29 subroutine fstr_newton( cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, &
30 restrt_step_num, sub_step, ctime, dtime )
32 integer,
intent(in) :: cstep
33 type (hecmwST_local_mesh) :: hecMESH
34 type (hecmwST_matrix) :: hecMAT
35 type (fstr_solid) :: fstrSOLID
36 integer,
intent(in) :: sub_step
37 real(kind=kreal),
intent(in) :: ctime
38 real(kind=kreal),
intent(in) :: dtime
39 type (fstr_param) :: fstrPARAM
40 type (hecmwST_matrix_lagrange) :: hecLagMAT
42 type (hecmwST_local_mesh),
pointer :: hecMESHmpc
43 type (hecmwST_matrix),
pointer :: hecMATmpc
44 integer(kind=kint) :: ndof
45 integer(kind=kint) :: i, iter
46 integer(kind=kint) :: stepcnt
47 integer(kind=kint) :: restrt_step_num
48 real(kind=kreal) :: tt0, tt, res, qnrm, rres, tincr, xnrm, dunrm, rxnrm
49 real(kind=kreal),
allocatable :: coord(:), p(:)
50 logical :: isLinear = .false.
51 integer(kind=kint) :: iterStatus
53 call hecmw_mpc_mat_init(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
55 if(.not.
fstrpr%nlgeom)
then
59 call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
61 allocate(p(hecmesh%n_node*ndof))
62 allocate(coord(hecmesh%n_node*ndof))
67 do iter=1,fstrsolid%step_ctrl(cstep)%max_iter
74 call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt)
75 call hecmw_mpc_mat_ass(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
76 call hecmw_mpc_trans_rhs(hecmesh, hecmat, hecmatmpc)
79 if( sub_step == restrt_step_num .and. iter == 1 ) hecmatmpc%Iarray(98) = 1
81 hecmatmpc%Iarray(97) = 2
83 hecmatmpc%Iarray(97) = 1
87 call solve_lineq(hecmeshmpc,hecmatmpc)
89 call hecmw_mpc_tback_sol(hecmesh, hecmat, hecmatmpc)
94 do i = 1, hecmesh%n_node*ndof
95 fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
109 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
110 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
114 do i=1,hecmesh%n_node*ndof
115 fstrsolid%unode(i) = fstrsolid%unode(i) + fstrsolid%dunode(i)
120 fstrsolid%CutBack_stat = 0
123 call hecmw_mpc_mat_finalize(hecmesh, hecmat, hecmeshmpc, hecmatmpc)
130 restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
135 integer,
intent(in) :: cstep
136 type (hecmwST_local_mesh) :: hecMESH
137 type (hecmwST_matrix) :: hecMAT
138 type (fstr_solid) :: fstrSOLID
139 integer,
intent(in) :: sub_step
140 real(kind=kreal),
intent(in) :: ctime
141 real(kind=kreal),
intent(in) :: dtime
142 type (fstr_param) :: fstrPARAM
143 type (fstr_info_contactChange) :: infoCTChange
144 type (hecmwST_matrix_lagrange) :: hecLagMAT
145 type (hecmwST_matrix) :: conMAT
147 integer(kind=kint) :: ndof
148 integer(kind=kint) :: ctAlgo
149 integer(kind=kint) :: i, iter
150 integer(kind=kint) :: al_step, n_al_step, stepcnt
151 real(kind=kreal) :: tt0, tt, res, res0, res1, maxv, relres, tincr
152 integer(kind=kint) :: restart_step_num, restart_substep_num
153 logical :: convg, ctchange
154 integer(kind=kint) :: n_node_global
155 integer(kind=kint) :: contact_changed_global
156 real(kind=kreal),
allocatable :: coord(:)
157 integer(kind=kint) :: istat
161 n_node_global = hecmesh%nn_internal
162 call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
164 ctalgo = fstrparam%contact_algo
166 call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
168 if( cstep == 1 .and. sub_step == restart_substep_num )
then
169 call fstr_save_originalmatrixstructure(hecmat)
170 if(hecmesh%my_rank==0)
write(*,*)
"---Scanning initial contact state---"
172 call hecmw_mat_copy_profile( hecmat, conmat )
174 call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat,
fstr_is_contact_active())
175 elseif( hecmat%Iarray(99)==4 )
then
176 write(*, *)
' This type of direct solver is not yet available in such case ! '
177 write(*, *)
' Please change the solver type to intel MKL direct solver !'
178 call hecmw_abort(hecmw_comm_get_comm())
186 allocate(coord(hecmesh%n_node*ndof))
188 call hecmw_mat_clear_b(conmat)
193 n_al_step = fstrsolid%step_ctrl(cstep)%max_contiter
196 do al_step = 1, n_al_step
198 if( hecmesh%my_rank == 0)
then
199 if( n_al_step > 1 )
then
200 print *,
"Augmentation step:", al_step
201 write(
imsg, *)
"Augmentation step:", al_step
210 do iter = 1,fstrsolid%step_ctrl(cstep)%max_iter
216 call hecmw_mat_clear( conmat )
220 if( al_step == 1 .and. stepcnt == 1 )
then
221 maxv = hecmw_mat_diag_max( hecmat, hecmesh )
229 call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
238 call hecmw_update_r (hecmesh, hecmat%X, hecmat%NP, hecmesh%n_dof)
243 do i = 1, hecmesh%n_node*ndof
244 fstrsolid%dunode(i) = fstrsolid%dunode(i)+hecmat%X(i)
250 if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 )
then
255 if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
256 call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
261 call hecmw_mat_clear_b( conmat )
271 res = sqrt(res)/n_node_global
272 if( iter == 1 ) res0 = res
273 if( res0 == 0.0d0 )
then
276 relres = dabs( res1-res )/res0
279 if( hecmesh%my_rank == 0 )
then
280 write(*,
'(a,i3,a,2e15.7)')
' - Residual(',iter,
') =', res, relres
284 if( res < fstrsolid%step_ctrl(cstep)%converg .or. &
285 relres < fstrsolid%step_ctrl(cstep)%converg_ddisp )
then
291 if( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res )
then
292 if( hecmesh%my_rank == 0)
then
293 write( *,
'(a,i5,a,i5)')
' ### Fail to Converge : at total_step=', cstep,
' sub_step=', sub_step
295 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
296 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
297 fstrsolid%NRstat_i(knstciter) = al_step
298 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
299 if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
300 if( res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) fstrsolid%NRstat_i(knstdresn) = 2
307 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
308 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
312 call fstr_scan_contact_state( cstep, sub_step, al_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
314 contact_changed_global = 0
316 call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat,
fstr_is_contact_active())
317 contact_changed_global = 1
319 call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
320 if (contact_changed_global > 0)
then
321 call hecmw_mat_clear_b( hecmat )
322 call hecmw_mat_clear_b( conmat )
329 if( al_step >= fstrsolid%step_ctrl(cstep)%max_contiter )
then
330 if( hecmesh%my_rank == 0)
then
331 write( *,
'(a,i5,a,i5)')
' ### Contact failed to Converge : at total_step=', cstep,
' sub_step=', sub_step
333 fstrsolid%NRstat_i(knstciter) = al_step
334 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
335 fstrsolid%NRstat_i(knstdresn) = 3
343 call hecmw_mat_clear_b( conmat )
352 do i=1,hecmesh%n_node*ndof
353 fstrsolid%unode(i) = fstrsolid%unode(i)+fstrsolid%dunode(i)
356 fstrsolid%NRstat_i(knstciter) = al_step-1
361 fstrsolid%CutBack_stat = 0
368 restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT )
374 integer,
intent(in) :: cstep
375 type (hecmwST_local_mesh) :: hecMESH
376 type (hecmwST_matrix) :: hecMAT
377 type (fstr_solid) :: fstrSOLID
378 integer,
intent(in) :: sub_step
379 real(kind=kreal),
intent(in) :: ctime
380 real(kind=kreal),
intent(in) :: dtime
381 type (fstr_param) :: fstrPARAM
382 type (fstr_info_contactChange) :: infoCTChange
383 type (hecmwST_matrix_lagrange) :: hecLagMAT
384 type (hecmwST_matrix) :: conMAT
386 integer(kind=kint) :: ndof
387 integer(kind=kint) :: ctAlgo
388 integer(kind=kint) :: i, iter, max_iter_contact
389 integer(kind=kint) :: stepcnt, count_step
390 real(kind=kreal) :: tt0, tt, res, res0, res1, relres, tincr, resx
391 integer(kind=kint) :: restart_step_num, restart_substep_num
392 logical :: is_mat_symmetric
393 integer(kind=kint) :: n_node_global
394 integer(kind=kint) :: contact_changed_global
395 integer(kint) :: nndof
396 real(kreal) :: q_residual,x_residual
397 real(kind=kreal),
allocatable :: coord(:)
398 integer(kind=kint) :: istat
400 ctalgo = fstrparam%contact_algo
403 n_node_global = hecmesh%nn_internal
404 call hecmw_allreduce_i1(hecmesh,n_node_global,hecmw_sum)
406 if( hecmat%Iarray(99) == 4 .and. .not. fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh) )
then
407 write(*, *)
' This type of direct solver is not yet available in such case ! '
408 write(*, *)
' Please use intel MKL direct solver !'
409 call hecmw_abort( hecmw_comm_get_comm() )
412 do i=1,fstrsolid%n_contacts
413 fstrsolid%contacts(i)%ctime = ctime + dtime
416 if( cstep==1 .and. sub_step==restart_substep_num )
then
417 call fstr_save_originalmatrixstructure(hecmat)
419 call hecmw_mat_copy_profile( hecmat, conmat )
421 call fstr_mat_con_contact(cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat,
fstr_is_contact_active())
422 elseif( hecmat%Iarray(99)==4 )
then
423 write(*, *)
' This type of direct solver is not yet available in such case ! '
424 write(*, *)
' Please change the solver type to intel MKL direct solver !'
425 call hecmw_abort(hecmw_comm_get_comm())
427 is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
431 call fstr_init_newton(hecmesh, hecmat, fstrsolid, ctime, tincr, iter, cstep, dtime, fstrparam, heclagmat, ndof)
433 call hecmw_mat_clear_b(conmat)
441 allocate(coord(hecmesh%n_node*ndof))
443 loopforcontactanalysis:
do while( .true. )
444 count_step = count_step+1
451 do iter = 1, fstrsolid%step_ctrl(cstep)%max_iter
452 call hecmw_barrier(hecmesh)
453 if(
myrank == 0 ) print *,
'-------------------------------------------------'
454 call hecmw_barrier(hecmesh)
460 call hecmw_mat_clear( conmat )
468 call fstr_addbc(cstep, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, stepcnt, conmat)
470 nndof = hecmat%N*hecmat%ndof
480 if( istat /= 0 )
then
481 if( hecmesh%my_rank == 0)
then
482 write( *,
'(a,i5,a,i5)')
' ### Fail to Converge : at total_step=', cstep,
' sub_step=', sub_step
484 fstrsolid%NRstat_i(knstdresn) = 4
485 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
491 call hecmw_innerproduct_r(hecmesh,ndof,hecmat%X,hecmat%X,resx)
492 resx = sqrt(resx)/n_node_global
494 if( hecmesh%my_rank==0 )
then
495 write(*,
'(a,i3,a,e15.7)')
' - ResidualX (',iter,
') =',resx
496 write(*,
'(a,i3,a,e15.7)')
' - ResidualX+LAG(',iter,
') =',sqrt(x_residual)/n_node_global
497 write(*,
'(a,i3,a,e15.7)')
' - ResidualQ (',iter,
') =',sqrt(q_residual)/n_node_global
501 do i = 1, hecmesh%n_node*ndof
502 fstrsolid%dunode(i) = fstrsolid%dunode(i) + hecmat%X(i)
507 do i = 1, heclagmat%num_lagrange
508 heclagmat%lagrange(i) = heclagmat%lagrange(i)+hecmat%X(hecmesh%n_node*ndof+i)
515 if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 )
then
520 if( fstrsolid%DLOAD_follow /= 0 .or. fstrsolid%CLOAD_ngrp_rot /= 0 ) &
521 call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
526 call hecmw_mat_clear_b( conmat )
532 res = sqrt(res)/n_node_global
533 if( iter == 1 ) res0 = res
534 if( res0 == 0.0d0 )
then
537 relres = dabs( res1-res )/res0
539 if( hecmesh%my_rank == 0 )
then
540 write(*,
'(a,i3,a,2e15.7)')
' - Residual(',iter,
') =',res,relres
544 if( res < fstrsolid%step_ctrl(cstep)%converg .or. &
545 relres < fstrsolid%step_ctrl(cstep)%converg_ddisp )
then
551 if( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res )
then
552 if( hecmesh%my_rank == 0)
then
553 write( *,
'(a,i5,a,i5)')
' ### Fail to Converge : at total_step=', cstep,
' sub_step=', sub_step
555 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
556 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
557 fstrsolid%NRstat_i(knstciter) = count_step
558 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
559 if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
560 if( res > fstrsolid%step_ctrl(cstep)%maxres .or. res /= res ) fstrsolid%NRstat_i(knstdresn) = 2
567 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
568 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
570 call fstr_scan_contact_state( cstep, sub_step, count_step, dtime, ctalgo, hecmesh, fstrsolid, infoctchange, hecmat%B )
573 write(*, *)
' This type of direct solver is not yet available in such case ! '
574 write(*, *)
' Please use intel MKL direct solver !'
575 call hecmw_abort( hecmw_comm_get_comm() )
578 is_mat_symmetric = fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh)
579 contact_changed_global = 0
581 call fstr_mat_con_contact( cstep, ctalgo, hecmat, fstrsolid, heclagmat, infoctchange, conmat,
fstr_is_contact_active())
582 contact_changed_global = 1
587 call hecmw_allreduce_i1(hecmesh, contact_changed_global, hecmw_max)
588 if (contact_changed_global > 0)
then
589 call hecmw_mat_clear_b( hecmat )
590 call hecmw_mat_clear_b( conmat )
595 if( count_step >= fstrsolid%step_ctrl(cstep)%max_contiter )
then
596 if( hecmesh%my_rank == 0)
then
597 write( *,
'(a,i5,a,i5)')
' ### Contact failed to Converge : at total_step=', cstep,
' sub_step=', sub_step
599 fstrsolid%NRstat_i(knstciter) = count_step
600 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
601 fstrsolid%NRstat_i(knstdresn) = 3
609 call hecmw_mat_clear_b( conmat )
613 enddo loopforcontactanalysis
615 fstrsolid%NRstat_i(knstciter) = count_step
619 do i = 1, hecmesh%n_node*ndof
620 fstrsolid%unode(i) = fstrsolid%unode(i)+fstrsolid%dunode(i)
631 fstrsolid%CutBack_stat = 0
634 subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof)
637 type (hecmwST_local_mesh) :: hecMESH
638 type (hecmwST_matrix) :: hecMAT
639 type (fstr_solid) :: fstrSOLID
640 real(kind=kreal),
intent(in) :: ctime
641 real(kind=kreal),
intent(in) :: dtime
642 type (fstr_param) :: fstrPARAM
643 real(kind=kreal),
intent(inout) :: tincr
644 integer(kind=kint) :: iter
645 integer,
intent(in) :: cstep
646 type (hecmwST_matrix_lagrange) :: hecLagMAT
647 integer(kind=kint),
intent(inout) :: ndof
649 hecmat%NDOF = hecmesh%n_dof
653 if( fstrsolid%step_ctrl(cstep)%solution == stepstatic ) tincr = 0.d0
655 fstrsolid%dunode(:) = 0.0d0
656 fstrsolid%NRstat_i(:) = 0
658 call fstr_ass_load(cstep, ctime+dtime, hecmesh, hecmat, fstrsolid, fstrparam)
659 if( fstrsolid%elemact%ELEMACT_egrp_tot > 0 )
then
669 integer(kind=kint) :: iterstatus
671 type (hecmwst_local_mesh) :: hecmesh
672 type (hecmwst_matrix) :: hecmat
674 integer(kind=kint),
intent(in) :: ndof
675 integer(kind=kint),
intent(in) :: iter
676 integer(kind=kint),
intent(in) :: sub_step, cstep
678 real(kind=kreal) :: res, qnrm, rres, xnrm, dunrm, rxnrm
682 call hecmw_innerproduct_r(hecmesh, ndof, hecmat%B, hecmat%B, res)
684 call hecmw_innerproduct_r(hecmesh, ndof, hecmat%X, hecmat%X, xnrm)
686 call hecmw_innerproduct_r(hecmesh, ndof, fstrsolid%QFORCE, fstrsolid%QFORCE, qnrm)
688 if (qnrm < 1.0d-8) qnrm = 1.0d0
692 call hecmw_innerproduct_r(hecmesh, ndof, fstrsolid%dunode, fstrsolid%dunode, dunrm)
698 if( hecmesh%my_rank == 0 )
then
699 if (qnrm == 1.0d0)
then
700 write(*,
"(a,i8,a,1pe11.4,a,1pe11.4)")
" iter:", iter,
", residual(abs):", rres,
", disp.corr.:", rxnrm
702 write(*,
"(a,i8,a,1pe11.4,a,1pe11.4)")
" iter:", iter,
", residual:", rres,
", disp.corr.:", rxnrm
705 if( hecmw_mat_get_flag_diverged(hecmat) ==
kno )
then
706 if( rres < fstrsolid%step_ctrl(cstep)%converg .or. &
707 rxnrm < fstrsolid%step_ctrl(cstep)%converg_ddisp )
then
714 if ( iter == fstrsolid%step_ctrl(cstep)%max_iter .or. rres > fstrsolid%step_ctrl(cstep)%maxres) &
718 if( hecmesh%my_rank == 0)
then
719 write(
ilog,
'(a,i5,a,i5)')
'### Fail to Converge : at total_step=', cstep,
' sub_step=', sub_step
720 write( *,
'(a,i5,a,i5)')
' ### Fail to Converge : at total_step=', cstep,
' sub_step=', sub_step
722 fstrsolid%NRstat_i(knstmaxit) = max(fstrsolid%NRstat_i(knstmaxit),iter)
723 fstrsolid%NRstat_i(knstsumit) = fstrsolid%NRstat_i(knstsumit) + iter
724 fstrsolid%CutBack_stat = fstrsolid%CutBack_stat + 1
725 if( iter == fstrsolid%step_ctrl(cstep)%max_iter ) fstrsolid%NRstat_i(knstdresn) = 1
726 if( rres > fstrsolid%step_ctrl(cstep)%maxres .or. rres /= rres ) fstrsolid%NRstat_i(knstdresn) = 2
This module provides a function to deal with prescribed displacement.
subroutine fstr_addbc(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, iter, conMAT, RHSvector)
Add Essential Boundary Conditions.
This module provides functions to take into account external load.
subroutine fstr_ass_load(cstep, ctime, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
This subroutine assmble following external force into fstrSOLIDGL and hecMATB afterwards.
This module provide a function to elemact elements.
subroutine fstr_update_elemact_solid_by_value(hecMESH, fstrSOLID, cstep, ctime)
This module provides functions on nonlinear analysis.
subroutine fstr_newton_contactslag(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, hecLagMAT, restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method....
subroutine fstr_newton_contactalag(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restart_step_num, restart_substep_num, sub_step, ctime, dtime, infoCTChange, conMAT)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method combined with Neste...
subroutine fstr_newton(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restrt_step_num, sub_step, ctime, dtime)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method.
integer(kind=kint) function fstr_check_iteration_converged(hecMESH, hecMAT, fstrSOLID, ndof, iter, sub_step, cstep)
\breaf This function check iteration status
subroutine fstr_init_newton(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, hecLagMAT, ndof)
This module provides function to calculate residual of nodal force.
subroutine, public fstr_update_ndforce_spc(cstep, hecMESH, fstrSOLID, B)
subroutine, public fstr_update_ndforce(cstep, hecMESH, hecMAT, fstrSOLID, conMAT)
real(kind=kreal) function, public fstr_get_norm_para_contact(hecMAT, hecLagMAT, conMAT, hecMESH)
subroutine fstr_calc_residual_vector(hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM)
\breaf This subroutine calculate residual vector
real(kind=kreal) function, public fstr_get_x_norm_contact(hecMAT, hecLagMAT, hecMESH)
This module provides functions to read in and write out restart files.
This module provides functions to deal with spring force.
subroutine fstr_addspring(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM)
This module provides function to calculate tangent stiffness matrix.
subroutine, public fstr_stiffmatrix(hecMESH, hecMAT, fstrSOLID, time, tincr)
This subroutine creates tangential stiffness matrix.
This module provides function to calculate to do updates.
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) myrank
PARALLEL EXECUTION.
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 ilog
FILE HANDLER.
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 kitrcontinue
iteration control
integer(kind=kint), parameter kno
integer(kind=kint), parameter kitrdiverged
logical paracontactflag
PARALLEL CONTACT FLAG.
This modules just summarizes all modules used in static analysis.
This module provides functions to output result.