FrontISTR  5.7.1
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
16  use m_fstr_ass_load
17  use m_fstr_addbc
18  use m_fstr_residual
19  use m_fstr_restart
20 
21  implicit none
22 
23  contains
24 
25 
28  subroutine fstr_newton( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, &
29  restrt_step_num, sub_step, ctime, dtime )
30 
31  integer, intent(in) :: cstep
32  type (hecmwST_local_mesh) :: hecMESH
33  type (hecmwST_matrix) :: hecMAT
34  type (fstr_solid) :: fstrSOLID
35  integer, intent(in) :: sub_step
36  real(kind=kreal), intent(in) :: ctime
37  real(kind=kreal), intent(in) :: dtime
38  type (fstr_param) :: fstrPARAM
39  type (hecmwST_matrix_lagrange) :: hecLagMAT
40 
41  type (hecmwST_local_mesh), pointer :: hecMESHmpc
42  type (hecmwST_matrix), pointer :: hecMATmpc
43  integer(kind=kint) :: ndof
44  integer(kind=kint) :: i, iter
45  integer(kind=kint) :: stepcnt
46  integer(kind=kint) :: restrt_step_num
47  real(kind=kreal) :: tt0, tt, res, qnrm, rres, tincr, xnrm, dunrm, rxnrm
48  real(kind=kreal), allocatable :: coord(:), p(:)
49  logical :: isLinear = .false.
50  integer(kind=kint) :: iterStatus
51 
52  call hecmw_mpc_mat_init(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
53 
54  if(.not. fstrpr%nlgeom)then
55  islinear = .true.
56  endif
57 
58  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
59 
60  allocate(p(hecmesh%n_node*ndof))
61  allocate(coord(hecmesh%n_node*ndof))
62  p = 0.0d0
63  stepcnt = 0
64 
65  ! ----- Inner Iteration, lagrange multiplier constant
66  do iter=1,fstrsolid%step_ctrl(cstep)%max_iter
67  stepcnt = stepcnt+1
68 
69  call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, ctime, tincr )
70  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
71 
72  ! ----- Set Boundary condition
73  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt)
74  call hecmw_mpc_mat_ass(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
75  call hecmw_mpc_trans_rhs(hecmesh, hecmat, hecmatmpc)
76 
77  !----- SOLVE [Kt]{du}={R}
78  if( sub_step == restrt_step_num .and. iter == 1 ) hecmatmpc%Iarray(98) = 1
79  if( iter == 1 ) then
80  hecmatmpc%Iarray(97) = 2 !Force numerical factorization
81  else
82  hecmatmpc%Iarray(97) = 1 !Need numerical factorization
83  endif
84  hecmatmpc%X = 0.0d0
85  call fstr_set_current_config_to_mesh(hecmeshmpc,fstrsolid,coord)
86  call solve_lineq(hecmeshmpc,hecmatmpc)
87  call fstr_recover_initial_config_to_mesh(hecmeshmpc,fstrsolid,coord)
88  call hecmw_mpc_tback_sol(hecmesh, hecmat, hecmatmpc)
89 
90  ! ----- update the small displacement and the displacement for 1step
91  ! \delta u^k => solver's solution
92  ! \Delta u_{n+1}^{k} = \Delta u_{n+1}^{k-1} + \delta u^k
93  do i = 1, hecmesh%n_node*ndof
94  fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
95  enddo
96 
97  call fstr_calc_residual_vector(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam)
98 
99  ! do i = 1, hecMESH%n_node*ndof
100  ! write(6,*) fstrSOLID%dunode(i), fstrSOLID%QFORCE(i), fstrSOLID%GL(i)
101  ! enddo
102 
103  if( islinear ) exit
104 
105  ! ----- check convergence
106  iterstatus = fstr_check_iteration_converged(hecmesh, hecmat, fstrsolid, ndof, iter, sub_step, cstep)
107  if (iterstatus == kitrconverged) exit
108  if (iterstatus == kitrdiverged .or. iterstatus==kitrfloatingerror) return
109  enddo
110  ! ----- end of inner loop
111 
112  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
113  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
114 
115  ! ----- update the total displacement
116  ! u_{n+1} = u_{n} + \Delta u_{n+1}
117  do i=1,hecmesh%n_node*ndof
118  fstrsolid%unode(i) = fstrsolid%unode(i) + fstrsolid%dunode(i)
119  enddo
120 
121  call fstr_updatestate( hecmesh, fstrsolid, tincr )
122 
123  fstrsolid%CutBack_stat = 0
124  deallocate(coord)
125  deallocate(p)
126  call hecmw_mpc_mat_finalize(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
127  end subroutine fstr_newton
128 
132  subroutine fstr_newton_contactalag( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, &
133  restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
134  use mcontact
137 
138  integer, intent(in) :: cstep
139  type (hecmwST_local_mesh) :: hecMESH
140  type (hecmwST_matrix) :: hecMAT
141  type (fstr_solid) :: fstrSOLID
142  integer, intent(in) :: sub_step
143  real(kind=kreal), intent(in) :: ctime
144  real(kind=kreal), intent(in) :: dtime
145  type (fstr_param) :: fstrPARAM
146  type (fstr_info_contactChange) :: infoCTChange
147  type (hecmwST_matrix_lagrange) :: hecLagMAT
148  type (hecmwST_matrix) :: conMAT
149 
150  integer(kind=kint) :: ndof
151  integer(kind=kint) :: ctAlgo
152  integer(kind=kint) :: i, iter
153  integer(kind=kint) :: al_step, n_al_step, stepcnt
154  real(kind=kreal) :: tt0, tt, res, res0, res1, maxv, relres, tincr
155  integer(kind=kint) :: restart_step_num, restart_substep_num
156  logical :: convg, ctchange
157  integer(kind=kint) :: n_node_global
158  integer(kind=kint) :: contact_changed_global
159  real(kind=kreal), allocatable :: coord(:)
160  integer(kind=kint) :: istat
161 
162 
163  ! sum of n_node among all subdomains (to be used to calc res)
164  n_node_global = hecmesh%nn_internal
165  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
166 
167  ctalgo = fstrparam%contact_algo
168 
169  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
170 
171  if( cstep == 1 .and. sub_step == restart_substep_num ) then
172  call fstr_save_originalmatrixstructure(hecmat)
173  if(hecmesh%my_rank==0) write(*,*) "---Scanning initial contact state---"
174  call fstr_scan_contact_state( cstep, sub_step, 0, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
175  call hecmw_mat_copy_profile( hecmat, conmat )
176  if ( fstr_is_contact_active() ) then
177  call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
178  elseif( hecmat%Iarray(99)==4 ) then
179  write(*, *) ' This type of direct solver is not yet available in such case ! '
180  write(*, *) ' Please change the solver type to intel MKL direct solver !'
181  call hecmw_abort(hecmw_comm_get_comm())
182  endif
183  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, .true.)
184  endif
185 
186  hecmat%X = 0.0d0
187 
188  stepcnt = 0
189  allocate(coord(hecmesh%n_node*ndof))
190 
191  call hecmw_mat_clear_b(conmat)
192 
193  if( fstr_is_contact_active() ) call fstr_ass_load_contactalag( hecmesh, fstrsolid, conmat%B )
194 
195  ! ----- Augmentation loop. In case of no contact, it is inactive
196  n_al_step = fstrsolid%step_ctrl(cstep)%max_contiter
197  if( .not. fstr_is_contact_active() ) n_al_step = 1
198 
199  do al_step = 1, n_al_step
200 
201  if( hecmesh%my_rank == 0) then
202  if( n_al_step > 1 ) then
203  print *, "Augmentation step:", al_step
204  write(imsg, *) "Augmentation step:", al_step
205  endif
206  end if
207 
208  ! ----- Inner Iteration, lagrange multiplier constant
209  res0 = 0.0d0
210  res1 = 0.0d0
211  relres = 1.0d0
212 
213  do iter = 1,fstrsolid%step_ctrl(cstep)%max_iter
214  stepcnt = stepcnt+1
215 
216  call fstr_stiffmatrix( hecmesh, hecmat, fstrsolid, ctime, tincr )
217  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
218 
219  call hecmw_mat_clear( conmat )
220  conmat%X = 0.0d0
221 
222  ! ----- Contact
223  if( al_step == 1 .and. stepcnt == 1 ) then
224  maxv = hecmw_mat_diag_max( hecmat, hecmesh )
225  call fstr_set_contact_penalty( maxv )
226  endif
227  if( fstr_is_contact_active() ) then
228  call fstr_contactbc( cstep, iter, hecmesh, conmat, fstrsolid )
229  endif
230 
231  ! ----- Set Boundary condition
232  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
233 
234  !----- SOLVE [Kt]{du}={R}
235  ! ---- For Parallel Contact with Multi-Partition Domains
236  hecmat%X = 0.0d0
237  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
238  call solve_lineq_contact(hecmesh, hecmat, heclagmat, conmat, istat, 1.0d0, fstr_is_contact_active())
239  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
240 
241  call hecmw_update_r (hecmesh, hecmat%X, hecmat%NP, hecmesh%n_dof)
242 
243  ! ----- update the small displacement and the displacement for 1step
244  ! \delta u^k => solver's solution
245  ! \Delta u_{n+1}^{k} = \Delta u_{n+1}^{k-1} + \delta u^k
246  do i = 1, hecmesh%n_node*ndof
247  fstrsolid%dunode(i) = fstrsolid%dunode(i)+hecmat%X(i)
248  enddo
249 
250  ! ----- update the strain, stress, and internal force
251  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter)
252 
253  ! ----- Set residual
254  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
255  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
256 
257  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid, conmat)
258 
259  if( fstr_is_contact_active() ) then
260  call hecmw_mat_clear_b( conmat )
261  call fstr_update_contact0(hecmesh, fstrsolid, conmat%B)
262  endif
263  ! Consider SPC condition
264  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
265  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
266 
267  !res = fstr_get_residual(hecMAT%B, hecMESH)
268  res = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
269  ! ----- Gather global residual
270  res = sqrt(res)/n_node_global
271  if( iter == 1 ) res0 = res
272  if( res0 == 0.0d0 ) then
273  res0 = 1.0d0
274  else
275  relres = dabs( res1-res )/res0
276  endif
277 
278  if( hecmesh%my_rank == 0 ) then
279  write(*, '(a,i3,a,2e15.7)') ' - Residual(',iter,') =', res, relres
280  endif
281 
282  ! ----- check convergence
283  if( res < fstrsolid%step_ctrl(cstep)%converg .or. &
284  relres < fstrsolid%step_ctrl(cstep)%converg_ddisp ) exit
285  res1 = res
286 
287  ! ----- check divergence and NaN
288  if( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) then
289  if( hecmesh%my_rank == 0) then
290  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
291  end if
292  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
293  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sumofiter)
294  fstrsolid%NRstat_i(knstciter) = al_step ! logging contact iteration
295  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
296  if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
297  if( res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) fstrsolid%NRstat_i(knstdresn) = 2
298  return
299  end if
300 
301  enddo
302  ! ----- end of inner loop
303 
304  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
305  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
306 
307  ! ----- deal with contact boundary
308  call fstr_update_contact_multiplier( hecmesh, fstrsolid, ctchange )
309  call fstr_scan_contact_state( cstep, sub_step, al_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
310 
311  contact_changed_global = 0
312  if( fstr_is_matrixstructure_changed(infoctchange) ) then
313  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
314  contact_changed_global = 1
315  endif
316  call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
317  if (contact_changed_global > 0) then
318  call hecmw_mat_clear_b( hecmat )
319  call hecmw_mat_clear_b( conmat )
320  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, .true.)
321  endif
322 
323  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) .and. .not. ctchange ) exit
324 
325  ! ----- check divergence
326  if( al_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
327  if( hecmesh%my_rank == 0) then
328  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', sub_step
329  end if
330  fstrsolid%NRstat_i(knstciter) = al_step ! logging contact iteration
331  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
332  fstrsolid%NRstat_i(knstdresn) = 3
333  return
334  end if
335 
336  ! ----- Set residual for next newton iteration
337  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
338 
339  if( fstr_is_contact_active() ) then
340  call hecmw_mat_clear_b( conmat )
341  call fstr_update_contact0(hecmesh, fstrsolid, conmat%B)
342  endif
343 
344  enddo
345  ! ----- end of augmentation loop
346 
347  ! ----- update the total displacement
348  ! u_{n+1} = u_{n} + \Delta u_{n+1}
349  do i=1,hecmesh%n_node*ndof
350  fstrsolid%unode(i) = fstrsolid%unode(i)+fstrsolid%dunode(i)
351  enddo
352 
353  fstrsolid%NRstat_i(knstciter) = al_step-1 ! logging contact iteration
354 
355  call fstr_updatestate( hecmesh, fstrsolid, tincr )
356 
357  deallocate(coord)
358  fstrsolid%CutBack_stat = 0
359  end subroutine fstr_newton_contactalag
360 
361 
364  subroutine fstr_newton_contactslag( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, &
365  restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
366 
367  use mcontact
370 
371  integer, intent(in) :: cstep
372  type (hecmwST_local_mesh) :: hecMESH
373  type (hecmwST_matrix) :: hecMAT
374  type (fstr_solid) :: fstrSOLID
375  integer, intent(in) :: sub_step
376  real(kind=kreal), intent(in) :: ctime
377  real(kind=kreal), intent(in) :: dtime
378  type (fstr_param) :: fstrPARAM
379  type (fstr_info_contactChange) :: infoCTChange
380  type (hecmwST_matrix_lagrange) :: hecLagMAT
381  type (hecmwST_matrix) :: conMAT
382 
383  integer(kind=kint) :: ndof
384  integer(kind=kint) :: ctAlgo
385  integer(kind=kint) :: i, iter, max_iter_contact
386  integer(kind=kint) :: stepcnt, count_step
387  real(kind=kreal) :: tt0, tt, res, res0, res1, relres, tincr, resx
388  integer(kind=kint) :: restart_step_num, restart_substep_num
389  logical :: is_mat_symmetric
390  integer(kind=kint) :: n_node_global
391  integer(kind=kint) :: contact_changed_global
392  integer(kint) :: nndof
393  real(kreal) :: q_residual,x_residual
394  real(kind=kreal), allocatable :: coord(:)
395  integer(kind=kint) :: istat
396 
397  ctalgo = fstrparam%contact_algo
398 
399  ! sum of n_node among all subdomains (to be used to calc res)
400  n_node_global = hecmesh%nn_internal
401  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
402 
403  if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh) ) then
404  write(*, *) ' This type of direct solver is not yet available in such case ! '
405  write(*, *) ' Please use intel MKL direct solver !'
406  call hecmw_abort( hecmw_comm_get_comm() )
407  endif
408 
409  do i=1,fstrsolid%n_contacts
410  fstrsolid%contacts(i)%ctime = ctime + dtime
411  enddo
412 
413  if( cstep==1 .and. sub_step==restart_substep_num ) then
414  call fstr_save_originalmatrixstructure(hecmat)
415  call fstr_scan_contact_state( cstep, sub_step, 0, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
416  call hecmw_mat_copy_profile( hecmat, conmat )
417  if ( fstr_is_contact_active() ) then
418  call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
419  elseif( hecmat%Iarray(99)==4 ) then
420  write(*, *) ' This type of direct solver is not yet available in such case ! '
421  write(*, *) ' Please change the solver type to intel MKL direct solver !'
422  call hecmw_abort(hecmw_comm_get_comm())
423  endif
424  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
425  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, is_mat_symmetric)
426  endif
427 
428  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
429 
430  call hecmw_mat_clear_b(conmat)
431 
432  if( fstr_is_contact_active() ) then
433  call fstr_ass_load_contact(cstep, hecmesh, conmat, fstrsolid, heclagmat)
434  endif
435 
436  stepcnt = 0
437  count_step = 0
438  allocate(coord(hecmesh%n_node*ndof))
439 
440  loopforcontactanalysis: do while( .true. )
441  count_step = count_step+1
442 
443  ! ----- Inner Iteration
444  res0 = 0.d0
445  res1 = 0.d0
446  relres = 1.d0
447 
448  do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
449  call hecmw_barrier(hecmesh)
450  if( myrank == 0 ) print *,'-------------------------------------------------'
451  call hecmw_barrier(hecmesh)
452  stepcnt = stepcnt+1
453 
454  call fstr_stiffmatrix(hecmesh, hecmat, fstrsolid, ctime, tincr)
455  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
456 
457  call hecmw_mat_clear( conmat )
458  conmat%X = 0.0d0
459 
460  if( fstr_is_contact_active() ) then
461  call fstr_addcontactstiffness(cstep,iter,conmat,heclagmat,fstrsolid)
462  endif
463 
464  ! ----- Set Boundary condition
465  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
466 
467  nndof = hecmat%N*hecmat%ndof
468 
469  !----- SOLVE [Kt]{du}={R}
470  ! ---- For Parallel Contact with Multi-Partition Domains
471  hecmat%X = 0.0d0
472  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
473  q_residual = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
474  call solve_lineq_contact(hecmesh, hecmat, heclagmat, conmat, istat, 1.0d0, fstr_is_contact_active())
475  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
476  ! ----- check matrix solver error
477  if( istat /= 0 ) then
478  if( hecmesh%my_rank == 0) then
479  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
480  end if
481  fstrsolid%NRstat_i(knstdresn) = 4
482  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
483  return
484  end if
485 
486  x_residual = fstr_get_x_norm_contact(hecmat,heclagmat,hecmesh)
487 
488  call hecmw_innerproduct_r(hecmesh,ndof,hecmat%X,hecmat%X,resx)
489  resx = sqrt(resx)/n_node_global
490 
491  if( hecmesh%my_rank==0 ) then
492  write(*,'(a,i3,a,e15.7)') ' - ResidualX (',iter,') =',resx
493  write(*,'(a,i3,a,e15.7)') ' - ResidualX+LAG(',iter,') =',sqrt(x_residual)/n_node_global
494  write(*,'(a,i3,a,e15.7)') ' - ResidualQ (',iter,') =',sqrt(q_residual)/n_node_global
495  endif
496 
497  ! ----- update the small displacement and the displacement for 1step
498  do i = 1, hecmesh%n_node*ndof
499  fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
500  enddo
501 
502  ! ----- update the Lagrange multipliers
503  if( fstr_is_contact_active() ) then
504  do i = 1, heclagmat%num_lagrange
505  heclagmat%lagrange(i) = heclagmat%lagrange(i)+hecmat%X(hecmesh%n_node*ndof+i)
506  enddo
507  endif
508 
509  ! ----- update the strain, stress, and internal force (only QFORCE)
510  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter)
511 
512  ! ----- Set residual
513  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
514  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
515 
516  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
517 
518  if( fstr_is_contact_active() ) then
519  call hecmw_mat_clear_b( conmat )
520  call fstr_update_ndforce_contact(cstep,hecmesh,hecmat,heclagmat,fstrsolid,conmat)
521  endif
522 
523  res = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
524 
525  res = sqrt(res)/n_node_global
526  if( iter == 1 ) res0 = res
527  if( res0 == 0.0d0 ) then
528  res0 =1.0d0
529  else
530  relres = dabs( res1-res )/res0
531  endif
532  if( hecmesh%my_rank == 0 ) then
533  write(*, '(a,i3,a,2e15.7)') ' - Residual(',iter,') =',res,relres
534  endif
535 
536  ! ----- check convergence
537  if( res < fstrsolid%step_ctrl(cstep)%converg .or. &
538  relres < fstrsolid%step_ctrl(cstep)%converg_ddisp ) exit
539  res1 = res
540 
541  ! ----- check divergence and NaN
542  if( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) then
543  if( hecmesh%my_rank == 0) then
544  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
545  end if
546  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
547  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sumofiter)
548  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
549  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
550  if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
551  if( res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) fstrsolid%NRstat_i(knstdresn) = 2
552  return
553  end if
554 
555  enddo
556  ! ----- end of inner loop
557 
558  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
559  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
560 
561  call fstr_scan_contact_state( cstep, sub_step, count_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
562 
563  if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_contact_active() ) then
564  write(*, *) ' This type of direct solver is not yet available in such case ! '
565  write(*, *) ' Please use intel MKL direct solver !'
566  call hecmw_abort( hecmw_comm_get_comm() )
567  endif
568 
569  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
570  contact_changed_global = 0
571  if( fstr_is_matrixstructure_changed(infoctchange) ) then
572  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
573  contact_changed_global = 1
574  endif
575 
576  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) ) exit loopforcontactanalysis
577 
578  call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
579  if (contact_changed_global > 0) then
580  call hecmw_mat_clear_b( hecmat )
581  call hecmw_mat_clear_b( conmat )
582  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, is_mat_symmetric)
583  endif
584 
585  ! ----- check divergence
586  if( count_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
587  if( hecmesh%my_rank == 0) then
588  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', sub_step
589  end if
590  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
591  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
592  fstrsolid%NRstat_i(knstdresn) = 3
593  return
594  end if
595 
596  ! ----- Set residual for next newton iteration
597  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
598 
599  if( fstr_is_contact_active() ) then
600  call hecmw_mat_clear_b( conmat )
601  call fstr_update_ndforce_contact(cstep,hecmesh,hecmat,heclagmat,fstrsolid,conmat)
602  endif
603 
604  enddo loopforcontactanalysis
605 
606  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
607 
608  ! ----- update the total displacement
609  ! u_{n+1} = u_{n} + \Delta u_{n+1}
610  do i = 1, hecmesh%n_node*ndof
611  fstrsolid%unode(i) = fstrsolid%unode(i)+fstrsolid%dunode(i)
612  enddo
613 
614  call fstr_updatestate(hecmesh, fstrsolid, tincr)
615  call fstr_update_contact_tangentforce( fstrsolid )
616  if( fstrsolid%n_embeds > 0 .and. paracontactflag ) then
617  call fstr_setup_parancon_contactvalue(hecmesh,ndof,fstrsolid%EMBED_NFORCE,1)
618  call fstr_update_ndforce_spc( cstep, hecmesh, fstrsolid, hecmat%B )
619  endif
620 
621  deallocate(coord)
622  fstrsolid%CutBack_stat = 0
623  end subroutine fstr_newton_contactslag
624 
625  subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof)
626  use m_fstr_update
627  implicit none
628  type (hecmwST_local_mesh) :: hecMESH
629  type (hecmwST_matrix) :: hecMAT
630  type (fstr_solid) :: fstrSOLID
631  real(kind=kreal), intent(in) :: ctime
632  real(kind=kreal), intent(in) :: dtime
633  type (fstr_param) :: fstrPARAM
634  real(kind=kreal), intent(inout) :: tincr
635  integer(kind=kint) :: iter
636  integer, intent(in) :: cstep
637  type (hecmwST_matrix_lagrange) :: hecLagMAT
638  integer(kind=kint), intent(inout) :: ndof
639 
640  hecmat%NDOF = hecmesh%n_dof
641  ndof = hecmat%ndof
642 
643  tincr = dtime
644  if( fstrsolid%step_ctrl(cstep)%solution == stepstatic ) tincr = 0.d0
645 
646  fstrsolid%dunode(:) = 0.0d0
647  fstrsolid%NRstat_i(:) = 0 ! logging newton iteration(init)
648 
649  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
650  end subroutine fstr_init_newton
651 
653  function fstr_check_iteration_converged(hecMESH, hecMAT, fstrSOLID, ndof, iter, sub_step, cstep) result(iterStatus)
654  implicit none
655  integer(kind=kint) :: iterstatus
656 
657  type (hecmwst_local_mesh) :: hecmesh
658  type (hecmwst_matrix) :: hecmat
659  type (fstr_solid) :: fstrsolid
660  integer(kind=kint), intent(in) :: ndof
661  integer(kind=kint), intent(in) :: iter
662  integer(kind=kint), intent(in) :: sub_step, cstep
663 
664  real(kind=kreal) :: res, qnrm, rres, xnrm, dunrm, rxnrm
665 
666  iterstatus = kitrcontinue
667 
668  call hecmw_innerproduct_r(hecmesh, ndof, hecmat%B, hecmat%B, res)
669  res = sqrt(res)
670  call hecmw_innerproduct_r(hecmesh, ndof, hecmat%X, hecmat%X, xnrm)
671  xnrm = sqrt(xnrm)
672  call hecmw_innerproduct_r(hecmesh, ndof, fstrsolid%QFORCE, fstrsolid%QFORCE, qnrm)
673  qnrm = sqrt(qnrm)
674  if (qnrm < 1.0d-8) qnrm = 1.0d0
675  if( iter == 1 ) then
676  dunrm = xnrm
677  else
678  call hecmw_innerproduct_r(hecmesh, ndof, fstrsolid%dunode, fstrsolid%dunode, dunrm)
679  dunrm = sqrt(dunrm)
680  endif
681  rres = res/qnrm
682  rxnrm = xnrm/dunrm
683 
684  if( hecmesh%my_rank == 0 ) then
685  if (qnrm == 1.0d0) then
686  write(*,"(a,i8,a,1pe11.4,a,1pe11.4)")" iter:", iter, ", residual(abs):", rres, ", disp.corr.:", rxnrm
687  else
688  write(*,"(a,i8,a,1pe11.4,a,1pe11.4)")" iter:", iter, ", residual:", rres, ", disp.corr.:", rxnrm
689  endif
690  endif
691  if( hecmw_mat_get_flag_diverged(hecmat) == kno ) then
692  if( rres < fstrsolid%step_ctrl(cstep)%converg .or. &
693  rxnrm < fstrsolid%step_ctrl(cstep)%converg_ddisp ) then
694  iterstatus=kitrconverged
695  return
696  endif
697  endif
698 
699  ! ----- check divergence and NaN
700  if ( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. rres > fstrsolid%step_ctrl(cstep)%maxres) &
701  iterstatus = kitrdiverged
702  if (rres /= rres ) iterstatus = kitrfloatingerror
703  if (iterstatus /= kitrcontinue) then
704  if( hecmesh%my_rank == 0) then
705  write(ilog,'(a,i5,a,i5)') '### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
706  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
707  end if
708  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
709  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sumofiter)
710  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
711  if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
712  if( rres > fstrsolid%step_ctrl(cstep)%maxres .or. rres /= rres ) fstrsolid%NRstat_i(knstdresn) = 2
713  return
714  end if
715  end function fstr_check_iteration_converged
716 
717 end module m_fstr_nonlinearmethod
This module provides functions: 1) obtain contact stiffness matrix of each contact pair and assemble ...
subroutine, public fstr_ass_load_contact(cstep, hecMESH, hecMAT, fstrSOLID, hecLagMAT)
This subroutine adds initial contact force vector to the right-hand side vector \at the beginning of ...
subroutine, public fstr_update_ndforce_contact(cstep, hecMESH, hecMAT, hecLagMAT, fstrSOLID, conMAT)
This subroutine obtains contact nodal force vector of each contact pair and assembles it into right-h...
subroutine, public fstr_addcontactstiffness(cstep, iter, hecMAT, hecLagMAT, fstrSOLID)
This subroutine obtains contact stiffness matrix of each contact pair and assembles it into global st...
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.
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_newton(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restrt_step_num, sub_step, ctime, dtime)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method.
integer(kind=kint) function fstr_check_iteration_converged(hecMESH, hecMAT, fstrSOLID, ndof, iter, sub_step, cstep)
\breaf This function check iteration status
subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof)
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
real(kind=kreal) function, public fstr_get_x_norm_contact(hecMAT, hecLagMAT, hecMESH)
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 tangent stiffness matrix.
subroutine, public fstr_stiffmatrix(hecMESH, hecMAT, fstrSOLID, time, tincr)
This subroutine creates tangential stiffness matrix.
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:96
subroutine fstr_recover_initial_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.f90:1112
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110
integer(kind=kint), parameter kitrfloatingerror
Definition: m_fstr.f90:92
integer(kind=kint), parameter kitrconverged
Definition: m_fstr.f90:90
integer(kind=kint), parameter ilog
FILE HANDLER.
Definition: m_fstr.f90:107
subroutine fstr_set_current_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.f90:1099
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.f90:209
integer(kind=kint), parameter kitrcontinue
iteration control
Definition: m_fstr.f90:89
integer(kind=kint), parameter kno
Definition: m_fstr.f90:31
integer(kind=kint), parameter kitrdiverged
Definition: m_fstr.f90:91
logical paracontactflag
PARALLEL CONTACT FLAG.
Definition: m_fstr.f90:100
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.
This module provides functions to calculate contact stiff matrix.
Definition: fstr_contact.f90:6
subroutine fstr_contactbc(cstep, iter, hecMESH, hecMAT, fstrSOLID)
Introduce contact stiff into global stiff matrix or mpc conditions into hecMESH.
subroutine fstr_update_contact_multiplier(hecMESH, fstrSOLID, ctchanged)
Update lagrangian multiplier.
logical function fstr_is_matrixstructure_changed(infoCTChange)
logical function fstr_is_contact_conv(ctAlgo, infoCTChange, hecMESH)
subroutine fstr_scan_contact_state(cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange, B)
Scanning contact state.
subroutine fstr_set_contact_penalty(maxv)
subroutine fstr_ass_load_contactalag(hecMESH, fstrSOLID, B)
Update lagrangian multiplier.
subroutine fstr_setup_parancon_contactvalue(hecMESH, ndof, vec, vtype)
subroutine fstr_update_contact_tangentforce(fstrSOLID)
Update tangent force.
subroutine fstr_update_contact0(hecMESH, fstrSOLID, B)
Update lagrangian multiplier.
logical function fstr_is_contact_active()