FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_solve_NonLinear.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
6 
8 
9  use m_fstr
10  use m_static_lib
11  use m_static_output
12 
13  use m_fstr_spring
15  use m_fstr_update
17  use m_fstr_ass_load
18  use m_fstr_addbc
19  use m_fstr_residual
20  use m_fstr_restart
21  use m_fstr_elemact
24  use mcontact
27 
28  implicit none
29 
30 contains
31 
32  subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, &
33  & ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof, ctAlgo, conMAT)
34  use m_fstr_update
35  implicit none
36  type (hecmwST_local_mesh) :: hecMESH
37  type (hecmwST_matrix) :: hecMAT
38  type (fstr_solid) :: fstrSOLID
39  real(kind=kreal), intent(in) :: ctime
40  real(kind=kreal), intent(in) :: dtime
41  type (fstr_param) :: fstrPARAM
42  real(kind=kreal), intent(inout) :: tincr
43  integer(kind=kint) :: iter
44  integer, intent(in) :: cstep
45  type (hecmwST_matrix_lagrange) :: hecLagMAT
46  integer(kind=kint), intent(inout) :: ndof
47  integer(kind=kint), intent(in) :: ctAlgo
48  type (hecmwST_matrix) :: conMAT
49 
50  hecmat%NDOF = hecmesh%n_dof
51  ndof = hecmat%ndof
52 
53  tincr = dtime
54  if( fstrsolid%step_ctrl(cstep)%solution == stepstatic ) tincr = 0.d0
55 
56  fstrsolid%dunode(:) = 0.0d0
57  fstrsolid%NRstat_i(:) = 0 ! logging newton iteration(init)
58 
59  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
60 
61  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
62  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, 0)
63  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid)
64  endif
65 
66  if( fstr_is_contact_active() ) then
67  call hecmw_mat_clear_b(conmat)
68  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
69  ! Consider SPC condition
70  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
71  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
72  endif
73 
74  end subroutine fstr_init_newton
75 
78  subroutine fstr_newton( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, &
79  restrt_step_num, sub_step, ctime, dtime )
80 
81  integer, intent(in) :: cstep
82  type (hecmwST_local_mesh) :: hecMESH
83  type (hecmwST_matrix) :: hecMAT
84  type (fstr_solid) :: fstrSOLID
85  integer, intent(in) :: sub_step
86  real(kind=kreal), intent(in) :: ctime
87  real(kind=kreal), intent(in) :: dtime
88  type (fstr_param) :: fstrPARAM
89  type (hecmwST_matrix_lagrange) :: hecLagMAT
90 
91  type (hecmwST_local_mesh), pointer :: hecMESHmpc
92  type (hecmwST_matrix), pointer :: hecMATmpc
93  integer(kind=kint) :: ndof
94  integer(kind=kint) :: i, iter
95  integer(kind=kint) :: stepcnt
96  integer(kind=kint) :: restrt_step_num
97  real(kind=kreal) :: tt0, tt, res, qnrm, rres, tincr, xnrm, dunrm, rxnrm
98  real(kind=kreal), allocatable :: coord(:), p(:)
99  logical :: isLinear = .false.
100  integer(kind=kint) :: iterStatus
101 
102  call hecmw_mpc_mat_init(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
103 
104  if(.not. fstrpr%nlgeom)then
105  islinear = .true.
106  endif
107 
108  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof, 0, hecmat)
109 
110  allocate(p(hecmesh%n_node*ndof))
111  allocate(coord(hecmesh%n_node*ndof))
112  p = 0.0d0
113  stepcnt = 0
114 
115  ! ----- Inner Iteration, lagrange multiplier constant
116  do iter=1,fstrsolid%step_ctrl(cstep)%max_iter
117  stepcnt = stepcnt+1
118 
119  call fstr_creatematrix_and_dampingforce( hecmesh, hecmat, fstrsolid, ctime, tincr )
120  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
121 
122  ! ----- Set Boundary condition
123  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt)
124  call hecmw_mpc_mat_ass(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
125  call hecmw_mpc_trans_rhs(hecmesh, hecmat, hecmatmpc)
126 
127  !----- SOLVE [Kt]{du}={R}
128  if( sub_step == restrt_step_num .and. iter == 1 ) hecmatmpc%Iarray(98) = 1
129  if( iter == 1 ) then
130  hecmatmpc%Iarray(97) = 2 !Force numerical factorization
131  else
132  hecmatmpc%Iarray(97) = 1 !Need numerical factorization
133  endif
134  hecmatmpc%X = 0.0d0
135  call fstr_set_current_config_to_mesh(hecmeshmpc,fstrsolid,coord)
136  call solve_lineq(hecmeshmpc,hecmatmpc)
137  call fstr_recover_initial_config_to_mesh(hecmeshmpc,fstrsolid,coord)
138  call hecmw_mpc_tback_sol(hecmesh, hecmat, hecmatmpc)
139 
140  ! ----- update the small displacement and the displacement for 1step
141  ! \delta u^k => solver's solution
142  ! \Delta u_{n+1}^{k} = \Delta u_{n+1}^{k-1} + \delta u^k
143  call fstr_apply_solution_increment( hecmesh, fstrsolid, ndof, hecmat%X )
144 
145  call fstr_calc_residual_vector(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam)
146 
147  if( islinear ) exit
148 
149  ! ----- check convergence
150  call fstr_check_convergence(hecmesh, hecmat, fstrsolid, fstrpr, &
151  ndof, iter, sub_step, cstep, &
152  hecmat%B, 0, &
153  res, res, &
154  0, &
155  iterstatus)
156  if (iterstatus == kitrconverged) exit
157  if (iterstatus == kitrdiverged .or. iterstatus==kitrfloatingerror) return
158  enddo
159  ! ----- end of inner loop
160 
161  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
162  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
163 
164  ! ----- update the total displacement
165  ! u_{n+1} = u_{n} + \Delta u_{n+1}
166  call fstr_commit_solution_increment( hecmesh, fstrsolid, ndof )
167 
168  call fstr_updatestate( hecmesh, fstrsolid, tincr )
169  ! Update REACTION using current QFORCE
170  call fstr_update_reaction_spc( cstep, hecmesh, fstrsolid )
171 
172  fstrsolid%CutBack_stat = 0
173  deallocate(coord)
174  deallocate(p)
175  call hecmw_mpc_mat_finalize(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
176  end subroutine fstr_newton
177 
181  subroutine fstr_newton_contactalag( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, &
182  restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
183 
184  integer, intent(in) :: cstep
185  type (hecmwST_local_mesh) :: hecMESH
186  type (hecmwST_matrix) :: hecMAT
187  type (fstr_solid) :: fstrSOLID
188  integer, intent(in) :: sub_step
189  real(kind=kreal), intent(in) :: ctime
190  real(kind=kreal), intent(in) :: dtime
191  type (fstr_param) :: fstrPARAM
192  type (fstr_info_contactChange) :: infoCTChange
193  type (hecmwST_matrix_lagrange) :: hecLagMAT
194  type (hecmwST_matrix) :: conMAT
195 
196  integer(kind=kint) :: ndof
197  integer(kind=kint) :: ctAlgo
198  integer(kind=kint) :: i, iter
199  integer(kind=kint) :: al_step, n_al_step, stepcnt, count_step
200  real(kind=kreal) :: tt0, tt, res, res0, res1, relres, tincr
201  integer(kind=kint) :: restart_step_num, restart_substep_num
202  logical :: convg, ctchange
203  integer(kind=kint) :: n_node_global
204  integer(kind=kint) :: contact_changed_global
205  real(kind=kreal), allocatable :: coord(:)
206  integer(kind=kint) :: istat
207  logical :: is_first_Stiffmatrixcall
208  integer(kind=kint) :: iterStatus, nresid
209  real(kind=kreal), allocatable :: resid_work(:)
210 
211 
212  ! sum of n_node among all subdomains (to be used to calc res)
213  n_node_global = hecmesh%nn_internal
214  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
215 
216  ctalgo = fstrparam%contact_algo
217 
218  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof, ctalgo, conmat)
219 
220  if( cstep == 1 .and. sub_step == restart_substep_num ) then
222  if(hecmesh%my_rank==0) write(*,*) "---Scanning initial contact state---"
223  call fstr_scan_contact_state( cstep, sub_step, 0, dtime, ctalgo, hecmesh, fstrsolid, infoctchange )
224  call hecmw_mat_copy_profile( hecmat, conmat )
225  if ( fstr_is_contact_active() ) then
226  call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
227  elseif( hecmat%Iarray(99)==4 ) then
228  write(*, *) ' This type of direct solver is not yet available in such case ! '
229  write(*, *) ' Please change the solver type to intel MKL direct solver !'
230  call hecmw_abort(hecmw_comm_get_comm())
231  endif
232  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, .true.)
233  endif
234 
235  hecmat%X = 0.0d0
236 
237  stepcnt = 0
238  allocate(coord(hecmesh%n_node*ndof))
239  allocate(resid_work(hecmesh%n_node*ndof + conmat%NP*ndof))
240 
241  ! ----- Augmentation loop. In case of no contact, it is inactive
242  n_al_step = fstrparam%augiter
243  count_step = 0
244  is_first_stiffmatrixcall = .true.
245 
246  loopforcontactanalysis: do while( .true. )
247  count_step = count_step + 1
248 
249  do al_step = 1, n_al_step
250 
251  if( hecmesh%my_rank == 0) then
252  write(*,*) "Contact iter: ", count_step, " Augmentation iter: ", al_step
253  end if
254 
255  ! ----- Inner Iteration, lagrange multiplier constant
256  res0 = 0.0d0
257  res1 = 0.0d0
258  relres = 1.0d0
259 
260  do iter = 1,fstrsolid%step_ctrl(cstep)%max_iter
261  stepcnt = stepcnt+1
262 
263  call fstr_creatematrix_and_dampingforce( hecmesh, hecmat, fstrsolid, ctime, tincr )
264  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
265 
266  call hecmw_mat_clear( conmat )
267  conmat%X = 0.0d0
268 
269  ! ----- Contact
270  if( is_first_stiffmatrixcall ) then
271  call fstr_calc_contact_refstiff(cstep, hecmesh, hecmat, fstrsolid)
272  is_first_stiffmatrixcall = .false.
273  endif
274  if( fstr_is_contact_active() ) then
275  call fstr_addcontactstiffness(cstep,ctalgo,iter,hecmesh,conmat,heclagmat,fstrsolid)
276  endif
277 
278  ! ----- Set Boundary condition
279  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
280 
281  !----- SOLVE [Kt]{du}={R}
282  ! ---- For Parallel Contact with Multi-Partition Domains
283  hecmat%X = 0.0d0
284  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
285  call solve_lineq_contact(hecmesh, hecmat, heclagmat, conmat, istat, 1.0d0, fstr_is_contact_active())
286  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
287 
288  call hecmw_update_r (hecmesh, hecmat%X, hecmat%NP, hecmesh%n_dof)
289 
290  ! ----- update the small displacement and the displacement for 1step
291  ! \delta u^k => solver's solution
292  ! \Delta u_{n+1}^{k} = \Delta u_{n+1}^{k-1} + \delta u^k
293  call fstr_apply_solution_increment( hecmesh, fstrsolid, ndof, hecmat%X )
294 
295  ! ----- update the strain, stress, and internal force
296  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter)
297 
298  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
299  call fstr_update_elemact_solid_by_value( hecmesh, fstrsolid, cstep, ctime )
300  endif
301 
302  ! ----- Set residual
303  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
304  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
305 
306  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid, conmat)
307  call hecmw_mat_clear_b( conmat )
308  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
309 
310  ! Consider SPC condition
311  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
312  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
313 
314  !res = fstr_get_residual(hecMAT%B, hecMESH)
315  call fstr_assemble_residual_contact(hecmat, heclagmat, conmat, hecmesh, resid_work, nresid)
316  call fstr_check_convergence(hecmesh, hecmat, fstrsolid, fstrpr, &
317  ndof, iter, sub_step, cstep, &
318  resid_work, nresid, &
319  res0, res1, &
320  n_node_global, &
321  iterstatus)
322  if (iterstatus == kitrconverged) exit
323  if (iterstatus == kitrdiverged .or. iterstatus == kitrfloatingerror) then
324  fstrsolid%NRstat_i(knstciter) = al_step
325  return
326  endif
327 
328  enddo
329  ! ----- end of inner loop
330 
331  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
332  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
333 
334  call fstr_update_contact_multiplier( cstep, ctalgo, hecmesh, heclagmat, fstrsolid, ctchange )
335 
336  ! ----- Set residual for next newton iteration
337  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid, conmat)
338  ! ----- deal with contact boundary
339  call hecmw_mat_clear_b( conmat )
340  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
341  ! Consider SPC condition
342  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
343  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
344 
345  enddo
346  ! ----- end of augmentation loop
347 
348  ! ----- compute CONT_NFORCE/CONT_FRIC for output
349  if( fstr_is_contact_active() ) &
350  call fstr_calc_contact_output_force(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
351 
352  ! ----- deal with contact boundary
353  call fstr_scan_contact_state( cstep, sub_step, count_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange )
354 
355  contact_changed_global = 0
356  if( fstr_is_matrixstructure_changed(infoctchange) ) then
357  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
358  contact_changed_global = 1
359  endif
360  call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
361  if (contact_changed_global > 0) then
362  call hecmw_mat_clear_b( hecmat )
363  call hecmw_mat_clear_b( conmat )
364  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, .true.)
365  endif
366 
367  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) .and. .not. ctchange ) exit loopforcontactanalysis
368 
369  ! ----- check divergence
370  if( count_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
371  if( hecmesh%my_rank == 0) then
372  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', sub_step
373  end if
374  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
375  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
376  fstrsolid%NRstat_i(knstdresn) = 3
377  return
378  end if
379 
380  ! ----- Set residual for next newton iteration
381  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
382 
383  if( fstr_is_contact_active() ) then
384  call hecmw_mat_clear_b( conmat )
385  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
386  ! Consider SPC condition
387  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
388  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
389  endif
390 
391  enddo loopforcontactanalysis
392 
393  ! ----- update the total displacement
394  ! u_{n+1} = u_{n} + \Delta u_{n+1}
395  call fstr_commit_solution_increment( hecmesh, fstrsolid, ndof )
396 
397  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
398 
399  call fstr_updatestate( hecmesh, fstrsolid, tincr )
400 
401  ! Update REACTION using current QFORCE
402  call fstr_update_reaction_spc( cstep, hecmesh, fstrsolid )
403 
404  deallocate(coord)
405  deallocate(resid_work)
406  fstrsolid%CutBack_stat = 0
407  end subroutine fstr_newton_contactalag
408 
411  subroutine fstr_newton_contactslag( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, &
412  restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
413 
414  integer, intent(in) :: cstep
415  type (hecmwST_local_mesh) :: hecMESH
416  type (hecmwST_matrix) :: hecMAT
417  type (fstr_solid) :: fstrSOLID
418  integer, intent(in) :: sub_step
419  real(kind=kreal), intent(in) :: ctime
420  real(kind=kreal), intent(in) :: dtime
421  type (fstr_param) :: fstrPARAM
422  type (fstr_info_contactChange) :: infoCTChange
423  type (hecmwST_matrix_lagrange) :: hecLagMAT
424  type (hecmwST_matrix) :: conMAT
425 
426  integer(kind=kint) :: ndof
427  integer(kind=kint) :: ctAlgo
428  integer(kind=kint) :: i, iter, max_iter_contact
429  integer(kind=kint) :: stepcnt, count_step
430  real(kind=kreal) :: tt0, tt, res, res0, res1, relres, tincr, resx
431  integer(kind=kint) :: restart_step_num, restart_substep_num
432  logical :: is_mat_symmetric
433  integer(kind=kint) :: n_node_global
434  integer(kind=kint) :: contact_changed_global
435  integer(kint) :: nndof
436  real(kreal) :: q_residual,x_residual
437  real(kind=kreal), allocatable :: coord(:)
438  integer(kind=kint) :: istat
439  integer(kind=kint) :: iterStatus, nresid
440  real(kind=kreal), allocatable :: resid_work(:)
441 
442  ctalgo = fstrparam%contact_algo
443 
444  ! sum of n_node among all subdomains (to be used to calc res)
445  n_node_global = hecmesh%nn_internal
446  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
447 
448  if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh) ) then
449  write(*, *) ' This type of direct solver is not yet available in such case ! '
450  write(*, *) ' Please use intel MKL direct solver !'
451  call hecmw_abort( hecmw_comm_get_comm() )
452  endif
453 
454  do i=1,fstrsolid%n_contacts
455  fstrsolid%contacts(i)%ctime = ctime + dtime
456  enddo
457 
458  if( cstep==1 .and. sub_step==restart_substep_num ) then
460  call fstr_scan_contact_state( cstep, sub_step, 0, dtime, ctalgo, hecmesh, fstrsolid, infoctchange )
461  call hecmw_mat_copy_profile( hecmat, conmat )
462  if ( fstr_is_contact_active() ) then
463  call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
464  elseif( hecmat%Iarray(99)==4 ) then
465  write(*, *) ' This type of direct solver is not yet available in such case ! '
466  write(*, *) ' Please change the solver type to intel MKL direct solver !'
467  call hecmw_abort(hecmw_comm_get_comm())
468  endif
469  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
470  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, is_mat_symmetric)
471  endif
472 
473  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof, ctalgo, conmat)
474 
475  stepcnt = 0
476  count_step = 0
477  allocate(coord(hecmesh%n_node*ndof))
478  allocate(resid_work(hecmesh%n_node*ndof + conmat%NP*ndof))
479 
480  loopforcontactanalysis: do while( .true. )
481  count_step = count_step+1
482 
483  ! ----- Inner Iteration
484  res0 = 0.d0
485  res1 = 0.d0
486  relres = 1.d0
487 
488  do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
489  call hecmw_barrier(hecmesh)
490  if( myrank == 0 ) print *,'-------------------------------------------------'
491  call hecmw_barrier(hecmesh)
492  stepcnt = stepcnt+1
493 
494  call fstr_creatematrix_and_dampingforce(hecmesh, hecmat, fstrsolid, ctime, tincr)
495  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
496 
497  call hecmw_mat_clear( conmat )
498  conmat%X = 0.0d0
499 
500  if( fstr_is_contact_active() ) then
501  call fstr_addcontactstiffness(cstep,ctalgo,iter,hecmesh,conmat,heclagmat,fstrsolid)
502  endif
503 
504  ! ----- Set Boundary condition
505  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
506 
507  nndof = hecmat%N*hecmat%ndof
508 
509  !----- SOLVE [Kt]{du}={R}
510  ! ---- For Parallel Contact with Multi-Partition Domains
511  hecmat%X = 0.0d0
512  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
513  q_residual = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
514  call solve_lineq_contact(hecmesh, hecmat, heclagmat, conmat, istat, 1.0d0, fstr_is_contact_active())
515  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
516  ! ----- check matrix solver error
517  if( istat /= 0 ) then
518  if( hecmesh%my_rank == 0) then
519  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
520  end if
521  fstrsolid%NRstat_i(knstdresn) = 4
522  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
523  return
524  end if
525 
526  x_residual = fstr_get_x_norm_contact(hecmat,heclagmat,hecmesh)
527 
528  call hecmw_innerproduct_r(hecmesh,ndof,hecmat%X,hecmat%X,resx)
529  resx = sqrt(resx)/n_node_global
530 
531  if( hecmesh%my_rank==0 ) then
532  write(*,'(a,i3,a,e15.7)') ' - ResidualX (',iter,') =',resx
533  write(*,'(a,i3,a,e15.7)') ' - ResidualX+LAG(',iter,') =',sqrt(x_residual)/n_node_global
534  write(*,'(a,i3,a,e15.7)') ' - ResidualQ (',iter,') =',sqrt(q_residual)/n_node_global
535  endif
536 
537  ! ----- update the small displacement and the displacement for 1step
538  call fstr_apply_solution_increment( hecmesh, fstrsolid, ndof, hecmat%X )
539 
540  ! ----- update the Lagrange multipliers
541  if( fstr_is_contact_active() ) then
542  do i = 1, heclagmat%num_lagrange
543  heclagmat%lagrange(i) = heclagmat%lagrange(i)+hecmat%X(hecmesh%n_node*ndof+i)
544  enddo
545  endif
546 
547  ! ----- update the strain, stress, and internal force (only QFORCE)
548  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter)
549 
550  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
551  call fstr_update_elemact_solid_by_value( hecmesh, fstrsolid, cstep, ctime )
552  endif
553 
554  ! ----- Set residual
555  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
556  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
557 
558  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
559 
560  if( fstr_is_contact_active() ) then
561  call hecmw_mat_clear_b( conmat )
562  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
563  ! Consider SPC condition
564  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
565  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
566  endif
567 
568  res = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
569 
570  call fstr_assemble_residual_contact(hecmat, heclagmat, conmat, hecmesh, resid_work, nresid)
571  call fstr_check_convergence(hecmesh, hecmat, fstrsolid, fstrpr, &
572  ndof, iter, sub_step, cstep, &
573  resid_work, nresid, &
574  res0, res1, &
575  n_node_global, &
576  iterstatus)
577  if (iterstatus == kitrconverged) then
578  exit
579  endif
580  if (iterstatus == kitrdiverged .or. iterstatus == kitrfloatingerror) then
581  fstrsolid%NRstat_i(knstciter) = count_step
582  return
583  endif
584 
585  enddo
586  ! ----- end of inner loop
587 
588  ! ----- compute CONT_NFORCE/CONT_FRIC for output
589  if( fstr_is_contact_active() ) &
590  call fstr_calc_contact_output_force(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
591 
592  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
593  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
594 
595  call fstr_scan_contact_state( cstep, sub_step, count_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange )
596 
597  if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_contact_active() ) then
598  write(*, *) ' This type of direct solver is not yet available in such case ! '
599  write(*, *) ' Please use intel MKL direct solver !'
600  call hecmw_abort( hecmw_comm_get_comm() )
601  endif
602 
603  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
604  contact_changed_global = 0
605  if( fstr_is_matrixstructure_changed(infoctchange) ) then
606  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
607  contact_changed_global = 1
608  endif
609 
610  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) ) exit loopforcontactanalysis
611 
612  call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
613  if (contact_changed_global > 0) then
614  call hecmw_mat_clear_b( hecmat )
615  call hecmw_mat_clear_b( conmat )
616  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, is_mat_symmetric)
617  endif
618 
619  ! ----- check divergence
620  if( count_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
621  if( hecmesh%my_rank == 0) then
622  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', sub_step
623  end if
624  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
625  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
626  fstrsolid%NRstat_i(knstdresn) = 3
627  return
628  end if
629 
630  ! ----- Set residual for next newton iteration
631  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
632 
633  if( fstr_is_contact_active() ) then
634  call hecmw_mat_clear_b( conmat )
635  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
636  ! Consider SPC condition
637  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
638  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
639  endif
640 
641  enddo loopforcontactanalysis
642 
643  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
644 
645  ! ----- update the total displacement
646  ! u_{n+1} = u_{n} + \Delta u_{n+1}
647  call fstr_commit_solution_increment( hecmesh, fstrsolid, ndof )
648 
649  call fstr_updatestate(hecmesh, fstrsolid, tincr)
650  call fstr_update_contact_tangentforce( cstep, fstrsolid )
651  if( fstrsolid%n_embeds > 0 .and. paracontactflag ) then
652  call fstr_setup_parancon_contactvalue(hecmesh,ndof,fstrsolid%EMBED_NFORCE,1)
653  call fstr_update_ndforce_spc( cstep, hecmesh, fstrsolid, hecmat%B )
654  endif
655  ! Update REACTION using current QFORCE
656  call fstr_update_reaction_spc( cstep, hecmesh, fstrsolid )
657 
658  deallocate(coord)
659  deallocate(resid_work)
660  fstrsolid%CutBack_stat = 0
661  end subroutine fstr_newton_contactslag
662 
663 end module m_fstr_nonlinearmethod
This module provides functions of reconstructing.
subroutine, public fstr_mat_con_contact(cstep, contact_algo, hecMAT, fstrSOLID, hecLagMAT, infoCTChange, conMAT, is_contact_active_flag)
this subroutine reconstructs node-based (stiffness) matrix structure \corresponding to contact state
subroutine, public fstr_save_originalmatrixstructure(hecMAT)
This subroutine saves original matrix structure constructed originally by hecMW_matrix.
logical function, public fstr_is_matrixstruct_symmetric(fstrSOLID, hecMESH)
this function judges whether sitiffness matrix is symmetric or not
This module provides a function to deal with prescribed displacement.
Definition: fstr_AddBC.f90:7
subroutine fstr_addbc(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, iter, conMAT, RHSvector)
Add Essential Boundary Conditions.
Definition: fstr_AddBC.f90:14
This module provides functions to take into account external load.
subroutine fstr_ass_load(cstep, ctime, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
This subroutine assmble following external force into fstrSOLIDGL and hecMATB afterwards.
Contact processing at assembly level (all pairs in one tContact object)
This module assembles the tangent stiffness matrix and, in the implicit dynamic case,...
subroutine, public fstr_creatematrix_and_dampingforce(hecMESH, hecMAT, fstrSOLID, time, tincr, fstrDYNAMIC, coef)
Assemble the system matrix and, optionally, the dynamic damping force.
This module provide a function to elemact elements.
subroutine fstr_update_elemact_solid_by_value(hecMESH, fstrSOLID, cstep, ctime)
This module provides a unified convergence check for Newton iteration.
subroutine, public fstr_check_convergence(hecMESH, hecMAT, fstrSOLID, fstrPR, ndof, iter, sub_step, cstep, residual_vec, nresid, resb, res_prev, n_node_global, iterStatus, maxDLag, converg_dlag)
Wrapper that calls fstr_check_convergence_main and applies the common divergence/NaN handling (status...
Finite-rotation nodal kinematics for NLGEOM.
subroutine, public fstr_apply_solution_increment(hecMESH, fstrSOLID, ndof, x)
Apply the linear-solver solution increment x to the step displacement dunode.
subroutine, public fstr_commit_solution_increment(hecMESH, fstrSOLID, ndof)
Commit the converged step increment dunode into the total displacement unode.
This module provides functions on nonlinear analysis.
subroutine fstr_newton_contactslag(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method....
subroutine fstr_newton_contactalag(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method combined with Neste...
subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof, ctAlgo, conMAT)
subroutine fstr_newton(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restrt_step_num, sub_step, ctime, dtime)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method.
This module provides function to calculate residual of nodal force.
subroutine, public fstr_update_ndforce_spc(cstep, hecMESH, fstrSOLID, B)
subroutine, public fstr_update_ndforce(cstep, hecMESH, hecMAT, fstrSOLID, conMAT)
real(kind=kreal) function, public fstr_get_norm_para_contact(hecMAT, hecLagMAT, conMAT, hecMESH)
subroutine fstr_calc_residual_vector(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM)
\breaf This subroutine calculate residual vector
subroutine, public fstr_assemble_residual_contact(hecMAT, hecLagMAT, conMAT, hecMESH, resid_vec, nresid)
Assemble contact residual vector (hecMATB + conMATB + Lagrange) into a single vector.
real(kind=kreal) function, public fstr_get_x_norm_contact(hecMAT, hecLagMAT, hecMESH)
subroutine, public fstr_update_reaction_spc(cstep, hecMESH, fstrSOLID)
Set fstrSOLIDREACTION at constrained DOFs using current fstrSOLIDQFORCE. Constrained DOFs are enumera...
This module provides functions to read in and write out restart files.
Definition: fstr_Restart.f90:8
This module provides functions to deal with spring force.
Definition: fstr_Spring.f90:7
subroutine fstr_addspring(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
Definition: fstr_Spring.f90:12
This module provides function to calculate to do updates.
Definition: fstr_Update.f90:6
subroutine fstr_updatestate(hecMESH, fstrSOLID, tincr)
Update elastiplastic status.
subroutine fstr_updatenewton(hecMESH, hecMAT, fstrSOLID, time, tincr, iter, strainEnergy)
Update displacement, stress, strain and internal forces.
Definition: fstr_Update.f90:26
This module defines common data and basic structures for analysis.
Definition: m_fstr.F90:15
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.F90:98
subroutine fstr_recover_initial_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.F90:1146
integer(kind=kint), parameter kitrfloatingerror
Definition: m_fstr.F90:94
integer(kind=kint), parameter kitrconverged
Definition: m_fstr.F90:92
subroutine fstr_set_current_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.F90:1133
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.F90:212
integer(kind=kint), parameter kitrdiverged
Definition: m_fstr.F90:93
logical paracontactflag
PARALLEL CONTACT FLAG.
Definition: m_fstr.F90:102
This module provides functions to solve sparse system of \linear equitions in the case of contact ana...
subroutine, public solve_lineq_contact_init(hecMESH, hecMAT, hecLagMAT, is_sym)
This subroutine.
subroutine, public solve_lineq_contact(hecMESH, hecMAT, hecLagMAT, conMAT, istat, rf, is_contact_active)
This subroutine.
This modules just summarizes all modules used in static analysis.
Definition: static_LIB.f90:6
This module provides functions to output result.
Top-level contact analysis module (System level)
subroutine fstr_update_ndforce_contact(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for residual vector (conMATB).
logical function fstr_is_matrixstructure_changed(infoCTChange)
logical function fstr_is_contact_conv(ctAlgo, infoCTChange, hecMESH)
subroutine fstr_addcontactstiffness(cstep, ctAlgo, iter, hecMESH, conMAT, hecLagMAT, fstrSOLID)
subroutine fstr_calc_contact_output_force(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for output (CONT_NFORCE/CONT_FRIC).
subroutine fstr_update_contact_multiplier(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, ctchanged)
subroutine fstr_calc_contact_refstiff(cstep, hecMESH, hecMAT, fstrSOLID)
Calculate reference stiffness for all contact pairs (System Level)
subroutine fstr_scan_contact_state(cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
logical function fstr_is_contact_active()
subroutine fstr_update_contact_tangentforce(cstep, fstrSOLID)
Update tangent force.