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
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  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
251  call fstr_update_elemact_solid_by_value( hecmesh, fstrsolid, cstep, ctime )
252  endif
253 
254  ! ----- Set residual
255  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
256  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
257 
258  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid, conmat)
259 
260  if( fstr_is_contact_active() ) then
261  call hecmw_mat_clear_b( conmat )
262  call fstr_update_contact0(hecmesh, fstrsolid, conmat%B)
263  endif
264  ! Consider SPC condition
265  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, hecmat%B)
266  call fstr_update_ndforce_spc(cstep, hecmesh, fstrsolid, conmat%B)
267 
268  !res = fstr_get_residual(hecMAT%B, hecMESH)
269  res = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
270  ! ----- Gather global residual
271  res = sqrt(res)/n_node_global
272  if( iter == 1 ) res0 = res
273  if( res0 == 0.0d0 ) then
274  res0 = 1.0d0
275  else
276  relres = dabs( res1-res )/res0
277  endif
278 
279  if( hecmesh%my_rank == 0 ) then
280  write(*, '(a,i3,a,2e15.7)') ' - Residual(',iter,') =', res, relres
281  endif
282 
283  ! ----- check convergence
284  if( res < fstrsolid%step_ctrl(cstep)%converg .or. &
285  relres < fstrsolid%step_ctrl(cstep)%converg_ddisp ) then
286  exit
287  endif
288  res1 = res
289 
290  ! ----- check divergence and NaN
291  if( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) then
292  if( hecmesh%my_rank == 0) then
293  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
294  end if
295  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
296  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sumofiter)
297  fstrsolid%NRstat_i(knstciter) = al_step ! logging contact iteration
298  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
299  if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
300  if( res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) fstrsolid%NRstat_i(knstdresn) = 2
301  return
302  end if
303 
304  enddo
305  ! ----- end of inner loop
306 
307  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
308  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
309 
310  ! ----- deal with contact boundary
311  call fstr_update_contact_multiplier( hecmesh, fstrsolid, ctchange )
312  call fstr_scan_contact_state( cstep, sub_step, al_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
313 
314  contact_changed_global = 0
315  if( fstr_is_matrixstructure_changed(infoctchange) ) then
316  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
317  contact_changed_global = 1
318  endif
319  call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
320  if (contact_changed_global > 0) then
321  call hecmw_mat_clear_b( hecmat )
322  call hecmw_mat_clear_b( conmat )
323  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, .true.)
324  endif
325 
326  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) .and. .not. ctchange ) exit
327 
328  ! ----- check divergence
329  if( al_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
330  if( hecmesh%my_rank == 0) then
331  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', sub_step
332  end if
333  fstrsolid%NRstat_i(knstciter) = al_step ! logging contact iteration
334  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
335  fstrsolid%NRstat_i(knstdresn) = 3
336  return
337  end if
338 
339  ! ----- Set residual for next newton iteration
340  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
341 
342  if( fstr_is_contact_active() ) then
343  call hecmw_mat_clear_b( conmat )
344  call fstr_update_contact0(hecmesh, fstrsolid, conmat%B)
345  endif
346 
347  enddo
348  ! ----- end of augmentation loop
349 
350  ! ----- update the total displacement
351  ! u_{n+1} = u_{n} + \Delta u_{n+1}
352  do i=1,hecmesh%n_node*ndof
353  fstrsolid%unode(i) = fstrsolid%unode(i)+fstrsolid%dunode(i)
354  enddo
355 
356  fstrsolid%NRstat_i(knstciter) = al_step-1 ! logging contact iteration
357 
358  call fstr_updatestate( hecmesh, fstrsolid, tincr )
359 
360  deallocate(coord)
361  fstrsolid%CutBack_stat = 0
362  end subroutine fstr_newton_contactalag
363 
364 
367  subroutine fstr_newton_contactslag( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, &
368  restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
369 
370  use mcontact
373 
374  integer, intent(in) :: cstep
375  type (hecmwST_local_mesh) :: hecMESH
376  type (hecmwST_matrix) :: hecMAT
377  type (fstr_solid) :: fstrSOLID
378  integer, intent(in) :: sub_step
379  real(kind=kreal), intent(in) :: ctime
380  real(kind=kreal), intent(in) :: dtime
381  type (fstr_param) :: fstrPARAM
382  type (fstr_info_contactChange) :: infoCTChange
383  type (hecmwST_matrix_lagrange) :: hecLagMAT
384  type (hecmwST_matrix) :: conMAT
385 
386  integer(kind=kint) :: ndof
387  integer(kind=kint) :: ctAlgo
388  integer(kind=kint) :: i, iter, max_iter_contact
389  integer(kind=kint) :: stepcnt, count_step
390  real(kind=kreal) :: tt0, tt, res, res0, res1, relres, tincr, resx
391  integer(kind=kint) :: restart_step_num, restart_substep_num
392  logical :: is_mat_symmetric
393  integer(kind=kint) :: n_node_global
394  integer(kind=kint) :: contact_changed_global
395  integer(kint) :: nndof
396  real(kreal) :: q_residual,x_residual
397  real(kind=kreal), allocatable :: coord(:)
398  integer(kind=kint) :: istat
399 
400  ctalgo = fstrparam%contact_algo
401 
402  ! sum of n_node among all subdomains (to be used to calc res)
403  n_node_global = hecmesh%nn_internal
404  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
405 
406  if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh) ) then
407  write(*, *) ' This type of direct solver is not yet available in such case ! '
408  write(*, *) ' Please use intel MKL direct solver !'
409  call hecmw_abort( hecmw_comm_get_comm() )
410  endif
411 
412  do i=1,fstrsolid%n_contacts
413  fstrsolid%contacts(i)%ctime = ctime + dtime
414  enddo
415 
416  if( cstep==1 .and. sub_step==restart_substep_num ) then
417  call fstr_save_originalmatrixstructure(hecmat)
418  call fstr_scan_contact_state( cstep, sub_step, 0, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
419  call hecmw_mat_copy_profile( hecmat, conmat )
420  if ( fstr_is_contact_active() ) then
421  call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
422  elseif( hecmat%Iarray(99)==4 ) then
423  write(*, *) ' This type of direct solver is not yet available in such case ! '
424  write(*, *) ' Please change the solver type to intel MKL direct solver !'
425  call hecmw_abort(hecmw_comm_get_comm())
426  endif
427  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
428  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, is_mat_symmetric)
429  endif
430 
431  call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
432 
433  call hecmw_mat_clear_b(conmat)
434 
435  if( fstr_is_contact_active() ) then
436  call fstr_ass_load_contact(cstep, hecmesh, conmat, fstrsolid, heclagmat)
437  endif
438 
439  stepcnt = 0
440  count_step = 0
441  allocate(coord(hecmesh%n_node*ndof))
442 
443  loopforcontactanalysis: do while( .true. )
444  count_step = count_step+1
445 
446  ! ----- Inner Iteration
447  res0 = 0.d0
448  res1 = 0.d0
449  relres = 1.d0
450 
451  do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
452  call hecmw_barrier(hecmesh)
453  if( myrank == 0 ) print *,'-------------------------------------------------'
454  call hecmw_barrier(hecmesh)
455  stepcnt = stepcnt+1
456 
457  call fstr_stiffmatrix(hecmesh, hecmat, fstrsolid, ctime, tincr)
458  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
459 
460  call hecmw_mat_clear( conmat )
461  conmat%X = 0.0d0
462 
463  if( fstr_is_contact_active() ) then
464  call fstr_addcontactstiffness(cstep,iter,conmat,heclagmat,fstrsolid)
465  endif
466 
467  ! ----- Set Boundary condition
468  call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
469 
470  nndof = hecmat%N*hecmat%ndof
471 
472  !----- SOLVE [Kt]{du}={R}
473  ! ---- For Parallel Contact with Multi-Partition Domains
474  hecmat%X = 0.0d0
475  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
476  q_residual = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
477  call solve_lineq_contact(hecmesh, hecmat, heclagmat, conmat, istat, 1.0d0, fstr_is_contact_active())
478  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
479  ! ----- check matrix solver error
480  if( istat /= 0 ) then
481  if( hecmesh%my_rank == 0) then
482  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
483  end if
484  fstrsolid%NRstat_i(knstdresn) = 4
485  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
486  return
487  end if
488 
489  x_residual = fstr_get_x_norm_contact(hecmat,heclagmat,hecmesh)
490 
491  call hecmw_innerproduct_r(hecmesh,ndof,hecmat%X,hecmat%X,resx)
492  resx = sqrt(resx)/n_node_global
493 
494  if( hecmesh%my_rank==0 ) then
495  write(*,'(a,i3,a,e15.7)') ' - ResidualX (',iter,') =',resx
496  write(*,'(a,i3,a,e15.7)') ' - ResidualX+LAG(',iter,') =',sqrt(x_residual)/n_node_global
497  write(*,'(a,i3,a,e15.7)') ' - ResidualQ (',iter,') =',sqrt(q_residual)/n_node_global
498  endif
499 
500  ! ----- update the small displacement and the displacement for 1step
501  do i = 1, hecmesh%n_node*ndof
502  fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
503  enddo
504 
505  ! ----- update the Lagrange multipliers
506  if( fstr_is_contact_active() ) then
507  do i = 1, heclagmat%num_lagrange
508  heclagmat%lagrange(i) = heclagmat%lagrange(i)+hecmat%X(hecmesh%n_node*ndof+i)
509  enddo
510  endif
511 
512  ! ----- update the strain, stress, and internal force (only QFORCE)
513  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter)
514 
515  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
516  call fstr_update_elemact_solid_by_value( hecmesh, fstrsolid, cstep, ctime )
517  endif
518 
519  ! ----- Set residual
520  if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
521  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
522 
523  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
524 
525  if( fstr_is_contact_active() ) then
526  call hecmw_mat_clear_b( conmat )
527  call fstr_update_ndforce_contact(cstep,hecmesh,hecmat,heclagmat,fstrsolid,conmat)
528  endif
529 
530  res = fstr_get_norm_para_contact(hecmat,heclagmat,conmat,hecmesh)
531 
532  res = sqrt(res)/n_node_global
533  if( iter == 1 ) res0 = res
534  if( res0 == 0.0d0 ) then
535  res0 =1.0d0
536  else
537  relres = dabs( res1-res )/res0
538  endif
539  if( hecmesh%my_rank == 0 ) then
540  write(*, '(a,i3,a,2e15.7)') ' - Residual(',iter,') =',res,relres
541  endif
542 
543  ! ----- check convergence
544  if( res < fstrsolid%step_ctrl(cstep)%converg .or. &
545  relres < fstrsolid%step_ctrl(cstep)%converg_ddisp ) then
546  exit
547  endif
548  res1 = res
549 
550  ! ----- check divergence and NaN
551  if( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) then
552  if( hecmesh%my_rank == 0) then
553  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
554  end if
555  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
556  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sumofiter)
557  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
558  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
559  if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
560  if( res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) fstrsolid%NRstat_i(knstdresn) = 2
561  return
562  end if
563 
564  enddo
565  ! ----- end of inner loop
566 
567  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
568  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
569 
570  call fstr_scan_contact_state( cstep, sub_step, count_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
571 
572  if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_contact_active() ) then
573  write(*, *) ' This type of direct solver is not yet available in such case ! '
574  write(*, *) ' Please use intel MKL direct solver !'
575  call hecmw_abort( hecmw_comm_get_comm() )
576  endif
577 
578  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
579  contact_changed_global = 0
580  if( fstr_is_matrixstructure_changed(infoctchange) ) then
581  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
582  contact_changed_global = 1
583  endif
584 
585  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) ) exit loopforcontactanalysis
586 
587  call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
588  if (contact_changed_global > 0) then
589  call hecmw_mat_clear_b( hecmat )
590  call hecmw_mat_clear_b( conmat )
591  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, is_mat_symmetric)
592  endif
593 
594  ! ----- check divergence
595  if( count_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
596  if( hecmesh%my_rank == 0) then
597  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', sub_step
598  end if
599  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
600  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
601  fstrsolid%NRstat_i(knstdresn) = 3
602  return
603  end if
604 
605  ! ----- Set residual for next newton iteration
606  call fstr_update_ndforce(cstep,hecmesh,hecmat,fstrsolid,conmat )
607 
608  if( fstr_is_contact_active() ) then
609  call hecmw_mat_clear_b( conmat )
610  call fstr_update_ndforce_contact(cstep,hecmesh,hecmat,heclagmat,fstrsolid,conmat)
611  endif
612 
613  enddo loopforcontactanalysis
614 
615  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
616 
617  ! ----- update the total displacement
618  ! u_{n+1} = u_{n} + \Delta u_{n+1}
619  do i = 1, hecmesh%n_node*ndof
620  fstrsolid%unode(i) = fstrsolid%unode(i)+fstrsolid%dunode(i)
621  enddo
622 
623  call fstr_updatestate(hecmesh, fstrsolid, tincr)
624  call fstr_update_contact_tangentforce( fstrsolid )
625  if( fstrsolid%n_embeds > 0 .and. paracontactflag ) then
626  call fstr_setup_parancon_contactvalue(hecmesh,ndof,fstrsolid%EMBED_NFORCE,1)
627  call fstr_update_ndforce_spc( cstep, hecmesh, fstrsolid, hecmat%B )
628  endif
629 
630  deallocate(coord)
631  fstrsolid%CutBack_stat = 0
632  end subroutine fstr_newton_contactslag
633 
634  subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof)
635  use m_fstr_update
636  implicit none
637  type (hecmwST_local_mesh) :: hecMESH
638  type (hecmwST_matrix) :: hecMAT
639  type (fstr_solid) :: fstrSOLID
640  real(kind=kreal), intent(in) :: ctime
641  real(kind=kreal), intent(in) :: dtime
642  type (fstr_param) :: fstrPARAM
643  real(kind=kreal), intent(inout) :: tincr
644  integer(kind=kint) :: iter
645  integer, intent(in) :: cstep
646  type (hecmwST_matrix_lagrange) :: hecLagMAT
647  integer(kind=kint), intent(inout) :: ndof
648 
649  hecmat%NDOF = hecmesh%n_dof
650  ndof = hecmat%ndof
651 
652  tincr = dtime
653  if( fstrsolid%step_ctrl(cstep)%solution == stepstatic ) tincr = 0.d0
654 
655  fstrsolid%dunode(:) = 0.0d0
656  fstrsolid%NRstat_i(:) = 0 ! logging newton iteration(init)
657 
658  call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
659  if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 ) then
660  call fstr_updatenewton(hecmesh, hecmat, fstrsolid, ctime, tincr, 0)
661  call fstr_update_ndforce(cstep, hecmesh, hecmat, fstrsolid)
662  endif
663 
664  end subroutine fstr_init_newton
665 
667  function fstr_check_iteration_converged(hecMESH, hecMAT, fstrSOLID, ndof, iter, sub_step, cstep) result(iterStatus)
668  implicit none
669  integer(kind=kint) :: iterstatus
670 
671  type (hecmwst_local_mesh) :: hecmesh
672  type (hecmwst_matrix) :: hecmat
673  type (fstr_solid) :: fstrsolid
674  integer(kind=kint), intent(in) :: ndof
675  integer(kind=kint), intent(in) :: iter
676  integer(kind=kint), intent(in) :: sub_step, cstep
677 
678  real(kind=kreal) :: res, qnrm, rres, xnrm, dunrm, rxnrm
679 
680  iterstatus = kitrcontinue
681 
682  call hecmw_innerproduct_r(hecmesh, ndof, hecmat%B, hecmat%B, res)
683  res = sqrt(res)
684  call hecmw_innerproduct_r(hecmesh, ndof, hecmat%X, hecmat%X, xnrm)
685  xnrm = sqrt(xnrm)
686  call hecmw_innerproduct_r(hecmesh, ndof, fstrsolid%QFORCE, fstrsolid%QFORCE, qnrm)
687  qnrm = sqrt(qnrm)
688  if (qnrm < 1.0d-8) qnrm = 1.0d0
689  if( iter == 1 ) then
690  dunrm = xnrm
691  else
692  call hecmw_innerproduct_r(hecmesh, ndof, fstrsolid%dunode, fstrsolid%dunode, dunrm)
693  dunrm = sqrt(dunrm)
694  endif
695  rres = res/qnrm
696  rxnrm = xnrm/dunrm
697 
698  if( hecmesh%my_rank == 0 ) then
699  if (qnrm == 1.0d0) then
700  write(*,"(a,i8,a,1pe11.4,a,1pe11.4)")" iter:", iter, ", residual(abs):", rres, ", disp.corr.:", rxnrm
701  else
702  write(*,"(a,i8,a,1pe11.4,a,1pe11.4)")" iter:", iter, ", residual:", rres, ", disp.corr.:", rxnrm
703  endif
704  endif
705  if( hecmw_mat_get_flag_diverged(hecmat) == kno ) then
706  if( rres < fstrsolid%step_ctrl(cstep)%converg .or. &
707  rxnrm < fstrsolid%step_ctrl(cstep)%converg_ddisp ) then
708  iterstatus=kitrconverged
709  return
710  endif
711  endif
712 
713  ! ----- check divergence and NaN
714  if ( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. rres > fstrsolid%step_ctrl(cstep)%maxres) &
715  iterstatus = kitrdiverged
716  if (rres /= rres ) iterstatus = kitrfloatingerror
717  if (iterstatus /= kitrcontinue) then
718  if( hecmesh%my_rank == 0) then
719  write(ilog,'(a,i5,a,i5)') '### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
720  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', sub_step
721  end if
722  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
723  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sumofiter)
724  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
725  if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
726  if( rres > fstrsolid%step_ctrl(cstep)%maxres .or. rres /= rres ) fstrsolid%NRstat_i(knstdresn) = 2
727  return
728  end if
729  end function fstr_check_iteration_converged
730 
731 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:1131
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:1118
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()