MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
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
10contains
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
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,wprim,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 if(crash) then
148 ! replace erroneous values with values at previous step
149 wnew=wct
150 if(stagger_grid) snew%ws=sct%ws
151 end if
152 endif
153
154 call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
155 dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim), &
156 ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
157
158 end associate
159
160 end subroutine fd
161
162 subroutine reconstructl(ixI^L,iL^L,idims,w,wLC)
164 use mod_mp5
165 use mod_limiter
166 use mod_comm_lib, only: mpistop
167
168 integer, intent(in) :: ixi^l, il^l, idims
169 double precision, intent(in) :: w(ixi^s,1:nw)
170
171 double precision, intent(out) :: wlc(ixi^s,1:nw)
172
173 double precision :: ldw(ixi^s), dwc(ixi^s)
174 integer :: jxr^l, ixc^l, jxc^l, kxc^l, iw
175
176 select case (type_limiter(block%level))
177 case (limiter_mp5)
178 call mp5limiterl(ixi^l,il^l,idims,w,wlc)
179 case (limiter_weno5)
180 call weno5limiterl(ixi^l,il^l,idims,w,wlc,1)
181 case (limiter_weno5nm)
182 call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,1)
183 case (limiter_wenoz5)
184 call weno5limiterl(ixi^l,il^l,idims,w,wlc,2)
185 case (limiter_wenoz5nm)
186 call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,2)
187 case (limiter_wenozp5)
188 call weno5limiterl(ixi^l,il^l,idims,w,wlc,3)
189 case (limiter_wenozp5nm)
190 call weno5nmlimiterl(ixi^l,il^l,idims,w,wlc,3)
191 case default
192
193 kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
194
195 wlc(kxc^s,iwstart:nwflux) = w(kxc^s,iwstart:nwflux)
196
197 jxr^l=il^l+kr(idims,^d);
198
199 ixcmax^d=jxrmax^d; ixcmin^d=ilmin^d-kr(idims,^d);
200 jxc^l=ixc^l+kr(idims,^d);
201
202 do iw=iwstart,nwflux
203 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
204
205 call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),ldw=ldw)
206
207 wlc(il^s,iw)=wlc(il^s,iw)+half*ldw(il^s)
208 end do
209 end select
210
211 end subroutine reconstructl
212
213 subroutine reconstructr(ixI^L,iL^L,idims,w,wRC)
215 use mod_mp5
216 use mod_limiter
217 use mod_comm_lib, only: mpistop
218
219 integer, intent(in) :: ixi^l, il^l, idims
220 double precision, intent(in) :: w(ixi^s,1:nw)
221
222 double precision, intent(out) :: wrc(ixi^s,1:nw)
223
224 double precision :: rdw(ixi^s), dwc(ixi^s)
225 integer :: jxr^l, ixc^l, jxc^l, kxc^l, kxr^l, iw
226
227 select case (type_limiter(block%level))
228 case (limiter_mp5)
229 call mp5limiterr(ixi^l,il^l,idims,w,wrc)
230 case (limiter_weno5)
231 call weno5limiterr(ixi^l,il^l,idims,w,wrc,1)
232 case (limiter_weno5nm)
233 call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,1)
234 case (limiter_wenoz5)
235 call weno5limiterr(ixi^l,il^l,idims,w,wrc,2)
236 case (limiter_wenoz5nm)
237 call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,2)
238 case (limiter_wenozp5)
239 call weno5limiterr(ixi^l,il^l,idims,w,wrc,3)
240 case (limiter_wenozp5nm)
241 call weno5nmlimiterr(ixi^l,il^l,idims,w,wrc,3)
242 case default
243
244 kxcmin^d=iximin^d; kxcmax^d=iximax^d-kr(idims,^d);
245 kxr^l=kxc^l+kr(idims,^d);
246
247 wrc(kxc^s,iwstart:nwflux)=w(kxr^s,iwstart:nwflux)
248
249 jxr^l=il^l+kr(idims,^d);
250 ixcmax^d=jxrmax^d; ixcmin^d=ilmin^d-kr(idims,^d);
251 jxc^l=ixc^l+kr(idims,^d);
252
253 do iw=iwstart,nwflux
254 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
255
256 call dwlimiter2(dwc,ixi^l,ixc^l,idims,type_limiter(block%level),rdw=rdw)
257
258 wrc(il^s,iw)=wrc(il^s,iw)-half*rdw(jxr^s)
259 end do
260 end select
261
262 end subroutine reconstructr
263
264 subroutine centdiff(method,qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,s,fC,fE,dxs,x)
265
266 ! Advance the flow variables from global_time to global_time+qdt within ixO^L by
267 ! fourth order centered differencing in space
268 ! for the dw/dt+dF_i(w)/dx_i=S type equation.
269 ! wCT contains the time centered variables at time qtC for flux and source.
270 ! w is the old value at qt on input and the new value at qt+qdt on output.
271 use mod_physics
274 use mod_source, only: addsource2
276 use mod_variables
277 use mod_comm_lib, only: mpistop
278
279 integer, intent(in) :: method
280 integer, intent(in) :: ixi^l, ixo^l, idims^lim
281 double precision, intent(in) :: qdt, dtfactor, qtc, qt, dxs(ndim)
282 type(state) :: sct, s
283 double precision, intent(in) :: x(ixi^s,1:ndim)
284 double precision :: fc(ixi^s,1:nwflux,1:ndim)
285 double precision :: fe(ixi^s,sdim:3)
286
287 double precision :: v(ixi^s,ndim), f(ixi^s, nwflux)
288
289 double precision, dimension(ixI^S,1:nw) :: wprim, wlc, wrc
290 ! left and right constructed status in primitive form, needed for better performance
291 double precision, dimension(ixI^S,1:nw) :: wlp, wrp
292 double precision, dimension(ixI^S) :: vlc, cmaxlc, cmaxrc
293 double precision, dimension(ixI^S,1:number_species) :: cmaxc
294 double precision, dimension(ixI^S,1:number_species) :: cminc
295 double precision, dimension(ixI^S) :: hspeed
296 double precision, dimension(ixO^S) :: inv_volume
297
298 double precision :: dxinv(1:ndim)
299 integer :: idims, iw, ix^l, hxo^l, ixc^l, jxc^l, hxc^l, kxc^l, kkxc^l, kkxr^l
300 type(ct_velocity) :: vcts
301 logical :: transport, new_cmax, patchw(ixi^s), active
302
303 associate(wct=>sct%w,w=>s%w)
304 ! two extra layers are needed in each direction for which fluxes are added.
305 ix^l=ixo^l;
306 do idims= idims^lim
307 ix^l=ix^l^ladd2*kr(idims,^d);
308 end do
309
310 if (ixi^l^ltix^l|.or.|.or.) then
311 call mpistop("Error in centdiff: Non-conforming input limits")
312 end if
313
314 fc=0.d0
315 wprim=wct
316 call phys_to_primitive(ixi^l,ixi^l,wprim,x)
317
318 ! get fluxes
319 do idims= idims^lim
320 b0i=idims
321
322 ix^l=ixo^l^ladd2*kr(idims,^d);
323 hxo^l=ixo^l-kr(idims,^d);
324
325 if(stagger_grid) then
326 ! ct needs all transverse cells
327 ixcmax^d=ixomax^d+nghostcells-nghostcells*kr(idims,^d);
328 ixcmin^d=hxomin^d-nghostcells+nghostcells*kr(idims,^d);
329 ixmax^d=ixmax^d+nghostcells-nghostcells*kr(idims,^d);
330 ixmin^d=ixmin^d-nghostcells+nghostcells*kr(idims,^d);
331 else
332 ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
333 ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
334 end if
335 hxc^l=ixc^l-kr(idims,^d);
336 jxc^l=ixc^l+kr(idims,^d);
337 kxc^l=ixc^l+2*kr(idims,^d);
338
339 kkxcmin^d=iximin^d; kkxcmax^d=iximax^d-kr(idims,^d);
340 kkxr^l=kkxc^l+kr(idims,^d);
341 wrp(kkxc^s,1:nwflux)=wprim(kkxr^s,1:nwflux)
342 wlp(kkxc^s,1:nwflux)=wprim(kkxc^s,1:nwflux)
343
344 ! apply limited reconstruction for left and right status at cell interfaces
345 call reconstruct_lr(ixi^l,ixc^l,ixc^l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
346
347 if(stagger_grid) then
348 if(h_correction) then
349 call phys_get_h_speed(wprim,x,ixi^l,ixo^l,idims,hspeed)
350 end if
351 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^l,ixc^l,idims,hspeed,cmaxc,cminc)
352 call phys_get_ct_velocity(vcts,wlp,wrp,ixi^l,ixc^l,idims,cmaxc,cminc)
353 end if
354
355 ! Calculate velocities from upwinded values
356 call phys_get_cmax(wlp,x,ixg^ll,ixc^l,idims,cmaxlc)
357 call phys_get_cmax(wrp,x,ixg^ll,ixc^l,idims,cmaxrc)
358 ! now take the maximum of left and right states
359 vlc(ixc^s)=max(cmaxrc(ixc^s),cmaxlc(ixc^s))
360
361 call phys_get_flux(wct,wprim,x,ixi^l,ix^l,idims,f)
362
363 ! Center flux to interface
364 if(method==fs_cd) then
365 fc(ixc^s,iwstart:nwflux,idims)=half*(f(ixc^s,iwstart:nwflux)+f(jxc^s,iwstart:nwflux))
366 else
367 ! f_i+1/2= (-f_(i+2) +7 f_(i+1) + 7 f_i - f_(i-1))/12
368 fc(ixc^s,iwstart:nwflux,idims)=(-f(kxc^s,iwstart:nwflux)+7.0d0*(f(jxc^s,iwstart:nwflux) + &
369 f(ixc^s,iwstart:nwflux))-f(hxc^s,iwstart:nwflux))/12.0d0
370
371 do iw=iwstart,nwflux
372 ! add rempel dissipative flux, only second order version for now
373 fc(ixc^s,iw,idims)=fc(ixc^s,iw,idims)-tvdlfeps*half*vlc(ixc^s) &
374 *(wrc(ixc^s,iw)-wlc(ixc^s,iw))
375 end do
376 end if
377
378 end do !next idims
379 b0i=0
380
381 if(stagger_grid) call phys_update_faces(ixi^l,ixo^l,qt,qdt,wprim,fc,fe,sct,s,vcts)
382
383 if(slab_uniform) then
384 dxinv=-qdt/dxs
385 do idims= idims^lim
386 hxo^l=ixo^l-kr(idims,^d);
387 if(local_timestep) then
388 do iw=iwstart,nwflux
389 fc(ixi^s,iw,idims)=-block%dt(ixi^s)*dtfactor/dxs(idims)*fc(ixi^s,iw,idims)
390 ! result: f_(i+1/2)-f_(i-1/2) = [-f_(i+2)+8(f_(i+1)-f_(i-1))+f_(i-2)]/12
391 w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
392 end do !next iw
393 else
394 do iw=iwstart,nwflux
395 fc(ixi^s,iw,idims)=dxinv(idims)*fc(ixi^s,iw,idims)
396 ! result: f_(i+1/2)-f_(i-1/2) = [-f_(i+2)+8(f_(i+1)-f_(i-1))+f_(i-2)]/12
397 w(ixo^s,iw)=w(ixo^s,iw)+(fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))
398 end do !next iw
399 endif
400 end do ! Next idims
401 else
402 inv_volume=1.d0/block%dvolume
403 do idims= idims^lim
404 hxo^l=ixo^l-kr(idims,^d);
405 if(local_timestep) then
406 do iw=iwstart,nwflux
407 fc(ixi^s,iw,idims)=-block%dt(ixi^s)*dtfactor*block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
408 w(ixo^s,iw)=w(ixo^s,iw)+ &
409 (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
410 end do !next iw
411 else
412 do iw=iwstart,nwflux
413 fc(ixi^s,iw,idims)=-qdt*block%surfaceC(ixi^s,idims)*fc(ixi^s,iw,idims)
414 w(ixo^s,iw)=w(ixo^s,iw)+ &
415 (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims))*inv_volume(ixo^s)
416 end do !next iw
417 endif
418 end do ! Next idims
419 end if
420
421 if (.not.slab.and.idimsmin==1) call phys_add_source_geom(qdt,dtfactor,ixi^l,ixo^l,wct,wprim,w,x)
422
423 if(stagger_grid) call phys_face_to_center(ixo^l,s)
424
425 ! check and optionally correct unphysical values
426 if(fix_small_values) then
427 call phys_handle_small_values(.false.,w,x,ixi^l,ixo^l,'centdiff')
428 endif
429
430 call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim), &
431 dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim), &
432 ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,w,x,.false.,active)
433
434 end associate
435 end subroutine centdiff
436
437end module mod_finite_difference
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with finite difference methods for fluxes.
subroutine, public fd(qdt, dtfactor, ixil, ixol, idimslim, qtc, sct, qt, snew, fc, fe, dxs, x)
subroutine, public centdiff(method, qdt, dtfactor, 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
double precision, dimension(:), allocatable, parameter d
logical local_timestep
each cell has its own timestep or not
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 fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
Module with slope/flux limiters.
Definition mod_limiter.t:2
subroutine dwlimiter2(dwc, ixil, ixcl, idims, typelim, ldw, rdw)
Limit the centered dwC differences within ixC for iw in direction idim. The limiter is chosen accordi...
integer, parameter limiter_wenozp5
Definition mod_limiter.t:33
integer, parameter limiter_weno5
Definition mod_limiter.t:29
integer, parameter limiter_wenoz5
Definition mod_limiter.t:31
integer, parameter limiter_wenoz5nm
Definition mod_limiter.t:32
integer, parameter limiter_wenozp5nm
Definition mod_limiter.t:34
integer, parameter limiter_mp5
Definition mod_limiter.t:26
integer, parameter limiter_weno5nm
Definition mod_limiter.t:30
Module containing the MP5 (fifth order) flux scheme.
Definition mod_mp5.t:2
subroutine, public mp5limiterl(ixil, ill, idims, w, wlc)
Definition mod_mp5.t:204
subroutine, public mp5limiterr(ixil, ill, idims, w, wrc)
Definition mod_mp5.t:301
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:75
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:54
procedure(sub_small_values), pointer phys_handle_small_values
Definition mod_physics.t:74
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:61
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:60
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:65
procedure(sub_update_faces), pointer phys_update_faces
Definition mod_physics.t:76
procedure(sub_face_to_center), pointer phys_face_to_center
Definition mod_physics.t:77
procedure(sub_get_h_speed), pointer phys_get_h_speed
Definition mod_physics.t:59
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:56
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.
integer nwflux
Number of flux variables.