32 real(kind=kreal),
parameter ::
pi = 3.14159265358979323846d0
39 ,fstrDYNAMIC,fstrRESULT,fstrPARAM &
40 ,fstrCPL,hecLagMAT,restart_step_num,restart_substep_num,infoCTChange &
41 ,conMAT,restart_step_count )
44 type(hecmwst_local_mesh) :: hecmesh
45 type(hecmwst_matrix) :: hecMAT
46 type(hecmwst_matrix),
pointer :: hecMAT0
49 type(hecmwst_result_data) :: fstrRESULT
53 type(hecmwst_matrix_lagrange) :: hecLagMAT
54 type(fstr_info_contactchange) :: infoCTChange
55 type(hecmwst_matrix) :: conMAT
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
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
69 logical :: is_outpoint
70 integer(kind=kint) :: n_node_global
71 logical :: is_mat_symmetric, is_interaction_active
75 is_interaction_active = (
associated( fstrsolid%contacts ) .or.
associated( fstrsolid%embeds ) )
80 n_node_global = hecmesh%nn_internal
81 call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
83 ctalgo = fstrparam%contact_algo
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())
91 hecmat%NDOF=hecmesh%n_dof
97 if(
associated( fstrsolid%contacts ) )
call initialize_contact_output_vectors(fstrsolid,hecmat)
100 time_1 = hecmw_wtime()
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'
107 call hecmw_abort( hecmw_comm_get_comm())
111 if(fstrdynamic%idx_mas == 1)
then
112 call setmass(fstrsolid,hecmesh,hecmat,fstreig)
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 !'
119 call hecmw_abort( hecmw_comm_get_comm())
122 hecmat%Iarray(98) = 1
123 hecmat%Iarray(97) = 1
126 if( restart_step_num == 1 .and. fstrdynamic%VarInitialize .and. abs(fstrdynamic%ray_m) > 1.0d-15 ) &
130 if( restart_step_num == 1 )
then
135 fstrdynamic%VEC3(:) =0.d0
140 ctalgo, hecmesh, fstrsolid, infoctchange)
142 call hecmw_mat_copy_profile( hecmat, conmat )
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())
155 fstrsolid%FACTOR = 0.0d0
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
165 sub_step = restart_substep_num
167 if (ndof == 4 .and. hecmesh%my_rank==0)
write(*,
'(a,i5)')
"iter: ",sub_step
170 & fstrsolid%NRstat_i, fstrsolid%NRstat_r, fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
176 fstrdynamic, fstrparam, fstrcpl, heclagmat, infoctchange, conmat, &
177 restart_step_num, hecmat0, sub_step, fstrdynamic%t_curr, fstrdynamic%t_delta)
181 & tot_step_print, sub_step, fstrsolid%NRstat_i, fstrsolid%NRstat_r, &
182 & fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
185 if( fstrsolid%CutBack_stat == 0 )
then
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
197 call hecmw_abort( hecmw_comm_get_comm() )
203 if( is_interaction_active .and. fstrparam%contact_algo ==
kcaslagrange )
then
209 if( hecmesh%my_rank == 0 )
write(*,*)
'### State has been restored at time =',
fstr_get_time()
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=', &
217 call hecmw_abort( hecmw_comm_get_comm())
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
227 if( fstrsolid%CutBack_stat > 0 ) stop
232 step_count = step_count + 1
238 call fstr_dynamic_output(tot_step, step_count, fstrdynamic%t_curr, hecmesh, fstrsolid, fstrdynamic, fstrparam, is_outpoint)
241 call dynamic_output_monit(tot_step, i, fstrdynamic%t_curr, hecmesh, fstrparam, fstrdynamic, fstreig, fstrsolid)
244 if( fstrdynamic%restart_nout > 0 )
then
245 if( mod(step_count,fstrdynamic%restart_nout).eq.0 )
then
247 .false.,infoctchange%contactNode_current,step_count)
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=', &
263 sub_step = sub_step + 1
267 if( fstrdynamic%restart_nout > 0 )
then
269 .true.,infoctchange%contactNode_current,step_count)
272 restart_substep_num = 1
276 if (
associated(hecmat0))
then
277 call hecmw_mat_finalize(hecmat0)
282 if( hecmesh%my_rank == 0 )
then
284 write(
imsg,
'("### FSTR_SOLVE_NLGEOM FINISHED!")')
285 write(*,
'("### FSTR_SOLVE_NLGEOM FINISHED!")')
291 fstrDYNAMIC, fstrPARAM, fstrCPL, hecLagMAT, infoCTChange, conMAT, &
292 restart_step_num, hecMAT0, istep, t_curr, t_delta)
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
305 type(hecmwst_matrix_lagrange) :: hecLagMAT
306 type(fstr_info_contactchange) :: infoCTChange
307 type(hecmwst_matrix) :: conMAT
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
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(:)
328 fstrsolid%NRstat_i(:) = 0
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
339 allocate(coord(hecmesh%n_node*ndof))
340 allocate(resid_work(hecmesh%n_node*ndof + conmat%NP*ndof))
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
351 coef(1) = a1; coef(2) = a2; coef(3) = a3
352 coef(4) = b1; coef(5) = b2; coef(6) = b3
354 if(hecmesh%my_rank==0)
then
355 write(*,
'(A)')
'-------------------------------------------------'
356 write(*,
'('' time step='',i10,'' time='',1pe13.4e3)') istep,t_curr
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)
369 fstrsolid%dunode(:) =0.d0
373 loopforcontactanalysis:
do while( .true. )
374 count_step = count_step + 1
381 do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
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)
397 & fstrparam, fstrdynamic, fstrcpl, restart_step_num, istep)
399 call hecmw_mat_clear( conmat )
400 call hecmw_mat_clear_b( conmat )
410 & fstrparam, heclagmat, t_curr+t_delta, stepcnt, conmat=conmat)
412 & fstrparam, heclagmat, t_curr+t_delta, stepcnt, conmat=conmat)
414 & fstrparam, heclagmat, t_curr+t_delta, stepcnt, conmat=conmat)
421 elseif( abs(maxdlag) < 1.0d-15)
then
424 call hecmw_allreduce_r1(hecmesh, maxdlag, hecmw_max)
427 ndof, iter, istep, cstep, &
428 resid_work, nresid, &
432 maxdlag, converg_dlag)
435 fstrsolid%NRstat_i(knstciter) = count_step
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
449 fstrsolid%NRstat_i(knstdresn) = 4
450 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
455 call hecmw_update_r (hecmesh, hecmat%X, hecmat%NP, hecmat%NDOF)
458 do j=1,hecmesh%n_node*ndof
459 fstrsolid%dunode(j) = fstrsolid%dunode(j)+hecmat%X(j)
462 & t_delta,iter, fstrdynamic%strainEnergy )
467 if(.not. fstrparam%nlgeom)
exit
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))
480 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
481 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
488 hecmesh, fstrsolid, infoctchange)
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())
497 contact_changed_global=0
499 exit loopforcontactanalysis
502 contact_changed_global=1
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 )
511 if( count_step > max_iter_contact )
exit loopforcontactanalysis
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
518 fstrsolid%NRstat_i(knstciter) = count_step
519 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
520 fstrsolid%NRstat_i(knstdresn) = 3
524 enddo loopforcontactanalysis
526 fstrsolid%NRstat_i(knstciter) = count_step
530 & fstrparam, fstrdynamic, fstrcpl, a1, a2, a3, b1, b2, b3, istep, is_cycle)
536 fstrdynamic%kineticEnergy = 0.0d0
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)
545 fstrsolid%unode(j) = fstrsolid%unode(j)+fstrsolid%dunode(j)
546 fstrdynamic%DISP(j,2) = fstrsolid%unode(j)
548 fstrdynamic%kineticEnergy = fstrdynamic%kineticEnergy + &
549 0.5d0*fstreig%mass(j)*fstrdynamic%VEL(j,2)*fstrdynamic%VEL(j,2)
555 deallocate(resid_work)
556 fstrsolid%CutBack_stat = 0
563 if( fstrparam%fg_couple == 1)
then
564 if( fstrparam%fg_couple_type==1 .or. &
565 fstrparam%fg_couple_type==3 .or. &
571 & fstrPARAM, fstrDYNAMIC, fstrCPL, restart_step_num, i)
573 type(hecmwst_local_mesh) :: hecMESH
574 type(hecmwst_matrix) :: hecMAT
579 integer(kint) :: kkk0, kkk1, j, kk, i, restart_step_num
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
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 )
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
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 )
609 & fstrPARAM, fstrDYNAMIC, fstrCPL, a1, a2, a3, b1, b2, b3, i, is_cycle)
611 type(hecmwst_local_mesh) :: hecMESH
612 type(hecmwst_matrix) :: hecMAT
617 integer(kint) :: kkk0, kkk1, j, i, revocap_flag
618 real(kreal) :: a1, a2, a3, b1, b2, b3
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
628 kkk1 = fstrcpl%coupled_node(j)*3
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 )
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)
648 kkk1 = fstrcpl%coupled_node(j)*2
650 fstrcpl%disp (kkk0-1) = fstrsolid%unode(kkk1-1) + fstrsolid%dunode(kkk1-1)
651 fstrcpl%disp (kkk0 ) = fstrsolid%unode(kkk1 ) + fstrsolid%dunode(kkk1 )
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)
666 select case ( fstrparam%fg_couple_type )
671 if( revocap_flag==0 ) is_cycle = .true.
674 if( revocap_flag==0 )
then
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 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.
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.
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.
subroutine fstr_update_ndforce_spring(cstep, hecMESH, fstrSOLID, B)
subroutine fstr_addspring(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
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.
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.
subroutine fstr_recover_initial_config_to_mesh(hecMESH, fstrSOLID, coord)
integer(kind=kint), parameter imsg
integer(kind=kint), parameter kitrfloatingerror
integer(kind=kint), parameter kitrconverged
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
subroutine fstr_set_current_config_to_mesh(hecMESH, fstrSOLID, coord)
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
integer(kind=kint), parameter kitrdiverged
Data for coupling analysis.
Data for DYNAMIC ANSLYSIS (fstrDYNAMIC)
Package of data used by Lanczos eigenvalue solver.
FSTR INNER CONTROL PARAMETERS (fstrPARAM)