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 double precision,
parameter :: qsmall=1.
d-12, qsmall2=2.
d-12
144 double precision,
parameter :: eps = sqrt(epsilon(1.0d0))
147 double precision,
parameter :: c_mcbeta=1.4d0
149 double precision,
parameter :: cadalfa=0.5d0, cadbeta=2.0d0, cadgamma=1.6d0
151 double precision :: rdelinv
152 double precision :: ldwA(ixI^S),ldwB(ixI^S),tmpeta(ixI^S)
153 double precision,
parameter :: cadepsilon=1.
d-14, invcadepsilon=1.d14,
cada3_radius=0.1d0
154 integer :: ixO^L, hxO^L, ix^D
158 ixomin^d=ixcmin^d+
kr(idims,^d); ixomax^d=ixcmax^d;
159 hxo^l=ixo^l-
kr(idims,^d);
171 select case (typelim)
174 tmp(ixo^s)=sign(one,dwc(ixo^s))
175 tmp(ixo^s)=tmp(ixo^s)* &
176 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)))
177 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
178 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
181 tmp(ixo^s)=sign(one,dwc(ixo^s))
182 tmp(ixo^s)=2*tmp(ixo^s)* &
183 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s),&
184 tmp(ixo^s)*quarter*(dwc(hxo^s)+dwc(ixo^s))))
185 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
186 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
189 tmp(ixo^s)=sign(one,dwc(ixo^s))
190 tmp(ixo^s)=tmp(ixo^s)* &
191 max(zero,min(c_mcbeta*abs(dwc(ixo^s)),c_mcbeta*tmp(ixo^s)*dwc(hxo^s),&
192 tmp(ixo^s)*half*(dwc(hxo^s)+dwc(ixo^s))))
193 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
194 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
197 tmp(ixo^s)=sign(one,dwc(ixo^s))
198 tmp(ixo^s)=tmp(ixo^s)* &
199 max(zero,min(2*abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)),&
200 min(abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s)))
201 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
202 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
205 tmp(ixo^s)=2*max(dwc(hxo^s)*dwc(ixo^s),zero) &
206 /(dwc(ixo^s)+dwc(hxo^s)+qsmall)
207 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
208 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
211 tmp(ixo^s)=(dwc(hxo^s)*(dwc(ixo^s)**2+qsmall)&
212 +dwc(ixo^s)*(dwc(hxo^s)**2+qsmall))&
213 /(dwc(ixo^s)**2+dwc(hxo^s)**2+qsmall2)
214 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
215 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
217 tmp(ixo^s)=sign(one,dwc(ixo^s))
218 tmp2(ixo^s)=min(2*abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s))
219 if (
present(ldw))
then
220 ldw(ixo^s)=tmp(ixo^s)* &
221 max(zero,min(tmp2(ixo^s),&
222 (dwc(hxo^s)*tmp(ixo^s)+2*abs(dwc(ixo^s)))*third))
224 if (
present(rdw))
then
225 rdw(ixo^s)=tmp(ixo^s)* &
226 max(zero,min(tmp2(ixo^s),&
227 (2*dwc(hxo^s)*tmp(ixo^s)+abs(dwc(ixo^s)))*third))
232 if (
present(ldw))
then
235 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
236 tmp2(ixo^s)=(2+tmp(ixo^s))*third
237 ldw(ixo^s)= max(zero,min(tmp2(ixo^s), &
238 max(-cadalfa*tmp(ixo^s), &
239 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
240 cadgamma)))) * dwc(ixo^s)
243 if (
present(rdw))
then
245 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
246 tmp2(ixo^s)=(2+tmp(ixo^s))*third
247 rdw(ixo^s)= max(zero,min(tmp2(ixo^s), &
248 max(-cadalfa*tmp(ixo^s), &
249 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
250 cadgamma)))) * dwc(hxo^s)
254 tmpeta(ixo^s)=(dwc(ixo^s)**2+dwc(hxo^s)**2)*rdelinv
255 if (
present(ldw))
then
256 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
257 ldwa(ixo^s)=(two+tmp(ixo^s))*third
258 where(tmpeta(ixo^s)<=one-cadepsilon)
259 ldw(ixo^s)=ldwa(ixo^s)
260 elsewhere(tmpeta(ixo^s)>=one+cadepsilon)
261 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
262 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
263 ldw(ixo^s)=ldwb(ixo^s)
265 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
266 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
267 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
268 ldw(ixo^s)=half*( (one-tmp2(ixo^s))*ldwa(ixo^s) &
269 +(one+tmp2(ixo^s))*ldwb(ixo^s))
271 ldw(ixo^s)=ldw(ixo^s) * dwc(ixo^s)
274 if (
present(rdw))
then
275 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
276 ldwa(ixo^s)=(two+tmp(ixo^s))*third
277 where(tmpeta(ixo^s)<=one-cadepsilon)
278 rdw(ixo^s)=ldwa(ixo^s)
279 elsewhere(tmpeta(ixo^s)>=one+cadepsilon)
280 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
281 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
282 rdw(ixo^s)=ldwb(ixo^s)
284 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
285 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
286 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
287 rdw(ixo^s)=half*( (one-tmp2(ixo^s))*ldwa(ixo^s) &
288 +(one+tmp2(ixo^s))*ldwb(ixo^s))
290 rdw(ixo^s)=rdw(ixo^s) * dwc(hxo^s)
293 tmpeta(ixo^s)=(sqrt(0.4d0*(dwc(ixo^s)**2+dwc(hxo^s)**2)))&
294 /((a2max+cadepsilon)*
dxlevel(idims)**2)
295 if(
present(ldw))
then
296 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s)+sign(eps,dwc(ixo^s)))
297 ldwa(ixo^s)=(two+tmp(ixo^s))*third
298 where(tmpeta(ixo^s)<=one-cadepsilon)
299 ldw(ixo^s)=ldwa(ixo^s)
300 else where(tmpeta(ixo^s)>=one+cadepsilon)
301 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
302 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
303 ldw(ixo^s)=ldwb(ixo^s)
305 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
306 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
307 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
308 ldw(ixo^s)=half*((one-tmp2(ixo^s))*ldwa(ixo^s)&
309 +(one+tmp2(ixo^s))*ldwb(ixo^s))
311 ldw(ixo^s)=ldw(ixo^s)*dwc(ixo^s)
313 if(
present(rdw))
then
314 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s)+sign(eps,dwc(hxo^s)))
315 ldwa(ixo^s)=(two+tmp(ixo^s))*third
316 where(tmpeta(ixo^s)<=one-cadepsilon)
317 rdw(ixo^s)=ldwa(ixo^s)
318 else where(tmpeta(ixo^s)>=one+cadepsilon)
319 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
320 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
321 rdw(ixo^s)=ldwb(ixo^s)
323 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
324 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), 1.5d0))))
325 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
326 rdw(ixo^s)=half*((one-tmp2(ixo^s))*ldwa(ixo^s)&
327 +(one+tmp2(ixo^s))*ldwb(ixo^s))
329 rdw(ixo^s)=rdw(ixo^s)*dwc(hxo^s)
332 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(^nd) dxlevel
store unstretched cell size of current level
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.