45 character(len=*),
intent(in) :: namelim
103 write(*,*)
'Unknown limiter: ', namelim
104 call mpistop(
"No such limiter")
109 integer,
intent(in) :: typelim
111 select case (typelim)
127 subroutine dwlimiter2(dwC,ixI^L,ixC^L,idims,typelim,ldw,rdw,a2max)
131 integer,
intent(in) :: ixI^L, ixC^L, idims
132 double precision,
intent(in) :: dwC(ixI^S)
133 integer,
intent(in) :: typelim
135 double precision,
intent(out),
optional :: ldw(ixI^S)
137 double precision,
intent(out),
optional :: rdw(ixI^S)
138 double precision,
intent(in),
optional :: a2max
140 double precision :: tmp(ixI^S), tmp2(ixI^S)
141 integer :: ixO^L, hxO^L
142 double precision,
parameter :: qsmall=1.
d-12, qsmall2=2.
d-12
143 double precision,
parameter :: eps = sqrt(epsilon(1.0d0))
146 double precision,
parameter :: c_mcbeta=1.4d0
148 double precision,
parameter :: cadalfa=0.5d0, cadbeta=2.0d0, cadgamma=1.6d0
150 double precision :: rdelinv
151 double precision :: ldwA(ixI^S),ldwB(ixI^S),tmpeta(ixI^S)
152 double precision,
parameter :: cadepsilon=1.
d-14, invcadepsilon=1.d14,
cada3_radius=0.1d0
157 ixomin^d=ixcmin^d+
kr(idims,^d); ixomax^d=ixcmax^d;
158 hxo^l=ixo^l-
kr(idims,^d);
170 select case (typelim)
173 tmp(ixo^s)=sign(one,dwc(ixo^s))
174 tmp(ixo^s)=tmp(ixo^s)* &
175 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)))
176 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
177 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
180 tmp(ixo^s)=sign(one,dwc(ixo^s))
181 tmp(ixo^s)=2*tmp(ixo^s)* &
182 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s),&
183 tmp(ixo^s)*quarter*(dwc(hxo^s)+dwc(ixo^s))))
184 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
185 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
188 tmp(ixo^s)=sign(one,dwc(ixo^s))
189 tmp(ixo^s)=tmp(ixo^s)* &
190 max(zero,min(c_mcbeta*abs(dwc(ixo^s)),c_mcbeta*tmp(ixo^s)*dwc(hxo^s),&
191 tmp(ixo^s)*half*(dwc(hxo^s)+dwc(ixo^s))))
192 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
193 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
196 tmp(ixo^s)=sign(one,dwc(ixo^s))
197 tmp(ixo^s)=tmp(ixo^s)* &
198 max(zero,min(2*abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)),&
199 min(abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s)))
200 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
201 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
204 tmp(ixo^s)=2*max(dwc(hxo^s)*dwc(ixo^s),zero) &
205 /(dwc(ixo^s)+dwc(hxo^s)+qsmall)
206 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
207 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
210 tmp(ixo^s)=(dwc(hxo^s)*(dwc(ixo^s)**2+qsmall)&
211 +dwc(ixo^s)*(dwc(hxo^s)**2+qsmall))&
212 /(dwc(ixo^s)**2+dwc(hxo^s)**2+qsmall2)
213 if (
present(ldw)) ldw(ixo^s) = tmp(ixo^s)
214 if (
present(rdw)) rdw(ixo^s) = tmp(ixo^s)
216 tmp(ixo^s)=sign(one,dwc(ixo^s))
217 tmp2(ixo^s)=min(2*abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s))
218 if (
present(ldw))
then
219 ldw(ixo^s)=tmp(ixo^s)* &
220 max(zero,min(tmp2(ixo^s),&
221 (dwc(hxo^s)*tmp(ixo^s)+2*abs(dwc(ixo^s)))*third))
223 if (
present(rdw))
then
224 rdw(ixo^s)=tmp(ixo^s)* &
225 max(zero,min(tmp2(ixo^s),&
226 (2*dwc(hxo^s)*tmp(ixo^s)+abs(dwc(ixo^s)))*third))
231 if (
present(ldw))
then
234 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
235 tmp2(ixo^s)=(2+tmp(ixo^s))*third
236 ldw(ixo^s)= max(zero,min(tmp2(ixo^s), &
237 max(-cadalfa*tmp(ixo^s), &
238 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
239 cadgamma)))) * dwc(ixo^s)
242 if (
present(rdw))
then
244 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
245 tmp2(ixo^s)=(2+tmp(ixo^s))*third
246 rdw(ixo^s)= max(zero,min(tmp2(ixo^s), &
247 max(-cadalfa*tmp(ixo^s), &
248 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
249 cadgamma)))) * dwc(hxo^s)
253 tmpeta(ixo^s)=(dwc(ixo^s)**2+dwc(hxo^s)**2)*rdelinv
254 if (
present(ldw))
then
255 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
256 ldwa(ixo^s)=(two+tmp(ixo^s))*third
257 where(tmpeta(ixo^s)<=one-cadepsilon)
258 ldw(ixo^s)=ldwa(ixo^s)
259 elsewhere(tmpeta(ixo^s)>=one+cadepsilon)
260 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
261 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
262 ldw(ixo^s)=ldwb(ixo^s)
264 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
265 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
266 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
267 ldw(ixo^s)=half*( (one-tmp2(ixo^s))*ldwa(ixo^s) &
268 +(one+tmp2(ixo^s))*ldwb(ixo^s))
270 ldw(ixo^s)=ldw(ixo^s) * dwc(ixo^s)
273 if (
present(rdw))
then
274 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
275 ldwa(ixo^s)=(two+tmp(ixo^s))*third
276 where(tmpeta(ixo^s)<=one-cadepsilon)
277 rdw(ixo^s)=ldwa(ixo^s)
278 elsewhere(tmpeta(ixo^s)>=one+cadepsilon)
279 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
280 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
281 rdw(ixo^s)=ldwb(ixo^s)
283 ldwb(ixo^s)= max(zero,min(ldwa(ixo^s), max(-cadalfa*tmp(ixo^s), &
284 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), cadgamma))))
285 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
286 rdw(ixo^s)=half*( (one-tmp2(ixo^s))*ldwa(ixo^s) &
287 +(one+tmp2(ixo^s))*ldwb(ixo^s))
289 rdw(ixo^s)=rdw(ixo^s) * dwc(hxo^s)
292 tmpeta(ixo^s)=(sqrt(0.4d0*(dwc(ixo^s)**2+dwc(hxo^s)**2)))&
293 /((a2max+cadepsilon)*
dxlevel(idims)**2)
294 if(
present(ldw))
then
295 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s)+sign(eps,dwc(ixo^s)))
296 ldwa(ixo^s)=(two+tmp(ixo^s))*third
297 where(tmpeta(ixo^s)<=one-cadepsilon)
298 ldw(ixo^s)=ldwa(ixo^s)
299 else where(tmpeta(ixo^s)>=one+cadepsilon)
300 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
301 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
302 ldw(ixo^s)=ldwb(ixo^s)
304 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
305 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
306 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
307 ldw(ixo^s)=half*((one-tmp2(ixo^s))*ldwa(ixo^s)&
308 +(one+tmp2(ixo^s))*ldwb(ixo^s))
310 ldw(ixo^s)=ldw(ixo^s)*dwc(ixo^s)
312 if(
present(rdw))
then
313 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s)+sign(eps,dwc(hxo^s)))
314 ldwa(ixo^s)=(two+tmp(ixo^s))*third
315 where(tmpeta(ixo^s)<=one-cadepsilon)
316 rdw(ixo^s)=ldwa(ixo^s)
317 else where(tmpeta(ixo^s)>=one+cadepsilon)
318 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
319 min(cadbeta*tmp(ixo^s),ldwa(ixo^s),1.5d0))))
320 rdw(ixo^s)=ldwb(ixo^s)
322 ldwb(ixo^s)=max(zero,min(ldwa(ixo^s),max(-tmp(ixo^s),&
323 min(cadbeta*tmp(ixo^s), ldwa(ixo^s), 1.5d0))))
324 tmp2(ixo^s)=(tmpeta(ixo^s)-one)*invcadepsilon
325 rdw(ixo^s)=half*((one-tmp2(ixo^s))*ldwa(ixo^s)&
326 +(one+tmp2(ixo^s))*ldwb(ixo^s))
328 rdw(ixo^s)=rdw(ixo^s)*dwc(hxo^s)
331 call mpistop(
"Error in dwLimiter: unknown limiter")
subroutine 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.
integer, dimension(:), allocatable, parameter d
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.