MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_ppm.t
Go to the documentation of this file.
1 module mod_ppm
2 
3  implicit none
4  private
5 
6  public :: ppmlimiter
7  public :: ppmlimitervar
8 
9 contains
10 
11  subroutine ppmlimitervar(ixI^L,ix^L,idims,q,qCT,qLC,qRC)
12 
13  ! references:
14  ! Mignone et al 2005, ApJS 160, 199,
15  ! Miller and Colella 2002, JCP 183, 26
16  ! Fryxell et al. 2000 ApJ, 131, 273 (Flash)
17  ! baciotti Phd (http://www.aei.mpg.de/~baiotti/Baiotti_PhD.pdf)
18  ! version : april 2009
19  ! author: zakaria.meliani@wis.kuleuven.be
20 
22 
23  integer, intent(in) :: ixi^l, ix^l, idims
24  double precision, intent(in) :: q(ixi^s),qct(ixi^s)
25 
26  double precision, intent(inout) :: qrc(ixg^t),qlc(ixg^t)
27 
28  double precision,dimension(ixG^T) :: dqc,d2qc,ldq
29  double precision,dimension(ixG^T) :: qmin,qmax,tmp
30 
31  integer :: lxc^l,lxr^l
32  integer :: ixll^l,ixl^l,ixo^l,ixr^l,ixrr^l
33  integer :: hxl^l,hxc^l,hxr^l
34  integer :: kxll^l,kxl^l,kxc^l,kxr^l,kxrr^l
35 
36  ixomin^d=ixmin^d-kr(idims,^d);ixomax^d=ixmax^d+kr(idims,^d);!ixO[ixMmin1-1,ixMmax1+1]
37  ixl^l=ixo^l-kr(idims,^d); !ixL[ixMmin1-2,ixMmax1]
38  ixll^l=ixl^l-kr(idims,^d); !ixLL[ixMmin1-3,ixMmax1-1]
39  ixr^l=ixo^l+kr(idims,^d); !ixR=[iMmin1,ixMmax+2]
40  ixrr^l=ixr^l+kr(idims,^d); !ixRR=[iMmin1+1,ixMmax+3]
41 
42  hxcmin^d=ixomin^d;hxcmax^d=ixmax^d; ! hxC = [ixMmin-1,ixMmax]
43  hxl^l=hxc^l-kr(idims,^d); ! hxL = [ixMmin-2,ixMmax-1]
44  hxr^l=hxc^l+kr(idims,^d); ! hxR = [ixMmin,ixMmax+1]
45 
46  kxcmin^d=ixllmin^d; kxcmax^d=ixrmax^d; ! kxC=[iMmin1-3,ixMmax1+2]
47  kxl^l=kxc^l-kr(idims,^d); ! kxL=[iMmin1-4,ixMmax1+1]
48  kxr^l=kxc^l+kr(idims,^d); ! kxR=[iMmin1-2,ixMmax1+3]
49 
50  lxcmin^d=ixllmin^d-kr(idims,^d);lxcmax^d=ixrrmax^d;! ixC=[iMmin1-4,ixMmax1+3]
51  lxr^l=lxc^l+kr(idims,^d); ! lxR=[iMmin1-3,ixMmax1+4]
52 
53 
54  dqc(lxc^s)=q(lxr^s)-q(lxc^s)
55  ! Eq. 64, Miller and Colella 2002, JCP 183, 26
56  d2qc(kxc^s)=half*(q(kxr^s)-q(kxl^s))
57  where(dqc(kxc^s)*dqc(kxl^s)>zero)
58  ! Store the sign of d2qC in qMin
59  qmin(kxc^s)= sign(one,d2qc(kxc^s))
60  ! Eq. 65, Miller and Colella 2002, JCP 183, 26
61  ldq(kxc^s)= qmin(kxc^s)*min(dabs(d2qc(kxc^s)),&
62  2.0d0*dabs(dqc(kxl^s)),&
63  2.0d0*dabs(dqc(kxc^s)))
64  elsewhere
65  ldq(kxc^s)=zero
66  end where
67 
68  ! Eq. 66, Miller and Colella 2002, JCP 183, 26
69  qlc(ixo^s)=qlc(ixo^s)+half*dqc(ixo^s)&
70  +(ldq(ixo^s)-ldq(ixr^s))/6.0d0
71  qrc(ixl^s)=qrc(ixl^s) -(half*dqc(ixl^s)&
72  -(ldq(ixl^s)-ldq(ixo^s))/6.0d0)
73 
74  ! make sure that min wCT(i)<wLC(i)<wCT(i+1) same for wRC(i)
75  call extremaq(ixi^l,kxc^l,qct,1,qmax,qmin)
76 
77  ! Eq. B8, page 217, Mignone et al 2005, ApJS
78  qrc(ixl^s)=max(qmin(ixo^s),min(qmax(ixo^s),qrc(ixl^s)))
79  qlc(ixo^s)=max(qmin(ixo^s),min(qmax(ixo^s),qlc(ixo^s)))
80 
81  ! Eq. B9, page 217, Mignone et al 2005, ApJS
82  where((qrc(ixl^s)-qct(ixo^s))*(qct(ixo^s)-qlc(ixo^s))<=zero)
83  qrc(ixl^s)=qct(ixo^s)
84  qlc(ixo^s)=qct(ixo^s)
85  end where
86 
87  qmax(ixo^s)=(qlc(ixo^s)-qrc(ixl^s))&
88  *(qct(ixo^s)-(qlc(ixo^s)+qrc(ixl^s))/2.0d0)
89  qmin(ixo^s)=(qlc(ixo^s)-qrc(ixl^s))**2.0d0/6.0d0
90  tmp(ixl^s)=qrc(ixl^s)
91 
92  ! Eq. B10, page 218, Mignone et al 2005, ApJS
93  where(qmax(hxr^s)>qmin(hxr^s))
94  qrc(hxc^s)= 3.0d0*qct(hxr^s)-2.0d0*qlc(hxr^s)
95  end where
96  ! Eq. B11, page 218, Mignone et al 2005, ApJS
97  where(qmax(hxc^s)<-qmin(hxc^s))
98  qlc(hxc^s)= 3.0d0*qct(hxc^s)-2.0d0*tmp(hxl^s)
99  end where
100 
101  end subroutine ppmlimitervar
102 
103  subroutine ppmlimiter(ixI^L,ix^L,idims,w,wCT,wLC,wRC)
104 
105  ! references:
106  ! Mignone et al 2005, ApJS 160, 199,
107  ! Miller and Colella 2002, JCP 183, 26
108  ! Fryxell et al. 2000 ApJ, 131, 273 (Flash)
109  ! baciotti Phd (http://www.aei.mpg.de/~baiotti/Baiotti_PhD.pdf)
110  ! version : april 2009
111  ! author: zakaria.meliani@wis.kuleuven.be
112 
114  use mod_physics, only: phys_ppm_flatcd, phys_ppm_flatsh
115 
116  integer, intent(in) :: ixi^l, ix^l, idims
117  double precision, intent(in) :: w(ixi^s,1:nw),wct(ixi^s,1:nw)
118 
119  double precision, intent(inout) :: wrc(ixg^t,1:nw),wlc(ixg^t,1:nw)
120 
121  double precision,dimension(ixG^T,1:nwflux) :: dwc,d2wc,ldw
122  double precision,dimension(ixG^T,1:nwflux) :: wmin,wmax,tmp
123  double precision,dimension(ixG^T) :: aa, ab, ac, dv
124  double precision,dimension(ixG^T,1:ndim) :: exi
125 
126  integer :: lxc^l,lxr^l
127  integer :: ixll^l,ixl^l,ixo^l,ixr^l,ixrr^l
128  integer :: hxl^l,hxc^l,hxr^l
129  integer :: kxll^l,kxl^l,kxc^l,kxr^l,kxrr^l
130  integer :: iw, idimss
131 
132  double precision, parameter :: betamin=0.75d0, betamax=0.85d0,&
133  zmin=0.25d0, zmax=0.75d0,&
134  eta1=20.0d0,eta2=0.05d0,eps=0.01d0,kappa=0.1d0
135 
136  ixomin^d=ixmin^d-kr(idims,^d);ixomax^d=ixmax^d+kr(idims,^d);!ixO[ixMmin1-1,ixMmax1+1]
137  ixl^l=ixo^l-kr(idims,^d); !ixL[ixMmin1-2,ixMmax1]
138  ixll^l=ixl^l-kr(idims,^d); !ixLL[ixMmin1-3,ixMmax1-1]
139  ixr^l=ixo^l+kr(idims,^d); !ixR=[iMmin1,ixMmax+2]
140  ixrr^l=ixr^l+kr(idims,^d); !ixRR=[iMmin1+1,ixMmax+3]
141 
142  hxcmin^d=ixomin^d;hxcmax^d=ixmax^d; ! hxC = [ixMmin-1,ixMmax]
143  hxl^l=hxc^l-kr(idims,^d); ! hxL = [ixMmin-2,ixMmax-1]
144  hxr^l=hxc^l+kr(idims,^d); ! hxR = [ixMmin,ixMmax+1]
145 
146  kxcmin^d=ixllmin^d; kxcmax^d=ixrmax^d; ! kxC=[iMmin1-3,ixMmax1+2]
147  kxl^l=kxc^l-kr(idims,^d); ! kxL=[iMmin1-4,ixMmax1+1]
148  kxr^l=kxc^l+kr(idims,^d); ! kxR=[iMmin1-2,ixMmax1+3]
149 
150  lxcmin^d=ixllmin^d-kr(idims,^d);lxcmax^d=ixrrmax^d;! ixC=[iMmin1-4,ixMmax1+3]
151  lxr^l=lxc^l+kr(idims,^d); ! lxR=[iMmin1-3,ixMmax1+4]
152 
153  dwc(lxc^s,1:nwflux)=w(lxr^s,1:nwflux)-w(lxc^s,1:nwflux)
154  ! Eq. 64, Miller and Colella 2002, JCP 183, 26
155  d2wc(kxc^s,1:nwflux)=half*(w(kxr^s,1:nwflux)-w(kxl^s,1:nwflux))
156  where(dwc(kxc^s,1:nwflux)*dwc(kxl^s,1:nwflux)>zero)
157  ! Store the sign of dwC in wMin
158  wmin(kxc^s,1:nwflux)= sign(one,d2wc(kxc^s,1:nwflux))
159  ! Eq. 65, Miller and Colella 2002, JCP 183, 26
160  ldw(kxc^s,1:nwflux)= wmin(kxc^s,1:nwflux)*min(dabs(d2wc(kxc^s,1:nwflux)),&
161  2.0d0*dabs(dwc(kxl^s,1:nwflux)),&
162  2.0d0*dabs(dwc(kxc^s,1:nwflux)))
163  elsewhere
164  ldw(kxc^s,1:nwflux)=zero
165  endwhere
166 
167  ! Eq. 66, Miller and Colella 2002, JCP 183, 26
168  wlc(ixo^s,1:nwflux)=wlc(ixo^s,1:nwflux)+half*dwc(ixo^s,1:nwflux)&
169  +(ldw(ixo^s,1:nwflux)-ldw(ixr^s,1:nwflux))/6.0d0
170 
171  wrc(ixl^s,1:nwflux)=wrc(ixl^s,1:nwflux)-(half*dwc(ixl^s,1:nwflux)&
172  -(ldw(ixl^s,1:nwflux)-ldw(ixo^s,1:nwflux))/6.0d0)
173 
174  ! make sure that min wCT(i)<wLC(i)<wCT(i+1) same for wRC(i)
175  call extremaw(ixi^l,kxc^l,wct,1,wmax,wmin)
176 
177  ! Eq. B8, page 217, Mignone et al 2005, ApJS
178  wrc(ixl^s,1:nwflux)=max(wmin(ixo^s,1:nwflux)&
179  ,min(wmax(ixo^s,1:nwflux),wrc(ixl^s,1:nwflux)))
180  wlc(ixo^s,1:nwflux)=max(wmin(ixo^s,1:nwflux)&
181  ,min(wmax(ixo^s,1:nwflux),wlc(ixo^s,1:nwflux)))
182 
183 
184  ! Eq. B9, page 217, Mignone et al 2005, ApJS
185  where((wrc(ixl^s,1:nwflux)-wct(ixo^s,1:nwflux))&
186  *(wct(ixo^s,1:nwflux)-wlc(ixo^s,1:nwflux))<=zero)
187  wrc(ixl^s,1:nwflux)=wct(ixo^s,1:nwflux)
188  wlc(ixo^s,1:nwflux)=wct(ixo^s,1:nwflux)
189  end where
190 
191  wmax(ixo^s,1:nwflux)=(wlc(ixo^s,1:nwflux)-wrc(ixl^s,1:nwflux))*&
192  (wct(ixo^s,1:nwflux)-&
193  (wlc(ixo^s,1:nwflux)+wrc(ixl^s,1:nwflux))/2.0d0)
194  wmin(ixo^s,1:nwflux)=(wlc(ixo^s,1:nwflux)-wrc(ixl^s,1:nwflux))**2.0d0/6.0d0
195  tmp(ixl^s,1:nwflux)=wrc(ixl^s,1:nwflux)
196  ! Eq. B10, page 218, Mignone et al 2005, ApJS
197  where(wmax(hxr^s,1:nwflux)>wmin(hxr^s,1:nwflux))
198  wrc(hxc^s,1:nwflux)= 3.0d0*wct(hxr^s,1:nwflux)&
199  -2.0d0*wlc(hxr^s,1:nwflux)
200  endwhere
201  ! Eq. B11, page 218, Mignone et al 2005, ApJS
202  where(wmax(hxc^s,1:nwflux)<-wmin(hxc^s,1:nwflux))
203  wlc(hxc^s,1:nwflux)= 3.0d0*wct(hxc^s,1:nwflux)&
204  -2.0d0*tmp(hxl^s,1:nwflux)
205  endwhere
206 
207  ! flattening at the contact discontinuities
208  if(flatcd)then
209  call phys_ppm_flatcd(ixi^l,kxc^l,kxl^l,kxr^l,wct,d2wc,aa,ab)
210  if(any(kappa*aa(kxc^s)>=ab(kxc^s)))then
211  do iw=1,nwflux
212  where(kappa*aa(kxc^s)>=ab(kxc^s).and. dabs(dwc(kxc^s,iw))>smalldouble)
213  wmax(kxc^s,iw) = wct(kxr^s,iw)-2.0d0*wct(kxc^s,iw)+wct(kxl^s,iw)
214  end where
215 
216  where(wmax(ixr^s,iw)*wmax(ixl^s,iw)<zero .and.&
217  dabs(wct(ixr^s,iw)-wct(ixl^s,iw))&
218  -eps*min(dabs(wct(ixr^s,iw)),dabs(wct(ixl^s,iw)))>zero &
219  .and. kappa*aa(ixo^s)>=ab(ixo^s)&
220  .and. dabs(dwc(ixo^s,iw))>smalldouble)
221 
222  ac(ixo^s)=(wct(ixll^s,iw)-wct(ixrr^s,iw)+4.0d0*dwc(ixo^s,iw))&
223  /(12.0d0*dwc(ixo^s,iw))
224  wmin(ixo^s,iw)=max(zero,min(eta1*(ac(ixo^s)-eta2),one))
225  elsewhere
226  wmin(ixo^s,iw)=zero
227  end where
228 
229  where(wmin(hxc^s,iw)>zero)
230  wlc(hxc^s,iw) = wlc(hxc^s,iw)*(one-wmin(hxc^s,iw))&
231  +(wct(hxc^s,iw)+half*ldw(hxc^s,iw))*wmin(hxc^s,iw)
232  end where
233  where(wmin(hxr^s,iw)>zero)
234  wrc(hxc^s,iw) = wrc(hxc^s,iw)*(one-wmin(hxr^s,iw))&
235  +(wct(hxr^s,iw)-half*ldw(hxr^s,iw))*wmin(hxr^s,iw)
236  end where
237  end do
238  endif
239  endif
240 
241  ! flattening at the shocks
242  if(flatsh)then
243  ! following MILLER and COLELLA 2002 JCP 183, 26
244  kxcmin^d=ixmin^d-2; kxcmax^d=ixmax^d+2; ! kxC=[ixMmin1-2,ixMmax1+2]
245  do idimss=1,ndim
246  kxl^l=kxc^l-kr(idimss,^d); ! kxL=[ixMmin1-3,ixMmax1+1]
247  kxr^l=kxc^l+kr(idimss,^d); ! kxR=[ixMmin1-1,ixMmax1+3]
248  kxll^l=kxl^l-kr(idimss,^d);! kxLL=[ixMmin-4,ixMmax]
249  kxrr^l=kxr^l+kr(idimss,^d);! kxRR=[ixMmin,ixMmax+4]
250 
251  call phys_ppm_flatsh(ixi^l,kxc^l,kxll^l,kxl^l,kxr^l,kxrr^l,idimss,wct,aa,ab,dv)
252 
253  ! eq. B17, page 218, Mignone et al 2005, ApJS (had(Xi1))
254  ac(kxc^s) = max(zero,min(one,(betamax-aa(kxc^s))/(betamax-betamin)))
255  ! eq. B18, page 218, Mignone et al 2005, ApJS (had(Xi1))
256  ! recycling aa(ixL^S)
257  where (dabs(dv(kxc^s))<smalldouble)
258  aa(kxc^s) = max(ac(kxc^s), min(one,(zmax-ab(kxc^s))/(zmax-zmin)))
259  elsewhere
260  aa(kxc^s) = one
261  endwhere
262 
263  {^ifoned call extremaa(ixi^l,ixo^l,aa,1,ab)}
264  {^nooned call extremaa(ixi^l,ixo^l,aa,1,exi(ixg^t,idimss))}
265  enddo
266  {^nooned ab(ixo^s)=min(^d&exi(ixo^s,^d))}
267  ! recycling wMax
268  do iw=1,nwflux
269  where(dabs(ab(ixo^s)-one)>smalldouble)
270  wmax(ixo^s,iw) = (one-ab(ixo^s))*wct(ixo^s,iw)
271  endwhere
272 
273  where(dabs(ab(hxc^s)-one)>smalldouble)
274  wlc(hxc^s,iw) = ab(hxc^s)*wlc(hxc^s,iw)+wmax(hxc^s,iw)
275  endwhere
276 
277  where(dabs(ab(hxr^s)-one)>smalldouble)
278  wrc(hxc^s,iw) = ab(hxr^s)*wrc(hxc^s,iw)+wmax(hxr^s,iw)
279  endwhere
280  enddo
281  endif
282 
283  end subroutine ppmlimiter
284 
285  subroutine extremaq(ixI^L,ixO^L,q,nshift,qMax,qMin)
286 
288 
289  integer,intent(in) :: ixI^L,ixO^L
290  double precision, intent(in) :: q(ixI^S)
291  integer,intent(in) :: nshift
292 
293  double precision, intent(out) :: qMax(ixI^S),qMin(ixI^S)
294 
295  integer :: ixs^L,ixsR^L,ixsL^L,idims,jdims,kdims,ishift,i,j
296 
297  do ishift=1,nshift
298  idims=1
299  ixsr^l=ixo^l+ishift*kr(idims,^d);
300  ixsl^l=ixo^l-ishift*kr(idims,^d);
301  if (ishift==1) then
302  qmax(ixo^s)=max(q(ixo^s),q(ixsr^s),q(ixsl^s))
303  qmin(ixo^s)=min(q(ixo^s),q(ixsr^s),q(ixsl^s))
304  else
305  qmax(ixo^s)=max(qmax(ixo^s),q(ixsr^s),q(ixsl^s))
306  qmin(ixo^s)=min(qmin(ixo^s),q(ixsr^s),q(ixsl^s))
307  end if
308  {^nooned
309  idims=1
310  jdims=idims+1
311  do i=-1,1
312  ixs^l=ixo^l+i*ishift*kr(idims,^d);
313  ixsr^l=ixs^l+ishift*kr(jdims,^d);
314  ixsl^l=ixs^l-ishift*kr(jdims,^d);
315  qmax(ixo^s)=max(qmax(ixo^s),q(ixsr^s),q(ixsl^s))
316  qmin(ixo^s)=min(qmin(ixo^s),q(ixsr^s),q(ixsl^s))
317  end do
318  }
319  {^ifthreed
320  idims=1
321  jdims=idims+1
322  kdims=jdims+1
323  do i=-1,1
324  ixs^l=ixo^l+i*ishift*kr(idims,^d);
325  do j=-1,1
326  ixs^l=ixo^l+j*ishift*kr(jdims,^d);
327  ixsr^l=ixs^l+ishift*kr(kdims,^d);
328  ixsl^l=ixs^l-ishift*kr(kdims,^d);
329  qmax(ixo^s)=max(qmax(ixo^s),q(ixsr^s),q(ixsl^s))
330  qmin(ixo^s)=min(qmin(ixo^s),q(ixsr^s),q(ixsl^s))
331  end do
332  end do
333  }
334  enddo
335 
336  end subroutine extremaq
337 
338  subroutine extremaa(ixI^L,ixO^L,a,nshift,aMin)
340 
341  integer,intent(in) :: ixI^L,ixO^L
342  double precision, intent(in) :: a(ixI^S)
343  integer,intent(in) :: nshift
344 
345  double precision, intent(out) :: aMin(ixI^S)
346 
347  integer :: ixs^L,ixsR^L,ixsL^L,idims,jdims,kdims,ishift,i,j
348 
349  do ishift=1,nshift
350  idims=1
351  ixsr^l=ixo^l+ishift*kr(idims,^d);
352  ixsl^l=ixo^l-ishift*kr(idims,^d);
353  amin(ixo^s)=min(a(ixsr^s),a(ixo^s),a(ixsl^s))
354  {^nooned
355  idims=1
356  jdims=idims+1
357  do i=-1,1
358  ixs^l=ixo^l+i*ishift*kr(idims,^d);
359  ixsr^l=ixs^l+ishift*kr(jdims,^d);
360  ixsl^l=ixs^l-ishift*kr(jdims,^d);
361  amin(ixo^s)=min(amin(ixo^s),a(ixsr^s),a(ixsl^s))
362  end do
363  }
364  {^ifthreed
365  idims=1
366  jdims=idims+1
367  kdims=jdims+1
368  do i=-1,1
369  ixs^l=ixo^l+i*ishift*kr(idims,^d);
370  do j=-1,1
371  ixs^l=ixo^l+j*ishift*kr(jdims,^d);
372  ixsr^l=ixs^l+ishift*kr(kdims,^d);
373  ixsl^l=ixs^l-ishift*kr(kdims,^d);
374  amin(ixo^s)=min(amin(ixo^s),a(ixsr^s),a(ixsl^s))
375  end do
376  end do
377  }
378  end do
379 
380  end subroutine extremaa
381 
382  subroutine extremaw(ixI^L,ixO^L,w,nshift,wMax,wMin)
384 
385  integer,intent(in) :: ixI^L,ixO^L
386  double precision, intent(in) :: w(ixI^S,1:nw)
387  integer,intent(in) :: nshift
388 
389  double precision, intent(out) :: wMax(ixI^S,1:nwflux),wMin(ixI^S,1:nwflux)
390 
391  integer :: ixs^L,ixsR^L,ixsL^L,idims,jdims,kdims,ishift,i,j
392 
393  do ishift=1,nshift
394  idims=1
395  ixsr^l=ixo^l+ishift*kr(idims,^d);
396  ixsl^l=ixo^l-ishift*kr(idims,^d);
397  if (ishift==1) then
398  wmax(ixo^s,1:nwflux)= &
399  max(w(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
400  wmin(ixo^s,1:nwflux)= &
401  min(w(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
402  else
403  wmax(ixo^s,1:nwflux)= &
404  max(wmax(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
405  wmin(ixo^s,1:nwflux)= &
406  min(wmin(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
407  end if
408  {^nooned
409  idims=1
410  jdims=idims+1
411  do i=-1,1
412  ixs^l=ixo^l+i*ishift*kr(idims,^d);
413  ixsr^l=ixs^l+ishift*kr(jdims,^d);
414  ixsl^l=ixs^l-ishift*kr(jdims,^d);
415  wmax(ixo^s,1:nwflux)= &
416  max(wmax(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
417  wmin(ixo^s,1:nwflux)= &
418  min(wmin(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
419  end do
420  }
421  {^ifthreed
422  idims=1
423  jdims=idims+1
424  kdims=jdims+1
425  do i=-1,1
426  ixs^l=ixo^l+i*ishift*kr(idims,^d);
427  do j=-1,1
428  ixs^l=ixo^l+j*ishift*kr(jdims,^d);
429  ixsr^l=ixs^l+ishift*kr(kdims,^d);
430  ixsl^l=ixs^l-ishift*kr(kdims,^d);
431  wmax(ixo^s,1:nwflux)= &
432  max(wmax(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
433  wmin(ixo^s,1:nwflux)= &
434  min(wmin(ixo^s,1:nwflux),w(ixsr^s,1:nwflux),w(ixsl^s,1:nwflux))
435  end do
436  end do
437  }
438  enddo
439 
440  end subroutine extremaw
441 
442 end module mod_ppm
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer, dimension(:), allocatable, parameter d
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
Definition: mod_ppm.t:1
subroutine, public ppmlimitervar(ixIL, ixL, idims, q, qCT, qLC, qRC)
Definition: mod_ppm.t:12
subroutine extremaq(ixIL, ixOL, q, nshift, qMax, qMin)
Definition: mod_ppm.t:286
subroutine extremaw(ixIL, ixOL, w, nshift, wMax, wMin)
Definition: mod_ppm.t:383
subroutine, public ppmlimiter(ixIL, ixL, idims, w, wCT, wLC, wRC)
Definition: mod_ppm.t:104
subroutine extremaa(ixIL, ixOL, a, nshift, aMin)
Definition: mod_ppm.t:339