FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_solve_NLGEOM.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
10  use m_static_output
13  use m_fstr_restart
15  use m_fstr_timeinc
16  use m_fstr_cutback
17  use mcontact
19 
20  implicit none
21 
22 contains
23 
24  !======================================================================!
30  subroutine fstr_solve_nlgeom(hecMESH,hecMAT,fstrSOLID,hecLagMAT,fstrPARAM,conMAT)
31  type (hecmwST_local_mesh) :: hecMESH
32  type (hecmwST_matrix ) :: hecMAT
33  type (fstr_param ) :: fstrPARAM
34  type (fstr_solid ) :: fstrSOLID
35  type (hecmwST_matrix_lagrange) :: hecLagMAT
36  type (fstr_info_contactChange) :: infoCTChange, infoCTChange_bak
37  type (hecmwST_matrix ) :: conMAT
38 
39  integer(kind=kint) :: ndof, nn
40  integer(kind=kint) :: j, i, tot_step, step_count, tot_step_print, CBbound
41  integer(kind=kint) :: sub_step
42  integer(kind=kint) :: restart_step_num, restart_substep_num
43  real(kind=kreal) :: ctime, dtime, endtime, factor
44  real(kind=kreal) :: time_1, time_2
45  logical :: ctchanged, is_OutPoint, is_interaction_active
46 
47  if(hecmesh%my_rank==0) call fstr_timeinc_printstatus_init
48 
49  hecmat%NDOF = hecmesh%n_dof
50 
51  ndof = hecmat%NDOF
52  nn = ndof*ndof
53 
54  is_interaction_active = ( associated( fstrsolid%contacts ) .or. associated( fstrsolid%embeds ) )
55 
56  if( fstrsolid%TEMP_ngrp_tot>0 .and. hecmesh%hecmw_flag_initcon==1 ) then
57  fstrsolid%last_temp = 0.0d0
58  fstrsolid%temperature = 0.0d0
59  do j=1, size(hecmesh%node_init_val_item)
60  i = hecmesh%node_init_val_index(j)
61  fstrsolid%last_temp(j) = hecmesh%node_init_val_item(i)
62  fstrsolid%temperature(j) = hecmesh%node_init_val_item(i)
63  end do
64  endif
65  if( fstrsolid%TEMP_ngrp_tot>0 .and. associated(g_initialcnd) ) then
66  fstrsolid%last_temp = 0.0d0
67  fstrsolid%temperature = 0.0d0
68  do j=1,size(g_initialcnd)
69  if( g_initialcnd(j)%cond_name=="temperature" ) then
70  if( .not. associated(fstrsolid%temperature) ) then
71  allocate( fstrsolid%temperature( hecmesh%n_node ) )
72  allocate( fstrsolid%temp_bak( hecmesh%n_node ) )
73  allocate( fstrsolid%last_temp( hecmesh%n_node ) )
74  endif
75  do i= 1, hecmesh%n_node
76  fstrsolid%last_temp(i) = g_initialcnd(j)%realval(i)
77  fstrsolid%temperature(i) = fstrsolid%last_temp(i)
78  enddo
79  endif
80  end do
81  endif
82 
83  if( associated( fstrsolid%contacts ) ) then
84  call initialize_contact_output_vectors(fstrsolid,hecmat)
85  call setup_contact_elesurf_for_area( 1, hecmesh, fstrsolid )
86  endif
87  if( fstrsolid%n_embeds > 0 ) call initialize_embed_vectors(fstrsolid,hecmat)
88 
89  restart_step_num = 1
90  restart_substep_num = 1
91  fstrsolid%unode = 0.0d0
92  step_count = 0 !**
93  infoctchange%contactNode_previous = 0
94  infoctchange%contactNode_current = 0
95  if( fstrsolid%restart_nout < 0 ) then
96  call fstr_read_restart(restart_step_num,restart_substep_num,step_count,ctime,dtime,hecmesh,fstrsolid, &
97  fstrparam,infoctchange%contactNode_previous)
98  hecmat%Iarray(98) = 1
99  call fstr_set_time( ctime )
100  call fstr_set_timeinc_base( dtime )
101  fstrsolid%restart_nout = - fstrsolid%restart_nout
102  else
103  call fstr_static_output( 1, 0, 0.d0, hecmesh, fstrsolid, fstrparam, fstrpr%solution_type, .true. )
104  endif
105 
106  fstrsolid%FACTOR = 0.0d0
107  call fstr_cutback_init( hecmesh, fstrsolid, fstrparam )
108  call fstr_cutback_save( fstrsolid, infoctchange, infoctchange_bak )
109 
110  do tot_step=1, fstrsolid%nstep_tot
111  tot_step_print = tot_step+restart_step_num-1
112  if(hecmesh%my_rank==0) write(*,*) ''
113  if(hecmesh%my_rank==0) write(*,'(a,i5)') ' loading step=',tot_step_print
114 
115  if( fstrsolid%TEMP_ngrp_tot>0 ) then
116  do j=1, hecmesh%n_node
117  fstrsolid%temp_bak(j) = fstrsolid%temperature(j)
118  end do
119  endif
120  call fstr_updatestate( hecmesh, fstrsolid, 0.0d0 )
121 
122  fstrsolid%unode_bak(:) = fstrsolid%unode(:)
123 
124  ! -------------------------------------------------------------------------
125  ! STEP LOOP
126  ! -------------------------------------------------------------------------
127  sub_step = restart_substep_num
128  do while(.true.)
129 
130  ! ----- time history of factor
131  call fstr_timeinc_settimeincrement( fstrsolid%step_ctrl(tot_step), fstrparam, sub_step, &
132  & fstrsolid%NRstat_i, fstrsolid%NRstat_r, fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
133  if( fstrsolid%TEMP_irres > 0 ) then
134  fstrsolid%FACTOR(1) = 0.d0
135  fstrsolid%FACTOR(2) = 1.d0
136  call table_nlsta(hecmesh,fstrsolid,tot_step,fstr_get_time()+fstr_get_timeinc(), factor)
137  fstrsolid%TEMP_FACTOR = factor
138  else
139  call table_nlsta(hecmesh,fstrsolid,tot_step,fstr_get_time(),factor)
140  fstrsolid%FACTOR(1) = factor
141  call table_nlsta(hecmesh,fstrsolid,tot_step,fstr_get_time()+fstr_get_timeinc(), factor)
142  fstrsolid%FACTOR(2) = factor
143  endif
144 
145  if(hecmesh%my_rank==0) then
146  write(*,'(A,I0,2(A,E12.4))') ' sub_step= ',sub_step,', &
147  & current_time=',fstr_get_time(), ', time_inc=',fstr_get_timeinc()
148  write(*,'(A,2f12.7)') ' loading_factor= ', fstrsolid%FACTOR
149  if( fstrsolid%TEMP_irres > 0 ) write(*,'(A,2f12.7)') ' readtemp_factor= ', fstrsolid%TEMP_FACTOR
150  endif
151 
152  time_1 = hecmw_wtime()
153 
154  ! analysis algorithm ( Newton-Rapshon Method )
155  if( .not. is_interaction_active ) then
156  if( fstrparam%nlsolver_method == knsmnewton ) then
157  call fstr_newton( tot_step, hecmesh, hecmat, fstrsolid, fstrparam, &
158  restart_step_num, sub_step, fstr_get_time(), fstr_get_timeinc() )
159  else if( fstrparam%nlsolver_method == knsmquasinewton ) then
160  call fstr_quasi_newton( tot_step, hecmesh, hecmat, fstrsolid, fstrparam, &
161  restart_step_num, sub_step, fstr_get_time(), fstr_get_timeinc() )
162  endif
163  else
164  if( fstrparam%contact_algo == kcaslagrange ) then
165  call fstr_newton_contactslag( tot_step, hecmesh, hecmat, fstrsolid, fstrparam, heclagmat, &
166  restart_step_num, restart_substep_num, sub_step, fstr_get_time(), fstr_get_timeinc(), infoctchange, conmat )
167  else if( fstrparam%contact_algo == kcaalagrange ) then
168  call fstr_newton_contactalag( tot_step, hecmesh, hecmat, fstrsolid, fstrparam, &
169  restart_step_num, restart_substep_num, sub_step, fstr_get_time(), fstr_get_timeinc(), infoctchange, conmat )
170  endif
171  endif
172 
173  ! Time Increment
174  if( hecmesh%my_rank == 0 ) call fstr_timeinc_printstatus( fstrsolid%step_ctrl(tot_step), fstrparam, &
175  & tot_step_print, sub_step, fstrsolid%NRstat_i, fstrsolid%NRstat_r, &
176  & fstrsolid%AutoINC_stat, fstrsolid%CutBack_stat )
177  if( fstr_cutback_active() ) then
178 
179  if( fstrsolid%CutBack_stat == 0 ) then ! converged
180  call fstr_cutback_save( fstrsolid, infoctchange, infoctchange_bak ) ! save analysis state
181  call fstr_proceed_time() ! current time += time increment
182 
183  else ! not converged
184  cbbound = fstrparam%ainc(fstrsolid%step_ctrl(tot_step)%AincParam_id)%CBbound
185  if( fstrsolid%CutBack_stat == cbbound ) then
186  if( hecmesh%my_rank == 0 ) then
187  write(*,*) 'Number of successive cutback reached max number: ',cbbound
188  call fstr_timeinc_printstatus_final(.false.)
189  endif
190  call hecmw_abort( hecmw_comm_get_comm() )
191  endif
192  call fstr_cutback_load( fstrsolid, infoctchange, infoctchange_bak ) ! load analysis state
193  call fstr_set_contact_active( infoctchange%contactNode_current > 0 )
194 
195  ! restore matrix structure for slagrange contact analysis
196  if( is_interaction_active .and. fstrparam%contact_algo == kcaslagrange ) then
197  call fstr_mat_con_contact( tot_step, fstrparam%contact_algo, hecmat, fstrsolid, heclagmat, &
198  & infoctchange, conmat, fstr_is_contact_active())
199  conmat%B(:) = 0.0d0
200  call solve_lineq_contact_init(hecmesh, hecmat, heclagmat, fstr_is_matrixstruct_symmetric(fstrsolid, hecmesh))
201  endif
202  if( hecmesh%my_rank == 0 ) write(*,*) '### State has been restored at time =',fstr_get_time()
203 
204  !stop if # of substeps reached upper bound.
205  if( sub_step == fstrsolid%step_ctrl(tot_step)%num_substep ) then
206  if( hecmesh%my_rank == 0 ) then
207  write(*,'(a,i5,a,f6.3)') '### Number of substeps reached max number: at total_step=', &
208  & tot_step_print, ' time=', fstr_get_time()
209  endif
210  call hecmw_abort( hecmw_comm_get_comm())
211  endif
212 
213  ! output time
214  time_2 = hecmw_wtime()
215  if( hecmesh%my_rank==0) write(imsg,'(a,",",2(I8,","),f10.2)') &
216  & 'step, substep, solve (sec) :', tot_step_print, sub_step, time_2 - time_1
217  cycle
218  endif
219  else
220  if( fstrsolid%CutBack_stat > 0 ) stop
221  call fstr_proceed_time() ! current time += time increment
222  endif
223 
224  step_count = step_count + 1
225 
226  ! ----- Restart
227  if( fstrsolid%restart_nout > 0) then
228  if( mod(step_count,fstrsolid%restart_nout) == 0 ) then
229  call fstr_write_restart(tot_step,tot_step_print,sub_step,step_count,fstr_get_time(), &
230  & fstr_get_timeinc_base(), hecmesh,fstrsolid,fstrparam,.false.,infoctchange%contactNode_current)
231  endif
232  endif
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  call fstr_static_output( tot_step, step_count, fstr_get_time(), hecmesh, fstrsolid, fstrparam, &
238  & fstrpr%solution_type, is_outpoint )
239 
240  time_2 = hecmw_wtime()
241  if( hecmesh%my_rank==0 ) then
242  write(imsg,'(A,",",2(I8,","),f10.2)') 'step, substep, solve (sec) :', tot_step_print, sub_step, time_2 - time_1
243  write(imsg,'(A,I0,",",1pE15.8)') '### stepcount (for output), time :', step_count, fstr_get_time()
244  endif
245 
246  !if time reached the end time of step, exit loop.
247  if( fstr_timeinc_isstepfinished( fstrsolid%step_ctrl(tot_step) ) ) exit
248 
249  if( sub_step == fstrsolid%step_ctrl(tot_step)%num_substep ) then
250  if( hecmesh%my_rank == 0 ) then
251  write(*,'(a,i5,a,f6.3)') '### Number of substeps reached max number: at total_step=', &
252  & tot_step_print, ' time=', fstr_get_time()
253  endif
254  if( hecmesh%my_rank == 0 ) call fstr_timeinc_printstatus_final(.false.)
255  stop !stop if # of substeps reached upper bound.
256  endif
257 
258  sub_step = sub_step + 1
259  enddo !--- end of substep loop
260 
261  ! ----- Restart at the end of step
262  if( fstrsolid%restart_nout > 0 ) then
263  call fstr_write_restart(tot_step,tot_step_print,sub_step,step_count,fstr_get_time(),fstr_get_timeinc_base(), &
264  & hecmesh,fstrsolid,fstrparam,.true.,infoctchange%contactNode_current)
265  endif
266  restart_substep_num = 1
267  if( fstrsolid%TEMP_irres > 0 ) exit
268  enddo !--- end of tot_step loop
269 
270  call fstr_cutback_finalize( fstrsolid )
271 
272  ! message
273  if(myrank == 0)then
274  call fstr_timeinc_printstatus_final(.true.)
275  write(imsg,'("### FSTR_SOLVE_NLGEOM FINISHED!")')
276  write(*,'("### FSTR_SOLVE_NLGEOM FINISHED!")')
277  endif
278 
279  end subroutine fstr_solve_nlgeom
280 
281  !C================================================================C
284  !C================================================================C
285  subroutine table_nlsta(hecMESH, fstrSOLID, cstep, time, f_t)
286  type ( hecmwST_local_mesh ), intent(in) :: hecMESH
287  type ( fstr_solid ), intent(in) :: fstrSOLID
288  integer(kind=kint), intent(in) :: cstep
289  real(kind=kreal), intent(in) :: time
290  real(kind=kreal), intent(out) :: f_t
291 
292  integer(kind=kint) :: jj_n_amp
293 
294  jj_n_amp = fstrsolid%step_ctrl( cstep )%amp_id
295 
296  if( jj_n_amp <= 0 ) then ! Amplitude not defined
297  f_t = (time-fstrsolid%step_ctrl(cstep)%starttime)/fstrsolid%step_ctrl(cstep)%elapsetime
298  if( f_t>1.d0 ) f_t=1.d0
299  else
300  call table_amp(hecmesh, fstrsolid, cstep, jj_n_amp, time, f_t)
301  endif
302 
303  end subroutine table_nlsta
304 
305 end module m_fstr_solve_nlgeom
m_fstr_cutback::fstr_cutback_init
subroutine fstr_cutback_init(hecMESH, fstrSOLID, fstrPARAM)
Initializer of cutback variables.
Definition: fstr_Cutback.f90:25
m_fstr_solve_nlgeom::fstr_solve_nlgeom
subroutine fstr_solve_nlgeom(hecMESH, hecMAT, fstrSOLID, hecLagMAT, fstrPARAM, conMAT)
This module provides main subroutine for nonlinear calculation.
Definition: fstr_solve_NLGEOM.f90:31
m_fstr_timeinc
This module provides functions to deal with time and increment of stress analysis.
Definition: fstr_Ctrl_TimeInc.f90:7
m_fstr_cutback::fstr_cutback_active
logical function fstr_cutback_active()
Definition: fstr_Cutback.f90:20
m_fstr_timeinc::fstr_set_timeinc_base
subroutine fstr_set_timeinc_base(dtime_base)
Definition: fstr_Ctrl_TimeInc.f90:45
m_fstr_timeinc::fstr_timeinc_printstatus
subroutine fstr_timeinc_printstatus(stepinfo, fstrPARAM, totstep, substep, NRstatI, NRstatR, AutoINC_stat, Cutback_stat)
Definition: fstr_Ctrl_TimeInc.f90:65
m_fstr_solve_nlgeom
This module provides main suboruitne for nonliear calculation.
Definition: fstr_solve_NLGEOM.f90:7
mcontact::initialize_embed_vectors
subroutine initialize_embed_vectors(fstrSOLID, hecMAT)
Definition: fstr_contact.f90:606
m_fstr_quasinewton::fstr_quasi_newton
subroutine fstr_quasi_newton(cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restrt_step_num, sub_step, ctime, dtime)
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method.
Definition: fstr_solve_QuasiNewton.f90:27
m_fstr_cutback::fstr_cutback_load
subroutine fstr_cutback_load(fstrSOLID, infoCTChange, infoCTChange_bak)
Load analysis status.
Definition: fstr_Cutback.f90:204
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_fstr_restart::fstr_read_restart
subroutine fstr_read_restart(cstep, substep, step_count, ctime, dtime, hecMESH, fstrSOLID, fstrPARAM, contactNode)
Read in restart file.
Definition: fstr_Restart.f90:18
m_solve_lineq_contact
This module provides functions to solve sparse system of \linear equitions in the case of contact ana...
Definition: hecmw_solver_contact.f90:9
mcontact::setup_contact_elesurf_for_area
subroutine setup_contact_elesurf_for_area(cstep, hecMESH, fstrSOLID)
Definition: fstr_contact.f90:616
m_fstr_nonlinearmethod::fstr_newton_contactslag
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....
Definition: fstr_solve_NonLinear.f90:366
m_fstr_timeinc::fstr_timeinc_isstepfinished
logical function fstr_timeinc_isstepfinished(stepinfo)
Definition: fstr_Ctrl_TimeInc.f90:121
m_fstr_timeinc::fstr_timeinc_settimeincrement
subroutine fstr_timeinc_settimeincrement(stepinfo, fstrPARAM, substep, NRstatI, NRstatR, AutoINC_stat, Cutback_stat)
Definition: fstr_Ctrl_TimeInc.f90:145
m_fstr_timeinc::fstr_get_time
real(kind=kreal) function fstr_get_time()
Definition: fstr_Ctrl_TimeInc.f90:22
fstr_matrix_con_contact::fstr_is_matrixstruct_symmetric
logical function, public fstr_is_matrixstruct_symmetric(fstrSOLID, hecMESH)
this function judges whether sitiffness matrix is symmetric or not
Definition: fstr_mat_con_contact.f90:240
m_static_output
This module provides functions to output result.
Definition: static_output.f90:6
m_fstr::knsmquasinewton
integer(kind=kint), parameter knsmquasinewton
Definition: m_fstr.f90:56
m_fstr::fstrpr
type(fstr_param), target fstrpr
GLOBAL VARIABLE INITIALIZED IN FSTR_SETUP.
Definition: m_fstr.f90:208
m_fstr::kcaalagrange
integer(kind=kint), parameter kcaalagrange
Definition: m_fstr.f90:60
m_fstr_restart::fstr_write_restart
subroutine fstr_write_restart(cstep, cstep_ext, substep, step_count, ctime, dtime, hecMESH, fstrSOLID, fstrPARAM, is_StepFinished, contactNode)
write out restart file
Definition: fstr_Restart.f90:121
m_fstr::kcaslagrange
integer(kind=kint), parameter kcaslagrange
contact analysis algorithm
Definition: m_fstr.f90:59
m_fstr_nonlinearmethod::fstr_newton_contactalag
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...
Definition: fstr_solve_NonLinear.f90:134
m_fstr::knsmnewton
integer(kind=kint), parameter knsmnewton
nonlinear solver method (nsm)
Definition: m_fstr.f90:55
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
mcontact::initialize_contact_output_vectors
subroutine initialize_contact_output_vectors(fstrSOLID, hecMAT)
Definition: fstr_contact.f90:565
m_fstr::g_initialcnd
type(tinitialcondition), dimension(:), pointer, save g_initialcnd
Definition: m_fstr.f90:151
m_fstr_cutback::fstr_cutback_save
subroutine fstr_cutback_save(fstrSOLID, infoCTChange, infoCTChange_bak)
Save analysis status.
Definition: fstr_Cutback.f90:141
m_fstr_cutback
This module provides functions to deal with cutback.
Definition: fstr_Cutback.f90:7
m_fstr_timeinc::fstr_get_timeinc_base
real(kind=kreal) function fstr_get_timeinc_base()
Definition: fstr_Ctrl_TimeInc.f90:30
mcontact::fstr_is_contact_active
logical function fstr_is_contact_active()
Definition: fstr_contact.f90:52
m_fstr_timeinc::fstr_timeinc_istimepoint
logical function fstr_timeinc_istimepoint(stepinfo, fstrPARAM)
Definition: fstr_Ctrl_TimeInc.f90:132
m_fstr_quasinewton
This module provides functions on nonlinear analysis.
Definition: fstr_solve_QuasiNewton.f90:7
m_fstr::table_amp
subroutine table_amp(hecMESH, fstrSOLID, cstep, jj_n_amp, time, f_t)
Definition: m_fstr.f90:1155
m_fstr_timeinc::fstr_set_time
subroutine fstr_set_time(time)
Definition: fstr_Ctrl_TimeInc.f90:35
fstr_matrix_con_contact
This module provides functions of reconstructing.
Definition: fstr_mat_con_contact.f90:9
m_solve_lineq_contact::solve_lineq_contact_init
subroutine, public solve_lineq_contact_init(hecMESH, hecMAT, hecLagMAT, is_sym)
This subroutine.
Definition: hecmw_solver_contact.f90:30
m_static_lib
This modules just summarizes all modules used in static analysis.
Definition: static_LIB.f90:6
m_fstr_nonlinearmethod
This module provides functions on nonlinear analysis.
Definition: fstr_solve_NonLinear.f90:7
m_fstr_nonlinearmethod::fstr_newton
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.
Definition: fstr_solve_NonLinear.f90:30
m_static_output::fstr_static_output
subroutine fstr_static_output(cstep, istep, time, hecMESH, fstrSOLID, fstrPARAM, flag, outflag)
Output result.
Definition: static_output.f90:13
m_fstr_timeinc::fstr_proceed_time
subroutine fstr_proceed_time()
Definition: fstr_Ctrl_TimeInc.f90:50
mcontact::fstr_set_contact_active
subroutine fstr_set_contact_active(a)
Definition: fstr_contact.f90:56
m_fstr_timeinc::fstr_timeinc_printstatus_final
subroutine fstr_timeinc_printstatus_final(success_flag)
Definition: fstr_Ctrl_TimeInc.f90:110
fstr_matrix_con_contact::fstr_mat_con_contact
subroutine, public fstr_mat_con_contact(cstep, contact_algo, hecMAT, fstrSOLID, hecLagMAT, infoCTChange, conMAT, is_contact_active)
this subroutine reconstructs node-based (stiffness) matrix structure \corresponding to contact state
Definition: fstr_mat_con_contact.f90:53
m_fstr_solve_nlgeom::table_nlsta
subroutine table_nlsta(hecMESH, fstrSOLID, cstep, time, f_t)
This subroutine decide the loading increment considering the amplitude definition.
Definition: fstr_solve_NLGEOM.f90:286
m_fstr_cutback::fstr_cutback_finalize
subroutine fstr_cutback_finalize(fstrSOLID)
Finalizer of cutback variables.
Definition: fstr_Cutback.f90:87
m_fstr_timeinc::fstr_get_timeinc
real(kind=kreal) function fstr_get_timeinc()
Definition: fstr_Ctrl_TimeInc.f90:26
m_fstr_restart
This module provides functions to read in and write out restart files.
Definition: fstr_Restart.f90:8
mcontact
This module provides functions to calculate contact stiff matrix.
Definition: fstr_contact.f90:6
m_fstr_timeinc::fstr_timeinc_printstatus_init
subroutine fstr_timeinc_printstatus_init
Definition: fstr_Ctrl_TimeInc.f90:55
m_fstr::imsg
integer(kind=kint), parameter imsg
Definition: m_fstr.f90:110