FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
fstr_Ctrl_TimeInc.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_step
10  use m_timepoint
11 
12  implicit none
13 
14  real(kind=kreal), private :: current_time = 0.d0
15  real(kind=kreal), private :: time_inc = 0.d0
16  real(kind=kreal), private :: time_inc_base = 0.d0
17 
18 contains
19 
20  !! basic functions to get time/timeinc
21  real(kind=kreal) function fstr_get_time()
22  fstr_get_time = current_time
23  end function
24 
25  real(kind=kreal) function fstr_get_timeinc()
26  fstr_get_timeinc = time_inc
27  end function
28 
29  real(kind=kreal) function fstr_get_timeinc_base()
30  fstr_get_timeinc_base = time_inc_base
31  end function
32 
33  !! basic functions to set/update time/timeinc
34  subroutine fstr_set_time( time )
35  real(kind=kreal), intent(in) :: time
36  current_time = time
37  end subroutine
38 
39  subroutine fstr_set_timeinc( dtime )
40  real(kind=kreal), intent(in) :: dtime
41  time_inc = dtime
42  end subroutine
43 
44  subroutine fstr_set_timeinc_base( dtime_base )
45  real(kind=kreal), intent(in) :: dtime_base
46  time_inc_base = dtime_base
47  end subroutine
48 
49  subroutine fstr_proceed_time()
50  current_time = current_time + time_inc
51  end subroutine
52 
53  !! functions to print status
55  write(ista,'(A10,"-+-",A60,"-+-",A40)') repeat("-",10),repeat("-",60),repeat("-",40)
56  write(ista,'(2A5," | ",2A5,2A7,3A12," | ",A40)') ' ', ' ', ' ', ' # of', 'MAX #', 'TOT #', '','', '', ''
57  write(ista,'(2A5," | ",2A5,2A7,3A12," | ",A7,A33)') &
58  & 'STEP', 'SUB', 'STAT', ' CONT', 'NEWTON', 'NEWTON', 'START', 'TIME','END', 'MESSAGE',''
59  write(ista,'(2A5," | ",2A5,2A7,3A12," | ",A40)') &
60  & '', 'STEP', ' ', 'ITER', 'ITER', 'ITER', 'TIME', 'INC', 'TIME',''
61  write(ista,'(A10,"-+-",A60,"-+-",A40)') repeat("-",10),repeat("-",60),repeat("-",40)
62  end subroutine
63 
64  subroutine fstr_timeinc_printstatus( stepinfo, fstrPARAM, totstep, substep, NRstatI, NRstatR, AutoINC_stat, Cutback_stat )
65  type(step_info), intent(in) :: stepinfo
66  type(fstr_param), intent(in) :: fstrPARAM
67  integer(kind=kint), intent(in) :: totstep
68  integer(kind=kint), intent(in) :: substep
69  integer(kind=kint), intent(in) :: NRstatI(:) !previous Newton Raphson Iteration Log
70  real(kind=kreal), intent(in) :: nrstatr(:) !previous Newton Raphson Iteration Log
71  integer(kind=kint), intent(inout) :: AutoINC_stat
72  integer(kind=kint), intent(in) :: Cutback_stat
73 
74  character(len=5) :: cstep, csstep, cstate, ccont
75  character(len=7) :: cmaxn, ctotn
76  character(len=12) :: ctime, cdtime, etime
77  character(len=40) :: message
78  type(tparamautoinc) :: pAinc
79 
80  painc = fstrparam%ainc(stepinfo%AincParam_id)
81 
82  write(cstep,'(I5)') totstep
83  write(csstep,'(I5)') substep
84  write(ccont,'(I5)') nrstati(knstciter)
85  write(cmaxn,'(I7)') nrstati(knstmaxit)
86  write(ctotn,'(I7)') nrstati(knstsumit)
87  write(ctime,'(1pE12.4)') current_time
88  write(cdtime,'(1pE12.4)') time_inc
89 
90  if( cutback_stat > 0 ) then
91  write(etime,'(1pE12.4)') current_time
92  write(cstate,'(I4,A)') cutback_stat,'F'
93  if( nrstati(knstdresn) == 1 ) write(message,'(A)') 'Failed to converge due to MAXITER.'
94  if( nrstati(knstdresn) == 2 ) write(message,'(A)') 'Failed to converge due to MAXRES.'
95  if( nrstati(knstdresn) == 3 ) write(message,'(A)') 'Failed to converge due to MAXCONTITER.'
96  if( nrstati(knstdresn) == 4 ) write(message,'(A)') 'Failed to converge due to MatSolveError.'
97  if( cutback_stat == painc%CBbound ) write(message,'(A)') '# of successive cutback reached max.'
98  else
99  write(etime,'(1pE12.4)') current_time+time_inc
100  write(cstate,'(A5)') 'S'
101  write(message,'(A)') ''
102  endif
103 
104  write(ista,'(2A5," | ",2A5,2A7,3A12," | ",A40)') cstep, csstep, cstate, ccont, &
105  & cmaxn, ctotn, ctime, cdtime, etime, message
106 
107  end subroutine
108 
109  subroutine fstr_timeinc_printstatus_final(success_flag)
110  logical, intent(in) :: success_flag
111  write(ista,'(A10,"-+-",A60,"-+-",A40)') repeat("-",10),repeat("-",60),repeat("-",40)
112  if(success_flag) then
113  write(ista,'(A)') 'FSTR_SOLVE_NLGEOM HAS COMPLETED SUCCESSFULLY'
114  else
115  write(ista,'(A)') 'FSTR_SOLVE_NLGEOM HAS NOT COMPLETED SUCCESSFULLY'
116  end if
117  end subroutine
118 
119  !! true if step finished
120  logical function fstr_timeinc_isstepfinished( stepinfo )
121  type(step_info), intent(in) :: stepinfo
122 
123  real(kind=kreal) :: endtime, remain_time
124 
125  endtime = stepinfo%starttime + stepinfo%elapsetime
126  remain_time = (endtime-current_time)/stepinfo%elapsetime
127  fstr_timeinc_isstepfinished = (remain_time < 1.d-12)
128  end function
129 
130  !! true if current time is time point
131  logical function fstr_timeinc_istimepoint( stepinfo, fstrPARAM )
132  type(step_info), intent(in) :: stepinfo
133  type(fstr_param), intent(in) :: fstrparam
134 
135  fstr_timeinc_istimepoint = .false.
136  if( stepinfo%inc_type == stepfixedinc ) return
137  if( stepinfo%timepoint_id == 0 ) return
138 
139  fstr_timeinc_istimepoint = is_at_timepoints(current_time,stepinfo%starttime, &
140  & fstrparam%timepoints(stepinfo%timepoint_id))
141  end function
142 
143  !! set time_inc from stepinfo and Newton Raphson iteration status
144  subroutine fstr_timeinc_settimeincrement( stepinfo, fstrPARAM, substep, NRstatI, NRstatR, AutoINC_stat, Cutback_stat )
145  type(step_info), intent(in) :: stepinfo
146  type(fstr_param), intent(in) :: fstrPARAM
147  integer(kind=kint), intent(in) :: substep
148  integer(kind=kint), intent(in) :: NRstatI(:) !previous Newton Raphson Iteration Log
149  real(kind=kreal), intent(in) :: nrstatr(:) !previous Newton Raphson Iteration Log
150  integer(kind=kint), intent(inout) :: AutoINC_stat
151  integer(kind=kint), intent(in) :: Cutback_stat
152 
153  real(kind=kreal) :: timeinc0
154  real(kind=kreal) :: endtime, remain
155  logical :: to_be_decreased, to_be_increased
156  type(tparamautoinc) :: pAinc
157 
158  painc = fstrparam%ainc(stepinfo%AincParam_id)
159 
160  if( stepinfo%inc_type == stepfixedinc ) then
161  timeinc0 = stepinfo%initdt
162  else ! INCTYPE==AUTO
163  if( cutback_stat > 0 ) then
164  timeinc0 = painc%ainc_Rc*time_inc
165  if(myrank == 0) write(*,'(2(A,E10.3))') 'time increment is decreased from ', time_inc_base, ' to ', timeinc0
166  autoinc_stat = -1
167  else
168  if( substep == 1 ) then
169  timeinc0 = stepinfo%initdt
170  autoinc_stat = 0
171  else
172  timeinc0 = time_inc_base
173 
174  !decrease condition
175  to_be_decreased = .false.
176  if( nrstati(knstmaxit) > painc%NRbound_s(knstmaxit) ) to_be_decreased = .true.
177  if( nrstati(knstsumit) > painc%NRbound_s(knstsumit) ) to_be_decreased = .true.
178  if( nrstati(knstciter) > painc%NRbound_s(knstciter) ) to_be_decreased = .true.
179 
180  !increase condition
181  to_be_increased = .true.
182  if( nrstati(knstmaxit) > painc%NRbound_l(knstmaxit) ) to_be_increased = .false.
183  if( nrstati(knstsumit) > painc%NRbound_l(knstsumit) ) to_be_increased = .false.
184  if( nrstati(knstciter) > painc%NRbound_l(knstciter) ) to_be_increased = .false.
185 
186  ! count # of times that increase/decrease condition has been satisfied
187  ! AutoINC_stat < 0 ... decrease condition has been satisfied -AutoINC_stat times
188  ! AutoINC_stat > 0 ... increase condition has been satisfied AutoINC_stat times
189  if(to_be_decreased) then
190  autoinc_stat = min(autoinc_stat,0) - 1
191  else if(to_be_increased) then
192  autoinc_stat = max(autoinc_stat,0) + 1
193  else
194  autoinc_stat = 0
195  end if
196 
197  if( autoinc_stat <= -painc%NRtimes_s ) then
198  timeinc0 = painc%ainc_Rs*timeinc0
199  if(myrank == 0) write(*,'(2(A,E10.3))') 'time increment is decreased from ', time_inc_base, ' to ', timeinc0
200  else if( autoinc_stat >= painc%NRtimes_l ) then
201  timeinc0 = dmin1(painc%ainc_Rl*timeinc0,stepinfo%maxdt)
202  if(myrank == 0) write(*,'(2(A,E10.3))') 'time increment is increased from ', time_inc_base, ' to ', timeinc0
203  end if
204  end if
205  end if
206  end if
207 
208  if( timeinc0 < stepinfo%mindt ) then
209  if(myrank == 0) then
210  write(*,'(2(A,E10.3))') 'Error: current time increment ',timeinc0,' is smaller than lower bound',stepinfo%mindt
211  call fstr_timeinc_printstatus_final(.false.)
212  endif
213  stop
214  end if
215 
216  endtime = stepinfo%starttime + stepinfo%elapsetime
217  time_inc_base = dmin1(timeinc0,endtime-current_time)
218  if( stepinfo%timepoint_id > 0 ) then
219  remain = get_remain_to_next_timepoints(current_time,stepinfo%starttime,fstrparam%timepoints(stepinfo%timepoint_id))
220  time_inc = dmin1(time_inc_base,remain)
221  else
222  time_inc = time_inc_base
223  end if
224 
225  end subroutine
226 end module m_fstr_timeinc
m_step::knstmaxit
integer(kind=kint), parameter knstmaxit
Definition: m_step.f90:18
m_fstr_timeinc
This module provides functions to deal with time and increment of stress analysis.
Definition: fstr_Ctrl_TimeInc.f90:7
m_fstr_timeinc::fstr_set_timeinc_base
subroutine fstr_set_timeinc_base(dtime_base)
Definition: fstr_Ctrl_TimeInc.f90:45
m_step::tparamautoinc
Definition: m_step.f90:51
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_step::knstciter
integer(kind=kint), parameter knstciter
Definition: m_step.f90:20
m_fstr::myrank
integer(kind=kint) myrank
PARALLEL EXECUTION.
Definition: m_fstr.f90:96
m_fstr_timeinc::fstr_timeinc_isstepfinished
logical function fstr_timeinc_isstepfinished(stepinfo)
Definition: fstr_Ctrl_TimeInc.f90:121
m_step::stepfixedinc
integer, parameter stepfixedinc
Definition: m_step.f90:14
m_step::knstsumit
integer(kind=kint), parameter knstsumit
Definition: m_step.f90:19
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
m_fstr::fstr_param
FSTR INNER CONTROL PARAMETERS (fstrPARAM)
Definition: m_fstr.f90:154
m_timepoint::is_at_timepoints
logical function is_at_timepoints(totaltime, starttime, tp)
Definition: m_timepoint.f90:45
m_timepoint
This module manages timepoint information.
Definition: m_timepoint.f90:6
m_fstr
This module defines common data and basic structures for analysis.
Definition: m_fstr.f90:15
m_fstr_timeinc::fstr_get_timeinc_base
real(kind=kreal) function fstr_get_timeinc_base()
Definition: fstr_Ctrl_TimeInc.f90:30
m_fstr_timeinc::fstr_timeinc_istimepoint
logical function fstr_timeinc_istimepoint(stepinfo, fstrPARAM)
Definition: fstr_Ctrl_TimeInc.f90:132
m_fstr_timeinc::fstr_set_time
subroutine fstr_set_time(time)
Definition: fstr_Ctrl_TimeInc.f90:35
m_step::step_info
Step control such as active boundary condition, convergent condition etc.
Definition: m_step.f90:24
m_fstr::ista
integer(kind=kint), parameter ista
Definition: m_fstr.f90:108
m_fstr_timeinc::fstr_proceed_time
subroutine fstr_proceed_time()
Definition: fstr_Ctrl_TimeInc.f90:50
m_timepoint::get_remain_to_next_timepoints
real(kind=kreal) function get_remain_to_next_timepoints(totaltime, starttime, tp)
Definition: m_timepoint.f90:65
m_fstr_timeinc::fstr_timeinc_printstatus_final
subroutine fstr_timeinc_printstatus_final(success_flag)
Definition: fstr_Ctrl_TimeInc.f90:110
m_step
This module manages step information.
Definition: m_step.f90:6
m_fstr_timeinc::fstr_get_timeinc
real(kind=kreal) function fstr_get_timeinc()
Definition: fstr_Ctrl_TimeInc.f90:26
m_fstr_timeinc::fstr_set_timeinc
subroutine fstr_set_timeinc(dtime)
Definition: fstr_Ctrl_TimeInc.f90:40
m_step::knstdresn
integer(kind=kint), parameter knstdresn
Definition: m_step.f90:21
m_fstr_timeinc::fstr_timeinc_printstatus_init
subroutine fstr_timeinc_printstatus_init
Definition: fstr_Ctrl_TimeInc.f90:55