MPI-AMRVAC  3.0
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,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, 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,7-2*ndim: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  do iw=iwstart,nwflux
109  fc(ixi^s,iw,idims) = dxinv(idims) * fc(ixi^s,iw,idims)
110  wnew(ixo^s,iw)=wnew(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
111  end do ! iw loop
112  end do ! Next idims
113  else
114  inv_volume=1.d0/block%dvolume(ixo^s)
115  do idims= idims^lim
116  hxo^l=ixo^l-kr(idims,^d);
117  do iw=iwstart,nwflux
118  fc(ixi^s,iw,idims)=-qdt*fc(ixi^s,iw,idims)*block%surfaceC(ixi^s,idims)
119  wnew(ixo^s,iw)=wnew(ixo^s,iw)+ &
120  (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
121  end do ! iw loop
122  end do ! Next idims
123  end if
124 
125  if (.not.slab.and.idimsmin==1) call phys_add_source_geom(qdt,ixi^l,ixo^l,wct,wnew,x)
126 
127  if(stagger_grid) call phys_face_to_center(ixo^l,snew)
128 
129  ! check and optionally correct unphysical values
130  if(fix_small_values) then
131  call phys_handle_small_values(.false.,wnew,x,ixi^l,ixo^l,'fd')
132  endif
133 
134  call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
135  ixi^l,ixo^l,1,nw,qtc,wct,qt,wnew,x,.false.,active,wprim)
136 
137  if(phys_solve_eaux.and.levmin==levmax) then
138  ! synchronize internal energy for uniform grid
139  call phys_energy_synchro(ixi^l,ixo^l,wnew,x)
140  endif
141  end associate
142 
143  end subroutine fd
144 
145  subroutine reconstructl(ixI^L,iL^L,idims,w,wLC)
147  use mod_mp5
148  use mod_limiter
149 
150  integer, intent(in) :: ixI^L, iL^L, idims
151  double precision, intent(in) :: w(ixI^S,1:nw)
152 
153  double precision, intent(out) :: wLC(ixI^S,1:nw)
154 
155  double precision :: ldw(ixI^S), dwC(ixI^S)
156  integer :: jxR^L, ixC^L, jxC^L, kxC^L, iw
157  double precision :: a2max
158 
159  select case (type_limiter(block%level))
160  case (limiter_mp5)
161  call mp5limiterl(ixi^l,il^l,idims,w,wlc)
162  case (limiter_weno5)
163  call weno5limiterl(ixi^l,il^l,idims,w,wlc,1)
164  case (limiter_weno5nm)
165  call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,1)
166  case (limiter_wenoz5)
167  call weno5limiterl(ixi^l,il^l,idims,w,wlc,2)
168  case (limiter_wenoz5nm)
169  call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,2)
170  case (limiter_wenozp5)
171  call weno5limiterl(ixi^l,il^l,idims,w,wlc,3)
172  case (limiter_wenozp5nm)
173  call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,3)
174  case default
175 
176  kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
177 
178  wlc(kxc^s,iwstart:nwflux) = w(kxc^s,iwstart:nwflux)
179 
180  jxr^l=il^l+kr(idims,^d);
181 
182  ixcmax^d=jxrmax^d; ixcmin^d=ilmin^d-kr(idims,^d);
183  jxc^l=ixc^l+kr(idims,^d);
184 
185  do iw=iwstart,nwflux
186  dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
187 
188  if(need_global_a2max) then
189  a2max=a2max_global(idims)
190  else
191  select case(idims)
192  case(1)
193  a2max=schmid_rad1
194  {^iftwod
195  case(2)
196  a2max=schmid_rad2}
197  {^ifthreed
198  case(2)
199  a2max=schmid_rad2
200  case(3)
201  a2max=schmid_rad3}
202  case default
203  call mpistop("idims is wrong in mod_limiter")
204  end select
205  end if
206 
207  call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),ldw=ldw,a2max=a2max)
208 
209  wlc(il^s,iw)=wlc(il^s,iw)+half*ldw(il^s)
210  end do
211  end select
212 
213  end subroutine reconstructl
214 
215  subroutine reconstructr(ixI^L,iL^L,idims,w,wRC)
217  use mod_mp5
218  use mod_limiter
219 
220  integer, intent(in) :: ixI^L, iL^L, idims
221  double precision, intent(in) :: w(ixI^S,1:nw)
222 
223  double precision, intent(out) :: wRC(ixI^S,1:nw)
224 
225  double precision :: rdw(ixI^S), dwC(ixI^S)
226  integer :: jxR^L, ixC^L, jxC^L, kxC^L, kxR^L, iw
227  double precision :: a2max
228 
229  select case (type_limiter(block%level))
230  case (limiter_mp5)
231  call mp5limiterr(ixi^l,il^l,idims,w,wrc)
232  case (limiter_weno5)
233  call weno5limiterr(ixi^l,il^l,idims,w,wrc,1)
234  case (limiter_weno5nm)
235  call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,1)
236  case (limiter_wenoz5)
237  call weno5limiterr(ixi^l,il^l,idims,w,wrc,2)
238  case (limiter_wenoz5nm)
239  call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,2)
240  case (limiter_wenozp5)
241  call weno5limiterr(ixi^l,il^l,idims,w,wrc,3)
242  case (limiter_wenozp5nm)
243  call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,3)
244  case default
245 
246  kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
247  kxr^l=kxc^l+kr(idims,^d);
248 
249  wrc(kxc^s,iwstart:nwflux)=w(kxr^s,iwstart:nwflux)
250 
251  jxr^l=il^l+kr(idims,^d);
252  ixcmax^d=jxrmax^d; ixcmin^d=ilmin^d-kr(idims,^d);
253  jxc^l=ixc^l+kr(idims,^d);
254 
255  do iw=iwstart,nwflux
256  dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
257 
258  if(need_global_a2max) then
259  a2max=a2max_global(idims)
260  else
261  select case(idims)
262  case(1)
263  a2max=schmid_rad1
264  {^iftwod
265  case(2)
266  a2max=schmid_rad2}
267  {^ifthreed
268  case(2)
269  a2max=schmid_rad2
270  case(3)
271  a2max=schmid_rad3}
272  case default
273  call mpistop("idims is wrong in mod_limiter")
274  end select
275  end if
276 
277  call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),rdw=rdw,a2max=a2max)
278 
279  wrc(il^s,iw)=wrc(il^s,iw)-half*rdw(jxr^s)
280  end do
281  end select
282 
283  end subroutine reconstructr
284 
285  subroutine centdiff(method,qdt,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,s,fC,fE,dxs,x)
286 
287  ! Advance the flow variables from global_time to global_time+qdt within ixO^L by
288  ! fourth order centered differencing in space
289  ! for the dw/dt+dF_i(w)/dx_i=S type equation.
290  ! wCT contains the time centered variables at time qtC for flux and source.
291  ! w is the old value at qt on input and the new value at qt+qdt on output.
292  use mod_physics
295  use mod_source, only: addsource2
296  use mod_usr_methods
297  use mod_variables
298 
299  integer, intent(in) :: method
300  integer, intent(in) :: ixi^l, ixo^l, idims^lim
301  double precision, intent(in) :: qdt, qtc, qt, dxs(ndim)
302  type(state) :: sct, s
303  double precision, intent(in) :: x(ixi^s,1:ndim)
304  double precision :: fc(ixi^s,1:nwflux,1:ndim)
305  double precision :: fe(ixi^s,7-2*ndim:3)
306 
307  double precision :: v(ixi^s,ndim), f(ixi^s, nwflux)
308 
309  double precision, dimension(ixI^S,1:nw) :: wprim, wlc, wrc
310  ! left and right constructed status in primitive form, needed for better performance
311  double precision, dimension(ixI^S,1:nw) :: wlp, wrp
312  double precision, dimension(ixI^S) :: vlc, cmaxlc, cmaxrc
313  double precision, dimension(ixI^S,1:number_species) :: cmaxc
314  double precision, dimension(ixI^S,1:number_species) :: cminc
315  double precision, dimension(ixI^S) :: hspeed
316  double precision, dimension(ixO^S) :: inv_volume
317 
318  double precision :: dxinv(1:ndim)
319  integer :: idims, iw, ix^l, hxo^l, ixc^l, jxc^l, hxc^l, kxc^l, kkxc^l, kkxr^l
320  type(ct_velocity) :: vcts
321  logical :: transport, new_cmax, patchw(ixi^s), active
322 
323  associate(wct=>sct%w,w=>s%w)
324  ! two extra layers are needed in each direction for which fluxes are added.
325  ix^l=ixo^l;
326  do idims= idims^lim
327  ix^l=ix^l^ladd2*kr(idims,^d);
328  end do
329 
330  if (ixi^l^ltix^l|.or.|.or.) then
331  call mpistop("Error in centdiff: Non-conforming input limits")
332  end if
333 
334  fc=0.d0
335  wprim=wct
336  call phys_to_primitive(ixi^l,ixi^l,wprim,x)
337 
338  ! get fluxes
339  do idims= idims^lim
340  b0i=idims
341 
342  ix^l=ixo^l^ladd2*kr(idims,^d);
343  hxo^l=ixo^l-kr(idims,^d);
344 
345  if(stagger_grid) then
346  ! ct needs all transverse cells
347  ixcmax^d=ixomax^d+nghostcells-nghostcells*kr(idims,^d);
348  ixcmin^d=hxomin^d-nghostcells+nghostcells*kr(idims,^d);
349  ixmax^d=ixmax^d+nghostcells-nghostcells*kr(idims,^d);
350  ixmin^d=ixmin^d-nghostcells+nghostcells*kr(idims,^d);
351  else
352  ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
353  ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
354  end if
355  hxc^l=ixc^l-kr(idims,^d);
356  jxc^l=ixc^l+kr(idims,^d);
357  kxc^l=ixc^l+2*kr(idims,^d);
358 
359  kkxcmin^d=iximin^d; kkxcmax^d=iximax^d-kr(idims,^d);
360  kkxr^l=kkxc^l+kr(idims,^d);
361  wrp(kkxc^s,1:nwflux)=wprim(kkxr^s,1:nwflux)
362  wlp(kkxc^s,1:nwflux)=wprim(kkxc^s,1:nwflux)
363 
364  ! apply limited reconstruction for left and right status at cell interfaces
365  call reconstruct_lr(ixi^l,ixc^l,ixc^l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
366 
367  if(stagger_grid) then
368  if(h_correction) then
369  call phys_get_h_speed(wprim,x,ixi^l,ixo^l,idims,hspeed)
370  end if
371  call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^l,ixc^l,idims,hspeed,cmaxc,cminc)
372  call phys_get_ct_velocity(vcts,wlp,wrp,ixi^l,ixc^l,idims,cmaxc,cminc)
373  end if
374 
375  ! Calculate velocities from upwinded values
376  call phys_get_cmax(wlc,x,ixg^ll,ixc^l,idims,cmaxlc)
377  call phys_get_cmax(wrc,x,ixg^ll,ixc^l,idims,cmaxrc)
378  ! now take the maximum of left and right states
379  vlc(ixc^s)=max(cmaxrc(ixc^s),cmaxlc(ixc^s))
380 
381  call phys_get_flux(wct,wprim,x,ixi^l,ix^l,idims,f)
382 
383  ! Center flux to interface
384  if(method==fs_cd) then
385  fc(ixc^s,iwstart:nwflux,idims)=half*(f(ixc^s,iwstart:nwflux)+f(jxc^s,iwstart:nwflux))
386  else
387  ! f_i+1/2= (-f_(i+2) +7 f_(i+1) + 7 f_i - f_(i-1))/12
388  fc(ixc^s,iwstart:nwflux,idims)=(-f(kxc^s,iwstart:nwflux)+7.0d0*(f(jxc^s,iwstart:nwflux) + &
389  f(ixc^s,iwstart:nwflux))-f(hxc^s,iwstart:nwflux))/12.0d0
390 
391  do iw=iwstart,nwflux
392  ! add rempel dissipative flux, only second order version for now
393  fc(ixc^s,iw,idims)=fc(ixc^s,iw,idims)-tvdlfeps*half*vlc(ixc^s) &
394  *(wrc(ixc^s,iw)-wlc(ixc^s,iw))
395  end do
396  end if
397 
398  end do !next idims
399  b0i=0
400 
401  if(stagger_grid) call phys_update_faces(ixi^l,ixo^l,qt,qdt,wprim,fc,fe,sct,s,vcts)
402 
403  if(slab_uniform) then
404  dxinv=-qdt/dxs
405  do idims= idims^lim
406  hxo^l=ixo^l-kr(idims,^d);
407  do iw=iwstart,nwflux
408  fc(ixi^s,iw,idims)=dxinv(idims)*fc(ixi^s,iw,idims)
409  ! result: f_(i+1/2)-f_(i-1/2) = [-f_(i+2)+8(f_(i+1)-f_(i-1))+f_(i-2)]/12
410  w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
411  end do !next iw
412  end do ! Next idims
413  else
414  inv_volume=1.d0/block%dvolume
415  do idims= idims^lim
416  hxo^l=ixo^l-kr(idims,^d);
417  do iw=iwstart,nwflux
418  fc(ixi^s,iw,idims)=-qdt*block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
419  w(ixo^s,iw)=w(ixo^s,iw)+ &
420  (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
421  end do !next iw
422  end do ! Next idims
423  end if
424 
425  if (.not.slab.and.idimsmin==1) call phys_add_source_geom(qdt,ixi^l,ixo^l,wct,w,x)
426 
427  if(stagger_grid) call phys_face_to_center(ixo^l,s)
428 
429  ! check and optionally correct unphysical values
430  if(fix_small_values) then
431  call phys_handle_small_values(.false.,w,x,ixi^l,ixo^l,'centdiff')
432  endif
433 
434  call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
435  ixi^l,ixo^l,1,nw,qtc,wct,qt,w,x,.false.,active,wprim)
436 
437  if(phys_solve_eaux.and.levmin==levmax) then
438  ! synchronize internal energy for uniform grid
439  call phys_energy_synchro(ixi^l,ixo^l,w,x)
440  endif
441  end associate
442  end subroutine centdiff
443 
444 end module mod_finite_difference
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:194
Module with finite difference methods for fluxes.
subroutine, public fd(qdt, ixIL, ixOL, idimsLIM, qtC, sCT, qt, snew, fC, fE, dxs, x)
subroutine reconstructr(ixIL, iLL, idims, w, wRC)
subroutine reconstructl(ixIL, iLL, idims, w, wLC)
subroutine, public centdiff(method, qdt, ixIL, ixOL, idimsLIM, qtC, sCT, qt, s, 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 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.
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:128
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:82
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:81
procedure(sub_get_flux), pointer phys_get_flux
Definition: mod_physics.t:66
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition: mod_physics.t:65
logical phys_solve_eaux
Solve internal energy and total energy equations.
Definition: mod_physics.t:39
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:83
procedure(sub_energy_synchro), pointer phys_energy_synchro
Definition: mod_physics.t:67
procedure(sub_face_to_center), pointer phys_face_to_center
Definition: mod_physics.t:84
procedure(sub_get_h_speed), pointer phys_get_h_speed
Definition: mod_physics.t:64
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, ixIL, ixOL, iwLIM, qtC, wCT, qt, w, x, qsourcesplit, src_active, wCTprim)
Add source within ixO for iws: w=w+qdt*S[wCT].
Definition: mod_source.t:121
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