30 ,fstrDYN,fstrRESULT,fstrPARAM,infoCTChange &
31 ,fstrCPL, restrt_step_num, restrt_step_count )
33 type(hecmwst_local_mesh) :: hecMESH
34 type(hecmwst_matrix) :: hecMAT
37 type(hecmwst_result_data) :: fstrRESULT
40 type(hecmwst_matrix_lagrange) :: hecLagMAT
41 type(fstr_info_contactchange) :: infoCTChange
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
60 call hecmw_mpc_mat_init_explicit(hecmesh, hecmat, hecmatmpc)
62 hecmat%NDOF=hecmesh%n_dof
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 )
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
76 call hecmw_abort( hecmw_comm_get_comm())
81 fstrsolid%dunode(:) =0.d0
84 & ndof, nnod, restrt_step_count )
86 if( restrt_step_count == 0 )
then
91 if(
associated( fstrsolid%contacts ) )
then
92 call initialize_contact_output_vectors(fstrsolid,hecmat)
94 & fstrdyn%DISP(:,2),fstrsolid%ddunode)
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
101 sub_step = restrt_step_num
104 & fstrsolid%NRstat_i, fstrsolid%NRstat_r, fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
109 step_count = step_count + 1
112 hecmesh, hecmat, hecmatmpc, fstrsolid, fstreig, fstrdyn, fstrparam, &
113 fstrcpl, infoctchange, &
114 restrt_step_num, ndof, nnod, prevb, &
122 call dynamic_output_monit(tot_step, sub_step, fstrdyn%t_curr, hecmesh, fstrparam, fstrdyn, fstreig, fstrsolid)
124 if( fstrdyn%restart_nout > 0 )
then
125 if ( mod(step_count,fstrdyn%restart_nout).eq.0 )
then
127 .false.,infoctchange%contactNode_current,step_count)
138 call fstr_dynamic_output(tot_step, step_count, fstrdyn%t_curr, hecmesh, fstrsolid, fstrdyn, fstrparam, is_outpoint)
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=', &
147 call hecmw_abort( hecmw_comm_get_comm())
150 sub_step = sub_step + 1
153 if( fstrdyn%restart_nout > 0 )
then
155 .true.,infoctchange%contactNode_current,step_count)
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
168 call hecmw_abort( hecmw_comm_get_comm())
173 call hecmw_mpc_mat_finalize_explicit(hecmesh, hecmat, hecmatmpc)
184 ndof, nnod, restrt_step_count )
186 type(hecmwst_local_mesh),
intent(inout) :: hecMESH
187 type(hecmwst_matrix),
intent(inout) :: hecMAT
191 integer(kind=kint),
intent(in) :: ndof
192 integer(kind=kint),
intent(in) :: nnod
193 integer(kind=kint),
intent(in) :: restrt_step_count
195 integer(kind=kint),
allocatable :: mark(:)
196 integer(kind=kint) :: j
197 real(kind=kreal) :: a1, a2
199 a1 = 1.d0/fstrdyn%t_delta**2
200 a2 = 1.d0/(2.d0*fstrdyn%t_delta)
202 call setmass(fstrsolid,hecmesh,hecmat,fstreig)
203 call hecmw_mpc_trans_mass(hecmesh, hecmat, fstreig%mass)
205 allocate(mark(hecmat%NP * hecmat%NDOF))
206 call hecmw_mpc_mark_slave(hecmesh, hecmat, mark)
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
216 call hecmw_abort( hecmw_comm_get_comm())
223 if( restrt_step_count == 0 )
then
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
234 hecMESH, hecMAT, hecMATmpc, fstrSOLID, fstrEIG, fstrDYN, fstrPARAM, &
235 fstrCPL, infoCTChange, &
236 restrt_step_num, ndof, nnod, prevB, &
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
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
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
264 a1 = 1.d0/fstrdyn%t_delta**2
265 a2 = 1.d0/(2.d0*fstrdyn%t_delta)
268 b1 = 0.d0; b2 = 0.d0; b3 = 0.d0
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)
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)
287 if( fstrparam%fg_couple == 1 )
then
288 if( fstrparam%fg_couple_type==1 .or. &
289 fstrparam%fg_couple_type==3 .or. &
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
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 )
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
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 )
316 call hecmw_mpc_trans_rhs(hecmesh, hecmat, hecmatmpc)
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)
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)
338 call hecmw_abort( hecmw_comm_get_comm())
341 call hecmw_mpc_tback_sol(hecmesh, hecmat, hecmatmpc)
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
350 kkk1 = fstrcpl%coupled_node(j)*3
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 )
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) )
370 kkk1 = fstrcpl%coupled_node(j)*2
372 fstrcpl%disp (kkk0-1) = hecmat%X(kkk1-1)
373 fstrcpl%disp (kkk0 ) = hecmat%X(kkk1 )
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) )
388 select case ( fstrparam%fg_couple_type )
393 if( revocap_flag==0 )
then
394 do j = 1, hecmat%NP * ndof
395 hecmat%B(j) = prevb(j)
401 if( revocap_flag==0 )
then
402 do j = 1, hecmat%NP * ndof
403 hecmat%B(j) = prevb(j)
419 fstrsolid%unode(j) = fstrdyn%DISP(j,1)
420 fstrsolid%dunode(j) = hecmat%X(j)-fstrdyn%DISP(j,1)
422 if(
associated( fstrsolid%contacts ) )
then
425 & fstrdyn%DISP(:,2),fstrsolid%ddunode)
427 hecmat%X(j) = hecmat%X(j) + fstrsolid%ddunode(j)
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)
444 call fstr_updatenewton( hecmesh, hecmat, fstrsolid, fstrdyn%t_curr, fstrdyn%t_delta, 0, fstrdyn%strainEnergy )
450 fstrsolid%unode(j) = fstrsolid%unode(j) + fstrsolid%dunode(j)
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
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)
470 if( .not. infoctchange%active )
return
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
485 fstrsolid%contacts(i)%states(j)%multiplier(:) =0.d0
486 fstrsolid%contacts(i)%states(j)%wkdist =0.d0
489 slave = fstrsolid%contacts(i)%slave(j)
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)
497 iss = fstrsolid%contacts(i)%master(sid)%nodes(k)
498 wkarray( iss ) = wkarray( iss ) + shapefunc(k) * fstrsolid%contacts(i)%states(j)%multiplier(1)
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 )
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 )
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 )
534 iss = fstrsolid%contacts(i)%master(sid)%nodes(k)
535 fdum = fdum + shapefunc(k)*shapefunc(k)/mmat( (iss-1)*ndof+1 )
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)
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(:)
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)
556 if( dsqrt(conv)<1.d-8 )
exit
560 do i=1,hecmesh%n_node*ndof
561 uc(i) = wkarray(i)/mmat(i)
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.
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.
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.
This module defines common data and basic structures for analysis.
integer(kind=kint), parameter imsg
integer(kind=kint), parameter idbg
This modules just summarizes all modules used in static analysis.
Data for coupling analysis.
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Package of data used by Lanczos eigenvalue solver.
FSTR INNER CONTROL PARAMETERS (fstrPARAM)