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