46 character(len=*),
intent(in) :: namelim
104 write(*,*)
'Unknown limiter: ', namelim
105 call mpistop(
"No such limiter")
110 integer,
intent(in) :: typelim
112 select case (typelim)
128 subroutine dwlimiter2(dwC,ixI^L,ixC^L,idims,typelim,ldw,rdw,a2max)
133 integer,
intent(in) :: ixI^L, ixC^L, idims
134 double precision,
intent(in) :: dwC(ixI^S)
135 integer,
intent(in) :: typelim
137 double precision,
intent(out),
optional :: ldw(ixI^S)
139 double precision,
intent(out),
optional :: rdw(ixI^S)
140 double precision,
intent(in),
optional :: a2max
142 double precision :: tmp(ixI^S), tmp2(ixI^S)
143 integer :: ixO^L, hxO^L
144 double precision,
parameter :: qsmall=1.
d-12, qsmall2=2.
d-12
145 double precision,
parameter :: eps = sqrt(epsilon(1.0d0))
148 double precision,
parameter :: c_mcbeta=1.4d0
150 double precision,
parameter :: cadalfa=0.5d0, cadbeta=2.0d0, cadgamma=1.6d0
152 double precision :: rdelinv
153 double precision :: ldwA(ixI^S),ldwB(ixI^S),tmpeta(ixI^S)
154 double precision,
parameter :: cadepsilon=1.
d-14, invcadepsilon=1.d14,
cada3_radius=0.1d0
159 ixomin^d=ixcmin^d+
kr(idims,^d); ixomax^d=ixcmax^d;
160 hxo^l=ixo^l-
kr(idims,^d);
172 select case (typelim)
175 tmp(ixo^s)=sign(one,dwc(ixo^s))
176 tmp(ixo^s)=tmp(ixo^s)* &
177 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)))
178 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
179 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
182 tmp(ixo^s)=sign(one,dwc(ixo^s))
183 tmp(ixo^s)=2*tmp(ixo^s)* &
184 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s),&
185 tmp(ixo^s)*quarter*(dwc(hxo^s)+dwc(ixo^s))))
186 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
187 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
190 tmp(ixo^s)=sign(one,dwc(ixo^s))
191 tmp(ixo^s)=tmp(ixo^s)* &
192 max(zero,min(c_mcbeta*abs(dwc(ixo^s)),c_mcbeta*tmp(ixo^s)*dwc(hxo^s),&
193 tmp(ixo^s)*half*(dwc(hxo^s)+dwc(ixo^s))))
194 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
195 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
198 tmp(ixo^s)=sign(one,dwc(ixo^s))
199 tmp(ixo^s)=tmp(ixo^s)* &
200 max(zero,min(2*abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)),&
201 min(abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s)))
202 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
203 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
206 tmp(ixo^s)=2*max(dwc(hxo^s)*dwc(ixo^s),zero) &
207 /(dwc(ixo^s)+dwc(hxo^s)+qsmall)
208 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
209 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
212 tmp(ixo^s)=(dwc(hxo^s)*(dwc(ixo^s)**2+qsmall)&
213 +dwc(ixo^s)*(dwc(hxo^s)**2+qsmall))&
214 /(dwc(ixo^s)**2+dwc(hxo^s)**2+qsmall2)
215 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
216 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
218 tmp(ixo^s)=sign(one,dwc(ixo^s))
219 tmp2(ixo^s)=min(2*abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s))
220 if (
present(ldw))
then
221 ldw(ixo^s)=tmp(ixo^s)* &
222 max(zero,min(tmp2(ixo^s),&
223 (dwc(hxo^s)*tmp(ixo^s)+2*abs(dwc(ixo^s)))*third))
225 if (
present(rdw))
then
226 rdw(ixo^s)=tmp(ixo^s)* &
227 max(zero,min(tmp2(ixo^s),&
228 (2*dwc(hxo^s)*tmp(ixo^s)+abs(dwc(ixo^s)))*third))
233 if (
present(ldw))
then
236 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
237 tmp2(ixo^s)=(2+tmp(ixo^s))*third
238 ldw(ixo^s)= max(zero,min(tmp2(ixo^s), &
239 max(-cadalfa*tmp(ixo^s), &
240 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
241 cadgamma)))) * dwc(ixo^s)
244 if (
present(rdw))
then
246 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
247 tmp2(ixo^s)=(2+tmp(ixo^s))*third
248 rdw(ixo^s)= max(zero,min(tmp2(ixo^s), &
249 max(-cadalfa*tmp(ixo^s), &
250 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
251 cadgamma)))) * dwc(hxo^s)
255 tmpeta(ixo^s)=(dwc(ixo^s)**2+dwc(hxo^s)**2)*rdelinv
256 if (
present(ldw))
then
257 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
258 ldwa(ixo^s)=(two+tmp(ixo^s))*third
259 where(tmpeta(ixo^s)<=one-cadepsilon)
260 ldw(ixo^s)=ldwa(ixo^s)
261 elsewhere(tmpeta(ixo^s)>=one+cadepsilon)
262 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
263 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
264 ldw(ixo^s)=ldwb(ixo^s)
266 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
267 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
268 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
269 ldw(ixo^s)=half*( (one-tmp2(ixo^s))*ldwa(ixo^s) &
270 +(one+tmp2(ixo^s))*ldwb(ixo^s))
272 ldw(ixo^s)=ldw(ixo^s) * dwc(ixo^s)
275 if (
present(rdw))
then
276 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
277 ldwa(ixo^s)=(two+tmp(ixo^s))*third
278 where(tmpeta(ixo^s)<=one-cadepsilon)
279 rdw(ixo^s)=ldwa(ixo^s)
280 elsewhere(tmpeta(ixo^s)>=one+cadepsilon)
281 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
282 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
283 rdw(ixo^s)=ldwb(ixo^s)
285 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
286 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
287 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
288 rdw(ixo^s)=half*( (one-tmp2(ixo^s))*ldwa(ixo^s) &
289 +(one+tmp2(ixo^s))*ldwb(ixo^s))
291 rdw(ixo^s)=rdw(ixo^s) * dwc(hxo^s)
294 tmpeta(ixo^s)=(sqrt(0.4d0*(dwc(ixo^s)**2+dwc(hxo^s)**2)))&
295 /((a2max+cadepsilon)*
dxlevel(idims)**2)
296 if(
present(ldw))
then
297 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s)+sign(eps,dwc(ixo^s)))
298 ldwa(ixo^s)=(two+tmp(ixo^s))*third
299 where(tmpeta(ixo^s)<=one-cadepsilon)
300 ldw(ixo^s)=ldwa(ixo^s)
301 else where(tmpeta(ixo^s)>=one+cadepsilon)
302 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
303 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
304 ldw(ixo^s)=ldwb(ixo^s)
306 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
307 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
308 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
309 ldw(ixo^s)=half*((one-tmp2(ixo^s))*ldwa(ixo^s)&
310 +(one+tmp2(ixo^s))*ldwb(ixo^s))
312 ldw(ixo^s)=ldw(ixo^s)*dwc(ixo^s)
314 if(
present(rdw))
then
315 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s)+sign(eps,dwc(hxo^s)))
316 ldwa(ixo^s)=(two+tmp(ixo^s))*third
317 where(tmpeta(ixo^s)<=one-cadepsilon)
318 rdw(ixo^s)=ldwa(ixo^s)
319 else where(tmpeta(ixo^s)>=one+cadepsilon)
320 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
321 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
322 rdw(ixo^s)=ldwb(ixo^s)
324 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
325 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), 1.5d0))))
326 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
327 rdw(ixo^s)=half*((one-tmp2(ixo^s))*ldwa(ixo^s)&
328 +(one+tmp2(ixo^s))*ldwb(ixo^s))
330 rdw(ixo^s)=rdw(ixo^s)*dwc(hxo^s)
333 call mpistop(
"Error in dwLimiter: unknown limiter")
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision, dimension(ndim) dxlevel
Module with slope/flux limiters.
integer, parameter limiter_mcbeta
integer, parameter limiter_weno5cu6
integer, parameter limiter_mpweno7
integer, parameter limiter_schmid
pure logical function limiter_symmetric(typelim)
integer, parameter limiter_teno5ad
integer, parameter limiter_weno3
integer, parameter limiter_vanleer
integer, parameter limiter_ppm
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...
integer function limiter_type(namelim)
integer, parameter limiter_cada3
integer, parameter limiter_wenozp5
integer, parameter limiter_woodward
integer, parameter limiter_superbee
integer, parameter limiter_weno5
double precision cada3_radius
radius of the asymptotic region [0.001, 10], larger means more accurate in smooth region but more ove...
integer, parameter limiter_wenoz5
integer, parameter limiter_wenoz5nm
integer, parameter limiter_koren
integer, parameter limiter_weno7
integer, parameter limiter_wenozp5nm
integer, parameter limiter_cada
integer, parameter limiter_mp5
double precision schmid_rad
integer, parameter limiter_venk
integer, parameter limiter_weno5nm
integer, parameter limiter_wenoyc3
integer, parameter limiter_albada
integer, parameter limiter_minmod
Module containing the MP5 (fifth order) flux scheme.