FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
m_step.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 module m_step
7  use hecmw
8  implicit none
9 
10  include 'fstr_ctrl_util_f.inc'
11 
12  integer, parameter :: stepstatic = 1
13  integer, parameter :: stepvisco = 2
14  integer, parameter :: stepfixedinc = 1
15  integer, parameter :: stepautoinc = 2
16 
17  ! statistics of newton iteration
18  integer(kind=kint),parameter :: knstmaxit = 1 ! maximum number of newton iteration
19  integer(kind=kint),parameter :: knstsumit = 2 ! total number of newton iteration
20  integer(kind=kint),parameter :: knstciter = 3 ! number of contact iteration
21  integer(kind=kint),parameter :: knstdresn = 4 ! reason of not to converged
22 
24  type step_info
25  integer :: solution
26  integer :: inc_type
27  character( len=80 ) :: control
28  character( len=80 ) :: convcontrol
30  real(kind=kreal) :: converg
31  real(kind=kreal) :: converg_lag
32  real(kind=kreal) :: converg_ddisp
33  real(kind=kreal) :: maxres
34 
35  integer :: num_substep
36  integer :: max_iter
37  integer :: max_contiter
38  integer :: amp_id
39  real(kind=kreal) :: initdt
40  real(kind=kreal) :: elapsetime
41  real(kind=kreal) :: mindt
42  real(kind=kreal) :: maxdt
43  real(kind=kreal) :: starttime
44  integer, pointer :: boundary(:)=>null()
45  integer, pointer :: load(:)=>null()
46  integer, pointer :: contact(:)=>null()
47  integer :: timepoint_id
48  integer :: aincparam_id
49  end type
50 
52  character(HECMW_NAME_LEN) :: name
53  real(kind=kreal) :: ainc_rs
54  real(kind=kreal) :: ainc_rl
55  integer( kind=kint ) :: nrbound_s(10)
56  integer( kind=kint ) :: nrbound_l(10)
57  integer( kind=kint ) :: nrtimes_s
58  integer( kind=kint ) :: nrtimes_l
59  real(kind=kreal) :: ainc_rc
60  integer( kind=kint ) :: cbbound
61  end type
62 
63 contains
64 
66  subroutine init_stepinfo( stepinfo )
67  type( step_info ), intent(out) :: stepinfo
68  stepinfo%solution = stepstatic
69  stepinfo%inc_type = stepfixedinc
70  stepinfo%num_substep = 1
71  stepinfo%max_iter = 50
72  stepinfo%max_contiter = 10
73  stepinfo%amp_id = -1
74  stepinfo%initdt = 1.d0
75  stepinfo%mindt = 1.d-4
76  stepinfo%maxdt = 1.d0
77  stepinfo%elapsetime = 1.d0
78  stepinfo%starttime = 0.d0
79  stepinfo%converg = 1.d-3
80  stepinfo%converg_lag = 1.d-4
81  stepinfo%converg_ddisp = 1.d-8
82  stepinfo%maxres = 1.d+10
83  stepinfo%timepoint_id = 0
84  stepinfo%AincParam_id = 0
85  end subroutine
86 
87  subroutine setup_stepinfo_starttime( stepinfos )
88  type( step_info ), pointer, intent(inout) :: stepinfos(:)
89 
90  integer :: i
91 
92  stepinfos(1)%starttime = 0.d0
93  do i=1,size(stepinfos)-1
94  stepinfos(i+1)%starttime = stepinfos(i)%starttime + stepinfos(i)%elapsetime
95  end do
96  end subroutine
97 
99  logical function isboundaryactive( bnd, stepinfo )
100  integer, intent(in) :: bnd
101  type( step_info ), intent(in) :: stepinfo
102  isboundaryactive = .false.
103  if( .not. associated( stepinfo%Boundary ) ) return
104  if( any( stepinfo%Boundary== bnd ) ) isboundaryactive = .true.
105  end function
106 
108  logical function isloadactive( bnd, stepinfo )
109  integer, intent(in) :: bnd
110  type( step_info ), intent(in) :: stepinfo
111  isloadactive = .false.
112  if( .not. associated( stepinfo%Load ) ) return
113  if( any( stepinfo%Load == bnd ) ) isloadactive = .true.
114  end function
115 
117  logical function iscontactactive( bnd, stepinfo )
118  integer, intent(in) :: bnd
119  type( step_info ), intent(in) :: stepinfo
120  iscontactactive = .false.
121  if( .not. associated( stepinfo%Contact ) ) return
122  if( any( stepinfo%Contact== bnd ) ) iscontactactive = .true.
123  end function
124 
126  subroutine free_stepinfo( step )
127  type(step_info), intent(inout) :: step
128  if( associated( step%Boundary ) ) deallocate( step%Boundary )
129  if( associated( step%Load ) ) deallocate( step%Load )
130  if( associated( step%Contact ) ) deallocate( step%Contact )
131  end subroutine
132 
134  subroutine fstr_print_steps( nfile, steps )
135  integer, intent(in) :: nfile
136  type(step_info), intent(in) :: steps(:)
137  integer :: i, j, nstep, nbc
138  nstep = size(steps)
139 
140  write( nfile, * ) "-----Information of steps:",nstep
141 
142  do i=1,nstep
143  write( nfile, * ) " -----Step:",i
144  write(nfile,*) steps(i)%solution, steps(i)%elapsetime, steps(i)%converg, &
145  steps(i)%num_substep, steps(i)%max_iter
146  if( associated( steps(i)%Boundary ) ) then
147  nbc = size( steps(i)%Boundary )
148  write(nfile,*) " Boundary conditions"
149  write(nfile,*) ( steps(i)%Boundary(j),j=1,nbc )
150  endif
151  if( associated( steps(i)%Load ) ) then
152  nbc = size( steps(i)%Load )
153  write(nfile,*) " External load conditions"
154  write(nfile,*) ( steps(i)%Load(j),j=1,nbc )
155  endif
156  if( associated( steps(i)%Contact ) ) then
157  nbc = size( steps(i)%Contact )
158  write(nfile,*) " Contact conditions"
159  write(nfile,*) ( steps(i)%Contact(j),j=1,nbc )
160  endif
161  enddo
162  end subroutine
163 
165  subroutine init_aincparam( aincparam )
166  type( tparamautoinc ), intent(out) :: aincparam
167 
168  aincparam%name = ''
169  aincparam%ainc_Rs = 0.25d0
170  aincparam%ainc_Rl = 1.25d0
171  aincparam%NRbound_s = 0
172  aincparam%NRbound_s(knstmaxit) = 10
173  aincparam%NRbound_s(knstsumit) = 50
174  aincparam%NRbound_s(knstciter) = 10
175  aincparam%NRbound_l = 0
176  aincparam%NRbound_l(knstmaxit) = 1
177  aincparam%NRbound_l(knstsumit) = 1
178  aincparam%NRbound_l(knstciter) = 1
179  aincparam%NRtimes_s = 1
180  aincparam%NRtimes_l = 2
181  aincparam%ainc_Rc = 0.25d0
182  aincparam%CBbound = 5
183  end subroutine
184 
185 end module
m_step::knstmaxit
integer(kind=kint), parameter knstmaxit
Definition: m_step.f90:18
m_step::iscontactactive
logical function iscontactactive(bnd, stepinfo)
Is contact condition in this step active.
Definition: m_step.f90:118
m_step::tparamautoinc
Definition: m_step.f90:51
m_step::stepautoinc
integer, parameter stepautoinc
Definition: m_step.f90:15
m_step::knstciter
integer(kind=kint), parameter knstciter
Definition: m_step.f90:20
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_step::fstr_print_steps
subroutine fstr_print_steps(nfile, steps)
Print out step control.
Definition: m_step.f90:135
m_step::stepvisco
integer, parameter stepvisco
Definition: m_step.f90:13
m_step::setup_stepinfo_starttime
subroutine setup_stepinfo_starttime(stepinfos)
Definition: m_step.f90:88
m_step::stepstatic
integer, parameter stepstatic
Definition: m_step.f90:12
m_step::isboundaryactive
logical function isboundaryactive(bnd, stepinfo)
Is boundary condition in this step active.
Definition: m_step.f90:100
hecmw
Definition: hecmw.f90:6
m_step::step_info
Step control such as active boundary condition, convergent condition etc.
Definition: m_step.f90:24
m_step::init_aincparam
subroutine init_aincparam(aincparam)
Initializer.
Definition: m_step.f90:166
m_step::init_stepinfo
subroutine init_stepinfo(stepinfo)
Initializer.
Definition: m_step.f90:67
m_step::isloadactive
logical function isloadactive(bnd, stepinfo)
Is external load in this step active.
Definition: m_step.f90:109
m_step
This module manages step information.
Definition: m_step.f90:6
m_step::knstdresn
integer(kind=kint), parameter knstdresn
Definition: m_step.f90:21
m_step::free_stepinfo
subroutine free_stepinfo(step)
Finalizer.
Definition: m_step.f90:127