MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
Loading...
Searching...
No Matches
mod_dt.t
Go to the documentation of this file.
1module mod_dt
2 implicit none
3 private
4 public :: setdt
5
6contains
7 !>setdt - set dt for all levels between levmin and levmax.
8 !> dtpar>0 --> use fixed dtpar for all level
9 !> dtpar<=0 --> determine CFL limited timestep
10 subroutine setdt()
12 use mod_physics
15 use mod_comm_lib, only: mpistop
16
17 integer :: iigrid, igrid, ncycle, ncycle2, ifile, idim
18 double precision :: dtnew, qdtnew, dtmin_mype, factor, dx^d, dxmin^d
19 double precision :: dtmax, dxmin, cmax_mype
20 double precision :: a2max_mype(ndim), tco_mype, tco_global, tmax_mype, t_peak
21 double precision :: trac_alfa, trac_dmax, trac_tau, t_bott
22 integer, parameter :: niter_print = 2000
23
24 if (dtpar<=zero) then
25 dtmin_mype = bigdouble
26 cmax_mype = zero
27 a2max_mype = zero
28 tco_mype = zero
29 tmax_mype = zero
30 !$OMP PARALLEL DO PRIVATE(igrid,qdtnew,dtnew,dx^D) REDUCTION(min:dtmin_mype) REDUCTION(max:cmax_mype,a2max_mype)
31 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
32 dtnew=bigdouble
33 dx^d=rnode(rpdx^d_,igrid);
34 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
35 block=>ps(igrid)
36 if(local_timestep) then
37 ps(igrid)%dt(ixm^t)=bigdouble
38 endif
39 call getdt_courant(ps(igrid)%w,ixg^ll,ixm^ll,qdtnew,dx^d,ps(igrid)%x,&
40 cmax_mype,a2max_mype)
41 dtnew=min(dtnew,qdtnew)
42
43 call phys_get_dt(ps(igrid)%w,ixg^ll,ixm^ll,qdtnew,dx^d,ps(igrid)%x)
44 dtnew=min(dtnew,qdtnew)
45
46 if (associated(usr_get_dt)) then
47 call usr_get_dt(ps(igrid)%w,ixg^ll,ixm^ll,qdtnew,dx^d,ps(igrid)%x)
48 dtnew = min(dtnew,qdtnew)
49 end if
50 dtmin_mype = min(dtmin_mype,dtnew)
51 end do
52 !$OMP END PARALLEL DO
53 else
54 dtmin_mype=dtpar
55 end if
56
57 if (dtmin_mype<dtmin) then
58 write(unitterm,*)"Error: Time step too small!", dtmin_mype
59 write(unitterm,*)"on processor:", mype, "at time:", global_time," step:", it
60 write(unitterm,*)"Lower limit of time step:", dtmin
61 crash=.true.
62 end if
63
64 if (slowsteps>it-it_init+1) then
65 factor=one-(one-dble(it-it_init+1)/dble(slowsteps))**2
66 dtmin_mype=dtmin_mype*factor
67 end if
68
69 if(final_dt_reduction)then
70 !if (dtmin_mype>time_max-global_time) then
71 ! write(unitterm,*)"WARNING final timestep artificially reduced!"
72 ! write(unitterm,*)"on processor:", mype, "at time:", global_time," step:", it
73 !endif
74 if(time_max-global_time<=dtmin) then
75 !write(unitterm,*)'Forcing to leave timeloop as time is reached!'
76 final_dt_exit=.true.
77 endif
78 dtmin_mype=min(dtmin_mype,time_max-global_time)
79 end if
80
81 if (dtpar<=zero) then
82 call mpi_allreduce(dtmin_mype,dt,1,mpi_double_precision,mpi_min, &
84 else
85 dt=dtmin_mype
86 end if
87
88 if(any(dtsave(1:nfile)<bigdouble).or.any(tsave(isavet(1:nfile),1:nfile)<bigdouble))then
89 dtmax = minval(ceiling(global_time/dtsave(1:nfile))*dtsave(1:nfile))-global_time
90 do ifile=1,nfile
91 dtmax = min(tsave(isavet(ifile),ifile)-global_time,dtmax)
92 end do
93 if(dtmax > smalldouble)then
94 dt=min(dt,dtmax)
95 else
96 ! dtmax=0 means dtsave is divisible by global_time
97 dt=min(dt,minval(dtsave(1:nfile)))
98 end if
99 end if
100
101 if(mype==0) then
102 if(any(dtsave(1:nfile)<dt)) then
103 write(unitterm,*) 'Warning: timesteps: ',dt,' exceeding output intervals ', dtsave(1:nfile)
104 endif
105 endif
106 if(is_sts_initialized()) then
108 qdtnew = 0.5d0 * dt
109 if (set_dt_sts_ncycles(qdtnew)) then
110 dt = 2.d0*qdtnew
111 !a quick way to print the reduction of time only every niter_print iterations
112 !Note that niter_print is a parameter variable hardcoded to the value of 200
113 if(mype==0 .and. mod(it-1, niter_print) .eq. 0) then
114 write(*,*) 'Max number of STS cycles exceeded, reducing dt to',dt
115 endif
116 endif
117 else
118 if(set_dt_sts_ncycles(dt))then
119 if(mype==0 .and. mod(it-1, niter_print) .eq. 0) then
120 write(*,*) 'Max number of STS cycles exceeded, reducing dt to',dt
121 endif
122 endif
123 endif
124 endif
125
126 ! global Lax-Friedrich finite difference flux splitting needs fastest wave-speed
127 ! so does GLM:
128 if(need_global_cmax) call mpi_allreduce(cmax_mype, cmax_global, 1,&
129 mpi_double_precision,mpi_max,icomm,ierrmpi)
130 if(need_global_a2max) call mpi_allreduce(a2max_mype, a2max_global, ndim,&
131 mpi_double_precision,mpi_max,icomm,ierrmpi)
132
133 ! transition region adaptive thermal conduction (Johnston 2019 ApJL, 873, L22)
134 ! transition region adaptive thermal conduction (Johnston 2020 A&A, 635, 168)
135 if(phys_trac) then
136 t_bott=2.d4/unit_temperature
137 call mpi_allreduce(tmax_mype,t_peak,1,mpi_double_precision,&
138 mpi_max,icomm,ierrmpi)
139 ! TODO trac stuff should not be here at all
140 if(phys_trac_type==1) then
141 !> 1D TRAC method
142 trac_dmax=0.1d0
143 trac_tau=1.d0/unit_time
144 trac_alfa=trac_dmax**(dt/trac_tau)
145 tco_global=zero
146 {^ifoned
147 call mpi_allreduce(tco_mype,tco_global,1,mpi_double_precision,&
148 mpi_max,icomm,ierrmpi)
149 }
150 endif
151 if(.not. associated(phys_trac_after_setdt)) call mpistop("phys_trac_after_setdt not set")
152 ! trac_alfa,tco_global are set only for phys_trac_type=1, should not be a problem when not initialized
153 ! side effect of modifying T_bott from mod_trac -> T_bott sent as param
154 call phys_trac_after_setdt(tco_global,trac_alfa,t_peak, t_bott)
155 end if
156
157 contains
158
159 !> compute CFL limited dt (for variable time stepping)
160 subroutine getdt_courant(w,ixI^L,ixO^L,dtnew,dx^D,x,cmax_mype,a2max_mype)
164
165 integer, intent(in) :: ixI^L, ixO^L
166 double precision, intent(in) :: x(ixI^S,1:ndim)
167 double precision, intent(in) :: dx^D
168 double precision, intent(inout) :: w(ixI^S,1:nw), dtnew, cmax_mype, a2max_mype(ndim)
169
170 double precision :: courantmax, dxinv(1:ndim), courantmaxtot, courantmaxtots
171 double precision :: cmax(ixI^S), cmaxtot(ixI^S), wprim(ixI^S,1:nw)
172 double precision :: a2max(ndim), tco_local, Tmax_local
173 integer :: idims
174 integer :: hxO^L
175
176 dtnew=bigdouble
177
178 ! local timestep dt has to be calculated in the
179 ! extended region because of the calculation from the
180 ! div fluxes in mod_finite_volume
181 if(local_timestep) then
182 hxomin^d=ixomin^d-1;
183 hxomax^d=ixomax^d;
184 else
185 hxomin^d=ixomin^d;
186 hxomax^d=ixomax^d;
187 end if
188
189 ! use primitive variables to get sound speed faster
190 wprim=w
191 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
192
193 if(need_global_a2max) then
194 call phys_get_a2max(w,x,ixi^l,ixo^l,a2max)
195 do idims=1,ndim
196 a2max_mype(idims) = max(a2max_mype(idims),a2max(idims))
197 end do
198 end if
199
200 if(phys_trac) then
201 call phys_get_tcutoff(ixi^l,ixo^l,wprim,x,tco_local,tmax_local)
202 {^ifoned tco_mype=max(tco_mype,tco_local) }
203 tmax_mype=max(tmax_mype,tmax_local)
204 end if
205
206 ! these are also calculated in hxO because of local timestep
207 if(nwaux>0) call phys_get_auxiliary(ixi^l,hxo^l,w,x)
208
209 select case (type_courant)
210 case (type_maxsum)
211 if(slab_uniform) then
212 ^d&dxinv(^d)=one/dx^d;
213 do idims=1,ndim
214 call phys_get_cmax(wprim,x,ixi^l,hxo^l,idims,cmax)
215 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
216 if(idims==1) then
217 cmaxtot(hxo^s)=cmax(hxo^s)*dxinv(idims)
218 else
219 cmaxtot(hxo^s)=cmaxtot(hxo^s)+cmax(hxo^s)*dxinv(idims)
220 end if
221 end do
222 else
223 do idims=1,ndim
224 call phys_get_cmax(wprim,x,ixi^l,hxo^l,idims,cmax)
225 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
226 if(idims==1) then
227 cmaxtot(hxo^s)=cmax(hxo^s)/block%ds(hxo^s,idims)
228 else
229 cmaxtot(hxo^s)=cmaxtot(hxo^s)+cmax(hxo^s)/block%ds(hxo^s,idims)
230 end if
231 end do
232 end if
233 ! courantmaxtots='max(summed c/dx)'
234 courantmaxtots=maxval(cmaxtot(ixo^s))
235 if(courantmaxtots>smalldouble) dtnew=min(dtnew,courantpar/courantmaxtots)
236 if(local_timestep) then
237 block%dt(hxo^s) = courantpar/cmaxtot(hxo^s)
238 end if
239
240 case (type_summax)
241 !TODO this should be mod_input_output?
242 if(local_timestep) then
243 call mpistop("Type courant summax incompatible with local_timestep")
244 end if
245 courantmax=zero
246 courantmaxtot=zero
247 if(slab_uniform) then
248 ^d&dxinv(^d)=one/dx^d;
249 do idims=1,ndim
250 call phys_get_cmax(wprim,x,ixi^l,ixo^l,idims,cmax)
251 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
252 courantmax=max(courantmax,maxval(cmax(ixo^s)*dxinv(idims)))
253 courantmaxtot=courantmaxtot+courantmax
254 end do
255 else
256 do idims=1,ndim
257 call phys_get_cmax(wprim,x,ixi^l,ixo^l,idims,cmax)
258 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
259 courantmax=max(courantmax,maxval(cmax(ixo^s)/block%ds(ixo^s,idims)))
260 courantmaxtot=courantmaxtot+courantmax
261 end do
262 end if
263 ! courantmaxtot='summed max(c/dx)'
264 if (courantmaxtot>smalldouble) dtnew=min(dtnew,courantpar/courantmaxtot)
265 case (type_minimum)
266 if(local_timestep) then
267 call mpistop("Type courant not implemented for local_timestep, use maxsum")
268 endif
269 courantmax=zero
270 if(slab_uniform) then
271 ^d&dxinv(^d)=one/dx^d;
272 do idims=1,ndim
273 call phys_get_cmax(wprim,x,ixi^l,ixo^l,idims,cmax)
274 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
275 courantmax=max(courantmax,maxval(cmax(ixo^s)*dxinv(idims)))
276 end do
277 else
278 do idims=1,ndim
279 call phys_get_cmax(wprim,x,ixi^l,ixo^l,idims,cmax)
280 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
281 courantmax=max(courantmax,maxval(cmax(ixo^s)/block%ds(ixo^s,idims)))
282 end do
283 end if
284 ! courantmax='max(c/dx)'
285 if (courantmax>smalldouble) dtnew=min(dtnew,courantpar/courantmax)
286 end select
287 end subroutine getdt_courant
288 end subroutine setdt
289end module mod_dt
subroutine getdt_courant(w, ixil, ixol, dtnew, dxd, x, cmax_mype, a2max_mype)
compute CFL limited dt (for variable time stepping)
Definition mod_dt.t:161
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition mod_dt.t:1
subroutine, public setdt()
setdt - set dt for all levels between levmin and levmax. dtpar>0 --> use fixed dtpar for all level dt...
Definition mod_dt.t:11
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
double precision unit_time
Physical scaling factor for time.
double precision global_time
The global simulation time.
double precision time_max
End time for the simulation.
integer it
Number of time steps taken.
integer it_init
initial iteration count
integer, parameter type_maxsum
integer switchers for type courant
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
double precision dtpar
If dtpar is positive, it sets the timestep dt, otherwise courantpar is used to limit the time step ba...
double precision, dimension(:), allocatable, parameter d
logical local_timestep
each cell has its own timestep or not
double precision dt
global time step
logical need_global_a2max
global value for schmid scheme
double precision courantpar
The Courant (CFL) number used for the simulation.
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
integer slowsteps
If > 1, then in the first slowsteps-1 time steps dt is reduced by a factor .
integer type_courant
How to compute the CFL-limited time step.
integer, parameter unitterm
Unit for standard output.
double precision, dimension(nfile) dtsave
Repeatedly save output of type N when dtsave(N) simulation time has passed.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
logical final_dt_reduction
If true, allow final dt reduction for matching time_max on output.
integer, parameter type_summax
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
logical phys_trac
Use TRAC for MHD or 1D HD.
double precision, dimension(nsavehi, nfile) tsave
Save output of type N on times tsave(:, N)
logical need_global_cmax
need global maximal wave speed
logical crash
Save a snapshot before crash a run met unphysical values.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision, dimension(^nd) a2max_global
global largest a2 for schmid scheme
double precision dtmin
Stop the simulation when the time step becomes smaller than this value.
integer, parameter nfile
Number of output methods.
logical final_dt_exit
Force timeloop exit when final dt < dtmin.
integer, parameter type_minimum
integer, dimension(nfile) isavet
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_get_a2max), pointer phys_get_a2max
Definition mod_physics.t:60
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:57
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:68
procedure(sub_get_tcutoff), pointer phys_get_tcutoff
Definition mod_physics.t:61
procedure(sub_get_auxiliary), pointer phys_get_auxiliary
Definition mod_physics.t:92
procedure(sub_trac_after_setdt), pointer phys_trac_after_setdt
Definition mod_physics.t:62
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:59
Generic supertimestepping method 1) in amrvac.par in sts_list set the following parameters which have...
pure logical function, public is_sts_initialized()
logical function, public set_dt_sts_ncycles(my_dt)
This sets the explicit dt and calculates the number of cycles for each of the terms implemented with ...
integer, parameter, public sourcetype_sts_split
Module with all the methods that users can customize in AMRVAC.
procedure(get_dt), pointer usr_get_dt