FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
fstr_dynamic_nlexplicit.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
9  use m_static_lib
16  use m_fstr_update
17  use m_fstr_residual
18  use m_fstr_restart
20  use m_fstr_rcap_io
21  use mcontact
22  use m_fstr_timeinc
23 
24 contains
25 
26  !C================================================================C
27  !C-- subroutine fstr_solve_LINEAR_DYNAMIC
28  !C================================================================C
29  subroutine fstr_solve_dynamic_nlexplicit(hecMESH,hecMAT,fstrSOLID,fstrEIG &
30  ,fstrDYN,fstrRESULT,fstrPARAM,infoCTChange &
31  ,fstrCPL, restrt_step_num, restrt_step_count )
32  implicit none
33  type(hecmwst_local_mesh) :: hecMESH
34  type(hecmwst_matrix) :: hecMAT
35  type(fstr_eigen) :: fstrEIG
36  type(fstr_solid) :: fstrSOLID
37  type(hecmwst_result_data) :: fstrRESULT
38  type(fstr_param) :: fstrPARAM
39  type(fstr_dynamic) :: fstrDYN
40  type(hecmwst_matrix_lagrange) :: hecLagMAT
41  type(fstr_info_contactchange) :: infoCTChange
42  type(fstr_couple) :: fstrCPL !for COUPLE
43  type(hecmwst_matrix), pointer :: hecMATmpc
44  integer(kind=kint), allocatable :: mark(:)
45  integer(kind=kint) :: nnod, ndof, nn, numnp
46  integer(kind=kint) :: i, j, ids, ide, kk
47  integer(kind=kint) :: kkk0, kkk1
48  integer(kind=kint) :: ierror
49  integer(kind=kint) :: iiii5, iexit
50  integer(kind=kint) :: revocap_flag
51  real(kind=kreal), allocatable :: prevb(:)
52  real(kind=kreal) :: bsize, res
53  real(kind=kreal) :: time_1, time_2
54  integer(kind=kint) :: restrt_step_num
55  integer(kind=kint) :: restrt_step_count
56  real(kind=kreal), parameter :: pi = 3.14159265358979323846d0
57  integer(kind=kint) :: tot_step, sub_step, step_count
58  logical :: is_OutPoint
59 
60  call hecmw_mpc_mat_init_explicit(hecmesh, hecmat, hecmatmpc)
61 
62  hecmat%NDOF=hecmesh%n_dof
63  nnod=hecmesh%n_node
64  ndof=hecmat%NDOF
65  nn=ndof*ndof
66 
67  if( fstrparam%fg_couple == 1) then
68  if( fstrparam%fg_couple_type==5 .or. &
69  fstrparam%fg_couple_type==6 ) then
70  allocate( prevb(hecmat%NP*ndof) ,stat=ierror )
71  prevb = 0.0d0
72  if( ierror /= 0 ) then
73  write(idbg,*) 'stop due to allocation error <fstr_solve_NONLINEAR_DYNAMIC, prevB>'
74  write(idbg,*) ' rank = ', hecmesh%my_rank,' ierror = ',ierror
75  call flush(idbg)
76  call hecmw_abort( hecmw_comm_get_comm())
77  endif
78  endif
79  endif
80 
81  fstrsolid%dunode(:) =0.d0
82 
83  call fstr_prepare_dynamic_explicit( hecmesh, hecmat, fstrsolid, fstreig, fstrdyn, &
84  & ndof, nnod, restrt_step_count )
85 
86  if( restrt_step_count == 0 ) then
87  call fstr_dynamic_output(1, 0, 0.d0, hecmesh, fstrsolid, fstrdyn, fstrparam, .true.)
88  call dynamic_output_monit(1, 0, 0.d0, hecmesh, fstrparam, fstrdyn, fstreig, fstrsolid)
89  end if
90 
91  if( associated( fstrsolid%contacts ) ) then
92  call initialize_contact_output_vectors(fstrsolid,hecmat)
93  call forward_increment_lagrange(1,ndof,fstrdyn%VEC1,hecmesh,fstrsolid,infoctchange,&
94  & fstrdyn%DISP(:,2),fstrsolid%ddunode)
95  endif
96 
97  step_count = restrt_step_count
98  do tot_step = 1, fstrsolid%nstep_tot
99  if(hecmesh%my_rank==0) write(*,'(a,i5)') ' loading step=',tot_step
100 
101  sub_step = restrt_step_num
102  do while(.true.)
103  call fstr_timeinc_settimeincrement( fstrsolid%step_ctrl(tot_step), fstrparam, sub_step, &
104  & fstrsolid%NRstat_i, fstrsolid%NRstat_r, fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
105 
106  fstrdyn%t_curr = fstr_get_time()
107  fstrdyn%t_delta = fstr_get_timeinc()
108 
109  step_count = step_count + 1
110 
111  call fstr_advance_dynamic_explicit( tot_step, sub_step, &
112  hecmesh, hecmat, hecmatmpc, fstrsolid, fstreig, fstrdyn, fstrparam, &
113  fstrcpl, infoctchange, &
114  restrt_step_num, ndof, nnod, prevb, &
115  fstr_timeinc_isstepfinished( fstrsolid%step_ctrl(tot_step) ) )
116 
117  ! ----- Result output (include visualize output)
118  ! Evaluate isTimePoint before time advance, then OR with isStepFinished after
119  is_outpoint = fstr_timeinc_istimepoint( fstrsolid%step_ctrl(tot_step), fstrparam )
120 
121  !C-- output result of monitoring node
122  call dynamic_output_monit(tot_step, sub_step, fstrdyn%t_curr, hecmesh, fstrparam, fstrdyn, fstreig, fstrsolid)
123 
124  if( fstrdyn%restart_nout > 0 ) then
125  if ( mod(step_count,fstrdyn%restart_nout).eq.0 ) then
126  call fstr_write_restart_dyna_nl(tot_step,sub_step,hecmesh,fstrsolid,fstrdyn,fstrparam,&
127  .false.,infoctchange%contactNode_current,step_count)
128  end if
129  end if
130 
131  call fstr_proceed_time()
132  fstrdyn%t_curr = fstr_get_time()
133 
134  ! isStepFinished must be evaluated after fstr_proceed_time
135  is_outpoint = is_outpoint .or. fstr_timeinc_isstepfinished( fstrsolid%step_ctrl(tot_step) )
136 
137  !C-- output new displacement, velocity and acceleration
138  call fstr_dynamic_output(tot_step, step_count, fstrdyn%t_curr, hecmesh, fstrsolid, fstrdyn, fstrparam, is_outpoint)
139 
140  if( fstr_timeinc_isstepfinished( fstrsolid%step_ctrl(tot_step) ) ) exit
141 
142  if( sub_step == fstrsolid%step_ctrl(tot_step)%num_substep ) then
143  if( hecmesh%my_rank == 0 ) then
144  write(*,'(a,i5,a,f6.3)') '### Number of substeps reached max number: at total_step=', &
145  & tot_step, ' time=', fstr_get_time()
146  endif
147  call hecmw_abort( hecmw_comm_get_comm())
148  endif
149 
150  sub_step = sub_step + 1
151  enddo
152 
153  if( fstrdyn%restart_nout > 0 ) then
154  call fstr_write_restart_dyna_nl(tot_step,sub_step,hecmesh,fstrsolid,fstrdyn,fstrparam,&
155  .true.,infoctchange%contactNode_current,step_count)
156  end if
157  restrt_step_num = 1
158  enddo
159 
160  if( fstrparam%fg_couple == 1) then
161  if( fstrparam%fg_couple_type==5 .or. &
162  fstrparam%fg_couple_type==6 ) then
163  deallocate( prevb ,stat=ierror )
164  if( ierror /= 0 ) then
165  write(idbg,*) 'stop due to deallocation error <fstr_solve_NONLINEAR_DYNAMIC, prevB>'
166  write(idbg,*) ' rank = ', hecmesh%my_rank,' ierror = ',ierror
167  call flush(idbg)
168  call hecmw_abort( hecmw_comm_get_comm())
169  endif
170  endif
171  endif
172 
173  call hecmw_mpc_mat_finalize_explicit(hecmesh, hecmat, hecmatmpc)
174 
175  end subroutine fstr_solve_dynamic_nlexplicit
176 
183  subroutine fstr_prepare_dynamic_explicit( hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYN, &
184  ndof, nnod, restrt_step_count )
185  implicit none
186  type(hecmwst_local_mesh), intent(inout) :: hecMESH
187  type(hecmwst_matrix), intent(inout) :: hecMAT
188  type(fstr_solid), intent(inout) :: fstrSOLID
189  type(fstr_eigen), intent(inout) :: fstrEIG
190  type(fstr_dynamic), intent(inout) :: fstrDYN
191  integer(kind=kint), intent(in) :: ndof
192  integer(kind=kint), intent(in) :: nnod
193  integer(kind=kint), intent(in) :: restrt_step_count
194 
195  integer(kind=kint), allocatable :: mark(:)
196  integer(kind=kint) :: j
197  real(kind=kreal) :: a1, a2
198 
199  a1 = 1.d0/fstrdyn%t_delta**2
200  a2 = 1.d0/(2.d0*fstrdyn%t_delta)
201 
202  call setmass(fstrsolid,hecmesh,hecmat,fstreig)
203  call hecmw_mpc_trans_mass(hecmesh, hecmat, fstreig%mass)
204 
205  allocate(mark(hecmat%NP * hecmat%NDOF))
206  call hecmw_mpc_mark_slave(hecmesh, hecmat, mark)
207 
208  do j = 1 ,ndof*nnod
209  fstrdyn%VEC1(j) = (a1 + a2 *fstrdyn%ray_m) * fstreig%mass(j)
210  if(mark(j) == 1) fstrdyn%VEC1(j) = 1.d0
211  if(dabs(fstrdyn%VEC1(j)) < 1.0e-20) then
212  if( hecmesh%my_rank == 0 ) then
213  write(*,*) 'stop due to fstrDYN%VEC(j) = 0 , j = ', j
214  write(imsg,*) 'stop due to fstrDYN%VEC(j) = 0 , j = ', j
215  end if
216  call hecmw_abort( hecmw_comm_get_comm())
217  endif
218  end do
219 
220  deallocate(mark)
221 
222  !C-- virtual past displacements for central-difference startup
223  if( restrt_step_count == 0 ) then
224  do j = 1 ,ndof*nnod
225  fstrdyn%DISP(j,3) = fstrdyn%DISP(j,1) - fstrdyn%VEL (j,1)/(2.d0*a2) + fstrdyn%ACC (j,1)/ (2.d0*a1)
226  fstrdyn%DISP(j,2) = fstrdyn%DISP(j,1) - fstrdyn%VEL (j,1)/ a2 + fstrdyn%ACC (j,1)/ (2.d0*a1) * 4.d0
227  end do
228  endif
229 
230  end subroutine fstr_prepare_dynamic_explicit
231 
233  subroutine fstr_advance_dynamic_explicit( cstep, istep, &
234  hecMESH, hecMAT, hecMATmpc, fstrSOLID, fstrEIG, fstrDYN, fstrPARAM, &
235  fstrCPL, infoCTChange, &
236  restrt_step_num, ndof, nnod, prevB, &
237  is_last_step )
238  implicit none
239  integer(kind=kint), intent(in) :: cstep
240  integer(kind=kint), intent(in) :: istep
241  type(hecmwst_local_mesh), intent(inout) :: hecMESH
242  type(hecmwst_matrix), intent(inout) :: hecMAT
243  type(hecmwst_matrix), pointer, intent(inout) :: hecMATmpc
244  type(fstr_solid), intent(inout) :: fstrSOLID
245  type(fstr_eigen), intent(inout) :: fstrEIG
246  type(fstr_dynamic), intent(inout) :: fstrDYN
247  type(fstr_param), intent(inout) :: fstrPARAM
248  type(fstr_couple), intent(inout) :: fstrCPL
249  type(fstr_info_contactchange), intent(inout) :: infoCTChange
250  integer(kind=kint), intent(in) :: restrt_step_num
251  integer(kind=kint), intent(in) :: ndof
252  integer(kind=kint), intent(in) :: nnod
253  real(kind=kreal), allocatable, intent(inout) :: prevb(:)
254  logical, intent(in) :: is_last_step
255 
256  integer(kind=kint) :: j, kk, kkk0, kkk1
257  integer(kind=kint) :: revocap_flag
258  real(kind=kreal) :: bsize
259  real(kind=kreal), parameter :: pi = 3.14159265358979323846d0
260  real(kind=kreal) :: a1, a2
261  real(kind=kreal) :: b1, b2, b3, a3
262 
263  !C-- central-difference coefficients (depend only on dt)
264  a1 = 1.d0/fstrdyn%t_delta**2
265  a2 = 1.d0/(2.d0*fstrdyn%t_delta)
266  !C-- coupling-only coefficients (kept zero for central difference)
267  a3 = 0.d0
268  b1 = 0.d0; b2 = 0.d0; b3 = 0.d0
269 
270  !C-- mechanical boundary condition
271  call dynamic_mat_ass_load (cstep, fstrdyn%t_curr, hecmesh, hecmat, fstrsolid, fstrdyn, fstrparam)
272  do j=1, hecmesh%n_node* hecmesh%n_dof
273  hecmat%B(j)=hecmat%B(j)-fstrsolid%QFORCE(j)
274  end do
275 
276  !C ********************************************************************************
277  !C for couple analysis
278  if( fstrparam%fg_couple == 1 ) then
279  if( fstrparam%fg_couple_type==5 .or. &
280  fstrparam%fg_couple_type==6 ) then
281  do j = 1, hecmat%NP * ndof
282  prevb(j) = hecmat%B(j)
283  enddo
284  endif
285  endif
286  do
287  if( fstrparam%fg_couple == 1 ) then
288  if( fstrparam%fg_couple_type==1 .or. &
289  fstrparam%fg_couple_type==3 .or. &
290  fstrparam%fg_couple_type==5 ) call fstr_rcap_get( fstrcpl )
291  if( fstrparam%fg_couple_first /= 0 ) then
292  bsize = dfloat( istep ) / dfloat( fstrparam%fg_couple_first )
293  if( bsize > 1.0 ) bsize = 1.0
294  do kkk0 = 1, fstrcpl%coupled_node_n
295  kkk1 = 3 * kkk0
296  fstrcpl%trac(kkk1-2) = bsize * fstrcpl%trac(kkk1-2)
297  fstrcpl%trac(kkk1-1) = bsize * fstrcpl%trac(kkk1-1)
298  fstrcpl%trac(kkk1 ) = bsize * fstrcpl%trac(kkk1 )
299  enddo
300  endif
301  if( fstrparam%fg_couple_window > 0 ) then
302  j = istep - restrt_step_num + 1
303  kk = fstrdyn%n_step - restrt_step_num + 1
304  bsize = 0.5*(1.0-cos(2.0*pi*dfloat(j)/dfloat(kk)))
305  do kkk0 = 1, fstrcpl%coupled_node_n
306  kkk1 = 3 * kkk0
307  fstrcpl%trac(kkk1-2) = bsize * fstrcpl%trac(kkk1-2)
308  fstrcpl%trac(kkk1-1) = bsize * fstrcpl%trac(kkk1-1)
309  fstrcpl%trac(kkk1 ) = bsize * fstrcpl%trac(kkk1 )
310  enddo
311  endif
312  call dynamic_mat_ass_couple( hecmesh, hecmat, fstrsolid, fstrcpl )
313  endif
314  !C ********************************************************************************
315 
316  call hecmw_mpc_trans_rhs(hecmesh, hecmat, hecmatmpc)
317 
318  do j = 1 ,ndof*nnod
319  hecmatmpc%B(j) = hecmatmpc%B(j) + 2.d0*a1* fstreig%mass(j) * fstrdyn%DISP(j,1) &
320  + (- a1 + a2 * fstrdyn%ray_m) * fstreig%mass(j) * fstrdyn%DISP(j,3)
321  end do
322 
323  !C
324  !C-- geometrical boundary condition
325 
326  call dynamic_explicit_ass_bc(hecmesh, hecmatmpc, fstrsolid, fstrdyn, fstrdyn%t_curr)
327  call dynamic_explicit_ass_vl(hecmesh, hecmatmpc, fstrsolid, fstrdyn, fstrdyn%t_curr)
328  call dynamic_explicit_ass_ac(hecmesh, hecmatmpc, fstrsolid, fstrdyn, fstrdyn%t_curr)
329 
330  ! Finish the calculation
331  do j = 1 ,ndof*nnod
332  hecmatmpc%X(j) = hecmatmpc%B(j) / fstrdyn%VEC1(j)
333  if(dabs(hecmatmpc%X(j)) > 1.0d+5) then
334  if( hecmesh%my_rank == 0 ) then
335  print *, 'Displacement increment too large, please adjust your step size!',istep,hecmatmpc%X(j)
336  write(imsg,*) 'Displacement increment too large, please adjust your step size!',istep,hecmatmpc%B(j),fstrdyn%VEC1(j)
337  end if
338  call hecmw_abort( hecmw_comm_get_comm())
339  end if
340  end do
341  call hecmw_mpc_tback_sol(hecmesh, hecmat, hecmatmpc)
342 
343  !C *****************************************************
344  !C for couple analysis
345  if( fstrparam%fg_couple == 1 ) then
346  if( fstrparam%fg_couple_type>1 ) then
347  do j=1, fstrcpl%coupled_node_n
348  if( fstrcpl%dof == 3 ) then
349  kkk0 = j*3
350  kkk1 = fstrcpl%coupled_node(j)*3
351 
352  fstrcpl%disp (kkk0-2) = hecmat%X(kkk1-2)
353  fstrcpl%disp (kkk0-1) = hecmat%X(kkk1-1)
354  fstrcpl%disp (kkk0 ) = hecmat%X(kkk1 )
355 
356  fstrcpl%velo (kkk0-2) = -b1*fstrdyn%ACC(kkk1-2,1) - b2*fstrdyn%VEL(kkk1-2,1) + &
357  b3*( hecmat%X(kkk1-2) - fstrdyn%DISP(kkk1-2,1) )
358  fstrcpl%velo (kkk0-1) = -b1*fstrdyn%ACC(kkk1-1,1) - b2*fstrdyn%VEL(kkk1-1,1) + &
359  b3*( hecmat%X(kkk1-1) - fstrdyn%DISP(kkk1-1,1) )
360  fstrcpl%velo (kkk0 ) = -b1*fstrdyn%ACC(kkk1,1) - b2*fstrdyn%VEL(kkk1,1) + &
361  b3*( hecmat%X(kkk1) - fstrdyn%DISP(kkk1,1) )
362  fstrcpl%accel(kkk0-2) = -a1*fstrdyn%ACC(kkk1-2,1) - a2*fstrdyn%VEL(kkk1-2,1) + &
363  a3*( hecmat%X(kkk1-2) - fstrdyn%DISP(kkk1-2,1) )
364  fstrcpl%accel(kkk0-1) = -a1*fstrdyn%ACC(kkk1-1,1) - a2*fstrdyn%VEL(kkk1-1,1) + &
365  a3*( hecmat%X(kkk1-1) - fstrdyn%DISP(kkk1-1,1) )
366  fstrcpl%accel(kkk0 ) = -a1*fstrdyn%ACC(kkk1,1) - a2*fstrdyn%VEL(kkk1,1) + &
367  a3*( hecmat%X(kkk1) - fstrdyn%DISP(kkk1,1) )
368  else
369  kkk0 = j*2
370  kkk1 = fstrcpl%coupled_node(j)*2
371 
372  fstrcpl%disp (kkk0-1) = hecmat%X(kkk1-1)
373  fstrcpl%disp (kkk0 ) = hecmat%X(kkk1 )
374 
375  fstrcpl%velo (kkk0-1) = -b1*fstrdyn%ACC(kkk1-1,1) - b2*fstrdyn%VEL(kkk1-1,1) + &
376  b3*( hecmat%X(kkk1-1) - fstrdyn%DISP(kkk1-1,1) )
377  fstrcpl%velo (kkk0 ) = -b1*fstrdyn%ACC(kkk1,1) - b2*fstrdyn%VEL(kkk1,1) + &
378  b3*( hecmat%X(kkk1) - fstrdyn%DISP(kkk1,1) )
379  fstrcpl%accel(kkk0-1) = -a1*fstrdyn%ACC(kkk1-1,1) - a2*fstrdyn%VEL(kkk1-1,1) + &
380  a3*( hecmat%X(kkk1-1) - fstrdyn%DISP(kkk1-1,1) )
381  fstrcpl%accel(kkk0 ) = -a1*fstrdyn%ACC(kkk1,1) - a2*fstrdyn%VEL(kkk1,1) + &
382  a3*( hecmat%X(kkk1) - fstrdyn%DISP(kkk1,1) )
383  endif
384  end do
385  call fstr_rcap_send( fstrcpl )
386  endif
387 
388  select case ( fstrparam%fg_couple_type )
389  case (4)
390  call fstr_rcap_get( fstrcpl )
391  case (5)
392  call fstr_get_convergence( revocap_flag )
393  if( revocap_flag==0 ) then
394  do j = 1, hecmat%NP * ndof
395  hecmat%B(j) = prevb(j)
396  enddo
397  cycle
398  endif
399  case (6)
400  call fstr_get_convergence( revocap_flag )
401  if( revocap_flag==0 ) then
402  do j = 1, hecmat%NP * ndof
403  hecmat%B(j) = prevb(j)
404  enddo
405  call fstr_rcap_get( fstrcpl )
406  cycle
407  else
408  if( .not. is_last_step ) call fstr_rcap_get( fstrcpl )
409  endif
410  end select
411  endif
412  exit
413  enddo
414 
415  !C *****************************************************
416  !C-- contact corrector
417  !C
418  do j = 1 ,ndof*nnod
419  fstrsolid%unode(j) = fstrdyn%DISP(j,1)
420  fstrsolid%dunode(j) = hecmat%X(j)-fstrdyn%DISP(j,1)
421  enddo
422  if( associated( fstrsolid%contacts ) ) then
423  !call fstr_scan_contact_state( cstep, fstrDYN%t_delta, kcaSLAGRANGE, hecMESH, fstrSOLID, infoCTChange )
424  call forward_increment_lagrange(cstep,ndof,fstrdyn%VEC1,hecmesh,fstrsolid,infoctchange,&
425  & fstrdyn%DISP(:,2),fstrsolid%ddunode)
426  do j = 1 ,ndof*nnod
427  hecmat%X(j) = hecmat%X(j) + fstrsolid%ddunode(j)
428  enddo
429  endif
430 
431  !C-- new displacement, velocity and acceleration
432  do j = 1 ,ndof*nnod
433  fstrdyn%ACC (j,1) = a1*(hecmat%X(j) - 2.d0*fstrdyn%DISP(j,1) + fstrdyn%DISP(j,3))
434  fstrdyn%VEL (j,1) = a2*(hecmat%X(j) - fstrdyn%DISP(j,3))
435  fstrsolid%unode(j) = fstrdyn%DISP(j,1)
436  fstrsolid%dunode(j) = hecmat%X(j)-fstrdyn%DISP(j,1)
437  fstrdyn%DISP(j,3) = fstrdyn%DISP(j,1)
438  fstrdyn%DISP(j,1) = hecmat%X(j)
439  hecmat%X(j) = fstrsolid%dunode(j)
440  fstrdyn%kineticEnergy = fstrdyn%kineticEnergy + 0.5d0*fstreig%mass(j)*fstrdyn%VEL(j,1)*fstrdyn%VEL(j,1)
441  end do
442 
443  ! ----- update strain, stress, and internal force
444  call fstr_updatenewton( hecmesh, hecmat, fstrsolid, fstrdyn%t_curr, fstrdyn%t_delta, 0, fstrdyn%strainEnergy )
445 
446  ! ----- update reaction force at constrained DOFs using converged QFORCE
447  call fstr_update_reaction_spc( cstep, hecmesh, fstrsolid )
448 
449  do j = 1 ,ndof*nnod
450  fstrsolid%unode(j) = fstrsolid%unode(j) + fstrsolid%dunode(j)
451  end do
452  call fstr_updatestate( hecmesh, fstrsolid, fstrdyn%t_delta )
453 
454  end subroutine fstr_advance_dynamic_explicit
455 
456 
457  subroutine forward_increment_lagrange(cstep,ndof,mmat,hecMESH,fstrSOLID,infoCTChange,wkarray,uc)
458  integer, intent(in) :: cstep
459  integer, intent(in) :: ndof
460  real(kind=kreal), intent(in) :: mmat(:)
461  type( hecmwst_local_mesh ), intent(in) :: hecmesh
462  type(fstr_solid), intent(inout) :: fstrSOLID
463  type(fstr_info_contactchange) :: infoCTChange
464  real(kind=kreal), intent(out) :: wkarray(:)
465  real(kind=kreal), intent(out) :: uc(:)
466  integer :: i, j, k, m, grpid, slave, nn, iSS, sid, etype, iter
467  real(kind=kreal) :: fdum, conv, dlambda, shapefunc(l_max_surface_node), lambda(3)
468 
469  call fstr_scan_contact_state_exp( cstep, hecmesh, fstrsolid, infoctchange )
470  if( .not. infoctchange%active ) return
471 
472  uc = 0.0d0
473 
474  iter = 0
475  do
476  wkarray = 0.0d0
477  do i=1,fstrsolid%n_contacts
478  do j= 1, size(fstrsolid%contacts(i)%slave)
479  if( .not. is_contact_active(fstrsolid%contacts(i)%states(j)%state) ) cycle
480  if( fstrsolid%contacts(i)%states(j)%distance>epsilon(1.d0) ) then
481  fstrsolid%contacts(i)%states(j)%state = contactfree
482  cycle
483  endif
484  if( iter==0 ) then
485  fstrsolid%contacts(i)%states(j)%multiplier(:) =0.d0
486  fstrsolid%contacts(i)%states(j)%wkdist =0.d0
487  cycle
488  endif
489  slave = fstrsolid%contacts(i)%slave(j)
490 
491  sid = fstrsolid%contacts(i)%states(j)%surface
492  nn = size( fstrsolid%contacts(i)%master(sid)%nodes )
493  etype = fstrsolid%contacts(i)%master(sid)%etype
494  call getshapefunc( etype, fstrsolid%contacts(i)%states(j)%lpos(:), shapefunc )
495  wkarray( slave ) = -fstrsolid%contacts(i)%states(j)%multiplier(1)
496  do k=1,nn
497  iss = fstrsolid%contacts(i)%master(sid)%nodes(k)
498  wkarray( iss ) = wkarray( iss ) + shapefunc(k) * fstrsolid%contacts(i)%states(j)%multiplier(1)
499  enddo
500  enddo
501  enddo
502 
503  if(iter > 0)then
504  do i=1,fstrsolid%n_contacts
505  do j= 1, size(fstrsolid%contacts(i)%slave)
506  if( .not. is_contact_active(fstrsolid%contacts(i)%states(j)%state) ) cycle
507  slave = fstrsolid%contacts(i)%slave(j)
508  sid = fstrsolid%contacts(i)%states(j)%surface
509  nn = size( fstrsolid%contacts(i)%master(sid)%nodes )
510  etype = fstrsolid%contacts(i)%master(sid)%etype
511  call getshapefunc( etype, fstrsolid%contacts(i)%states(j)%lpos(:), shapefunc )
512  fstrsolid%contacts(i)%states(j)%wkdist = -wkarray( slave )/mmat( (slave-1)*ndof+1 )
513  do k=1,nn
514  iss = fstrsolid%contacts(i)%master(sid)%nodes(k)
515  fstrsolid%contacts(i)%states(j)%wkdist = fstrsolid%contacts(i)%states(j)%wkdist &
516  + shapefunc(k) * wkarray(iss) / mmat( (iss-1)*ndof+1 )
517  enddo
518  enddo
519  enddo
520  endif
521 
522  conv = 0.d0
523  wkarray = 0.d0
524  do i=1,fstrsolid%n_contacts
525  do j= 1, size(fstrsolid%contacts(i)%slave)
526  if( .not. is_contact_active(fstrsolid%contacts(i)%states(j)%state) ) cycle
527  slave = fstrsolid%contacts(i)%slave(j)
528  sid = fstrsolid%contacts(i)%states(j)%surface
529  nn = size( fstrsolid%contacts(i)%master(sid)%nodes )
530  etype = fstrsolid%contacts(i)%master(sid)%etype
531  call getshapefunc( etype, fstrsolid%contacts(i)%states(j)%lpos(:), shapefunc )
532  fdum = 1.d0/mmat( (slave-1)*ndof+1 )
533  do k=1,nn
534  iss = fstrsolid%contacts(i)%master(sid)%nodes(k)
535  fdum = fdum + shapefunc(k)*shapefunc(k)/mmat( (iss-1)*ndof+1 )
536  enddo
537  dlambda= (fstrsolid%contacts(i)%states(j)%distance-fstrsolid%contacts(i)%states(j)%wkdist) /fdum
538  conv = conv + dlambda*dlambda;
539  fstrsolid%contacts(i)%states(j)%multiplier(1) = fstrsolid%contacts(i)%states(j)%multiplier(1) + dlambda
540  if( fstrsolid%contacts(i)%fcoeff>0.d0 ) then
541  if( fstrsolid%contacts(i)%states(j)%state == contactslip ) then
542  fstrsolid%contacts(i)%states(j)%multiplier(2) = &
543  fstrsolid%contacts(i)%fcoeff * fstrsolid%contacts(i)%states(j)%multiplier(1)
544  else ! stick
545  ! fstrSOLID%contacts(i)%states(j)%multiplier(2) =
546  endif
547  endif
548  lambda = fstrsolid%contacts(i)%states(j)%multiplier(1)* fstrsolid%contacts(i)%states(j)%direction
549  wkarray((slave-1)*ndof+1:(slave-1)*ndof+3) = lambda(:)
550  do k=1,nn
551  iss = fstrsolid%contacts(i)%master(sid)%nodes(k)
552  wkarray((iss-1)*ndof+1:(iss-1)*ndof+3) = wkarray((iss-1)*ndof+1:(iss-1)*ndof+3) -lambda(:)*shapefunc(k)
553  enddo
554  enddo
555  enddo
556  if( dsqrt(conv)<1.d-8 ) exit
557  iter = iter+1
558  enddo
559 
560  do i=1,hecmesh%n_node*ndof
561  uc(i) = wkarray(i)/mmat(i)
562  enddo
563  end subroutine forward_increment_lagrange
564 
565 end module fstr_dynamic_nlexplicit
This module contains subroutines for nonlinear explicit dynamic analysis.
subroutine fstr_solve_dynamic_nlexplicit(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYN, fstrRESULT, fstrPARAM, infoCTChange, fstrCPL, restrt_step_num, restrt_step_count)
subroutine fstr_advance_dynamic_explicit(cstep, istep, hecMESH, hecMAT, hecMATmpc, fstrSOLID, fstrEIG, fstrDYN, fstrPARAM, fstrCPL, infoCTChange, restrt_step_num, ndof, nnod, prevB, is_last_step)
Advance one time step of explicit dynamic analysis.
subroutine forward_increment_lagrange(cstep, ndof, mmat, hecMESH, fstrSOLID, infoCTChange, wkarray, uc)
subroutine fstr_prepare_dynamic_explicit(hecMESH, hecMAT, fstrSOLID, fstrEIG, fstrDYN, ndof, nnod, restrt_step_count)
Prepare initial state for explicit dynamic analysis (central difference).
This module contains functions to set acceleration boundary condition in dynamic analysis.
subroutine dynamic_explicit_ass_ac(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, t_curr, iter)
This module contains functions to set velocity boundary condition in dynamic analysis.
subroutine dynamic_explicit_ass_vl(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, t_curr, iter)
This module contains functions to set displacement boundary condition in dynamic analysis.
subroutine dynamic_explicit_ass_bc(hecMESH, hecMAT, fstrSOLID, fstrDYNAMIC, t_curr, iter)
This subroutine setup disp boundary 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)
Set up lumped mass matrix.
subroutine setmass(fstrSOLID, hecMESH, hecMAT, fstrEIG)
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_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 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()
logical function fstr_timeinc_isstepfinished(stepinfo)
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), parameter imsg
Definition: m_fstr.F90:112
integer(kind=kint), parameter idbg
Definition: m_fstr.F90:113
This modules just summarizes all modules used in static analysis.
Definition: static_LIB.f90:6
Top-level contact analysis module (System level)
subroutine fstr_scan_contact_state_exp(cstep, hecMESH, fstrSOLID, infoCTChange)
Scanning contact state.
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