MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
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, cs2_mype
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 cs2_mype = zero
31 !$OMP PARALLEL DO PRIVATE(igrid,qdtnew,dtnew,dx^D) REDUCTION(min:dtmin_mype) REDUCTION(max:cmax_mype,a2max_mype,cs2_mype)
32 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
33 dtnew=bigdouble
34 dx^d=rnode(rpdx^d_,igrid);
35 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
36 block=>ps(igrid)
37 if(local_timestep) then
38 ps(igrid)%dt(ixm^t)=bigdouble
39 endif
40 call getdt_courant(ps(igrid)%w,ixg^ll,ixm^ll,qdtnew,dx^d,ps(igrid)%x,&
41 cmax_mype,a2max_mype,cs2_mype)
42 dtnew=min(dtnew,qdtnew)
43
44 call phys_get_dt(ps(igrid)%w,ixg^ll,ixm^ll,qdtnew,dx^d,ps(igrid)%x)
45 dtnew=min(dtnew,qdtnew)
46
47 if (associated(usr_get_dt)) then
48 call usr_get_dt(ps(igrid)%w,ixg^ll,ixm^ll,qdtnew,dx^d,ps(igrid)%x)
49 dtnew = min(dtnew,qdtnew)
50 end if
51 dtmin_mype = min(dtmin_mype,dtnew)
52 end do
53 !$OMP END PARALLEL DO
54 else
55 dtmin_mype=dtpar
56 end if
57
58 if (dtmin_mype<dtmin) then
59 write(unitterm,*)"Error: Time step too small!", dtmin_mype
60 write(unitterm,*)"on processor:", mype, "at time:", global_time," step:", it
61 write(unitterm,*)"Lower limit of time step:", dtmin
62 crash=.true.
63 end if
64
65 if (slowsteps>it-it_init+1) then
66 factor=one-(one-dble(it-it_init+1)/dble(slowsteps))**2
67 dtmin_mype=dtmin_mype*factor
68 end if
69
70 if(final_dt_reduction)then
71 !if (dtmin_mype>time_max-global_time) then
72 ! write(unitterm,*)"WARNING final timestep artificially reduced!"
73 ! write(unitterm,*)"on processor:", mype, "at time:", global_time," step:", it
74 !endif
75 if(time_max-global_time<=dtmin) then
76 !write(unitterm,*)'Forcing to leave timeloop as time is reached!'
77 final_dt_exit=.true.
78 endif
79 dtmin_mype=min(dtmin_mype,time_max-global_time)
80 end if
81
82 if (dtpar<=zero) then
83 call mpi_allreduce(dtmin_mype,dt,1,mpi_double_precision,mpi_min, &
85 else
86 dt=dtmin_mype
87 end if
88
89 if(any(dtsave(1:nfile)<bigdouble).or.any(tsave(isavet(1:nfile),1:nfile)<bigdouble))then
90 dtmax = minval(ceiling(global_time/dtsave(1:nfile))*dtsave(1:nfile))-global_time
91 do ifile=1,nfile
92 dtmax = min(tsave(isavet(ifile),ifile)-global_time,dtmax)
93 end do
94 if(dtmax > smalldouble)then
95 dt=min(dt,dtmax)
96 else
97 ! dtmax=0 means dtsave is divisible by global_time
98 dt=min(dt,minval(dtsave(1:nfile)))
99 end if
100 end if
101
102 if(mype==0) then
103 if(any(dtsave(1:nfile)<dt)) then
104 write(unitterm,*) 'Warning: timesteps: ',dt,' exceeding output intervals ', dtsave(1:nfile)
105 endif
106 endif
107 if(is_sts_initialized()) then
109 qdtnew = 0.5d0 * dt
110 if (set_dt_sts_ncycles(qdtnew)) then
111 dt = 2.d0*qdtnew
112 !a quick way to print the reduction of time only every niter_print iterations
113 !Note that niter_print is a parameter variable hardcoded to the value of 200
114 if(mype==0 .and. mod(it-1, niter_print) .eq. 0) then
115 write(*,*) 'Max number of STS cycles exceeded, reducing dt to',dt
116 endif
117 endif
118 else
119 if(set_dt_sts_ncycles(dt))then
120 if(mype==0 .and. mod(it-1, niter_print) .eq. 0) then
121 write(*,*) 'Max number of STS cycles exceeded, reducing dt to',dt
122 endif
123 endif
124 endif
125 endif
126
127 ! global Lax-Friedrich finite difference flux splitting needs fastest wave-speed
128 ! so does GLM:
129 if(need_global_cmax) call mpi_allreduce(cmax_mype, cmax_global, 1,&
130 mpi_double_precision,mpi_max,icomm,ierrmpi)
131 if(need_global_a2max) call mpi_allreduce(a2max_mype, a2max_global, ndim,&
132 mpi_double_precision,mpi_max,icomm,ierrmpi)
133 if(need_global_cs2) call mpi_allreduce(cs2_mype, cs2_global, 1,&
134 mpi_double_precision,mpi_max,icomm,ierrmpi)
135
136 ! transition region adaptive thermal conduction (Johnston 2019 ApJL, 873, L22)
137 ! transition region adaptive thermal conduction (Johnston 2020 A&A, 635, 168)
138 if(phys_trac) then
139 t_bott=2.d4/unit_temperature
140 call mpi_allreduce(tmax_mype,t_peak,1,mpi_double_precision,&
141 mpi_max,icomm,ierrmpi)
142 if(phys_trac_type==1) then
143 !> 1D TRAC method
144 trac_dmax=0.1d0
145 trac_tau=1.d0/unit_time
146 trac_alfa=trac_dmax**(dt/trac_tau)
147 tco_global=zero
148 {^ifoned
149 call mpi_allreduce(tco_mype,tco_global,1,mpi_double_precision,&
150 mpi_max,icomm,ierrmpi)
151 }
152 endif
153 if(.not. associated(phys_trac_after_setdt)) call mpistop("phys_trac_after_setdt not set")
154 ! trac_alfa,tco_global are set only for phys_trac_type=1, should not be a problem when not initialized
155 ! side effect of modifying T_bott from mod_trac -> T_bott sent as param
156 call phys_trac_after_setdt(tco_global,trac_alfa,t_peak, t_bott)
157 end if
158
159 contains
160
161 !> compute CFL limited dt (for variable time stepping)
162 subroutine getdt_courant(w,ixI^L,ixO^L,dtnew,dx^D,x,cmax_mype,a2max_mype,cs2_mype)
166
167 integer, intent(in) :: ixI^L, ixO^L
168 double precision, intent(in) :: x(ixI^S,1:ndim)
169 double precision, intent(in) :: dx^D
170 double precision, intent(inout) :: w(ixI^S,1:nw), dtnew, cmax_mype, a2max_mype(ndim), cs2_mype
171
172 double precision :: courantmax, dxinv(1:ndim), courantmaxtot, courantmaxtots
173 double precision :: cmax(ixI^S), cmaxtot(ixI^S), wprim(ixI^S,1:nw)
174 double precision :: a2max(ndim), tco_local, Tmax_local
175 integer :: idims
176 integer :: hxO^L
177
178 dtnew=bigdouble
179
180 ! local timestep dt has to be calculated in the
181 ! extended region because of the calculation from the
182 ! div fluxes in mod_finite_volume
183 if(local_timestep) then
184 hxomin^d=ixomin^d-1;
185 hxomax^d=ixomax^d;
186 else
187 hxomin^d=ixomin^d;
188 hxomax^d=ixomax^d;
189 end if
190
191 ! use primitive variables to get sound speed faster
192 wprim=w
193 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
194
195 if(need_global_a2max) then
196 call phys_get_a2max(w,x,ixi^l,ixo^l,a2max)
197 do idims=1,ndim
198 a2max_mype(idims) = max(a2max_mype(idims),a2max(idims))
199 end do
200 end if
201
202 if(need_global_cs2) then
203 call phys_get_cs2(wprim,x,ixi^l,ixo^l,cmax)
204 cs2_mype=max(cs2_mype,maxval(cmax(ixo^s)))
205 end if
206
207 if(phys_trac) then
208 call phys_get_tcutoff(ixi^l,ixo^l,wprim,x,tco_local,tmax_local)
209 {^ifoned tco_mype=max(tco_mype,tco_local) }
210 tmax_mype=max(tmax_mype,tmax_local)
211 end if
212
213 ! these are also calculated in hxO because of local timestep
214 if(nwaux>0) call phys_get_auxiliary(ixi^l,hxo^l,w,x)
215
216 select case (type_courant)
217 case (type_maxsum)
218 if(slab_uniform) then
219 ^d&dxinv(^d)=one/dx^d;
220 do idims=1,ndim
221 call phys_get_cmax(wprim,x,ixi^l,hxo^l,idims,cmax)
222 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
223 if(idims==1) then
224 cmaxtot(hxo^s)=cmax(hxo^s)*dxinv(idims)
225 else
226 cmaxtot(hxo^s)=cmaxtot(hxo^s)+cmax(hxo^s)*dxinv(idims)
227 end if
228 end do
229 else
230 do idims=1,ndim
231 call phys_get_cmax(wprim,x,ixi^l,hxo^l,idims,cmax)
232 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
233 if(idims==1) then
234 cmaxtot(hxo^s)=cmax(hxo^s)/block%ds(hxo^s,idims)
235 else
236 cmaxtot(hxo^s)=cmaxtot(hxo^s)+cmax(hxo^s)/block%ds(hxo^s,idims)
237 end if
238 end do
239 end if
240 ! courantmaxtots='max(summed c/dx)'
241 courantmaxtots=maxval(cmaxtot(ixo^s))
242 if(courantmaxtots>smalldouble) dtnew=min(dtnew,courantpar/courantmaxtots)
243 if(local_timestep) then
244 block%dt(hxo^s) = courantpar/cmaxtot(hxo^s)
245 end if
246
247 case (type_summax)
248 courantmax=zero
249 courantmaxtot=zero
250 if(slab_uniform) then
251 ^d&dxinv(^d)=one/dx^d;
252 do idims=1,ndim
253 call phys_get_cmax(wprim,x,ixi^l,ixo^l,idims,cmax)
254 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
255 courantmax=max(courantmax,maxval(cmax(ixo^s)*dxinv(idims)))
256 courantmaxtot=courantmaxtot+courantmax
257 end do
258 else
259 do idims=1,ndim
260 call phys_get_cmax(wprim,x,ixi^l,ixo^l,idims,cmax)
261 if(need_global_cmax) cmax_mype = max(cmax_mype,maxval(cmax(ixo^s)))
262 courantmax=max(courantmax,maxval(cmax(ixo^s)/block%ds(ixo^s,idims)))
263 courantmaxtot=courantmaxtot+courantmax
264 end do
265 end if
266 ! courantmaxtot='summed max(c/dx)'
267 if (courantmaxtot>smalldouble) dtnew=min(dtnew,courantpar/courantmaxtot)
268 case (type_minimum)
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, cs2_mype)
compute CFL limited dt (for variable time stepping)
Definition mod_dt.t:163
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.
logical need_global_cs2
need global squared sound speed
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 cs2_global
global sound speed squared for hyperbolic TC
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:58
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:54
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:66
procedure(sub_get_tcutoff), pointer phys_get_tcutoff
Definition mod_physics.t:59
procedure(sub_get_auxiliary), pointer phys_get_auxiliary
Definition mod_physics.t:90
procedure(sub_trac_after_setdt), pointer phys_trac_after_setdt
Definition mod_physics.t:60
procedure(sub_get_cs2), pointer phys_get_cs2
Definition mod_physics.t:57
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:56
Generic supertimestepping method which can be used for multiple source terms in the governing equatio...
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