MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_finite_difference.t
Go to the documentation of this file.
1 !> Module with finite difference methods for fluxes
3 
4  implicit none
5  private
6 
7  public :: fd
8  public :: centdiff
9 
10 contains
11 
12  subroutine fd(qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,snew,fC,fE,dxs,x)
13  use mod_physics
14  use mod_source, only: addsource2
17  use mod_usr_methods
18 
19  double precision, intent(in) :: qdt, dtfactor, qtc, qt, dxs(ndim)
20  integer, intent(in) :: ixi^l, ixo^l, idims^lim
21  double precision, dimension(ixI^S,1:ndim), intent(in) :: x
22 
23  type(state) :: sct, snew
24  double precision, dimension(ixI^S,1:nwflux,1:ndim), intent(out) :: fc
25  double precision, dimension(ixI^S,sdim:3) :: fe
26 
27  double precision, dimension(ixI^S,1:nwflux) :: fct
28  double precision, dimension(ixI^S,1:nw) :: fm, fp, fmr, fpl, wprim
29  ! left and right constructed status in conservative form
30  double precision, dimension(ixI^S,1:nw) :: wlc, wrc
31  ! left and right constructed status in primitive form, needed for better performance
32  double precision, dimension(ixI^S,1:nw) :: wlp, wrp
33  double precision, dimension(ixI^S) :: cmaxc, cminc
34  double precision, dimension(ixI^S) :: hspeed
35  double precision, dimension(ixO^S) :: inv_volume
36  double precision, dimension(1:ndim) :: dxinv
37  logical :: transport, active
38  integer :: idims, iw, ixc^l, ix^l, hxo^l, kxc^l, kxr^l
39  type(ct_velocity) :: vcts
40 
41  associate(wct=>sct%w,wnew=>snew%w)
42 
43  fc=0.d0
44  wprim=wct
45  call phys_to_primitive(ixi^l,ixi^l,wprim,x)
46 
47  b0i = 0 ! fd uses centered values in phys_get_flux
48  do idims= idims^lim
49 
50  !b0i=idims
51 
52  ! Get fluxes for the whole grid (mesh+nghostcells)
53  {^d& ixmin^d = ixomin^d - nghostcells * kr(idims,^d)\}
54  {^d& ixmax^d = ixomax^d + nghostcells * kr(idims,^d)\}
55 
56  hxo^l=ixo^l-kr(idims,^d);
57 
58  if(stagger_grid) then
59  ! ct needs all transverse cells
60  ixcmax^d=ixomax^d+nghostcells-nghostcells*kr(idims,^d); ixcmin^d=hxomin^d-nghostcells+nghostcells*kr(idims,^d);
61  ixmax^d=ixmax^d+nghostcells-nghostcells*kr(idims,^d); ixmin^d=ixmin^d-nghostcells+nghostcells*kr(idims,^d);
62  kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
63  kxr^l=kxc^l+kr(idims,^d);
64  ! wRp and wLp are defined at the same locations, and will correspond to
65  ! the left and right reconstructed values at a cell face. Their indexing
66  ! is similar to cell-centered values, but in direction idims they are
67  ! shifted half a cell towards the 'lower' direction.
68  wrp(kxc^s,1:nw)=wprim(kxr^s,1:nw)
69  wlp(kxc^s,1:nw)=wprim(kxc^s,1:nw)
70  else
71  ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
72  ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
73  end if
74 
75  call phys_get_flux(wct,wprim,x,ixg^ll,ix^l,idims,fct)
76 
77  do iw=iwstart,nwflux
78  ! Lax-Friedrich splitting:
79  fp(ix^s,iw) = half * (fct(ix^s,iw) + tvdlfeps * cmax_global * wct(ix^s,iw))
80  fm(ix^s,iw) = half * (fct(ix^s,iw) - tvdlfeps * cmax_global * wct(ix^s,iw))
81  end do ! iw loop
82 
83  ! now do the reconstruction of fp and fm:
84  call reconstructl(ixi^l,ixc^l,idims,fp,fpl)
85  call reconstructr(ixi^l,ixc^l,idims,fm,fmr)
86 
87  fc(ixc^s,iwstart:nwflux,idims) = fpl(ixc^s,iwstart:nwflux) + fmr(ixc^s,iwstart:nwflux)
88 
89  if(stagger_grid) then
90  ! apply limited reconstruction for left and right status at cell interfaces
91  call reconstruct_lr(ixi^l,ixc^l,ixc^l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
92  if(h_correction) then
93  call phys_get_h_speed(wprim,x,ixi^l,ixo^l,idims,hspeed)
94  end if
95  call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^l,ixc^l,idims,hspeed,cmaxc,cminc)
96  call phys_get_ct_velocity(vcts,wlp,wrp,ixi^l,ixc^l,idims,cmaxc,cminc)
97  end if
98 
99  end do !idims loop
100  !b0i=0
101 
102  if(stagger_grid) call phys_update_faces(ixi^l,ixo^l,qt,qdt,wprim,fc,fe,sct,snew,vcts)
103 
104  if(slab_uniform) then
105  dxinv=-qdt/dxs
106  do idims= idims^lim
107  hxo^l=ixo^l-kr(idims,^d);
108  if(local_timestep) then
109  do iw=iwstart,nwflux
110  fc(ixi^s,iw,idims) = -block%dt(ixi^s)*dtfactor/dxs(idims) * fc(ixi^s,iw,idims)
111  wnew(ixo^s,iw)=wnew(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
112  end do ! iw loop
113  else
114  do iw=iwstart,nwflux
115  fc(ixi^s,iw,idims) = dxinv(idims) * fc(ixi^s,iw,idims)
116  wnew(ixo^s,iw)=wnew(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
117  end do ! iw loop
118  endif
119  end do ! Next idims
120  else
121  inv_volume=1.d0/block%dvolume(ixo^s)
122  do idims= idims^lim
123  hxo^l=ixo^l-kr(idims,^d);
124  if(local_timestep) then
125  do iw=iwstart,nwflux
126  fc(ixi^s,iw,idims)=-block%dt(ixi^s)*dtfactor*fc(ixi^s,iw,idims)*block%surfaceC(ixi^s,idims)
127  wnew(ixo^s,iw)=wnew(ixo^s,iw)+ &
128  (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
129  end do ! iw loop
130  else
131  do iw=iwstart,nwflux
132  fc(ixi^s,iw,idims)=-qdt*fc(ixi^s,iw,idims)*block%surfaceC(ixi^s,idims)
133  wnew(ixo^s,iw)=wnew(ixo^s,iw)+ &
134  (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
135  end do ! iw loop
136  endif ! local_timestep
137  end do ! Next idims
138  end if
139 
140  if (.not.slab.and.idimsmin==1) call phys_add_source_geom(qdt,dtfactor,ixi^l,ixo^l,wct,wnew,x)
141 
142  if(stagger_grid) call phys_face_to_center(ixo^l,snew)
143 
144  ! check and optionally correct unphysical values
145  if(fix_small_values) then
146  call phys_handle_small_values(.false.,wnew,x,ixi^l,ixo^l,'fd')
147  endif
148 
149  call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
150  dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim), &
151  ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
152 
153  end associate
154 
155  end subroutine fd
156 
157  subroutine reconstructl(ixI^L,iL^L,idims,w,wLC)
159  use mod_mp5
160  use mod_limiter
161  use mod_comm_lib, only: mpistop
162 
163  integer, intent(in) :: ixI^L, iL^L, idims
164  double precision, intent(in) :: w(ixI^S,1:nw)
165 
166  double precision, intent(out) :: wLC(ixI^S,1:nw)
167 
168  double precision :: ldw(ixI^S), dwC(ixI^S)
169  integer :: jxR^L, ixC^L, jxC^L, kxC^L, iw
170  double precision :: a2max
171 
172  select case (type_limiter(block%level))
173  case (limiter_mp5)
174  call mp5limiterl(ixi^l,il^l,idims,w,wlc)
175  case (limiter_weno5)
176  call weno5limiterl(ixi^l,il^l,idims,w,wlc,1)
177  case (limiter_weno5nm)
178  call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,1)
179  case (limiter_wenoz5)
180  call weno5limiterl(ixi^l,il^l,idims,w,wlc,2)
181  case (limiter_wenoz5nm)
182  call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,2)
183  case (limiter_wenozp5)
184  call weno5limiterl(ixi^l,il^l,idims,w,wlc,3)
185  case (limiter_wenozp5nm)
186  call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,3)
187  case default
188 
189  kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
190 
191  wlc(kxc^s,iwstart:nwflux) = w(kxc^s,iwstart:nwflux)
192 
193  jxr^l=il^l+kr(idims,^d);
194 
195  ixcmax^d=jxrmax^d; ixcmin^d=ilmin^d-kr(idims,^d);
196  jxc^l=ixc^l+kr(idims,^d);
197 
198  do iw=iwstart,nwflux
199  dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
200 
201  if(need_global_a2max) then
202  a2max=a2max_global(idims)
203  else
204  select case(idims)
205  case(1)
206  a2max=schmid_rad1
207  {^iftwod
208  case(2)
209  a2max=schmid_rad2}
210  {^ifthreed
211  case(2)
212  a2max=schmid_rad2
213  case(3)
214  a2max=schmid_rad3}
215  case default
216  call mpistop("idims is wrong in mod_limiter")
217  end select
218  end if
219 
220  call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),ldw=ldw,a2max=a2max)
221 
222  wlc(il^s,iw)=wlc(il^s,iw)+half*ldw(il^s)
223  end do
224  end select
225 
226  end subroutine reconstructl
227 
228  subroutine reconstructr(ixI^L,iL^L,idims,w,wRC)
230  use mod_mp5
231  use mod_limiter
232  use mod_comm_lib, only: mpistop
233 
234  integer, intent(in) :: ixI^L, iL^L, idims
235  double precision, intent(in) :: w(ixI^S,1:nw)
236 
237  double precision, intent(out) :: wRC(ixI^S,1:nw)
238 
239  double precision :: rdw(ixI^S), dwC(ixI^S)
240  integer :: jxR^L, ixC^L, jxC^L, kxC^L, kxR^L, iw
241  double precision :: a2max
242 
243  select case (type_limiter(block%level))
244  case (limiter_mp5)
245  call mp5limiterr(ixi^l,il^l,idims,w,wrc)
246  case (limiter_weno5)
247  call weno5limiterr(ixi^l,il^l,idims,w,wrc,1)
248  case (limiter_weno5nm)
249  call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,1)
250  case (limiter_wenoz5)
251  call weno5limiterr(ixi^l,il^l,idims,w,wrc,2)
252  case (limiter_wenoz5nm)
253  call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,2)
254  case (limiter_wenozp5)
255  call weno5limiterr(ixi^l,il^l,idims,w,wrc,3)
256  case (limiter_wenozp5nm)
257  call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,3)
258  case default
259 
260  kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
261  kxr^l=kxc^l+kr(idims,^d);
262 
263  wrc(kxc^s,iwstart:nwflux)=w(kxr^s,iwstart:nwflux)
264 
265  jxr^l=il^l+kr(idims,^d);
266  ixcmax^d=jxrmax^d; ixcmin^d=ilmin^d-kr(idims,^d);
267  jxc^l=ixc^l+kr(idims,^d);
268 
269  do iw=iwstart,nwflux
270  dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
271 
272  if(need_global_a2max) then
273  a2max=a2max_global(idims)
274  else
275  select case(idims)
276  case(1)
277  a2max=schmid_rad1
278  {^iftwod
279  case(2)
280  a2max=schmid_rad2}
281  {^ifthreed
282  case(2)
283  a2max=schmid_rad2
284  case(3)
285  a2max=schmid_rad3}
286  case default
287  call mpistop("idims is wrong in mod_limiter")
288  end select
289  end if
290 
291  call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),rdw=rdw,a2max=a2max)
292 
293  wrc(il^s,iw)=wrc(il^s,iw)-half*rdw(jxr^s)
294  end do
295  end select
296 
297  end subroutine reconstructr
298 
299  subroutine centdiff(method,qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,s,fC,fE,dxs,x)
300 
301  ! Advance the flow variables from global_time to global_time+qdt within ixO^L by
302  ! fourth order centered differencing in space
303  ! for the dw/dt+dF_i(w)/dx_i=S type equation.
304  ! wCT contains the time centered variables at time qtC for flux and source.
305  ! w is the old value at qt on input and the new value at qt+qdt on output.
306  use mod_physics
309  use mod_source, only: addsource2
310  use mod_usr_methods
311  use mod_variables
312  use mod_comm_lib, only: mpistop
313 
314  integer, intent(in) :: method
315  integer, intent(in) :: ixi^l, ixo^l, idims^lim
316  double precision, intent(in) :: qdt, dtfactor, qtc, qt, dxs(ndim)
317  type(state) :: sct, s
318  double precision, intent(in) :: x(ixi^s,1:ndim)
319  double precision :: fc(ixi^s,1:nwflux,1:ndim)
320  double precision :: fe(ixi^s,sdim:3)
321 
322  double precision :: v(ixi^s,ndim), f(ixi^s, nwflux)
323 
324  double precision, dimension(ixI^S,1:nw) :: wprim, wlc, wrc
325  ! left and right constructed status in primitive form, needed for better performance
326  double precision, dimension(ixI^S,1:nw) :: wlp, wrp
327  double precision, dimension(ixI^S) :: vlc, cmaxlc, cmaxrc
328  double precision, dimension(ixI^S,1:number_species) :: cmaxc
329  double precision, dimension(ixI^S,1:number_species) :: cminc
330  double precision, dimension(ixI^S) :: hspeed
331  double precision, dimension(ixO^S) :: inv_volume
332 
333  double precision :: dxinv(1:ndim)
334  integer :: idims, iw, ix^l, hxo^l, ixc^l, jxc^l, hxc^l, kxc^l, kkxc^l, kkxr^l
335  type(ct_velocity) :: vcts
336  logical :: transport, new_cmax, patchw(ixi^s), active
337 
338  associate(wct=>sct%w,w=>s%w)
339  ! two extra layers are needed in each direction for which fluxes are added.
340  ix^l=ixo^l;
341  do idims= idims^lim
342  ix^l=ix^l^ladd2*kr(idims,^d);
343  end do
344 
345  if (ixi^l^ltix^l|.or.|.or.) then
346  call mpistop("Error in centdiff: Non-conforming input limits")
347  end if
348 
349  fc=0.d0
350  wprim=wct
351  call phys_to_primitive(ixi^l,ixi^l,wprim,x)
352 
353  ! get fluxes
354  do idims= idims^lim
355  b0i=idims
356 
357  ix^l=ixo^l^ladd2*kr(idims,^d);
358  hxo^l=ixo^l-kr(idims,^d);
359 
360  if(stagger_grid) then
361  ! ct needs all transverse cells
362  ixcmax^d=ixomax^d+nghostcells-nghostcells*kr(idims,^d);
363  ixcmin^d=hxomin^d-nghostcells+nghostcells*kr(idims,^d);
364  ixmax^d=ixmax^d+nghostcells-nghostcells*kr(idims,^d);
365  ixmin^d=ixmin^d-nghostcells+nghostcells*kr(idims,^d);
366  else
367  ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
368  ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
369  end if
370  hxc^l=ixc^l-kr(idims,^d);
371  jxc^l=ixc^l+kr(idims,^d);
372  kxc^l=ixc^l+2*kr(idims,^d);
373 
374  kkxcmin^d=iximin^d; kkxcmax^d=iximax^d-kr(idims,^d);
375  kkxr^l=kkxc^l+kr(idims,^d);
376  wrp(kkxc^s,1:nwflux)=wprim(kkxr^s,1:nwflux)
377  wlp(kkxc^s,1:nwflux)=wprim(kkxc^s,1:nwflux)
378 
379  ! apply limited reconstruction for left and right status at cell interfaces
380  call reconstruct_lr(ixi^l,ixc^l,ixc^l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
381 
382  if(stagger_grid) then
383  if(h_correction) then
384  call phys_get_h_speed(wprim,x,ixi^l,ixo^l,idims,hspeed)
385  end if
386  call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^l,ixc^l,idims,hspeed,cmaxc,cminc)
387  call phys_get_ct_velocity(vcts,wlp,wrp,ixi^l,ixc^l,idims,cmaxc,cminc)
388  end if
389 
390  ! Calculate velocities from upwinded values
391  call phys_get_cmax(wlc,x,ixg^ll,ixc^l,idims,cmaxlc)
392  call phys_get_cmax(wrc,x,ixg^ll,ixc^l,idims,cmaxrc)
393  ! now take the maximum of left and right states
394  vlc(ixc^s)=max(cmaxrc(ixc^s),cmaxlc(ixc^s))
395 
396  call phys_get_flux(wct,wprim,x,ixi^l,ix^l,idims,f)
397 
398  ! Center flux to interface
399  if(method==fs_cd) then
400  fc(ixc^s,iwstart:nwflux,idims)=half*(f(ixc^s,iwstart:nwflux)+f(jxc^s,iwstart:nwflux))
401  else
402  ! f_i+1/2= (-f_(i+2) +7 f_(i+1) + 7 f_i - f_(i-1))/12
403  fc(ixc^s,iwstart:nwflux,idims)=(-f(kxc^s,iwstart:nwflux)+7.0d0*(f(jxc^s,iwstart:nwflux) + &
404  f(ixc^s,iwstart:nwflux))-f(hxc^s,iwstart:nwflux))/12.0d0
405 
406  do iw=iwstart,nwflux
407  ! add rempel dissipative flux, only second order version for now
408  fc(ixc^s,iw,idims)=fc(ixc^s,iw,idims)-tvdlfeps*half*vlc(ixc^s) &
409  *(wrc(ixc^s,iw)-wlc(ixc^s,iw))
410  end do
411  end if
412 
413  end do !next idims
414  b0i=0
415 
416  if(stagger_grid) call phys_update_faces(ixi^l,ixo^l,qt,qdt,wprim,fc,fe,sct,s,vcts)
417 
418  if(slab_uniform) then
419  dxinv=-qdt/dxs
420  do idims= idims^lim
421  hxo^l=ixo^l-kr(idims,^d);
422  if(local_timestep) then
423  do iw=iwstart,nwflux
424  fc(ixi^s,iw,idims)=-block%dt(ixi^s)*dtfactor/dxs(idims)*fc(ixi^s,iw,idims)
425  ! result: f_(i+1/2)-f_(i-1/2) = [-f_(i+2)+8(f_(i+1)-f_(i-1))+f_(i-2)]/12
426  w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
427  end do !next iw
428  else
429  do iw=iwstart,nwflux
430  fc(ixi^s,iw,idims)=dxinv(idims)*fc(ixi^s,iw,idims)
431  ! result: f_(i+1/2)-f_(i-1/2) = [-f_(i+2)+8(f_(i+1)-f_(i-1))+f_(i-2)]/12
432  w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
433  end do !next iw
434  endif
435  end do ! Next idims
436  else
437  inv_volume=1.d0/block%dvolume
438  do idims= idims^lim
439  hxo^l=ixo^l-kr(idims,^d);
440  if(local_timestep) then
441  do iw=iwstart,nwflux
442  fc(ixi^s,iw,idims)=-block%dt(ixi^s)*dtfactor*block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
443  w(ixo^s,iw)=w(ixo^s,iw)+ &
444  (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
445  end do !next iw
446  else
447  do iw=iwstart,nwflux
448  fc(ixi^s,iw,idims)=-qdt*block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
449  w(ixo^s,iw)=w(ixo^s,iw)+ &
450  (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
451  end do !next iw
452  endif
453  end do ! Next idims
454  end if
455 
456  if (.not.slab.and.idimsmin==1) call phys_add_source_geom(qdt,dtfactor,ixi^l,ixo^l,wct,w,x)
457 
458  if(stagger_grid) call phys_face_to_center(ixo^l,s)
459 
460  ! check and optionally correct unphysical values
461  if(fix_small_values) then
462  call phys_handle_small_values(.false.,w,x,ixi^l,ixo^l,'centdiff')
463  endif
464 
465  call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
466  dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim), &
467  ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,w,x,.false.,active)
468 
469  end associate
470  end subroutine centdiff
471 
472 end module mod_finite_difference
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module with finite difference methods for fluxes.
subroutine, public centdiff(method, qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, s, fC, fE, dxs, x)
subroutine reconstructr(ixIL, iLL, idims, w, wRC)
subroutine reconstructl(ixIL, iLL, idims, w, wLC)
subroutine, public fd(qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, snew, fC, fE, dxs, x)
Module with finite volume methods for fluxes.
subroutine, public reconstruct_lr(ixIL, ixLL, ixRL, idims, w, wLC, wRC, wLp, wRp, x, dxdim)
Determine the upwinded wLC(ixL) and wRC(ixR) from w. the wCT is only used when PPM is exploited.
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.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
integer b0i
background magnetic field location indicator
integer, dimension(:), allocatable, parameter d
logical local_timestep
each cell has its own timestep or not
logical need_global_a2max
global value for schmid scheme
logical slab
Cartesian geometry or not.
integer, dimension(:), allocatable type_limiter
Type of slope limiter used for reconstructing variables on cell edges.
integer nghostcells
Number of ghost cells surrounding a grid.
integer, parameter sdim
starting dimension for electric field
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision, dimension(ndim) a2max_global
global largest a2 for schmid scheme
Module with slope/flux limiters.
Definition: mod_limiter.t:2
subroutine dwlimiter2(dwC, ixIL, ixCL, idims, typelim, ldw, rdw, a2max)
Limit the centered dwC differences within ixC for iw in direction idim. The limiter is chosen accordi...
Definition: mod_limiter.t:129
double precision d
Definition: mod_limiter.t:14
integer, parameter limiter_wenozp5
Definition: mod_limiter.t:35
integer, parameter limiter_weno5
Definition: mod_limiter.t:31
integer, parameter limiter_wenoz5
Definition: mod_limiter.t:33
integer, parameter limiter_wenoz5nm
Definition: mod_limiter.t:34
integer, parameter limiter_wenozp5nm
Definition: mod_limiter.t:36
integer, parameter limiter_mp5
Definition: mod_limiter.t:28
integer, parameter limiter_weno5nm
Definition: mod_limiter.t:32
Module containing the MP5 (fifth order) flux scheme.
Definition: mod_mp5.t:2
subroutine, public mp5limiterr(ixIL, iLL, idims, w, wRC)
Definition: mod_mp5.t:301
subroutine, public mp5limiterl(ixIL, iLL, idims, w, wLC)
Definition: mod_mp5.t:204
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_get_ct_velocity), pointer phys_get_ct_velocity
Definition: mod_physics.t:81
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59
procedure(sub_small_values), pointer phys_handle_small_values
Definition: mod_physics.t:80
procedure(sub_get_flux), pointer phys_get_flux
Definition: mod_physics.t:67
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition: mod_physics.t:66
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition: mod_physics.t:70
procedure(sub_update_faces), pointer phys_update_faces
Definition: mod_physics.t:82
procedure(sub_face_to_center), pointer phys_face_to_center
Definition: mod_physics.t:83
procedure(sub_get_h_speed), pointer phys_get_h_speed
Definition: mod_physics.t:65
procedure(sub_get_cmax), pointer phys_get_cmax
Definition: mod_physics.t:61
Module for handling split source terms (split from the fluxes)
Definition: mod_source.t:2
subroutine, public addsource2(qdt, dtfactor, ixIL, ixOL, iwLIM, qtC, wCT, wCTprim, qt, w, x, qsourcesplit, src_active)
Add source within ixO for iws: w=w+qdt*S[wCT].
Definition: mod_source.t:132
Module with all the methods that users can customize in AMRVAC.
integer nw
Total number of variables.
Definition: mod_variables.t:23
integer iwstart
Definition: mod_variables.t:38
integer nwflux
Number of flux variables.
Definition: mod_variables.t:8