FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_dynamic_nlimplicit.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  use m_fstr
16  use m_fstr_update
17  use m_fstr_restart
19  use m_fstr_residual
21  use mcontact
24  use m_fstr_timeinc
25  use m_fstr_cutback
26  use m_fstr_spring
27 
28  !-------- for couple -------
30  use m_fstr_rcap_io
31 
32  real(kind=kreal), parameter :: pi = 3.14159265358979323846d0
33 
34 contains
35 
38  subroutine fstr_solve_nlgeom_dynamic_implicit_contactslag(hecMESH,hecMAT,fstrSOLID,fstrEIG &
39  ,fstrDYNAMIC,fstrRESULT,fstrPARAM &
40  ,fstrCPL,hecLagMAT,restart_step_num,restart_substep_num,infoCTChange &
41  ,conMAT,restart_step_count )
42  implicit none
43  !C-- global variable
44  type(hecmwst_local_mesh) :: hecmesh
45  type(hecmwst_matrix) :: hecMAT
46  type(hecmwst_matrix), pointer :: hecMAT0
47  type(fstr_eigen) :: fstrEIG
48  type(fstr_solid) :: fstrSOLID
49  type(hecmwst_result_data) :: fstrRESULT ! unused - kept for interface compatibility
50  type(fstr_param) :: fstrPARAM
51  type(fstr_dynamic) :: fstrDYNAMIC
52  type(fstr_couple) :: fstrCPL !for COUPLE
53  type(hecmwst_matrix_lagrange) :: hecLagMAT
54  type(fstr_info_contactchange) :: infoCTChange
55  type(hecmwst_matrix) :: conMAT
56 
57  !C-- local variable
58  integer(kind=kint) :: nnod, ndof, nn
59  integer(kind=kint) :: i, tot_step_print, CBbound
60  real(kind=kreal) :: time_1, time_2, factor
61  integer(kind=kint) :: sub_step
62 
63  integer(kind=kint) :: restart_step_num, restart_substep_num, restart_step_count, tot_step, step_count
64  integer(kind=kint) :: ctAlgo
65  integer(kind=kint) :: max_iter_contact
66  real(kind=kreal) :: converg_dlag
67  type(fstr_info_contactchange) :: infoctchange_bak
68 
69  logical :: is_outpoint
70  integer(kind=kint) :: n_node_global
71  logical :: is_mat_symmetric, is_interaction_active
72 
73  if(hecmesh%my_rank==0) call fstr_timeinc_printstatus_init
74 
75  is_interaction_active = ( associated( fstrsolid%contacts ) .or. associated( fstrsolid%embeds ) )
76 
77  nullify(hecmat0)
78 
79  ! sum of n_node among all subdomains (to be used to calc res)
80  n_node_global = hecmesh%nn_internal
81  call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
82 
83  ctalgo = fstrparam%contact_algo
84 
85  if( hecmat%Iarray(99)==4 .and. .not.fstr_is_matrixstruct_symmetric(fstrsolid,hecmesh) ) then
86  write(*,*) ' This type of direct solver is not yet available in such case ! '
87  write(*,*) ' Please use intel MKL direct solver !'
88  call hecmw_abort(hecmw_comm_get_comm())
89  endif
90 
91  hecmat%NDOF=hecmesh%n_dof
92 
93  nnod=hecmesh%n_node
94  ndof=hecmat%NDOF
95  nn=ndof*ndof
96 
97  if( associated( fstrsolid%contacts ) ) call initialize_contact_output_vectors(fstrsolid,hecmat)
98 
99  !!-- initial value
100  time_1 = hecmw_wtime()
101 
102  !C-- check parameters
103  if(dabs(fstrdynamic%beta) < 1.0e-20) then
104  if( hecmesh%my_rank == 0 ) then
105  write(imsg,*) 'stop due to Newmark-beta = 0'
106  endif
107  call hecmw_abort( hecmw_comm_get_comm())
108  endif
109 
110  !C-- matrix [M] lumped mass matrix
111  if(fstrdynamic%idx_mas == 1) then
112  call setmass(fstrsolid,hecmesh,hecmat,fstreig)
113 
114  !C-- consistent mass matrix
115  else if(fstrdynamic%idx_mas == 2) then
116  if( hecmesh%my_rank .eq. 0 ) then
117  write(imsg,*) 'stop: consistent mass matrix is not yet available !'
118  endif
119  call hecmw_abort( hecmw_comm_get_comm())
120  endif
121 
122  hecmat%Iarray(98) = 1 !Assembly complete
123  hecmat%Iarray(97) = 1 !Need numerical factorization
124 
125  !C-- initialize variables
126  if( restart_step_num == 1 .and. fstrdynamic%VarInitialize .and. abs(fstrdynamic%ray_m) > 1.0d-15 ) &
127  call dynamic_init_varibles( hecmesh, hecmat, fstrsolid, fstreig, fstrdynamic, fstrparam )
128 
129  !C-- output of initial state
130  if( restart_step_num == 1 ) then
131  call fstr_dynamic_output(1, 0, 0.d0, hecmesh, fstrsolid, fstrdynamic, fstrparam, .true.)
132  call dynamic_output_monit(1, 0, 0.d0, hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
133  endif
134 
135  fstrdynamic%VEC3(:) =0.d0
136  hecmat%X(:) =0.d0
137 
139  call fstr_scan_contact_state(restart_step_num, restart_step_num, 0, fstrdynamic%t_delta, &
140  ctalgo, hecmesh, fstrsolid, infoctchange)
141 
142  call hecmw_mat_copy_profile( hecmat, conmat )
143 
144  if ( fstr_is_contact_active() ) then
145  call fstr_mat_con_contact( restart_step_num, ctalgo, hecmat, fstrsolid, heclagmat, &
146  infoctchange, conmat, fstr_is_contact_active())
147  elseif( hecmat%Iarray(99)==4 ) then
148  write(*,*) ' This type of direct solver is not yet available in such case ! '
149  write(*,*) ' Please change solver type to intel MKL direct solver !'
150  call hecmw_abort(hecmw_comm_get_comm())
151  endif
152  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid,hecmesh)
153  call solve_lineq_contact_init(hecmesh,hecmat,heclagmat,is_mat_symmetric)
154 
155  fstrsolid%FACTOR = 0.0d0
156  call fstr_cutback_init( hecmesh, fstrsolid, fstrparam )
157  call fstr_cutback_save( fstrsolid, infoctchange, infoctchange_bak )
158 
159  step_count = restart_step_count
160  do tot_step=1, fstrsolid%nstep_tot
161  tot_step_print = tot_step+restart_step_num-1
162  if(hecmesh%my_rank==0) write(*,*) ''
163  if(hecmesh%my_rank==0) write(*,'(a,i5)') ' loading step=',tot_step_print
164 
165  sub_step = restart_substep_num
166  do while(.true.)
167  if (ndof == 4 .and. hecmesh%my_rank==0) write(*,'(a,i5)')"iter: ",sub_step
168  ! ----- time history of factor
169  call fstr_timeinc_settimeincrement( fstrsolid%step_ctrl(tot_step), fstrparam, sub_step, &
170  & fstrsolid%NRstat_i, fstrsolid%NRstat_r, fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
171 
172  fstrdynamic%t_curr = fstr_get_time()
173  fstrdynamic%t_delta = fstr_get_timeinc()
174 
175  call fstr_newton_dynamic_contactslag(tot_step, hecmesh, hecmat, fstrsolid, fstreig, &
176  fstrdynamic, fstrparam, fstrcpl, heclagmat, infoctchange, conmat, &
177  restart_step_num, hecmat0, sub_step, fstrdynamic%t_curr, fstrdynamic%t_delta)
178 
179  ! Time Increment
180  if( hecmesh%my_rank == 0 ) call fstr_timeinc_printstatus( fstrsolid%step_ctrl(tot_step), fstrparam, &
181  & tot_step_print, sub_step, fstrsolid%NRstat_i, fstrsolid%NRstat_r, &
182  & fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
183  if( fstr_cutback_active() ) then
184 
185  if( fstrsolid%CutBack_stat == 0 ) then ! converged
186  call fstr_cutback_save( fstrsolid, infoctchange, infoctchange_bak ) ! save analysis state
187  call fstr_proceed_time() ! current time += time increment
188  fstrdynamic%t_curr = fstr_get_time()
189 
190  else ! not converged
191  cbbound = fstrparam%ainc(fstrsolid%step_ctrl(tot_step)%AincParam_id)%CBbound
192  if( fstrsolid%CutBack_stat == cbbound ) then
193  if( hecmesh%my_rank == 0 ) then
194  write(*,*) 'Number of successive cutback reached max number: ',cbbound
195  call fstr_timeinc_printstatus_final(.false.)
196  endif
197  call hecmw_abort( hecmw_comm_get_comm() )
198  endif
199  call fstr_cutback_load( fstrsolid, infoctchange, infoctchange_bak ) ! load analysis state
200  call fstr_set_contact_active( infoctchange%contactNode_current > 0 )
201 
202  ! restore matrix structure for slagrange contact analysis
203  if( is_interaction_active .and. fstrparam%contact_algo == kcaslagrange ) then
204  call fstr_mat_con_contact( tot_step, fstrparam%contact_algo, hecmat, fstrsolid, heclagmat, &
205  & infoctchange, conmat, fstr_is_contact_active())
206  conmat%B(:) = 0.0d0
207  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh))
208  endif
209  if( hecmesh%my_rank == 0 ) write(*,*) '### State has been restored at time =',fstr_get_time()
210 
211  !stop if # of substeps reached upper bound.
212  if( sub_step == fstrsolid%step_ctrl(tot_step)%num_substep ) then
213  if( hecmesh%my_rank == 0 ) then
214  write(*,'(a,i5,a,f6.3)') '### Number of substeps reached max number: at total_step=', &
215  & tot_step_print, ' time=', fstr_get_time()
216  endif
217  call hecmw_abort( hecmw_comm_get_comm())
218  endif
219 
220  ! output time
221  time_2 = hecmw_wtime()
222  if( hecmesh%my_rank==0) write(imsg,'(a,",",2(I8,","),f10.2)') &
223  & 'step, substep, solve (sec) :', tot_step_print, sub_step, time_2 - time_1
224  cycle
225  endif
226  else
227  if( fstrsolid%CutBack_stat > 0 ) stop
228  call fstr_proceed_time() ! current time += time increment
229  fstrdynamic%t_curr = fstr_get_time()
230  endif
231 
232  step_count = step_count + 1
233 
234  ! ----- Result output (include visualize output)
235  is_outpoint = fstr_timeinc_istimepoint( fstrsolid%step_ctrl(tot_step), fstrparam ) &
236  & .or. fstr_timeinc_isstepfinished( fstrsolid%step_ctrl(tot_step) )
237  !C-- output new displacement, velocity and acceleration
238  call fstr_dynamic_output(tot_step, step_count, fstrdynamic%t_curr, hecmesh, fstrsolid, fstrdynamic, fstrparam, is_outpoint)
239 
240  !C-- output result of monitoring node
241  call dynamic_output_monit(tot_step, i, fstrdynamic%t_curr, hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
242 
243  !--- Restart info
244  if( fstrdynamic%restart_nout > 0 ) then
245  if( mod(step_count,fstrdynamic%restart_nout).eq.0 ) then
246  call fstr_write_restart_dyna_nl(tot_step,sub_step,hecmesh,fstrsolid,fstrdynamic,fstrparam,&
247  .false.,infoctchange%contactNode_current,step_count)
248  endif
249  endif
250 
251  !if time reached the end time of step, exit loop.
252  if( fstr_timeinc_isstepfinished( fstrsolid%step_ctrl(tot_step) ) ) exit
253 
254  if( sub_step == fstrsolid%step_ctrl(tot_step)%num_substep ) then
255  if( hecmesh%my_rank == 0 ) then
256  write(*,'(a,i5,a,f6.3)') '### Number of substeps reached max number: at total_step=', &
257  & tot_step_print, ' time=', fstr_get_time()
258  endif
259  if( hecmesh%my_rank == 0 ) call fstr_timeinc_printstatus_final(.false.)
260  stop !stop if # of substeps reached upper bound.
261  endif
262 
263  sub_step = sub_step + 1
264  enddo
265 
266  !--- Restart at the end of step
267  if( fstrdynamic%restart_nout > 0 ) then
268  call fstr_write_restart_dyna_nl(tot_step,sub_step,hecmesh,fstrsolid,fstrdynamic,fstrparam,&
269  .true.,infoctchange%contactNode_current,step_count)
270  endif
271 
272  restart_substep_num = 1
273  !C-- end of time step loop
274  enddo
275 
276  if (associated(hecmat0)) then
277  call hecmw_mat_finalize(hecmat0)
278  deallocate(hecmat0)
279  endif
280 
281  ! message
282  if( hecmesh%my_rank == 0 ) then
283  call fstr_timeinc_printstatus_final(.true.)
284  write(imsg,'("### FSTR_SOLVE_NLGEOM FINISHED!")')
285  write(*,'("### FSTR_SOLVE_NLGEOM FINISHED!")')
286  endif
287 
289 
290  subroutine fstr_newton_dynamic_contactslag(cstep, hecMESH, hecMAT, fstrSOLID, fstrEIG, &
291  fstrDYNAMIC, fstrPARAM, fstrCPL, hecLagMAT, infoCTChange, conMAT, &
292  restart_step_num, hecMAT0, istep, t_curr, t_delta)
293  implicit none
294  !C-- arguments
295  integer(kind=kint), intent(in) :: cstep, restart_step_num, istep
296  real(kind=kreal), intent(in) :: t_curr, t_delta
297  type(hecmwst_local_mesh) :: hecmesh
298  type(hecmwst_matrix) :: hecMAT
299  type(hecmwst_matrix), pointer :: hecMAT0
300  type(fstr_solid) :: fstrSOLID
301  type(fstr_eigen) :: fstrEIG
302  type(fstr_dynamic) :: fstrDYNAMIC
303  type(fstr_param) :: fstrPARAM
304  type(fstr_couple) :: fstrCPL
305  type(hecmwst_matrix_lagrange) :: hecLagMAT
306  type(fstr_info_contactchange) :: infoCTChange
307  type(hecmwst_matrix) :: conMAT
308 
309  !C-- local variables
310  integer(kind=kint) :: j, kk, idm, imm
311  integer(kind=kint) :: iter
312  real(kind=kreal) :: a1, a2, a3, b1, b2, b3, c1, c2
313  real(kind=kreal) :: coef(6)
314  real(kind=kreal) :: res, res1, res0, relres
315  integer(kind=kint) :: count_step, stepcnt
316  real(kind=kreal) :: maxdlag
317  integer(kind=kint) :: contact_changed_global
318  logical :: is_mat_symmetric
319  integer :: istat
320  logical :: is_cycle
321  integer(kind=kint) :: ctAlgo, max_iter_contact
322  integer(kind=kint) :: nnod, ndof, nn
323  real(kind=kreal) :: converg_dlag
324  real(kind=kreal), allocatable :: coord(:)
325  integer(kind=kint) :: iterStatus, nresid, n_node_global
326  real(kind=kreal), allocatable :: resid_work(:)
327 
328  fstrsolid%NRstat_i(:) = 0 ! logging newton iteration(init)
329 
330  !C-- initialize local variables
331  n_node_global = hecmesh%nn_internal
332  call hecmw_allreduce_i1(hecmesh, n_node_global, hecmw_sum)
333  ctalgo = fstrparam%contact_algo
334  max_iter_contact = fstrsolid%step_ctrl(cstep)%max_contiter
335  converg_dlag = fstrsolid%step_ctrl(cstep)%converg_lag
336  nnod = hecmesh%n_node
337  ndof = hecmat%NDOF
338  nn = ndof*ndof
339  allocate(coord(hecmesh%n_node*ndof))
340  allocate(resid_work(hecmesh%n_node*ndof + conmat%NP*ndof))
341 
342  a1 = .5d0/fstrdynamic%beta - 1.d0
343  a2 = 1.d0/(fstrdynamic%beta*t_delta)
344  a3 = 1.d0/(fstrdynamic%beta*t_delta*t_delta)
345  b1 = ( .5d0*fstrdynamic%gamma/fstrdynamic%beta - 1.d0 )*t_delta
346  b2 = fstrdynamic%gamma/fstrdynamic%beta - 1.d0
347  b3 = fstrdynamic%gamma/(fstrdynamic%beta*t_delta)
348  c1 = 1.d0 + fstrdynamic%ray_k*b3
349  c2 = a3 + fstrdynamic%ray_m*b3
350 
351  coef(1) = a1; coef(2) = a2; coef(3) = a3
352  coef(4) = b1; coef(5) = b2; coef(6) = b3
353 
354  if(hecmesh%my_rank==0) then
355  write(*,'(A)')'-------------------------------------------------'
356  write(*,'('' time step='',i10,'' time='',1pe13.4e3)') istep,t_curr
357  endif
358 
359  do j = 1 ,ndof*nnod
360  fstrdynamic%VEC1(j) = a1*fstrdynamic%ACC(j,1) + a2*fstrdynamic%VEL(j,1)
361  fstrdynamic%VEC2(j) = b1*fstrdynamic%ACC(j,1) + b2*fstrdynamic%VEL(j,1)
362  enddo
363 
364  count_step = 0
365  stepcnt = 0
366 
367  !C for couple analysis
368  do
369  fstrsolid%dunode(:) =0.d0
370  ! call fstr_UpdateEPState( hecMESH, fstrSOLID )
371  call fstr_solve_dynamic_nlimplicit_couple_init(fstrparam, fstrcpl)
372 
373  loopforcontactanalysis: do while( .true. )
374  count_step = count_step + 1
375 
376  ! ----- Inner Iteration
377  res0 = 0.d0
378  res1 = 0.d0
379  relres = 1.d0
380 
381  do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
382  stepcnt=stepcnt+1
383  call fstr_creatematrix_and_dampingforce( hecmesh, hecmat, fstrsolid, t_curr, t_delta, fstrdynamic, coef )
384 
385  call fstr_addspring(cstep, hecmesh, hecmat, fstrsolid, fstrparam)
386 
387  !C-- mechanical boundary condition
388  call dynamic_mat_ass_load (cstep, t_curr+t_delta, hecmesh, hecmat, fstrsolid, fstrdynamic, fstrparam, iter )
389  do j=1, hecmesh%n_node* hecmesh%n_dof
390  hecmat%B(j) = hecmat%B(j) - fstrsolid%QFORCE(j) + fstrsolid%DFORCE(j)
391  enddo
392 
393  call fstr_update_ndforce_spring( cstep, hecmesh, fstrsolid, hecmat%B )
394 
395  !C for couple analysis
396  call fstr_solve_dynamic_nlimplicit_couple_pre(hecmesh, hecmat, fstrsolid, &
397  & fstrparam, fstrdynamic, fstrcpl, restart_step_num, istep)
398 
399  call hecmw_mat_clear( conmat )
400  call hecmw_mat_clear_b( conmat )
401  conmat%X = 0.0d0
402 
403  if( fstr_is_contact_active() ) then
404  call fstr_update_ndforce_contact(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
405  call fstr_addcontactstiffness(cstep,ctalgo,iter,hecmesh,conmat,heclagmat,fstrsolid)
406  endif
407 
408  !C-- geometrical boundary condition
409  call dynamic_mat_ass_bc (cstep, hecmesh, hecmat, fstrsolid, fstrdynamic, &
410  & fstrparam, heclagmat, t_curr+t_delta, stepcnt, conmat=conmat)
411  call dynamic_mat_ass_bc_vl(cstep, hecmesh, hecmat, fstrsolid, fstrdynamic, &
412  & fstrparam, heclagmat, t_curr+t_delta, stepcnt, conmat=conmat)
413  call dynamic_mat_ass_bc_ac(cstep, hecmesh, hecmat, fstrsolid, fstrdynamic, &
414  & fstrparam, heclagmat, t_curr+t_delta, stepcnt, conmat=conmat)
415 
416  ! ----- check convergence
417  call fstr_assemble_residual_contact(hecmat, heclagmat, conmat, hecmesh, resid_work, nresid)
418 
419  if( .not.fstr_is_contact_active() ) then
420  maxdlag = 0.0d0
421  elseif( abs(maxdlag) < 1.0d-15) then
422  maxdlag = 1.0d0
423  endif
424  call hecmw_allreduce_r1(hecmesh, maxdlag, hecmw_max)
425 
426  call fstr_check_convergence(hecmesh, hecmat, fstrsolid, fstrpr, &
427  ndof, iter, istep, cstep, &
428  resid_work, nresid, &
429  res0, res, &
430  n_node_global, &
431  iterstatus, &
432  maxdlag, converg_dlag)
433  if (iterstatus == kitrconverged) exit
434  if (iterstatus == kitrdiverged .or. iterstatus == kitrfloatingerror) then
435  fstrsolid%NRstat_i(knstciter) = count_step
436  return
437  endif
438 
439  ! ---- For Parallel Contact with Multi-Partition Domains
440  hecmat%X = 0.0d0
441  call fstr_set_current_config_to_mesh(hecmesh,fstrsolid,coord)
442  call solve_lineq_contact(hecmesh,hecmat,heclagmat,conmat,istat,1.0d0,fstr_is_contact_active())
443  call fstr_recover_initial_config_to_mesh(hecmesh,fstrsolid,coord)
444  ! ----- check matrix solver error
445  if( istat /= 0 ) then
446  if( hecmesh%my_rank == 0) then
447  write( *,'(a,i5,a,i5)') ' ### Fail to Converge : at total_step=', cstep, ' sub_step=', istep
448  end if
449  fstrsolid%NRstat_i(knstdresn) = 4
450  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
451  return
452  end if
453 
454  ! ----- update external nodal displacement increments
455  call hecmw_update_r (hecmesh, hecmat%X, hecmat%NP, hecmat%NDOF)
456 
457  ! ----- update the strain, stress, and internal force
458  do j=1,hecmesh%n_node*ndof
459  fstrsolid%dunode(j) = fstrsolid%dunode(j)+hecmat%X(j)
460  enddo
461  call fstr_updatenewton( hecmesh, hecmat, fstrsolid, t_curr, &
462  & t_delta,iter, fstrdynamic%strainEnergy )
463 
464  ! ----- update reaction force at constrained DOFs using converged QFORCE
465  call fstr_update_reaction_spc( cstep, hecmesh, fstrsolid )
466 
467  if(.not. fstrparam%nlgeom) exit
468 
469  ! ----- update the Lagrange multipliers
470  if( fstr_is_contact_active() ) then
471  maxdlag = 0.0d0
472  do j=1,heclagmat%num_lagrange
473  heclagmat%lagrange(j) = heclagmat%lagrange(j) + hecmat%X(hecmesh%n_node*ndof+j)
474  if(dabs(hecmat%X(hecmesh%n_node*ndof+j))>maxdlag) maxdlag=dabs(hecmat%X(hecmesh%n_node*ndof+j))
475  ! write(*,*)'Lagrange:', j,hecLagMAT%lagrange(j),hecMAT%X(hecMESH%n_node*ndof+j)
476  enddo
477  endif
478  enddo
479 
480  fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter) ! logging newton iteration(maxtier)
481  fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter ! logging newton iteration(sum of iter)
482 
483  ! ----- compute CONT_NFORCE/CONT_FRIC for output
484  if( fstr_is_contact_active() ) &
485  call fstr_calc_contact_output_force(cstep,ctalgo,hecmesh,heclagmat,fstrsolid,conmat)
486 
487  call fstr_scan_contact_state(cstep, istep, count_step, t_delta, ctalgo, &
488  hecmesh, fstrsolid, infoctchange)
489 
490  if( hecmat%Iarray(99)==4 .and. .not. fstr_is_contact_active() ) then
491  write(*,*) ' This type of direct solver is not yet available in such case ! '
492  write(*,*) ' Please use intel MKL direct solver !'
493  call hecmw_abort(hecmw_comm_get_comm())
494  endif
495 
496  is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid,hecmesh)
497  contact_changed_global=0
498  if( fstr_is_contact_conv(ctalgo,infoctchange,hecmesh) ) then
499  exit loopforcontactanalysis
500  elseif( fstr_is_matrixstructure_changed(infoctchange) ) then
501  call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat, fstr_is_contact_active())
502  contact_changed_global=1
503  endif
504  call hecmw_allreduce_i1(hecmesh,contact_changed_global,hecmw_max)
505  if (contact_changed_global > 0) then
506  call hecmw_mat_clear_b( hecmat )
507  call hecmw_mat_clear_b( conmat )
508  call solve_lineq_contact_init(hecmesh,hecmat,heclagmat,is_mat_symmetric)
509  endif
510 
511  if( count_step > max_iter_contact ) exit loopforcontactanalysis
512 
513  ! ----- check divergence
514  if( count_step >= fstrsolid%step_ctrl(cstep)%max_contiter ) then
515  if( hecmesh%my_rank == 0) then
516  write( *,'(a,i5,a,i5)') ' ### Contact failed to Converge : at total_step=', cstep, ' sub_step=', istep
517  end if
518  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
519  fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
520  fstrsolid%NRstat_i(knstdresn) = 3
521  return
522  end if
523 
524  enddo loopforcontactanalysis
525 
526  fstrsolid%NRstat_i(knstciter) = count_step ! logging contact iteration
527 
528  !C for couple analysis
529  call fstr_solve_dynamic_nlimplicit_couple_post(hecmesh, hecmat, fstrsolid, &
530  & fstrparam, fstrdynamic, fstrcpl, a1, a2, a3, b1, b2, b3, istep, is_cycle)
531  if(is_cycle) cycle
532  exit
533  enddo
534 
535  !C-- new displacement, velocity and acceleration
536  fstrdynamic%kineticEnergy = 0.0d0
537  do j = 1 ,ndof*nnod
538  fstrdynamic%ACC (j,2) = -a1*fstrdynamic%ACC(j,1) - a2*fstrdynamic%VEL(j,1) + &
539  a3*fstrsolid%dunode(j)
540  fstrdynamic%VEL (j,2) = -b1*fstrdynamic%ACC(j,1) - b2*fstrdynamic%VEL(j,1) + &
541  b3*fstrsolid%dunode(j)
542  fstrdynamic%ACC (j,1) = fstrdynamic%ACC (j,2)
543  fstrdynamic%VEL (j,1) = fstrdynamic%VEL (j,2)
544 
545  fstrsolid%unode(j) = fstrsolid%unode(j)+fstrsolid%dunode(j)
546  fstrdynamic%DISP(j,2) = fstrsolid%unode(j)
547 
548  fstrdynamic%kineticEnergy = fstrdynamic%kineticEnergy + &
549  0.5d0*fstreig%mass(j)*fstrdynamic%VEL(j,2)*fstrdynamic%VEL(j,2)
550  enddo
551 
552  call fstr_updatestate( hecmesh, fstrsolid, t_delta )
553 
554  deallocate(coord)
555  deallocate(resid_work)
556  fstrsolid%CutBack_stat = 0
557  end subroutine fstr_newton_dynamic_contactslag
558 
559  subroutine fstr_solve_dynamic_nlimplicit_couple_init(fstrPARAM, fstrCPL)
560  implicit none
561  type(fstr_param) :: fstrparam
562  type(fstr_couple) :: fstrCPL
563  if( fstrparam%fg_couple == 1) then
564  if( fstrparam%fg_couple_type==1 .or. &
565  fstrparam%fg_couple_type==3 .or. &
566  fstrparam%fg_couple_type==5 ) call fstr_rcap_get( fstrcpl )
567  endif
569 
570  subroutine fstr_solve_dynamic_nlimplicit_couple_pre(hecMESH, hecMAT, fstrSOLID, &
571  & fstrPARAM, fstrDYNAMIC, fstrCPL, restart_step_num, i)
572  implicit none
573  type(hecmwst_local_mesh) :: hecMESH
574  type(hecmwst_matrix) :: hecMAT
575  type(fstr_solid) :: fstrSOLID
576  type(fstr_param) :: fstrPARAM
577  type(fstr_dynamic) :: fstrDYNAMIC
578  type(fstr_couple) :: fstrCPL
579  integer(kint) :: kkk0, kkk1, j, kk, i, restart_step_num
580  real(kreal) :: bsize
581 
582  if( fstrparam%fg_couple == 1) then
583  if( fstrparam%fg_couple_first /= 0 ) then
584  bsize = dfloat( i ) / dfloat( fstrparam%fg_couple_first )
585  if( bsize > 1.0 ) bsize = 1.0
586  do kkk0 = 1, fstrcpl%coupled_node_n
587  kkk1 = 3 * kkk0
588  fstrcpl%trac(kkk1-2) = bsize * fstrcpl%trac(kkk1-2)
589  fstrcpl%trac(kkk1-1) = bsize * fstrcpl%trac(kkk1-1)
590  fstrcpl%trac(kkk1 ) = bsize * fstrcpl%trac(kkk1 )
591  enddo
592  endif
593  if( fstrparam%fg_couple_window > 0 ) then
594  j = i - restart_step_num + 1
595  kk = fstrdynamic%n_step - restart_step_num + 1
596  bsize = 0.5*(1.0-cos(2.0*pi*dfloat(j)/dfloat(kk)))
597  do kkk0 = 1, fstrcpl%coupled_node_n
598  kkk1 = 3 * kkk0
599  fstrcpl%trac(kkk1-2) = bsize * fstrcpl%trac(kkk1-2)
600  fstrcpl%trac(kkk1-1) = bsize * fstrcpl%trac(kkk1-1)
601  fstrcpl%trac(kkk1 ) = bsize * fstrcpl%trac(kkk1 )
602  enddo
603  endif
604  call dynamic_mat_ass_couple( hecmesh, hecmat, fstrsolid, fstrcpl )
605  endif
607 
608  subroutine fstr_solve_dynamic_nlimplicit_couple_post(hecMESH, hecMAT, fstrSOLID, &
609  & fstrPARAM, fstrDYNAMIC, fstrCPL, a1, a2, a3, b1, b2, b3, i, is_cycle)
610  implicit none
611  type(hecmwst_local_mesh) :: hecMESH ! unused - kept for interface compatibility
612  type(hecmwst_matrix) :: hecMAT ! unused - kept for interface compatibility
613  type(fstr_solid) :: fstrSOLID
614  type(fstr_param) :: fstrPARAM
615  type(fstr_dynamic) :: fstrDYNAMIC
616  type(fstr_couple) :: fstrCPL
617  integer(kint) :: kkk0, kkk1, j, i, revocap_flag
618  real(kreal) :: a1, a2, a3, b1, b2, b3
619  logical :: is_cycle
620 
621  is_cycle = .false.
622 
623  if( fstrparam%fg_couple == 1 ) then
624  if( fstrparam%fg_couple_type>1 ) then
625  do j=1, fstrcpl%coupled_node_n
626  if( fstrcpl%dof == 3 ) then
627  kkk0 = j*3
628  kkk1 = fstrcpl%coupled_node(j)*3
629 
630  fstrcpl%disp (kkk0-2) = fstrsolid%unode(kkk1-2) + fstrsolid%dunode(kkk1-2)
631  fstrcpl%disp (kkk0-1) = fstrsolid%unode(kkk1-1) + fstrsolid%dunode(kkk1-1)
632  fstrcpl%disp (kkk0 ) = fstrsolid%unode(kkk1 ) + fstrsolid%dunode(kkk1 )
633 
634  fstrcpl%velo (kkk0-2) = -b1*fstrdynamic%ACC(kkk1-2,1) - b2*fstrdynamic%VEL(kkk1-2,1) + &
635  b3*fstrsolid%dunode(kkk1-2)
636  fstrcpl%velo (kkk0-1) = -b1*fstrdynamic%ACC(kkk1-1,1) - b2*fstrdynamic%VEL(kkk1-1,1) + &
637  b3*fstrsolid%dunode(kkk1-1)
638  fstrcpl%velo (kkk0 ) = -b1*fstrdynamic%ACC(kkk1,1) - b2*fstrdynamic%VEL(kkk1,1) + &
639  b3*fstrsolid%dunode(kkk1)
640  fstrcpl%accel(kkk0-2) = -a1*fstrdynamic%ACC(kkk1-2,1) - a2*fstrdynamic%VEL(kkk1-2,1) + &
641  a3*fstrsolid%dunode(kkk1-2)
642  fstrcpl%accel(kkk0-1) = -a1*fstrdynamic%ACC(kkk1-1,1) - a2*fstrdynamic%VEL(kkk1-1,1) + &
643  a3*fstrsolid%dunode(kkk1-1)
644  fstrcpl%accel(kkk0 ) = -a1*fstrdynamic%ACC(kkk1,1) - a2*fstrdynamic%VEL(kkk1,1) + &
645  a3*fstrsolid%dunode(kkk1)
646  else
647  kkk0 = j*2
648  kkk1 = fstrcpl%coupled_node(j)*2
649 
650  fstrcpl%disp (kkk0-1) = fstrsolid%unode(kkk1-1) + fstrsolid%dunode(kkk1-1)
651  fstrcpl%disp (kkk0 ) = fstrsolid%unode(kkk1 ) + fstrsolid%dunode(kkk1 )
652 
653  fstrcpl%velo (kkk0-1) = -b1*fstrdynamic%ACC(kkk1-1,1) - b2*fstrdynamic%VEL(kkk1-1,1) + &
654  b3*fstrsolid%dunode(kkk1-1)
655  fstrcpl%velo (kkk0 ) = -b1*fstrdynamic%ACC(kkk1,1) - b2*fstrdynamic%VEL(kkk1,1) + &
656  b3*fstrsolid%dunode(kkk1)
657  fstrcpl%accel(kkk0-1) = -a1*fstrdynamic%ACC(kkk1-1,1) - a2*fstrdynamic%VEL(kkk1-1,1) + &
658  a3*fstrsolid%dunode(kkk1-1)
659  fstrcpl%accel(kkk0 ) = -a1*fstrdynamic%ACC(kkk1,1) - a2*fstrdynamic%VEL(kkk1,1) + &
660  a3*fstrsolid%dunode(kkk1)
661  endif
662  enddo
663  call fstr_rcap_send( fstrcpl )
664  endif
665 
666  select case ( fstrparam%fg_couple_type )
667  case (4)
668  call fstr_rcap_get( fstrcpl )
669  case (5)
670  call fstr_get_convergence( revocap_flag )
671  if( revocap_flag==0 ) is_cycle = .true.
672  case (6)
673  call fstr_get_convergence( revocap_flag )
674  if( revocap_flag==0 ) then
675  call fstr_rcap_get( fstrcpl )
676  is_cycle = .true.
677  else
678  if( i /= fstrdynamic%n_step ) call fstr_rcap_get( fstrcpl )
679  endif
680  end select
681  endif
683 
684 end module fstr_dynamic_nlimplicit
This module contains subroutines for nonlinear implicit dynamic analysis.
subroutine fstr_solve_nlgeom_dynamic_implicit_contactslag(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrRESULT, fstrPARAM, fstrCPL, hecLagMAT, restart_step_num, restart_substep_num, infoCTChange, conMAT, restart_step_count)
This subroutine provides function of nonlinear implicit dynamic analysis using the Newmark method....
subroutine fstr_solve_dynamic_nlimplicit_couple_init(fstrPARAM, fstrCPL)
subroutine fstr_solve_dynamic_nlimplicit_couple_pre(hecMESH, hecMAT, fstrSOLID, fstrPARAM, fstrDYNAMIC, fstrCPL, restart_step_num, i)
subroutine fstr_solve_dynamic_nlimplicit_couple_post(hecMESH, hecMAT, fstrSOLID, fstrPARAM, fstrDYNAMIC, fstrCPL, a1, a2, a3, b1, b2, b3, i, is_cycle)
real(kind=kreal), parameter pi
subroutine fstr_newton_dynamic_contactslag(cstep, hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrPARAM, fstrCPL, hecLagMAT, infoCTChange, conMAT, restart_step_num, hecMAT0, istep, t_curr, t_delta)
This module provides functions of reconstructing.
subroutine, public fstr_mat_con_contact(cstep, contact_algo, hecMAT, fstrSOLID, hecLagMAT, infoCTChange, conMAT, is_contact_active_flag)
this subroutine reconstructs node-based (stiffness) matrix structure \corresponding to contact state
subroutine, public fstr_save_originalmatrixstructure(hecMAT)
This subroutine saves original matrix structure constructed originally by hecMW_matrix.
logical function, public fstr_is_matrixstruct_symmetric(fstrSOLID, hecMESH)
this function judges whether sitiffness matrix is symmetric or not
This module provides functions to initialize variables when initial velocity or acceleration boundary...
subroutine dynamic_init_varibles(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYNAMIC, fstrPARAM)
This module contains functions to set acceleration boundary condition in dynamic analysis.
subroutine dynamic_mat_ass_bc_ac(cstep, hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, t_curr, iter, conMAT)
This subrouitne set acceleration boundary condition in dynamic analysis.
This module contains functions to set velocity boundary condition in dynamic analysis.
subroutine dynamic_mat_ass_bc_vl(cstep, hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, t_curr, iter, conMAT)
This subrouitne set velocity boundary condition in dynamic analysis.
This module contains functions to set displacement boundary condition in dynamic analysis.
subroutine dynamic_mat_ass_bc(cstep, hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, hecLagMAT, t_curr, iter, conMAT)
This subroutine setup disp bundary condition.
This module contains functions relates to coupling analysis.
subroutine dynamic_mat_ass_couple(hecMESH, hecMAT, fstrSOLID, fstrCPL)
This module contains function to set boundary condition of external load in dynamic analysis.
subroutine dynamic_mat_ass_load(cstep, t_curr, hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, fstrPARAM, iter)
This function sets boundary condition of external load.
This module provides functions to output result.
subroutine fstr_dynamic_output(cstep, istep, t_curr, hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM, outflag)
Output result.
subroutine dynamic_output_monit(cstep, istep, t_curr, hecMESH, fstrPARAM, fstrDYNAMIC, fstrEIG, fstrSOLID)
This module assembles the tangent stiffness matrix and, in the implicit dynamic case,...
subroutine, public fstr_creatematrix_and_dampingforce(hecMESH, hecMAT, fstrSOLID, time, tincr, fstrDYNAMIC, coef)
Assemble the system matrix and, optionally, the dynamic damping force.
This module provides functions to deal with cutback.
Definition: fstr_Cutback.f90:7
subroutine fstr_cutback_save(fstrSOLID, infoCTChange, infoCTChange_bak)
Save analysis status.
subroutine fstr_cutback_load(fstrSOLID, infoCTChange, infoCTChange_bak)
Load analysis status.
subroutine fstr_cutback_init(hecMESH, fstrSOLID, fstrPARAM)
Initializer of cutback variables.
logical function fstr_cutback_active()
Set up lumped mass matrix.
subroutine setmass(fstrSOLID, hecMESH, hecMAT, fstrEIG)
This module provides a unified convergence check for Newton iteration.
subroutine, public fstr_check_convergence(hecMESH, hecMAT, fstrSOLID, fstrPR, ndof, iter, sub_step, cstep, residual_vec, nresid, resb, res_prev, n_node_global, iterStatus, maxDLag, converg_dlag)
Wrapper that calls fstr_check_convergence_main and applies the common divergence/NaN handling (status...
subroutine, public fstr_rcap_send(fstrCPL)
subroutine, public fstr_rcap_get(fstrCPL)
subroutine fstr_get_convergence(revocap_flag)
This module provides function to calculate residual of nodal force.
subroutine, public fstr_assemble_residual_contact(hecMAT, hecLagMAT, conMAT, hecMESH, resid_vec, nresid)
Assemble contact residual vector (hecMATB + conMATB + Lagrange) into a single vector.
subroutine, public fstr_update_reaction_spc(cstep, hecMESH, fstrSOLID)
Set fstrSOLIDREACTION at constrained DOFs using current fstrSOLIDQFORCE. Constrained DOFs are enumera...
This module provides functions to read in and write out restart files.
Definition: fstr_Restart.f90:8
subroutine fstr_write_restart_dyna_nl(cstep, substep, hecMESH, fstrSOLID, fstrDYNAMIC, fstrPARAM, is_StepFinished, contactNode, step_count)
write out restart file for nonlinear dynamic analysis
This module provides functions to deal with spring force.
Definition: fstr_Spring.f90:7
subroutine fstr_update_ndforce_spring(cstep, hecMESH, fstrSOLID, B)
Definition: fstr_Spring.f90:46
subroutine fstr_addspring(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
Definition: fstr_Spring.f90:12
This module provides functions to deal with time and increment of stress analysis.
real(kind=kreal) function fstr_get_timeinc()
logical function fstr_timeinc_istimepoint(stepinfo, fstrPARAM)
subroutine fstr_timeinc_settimeincrement(stepinfo, fstrPARAM, substep, NRstatI, NRstatR, AutoINC_stat, Cutback_stat)
real(kind=kreal) function fstr_get_time()
subroutine fstr_proceed_time()
subroutine fstr_timeinc_printstatus_final(success_flag)
subroutine fstr_timeinc_printstatus_init
logical function fstr_timeinc_isstepfinished(stepinfo)
subroutine fstr_timeinc_printstatus(stepinfo, fstrPARAM, totstep, substep, NRstatI, NRstatR, AutoINC_stat, Cutback_stat)
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
subroutine fstr_recover_initial_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.F90:1146
integer(kind=kint), parameter imsg
Definition: m_fstr.F90:112
integer(kind=kint), parameter kitrfloatingerror
Definition: m_fstr.F90:94
integer(kind=kint), parameter kitrconverged
Definition: m_fstr.F90:92
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.F90:61
subroutine fstr_set_current_config_to_mesh(hecMESH, fstrSOLID, coord)
Definition: m_fstr.F90:1133
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.F90:212
integer(kind=kint), parameter kitrdiverged
Definition: m_fstr.F90:93
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.
Top-level contact analysis module (System level)
subroutine fstr_set_contact_active(a)
subroutine fstr_update_ndforce_contact(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for residual vector (conMATB).
logical function fstr_is_matrixstructure_changed(infoCTChange)
logical function fstr_is_contact_conv(ctAlgo, infoCTChange, hecMESH)
subroutine fstr_addcontactstiffness(cstep, ctAlgo, iter, hecMESH, conMAT, hecLagMAT, fstrSOLID)
subroutine fstr_calc_contact_output_force(cstep, ctAlgo, hecMESH, hecLagMAT, fstrSOLID, conMAT)
Compute contact forces for output (CONT_NFORCE/CONT_FRIC).
subroutine fstr_scan_contact_state(cstep, sub_step, cont_step, dt, ctAlgo, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
logical function fstr_is_contact_active()
Data for coupling analysis.
Definition: m_fstr.F90:636
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Definition: m_fstr.F90:530
Package of data used by Lanczos eigenvalue solver.
Definition: m_fstr.F90:618
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.F90:156