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