MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_limiter.t
Go to the documentation of this file.
1!> Module with slope/flux limiters
3 use mod_ppm
4 use mod_mp5
5 use mod_weno
6 use mod_venk
7
8 implicit none
9 public
10
11 !> radius of the asymptotic region [0.001, 10], larger means more accurate in smooth
12 !> region but more overshooting at discontinuities
13 double precision :: cada3_radius
14 integer, parameter :: limiter_minmod = 1
15 integer, parameter :: limiter_woodward = 2
16 integer, parameter :: limiter_mcbeta = 3
17 integer, parameter :: limiter_superbee = 4
18 integer, parameter :: limiter_vanleer = 5
19 integer, parameter :: limiter_albada = 6
20 integer, parameter :: limiter_koren = 7
21 integer, parameter :: limiter_cada = 8
22 integer, parameter :: limiter_cada3 = 9
23 integer, parameter :: limiter_venk = 11
24 ! Special cases
25 integer, parameter :: limiter_ppm = 12
26 integer, parameter :: limiter_mp5 = 13
27 integer, parameter :: limiter_weno3 = 14
28 integer, parameter :: limiter_wenoyc3 = 15
29 integer, parameter :: limiter_weno5 = 16
30 integer, parameter :: limiter_weno5nm = 17
31 integer, parameter :: limiter_wenoz5 = 18
32 integer, parameter :: limiter_wenoz5nm = 19
33 integer, parameter :: limiter_wenozp5 = 20
34 integer, parameter :: limiter_wenozp5nm = 21
35 integer, parameter :: limiter_weno5cu6 = 22
36 integer, parameter :: limiter_teno5ad = 23
37 integer, parameter :: limiter_weno7 = 24
38 integer, parameter :: limiter_mpweno7 = 25
39
40contains
41
42 integer function limiter_type(namelim)
43 use mod_comm_lib, only: mpistop
44 character(len=*), intent(in) :: namelim
45
46 select case (namelim)
47 case ('minmod')
49 case ('woodward')
51 case ('mcbeta')
53 case ('superbee')
55 case ('vanleer')
57 case ('albada')
59 case ('koren')
61 case ('cada')
63 case ('cada3')
65 case('venk')
67 case ('ppm')
69 case ('mp5')
71 case ('weno3')
73 case ('wenoyc3')
75 case ('weno5')
77 case ('weno5nm')
79 case ('wenoz5')
81 case ('wenoz5nm')
83 case ('wenozp5')
85 case ('wenozp5nm')
87 case ('weno5cu6')
89 case ('teno5ad')
91 case ('weno7')
93 case ('mpweno7')
95
96 case default
97 limiter_type = -1
98 write(*,*) 'Unknown limiter: ', namelim
99 call mpistop("No such limiter")
100 end select
101 end function limiter_type
102
103 pure logical function limiter_symmetric(typelim)
104 integer, intent(in) :: typelim
105
106 select case (typelim)
108 limiter_symmetric = .false.
109 case default
110 limiter_symmetric = .true.
111 end select
112 end function limiter_symmetric
113
114 !> Limit the centered dwC differences within ixC for iw in direction idim.
115 !> The limiter is chosen according to typelim.
116 !>
117 !> Note that this subroutine is called from upwindLR (hence from methods
118 !> like tvdlf, hancock, hll(c) etc) or directly from tvd.t,
119 !> but also from the gradientL and divvectorS subroutines in geometry.t
120 !> Accordingly, the typelim here corresponds to one of limiter
121 !> or one of gradient_limiter.
122 subroutine dwlimiter2(dwC,ixI^L,ixC^L,idims,typelim,ldw,rdw)
123
125 use mod_comm_lib, only: mpistop
126
127 integer, intent(in) :: ixI^L, ixC^L, idims
128 double precision, intent(in) :: dwC(ixI^S)
129 integer, intent(in) :: typelim
130 !> Result using left-limiter (same as right for symmetric)
131 double precision, intent(out), optional :: ldw(ixI^S)
132 !> Result using right-limiter (same as left for symmetric)
133 double precision, intent(out), optional :: rdw(ixI^S)
134
135 double precision :: tmp(ixI^S), tmp2(ixI^S)
136 double precision, parameter :: qsmall=1.d-12, qsmall2=2.d-12
137 double precision, parameter :: eps = sqrt(epsilon(1.0d0))
138
139 ! mcbeta limiter parameter value
140 double precision, parameter :: c_mcbeta=1.4d0
141 ! cada limiter parameter values
142 double precision, parameter :: cadalfa=0.5d0, cadbeta=2.0d0, cadgamma=1.6d0
143 ! full third order cada limiter
144 double precision :: rdelinv
145 double precision :: ldwA(ixI^S),ldwB(ixI^S),tmpeta(ixI^S)
146 double precision, parameter :: cadepsilon=1.d-14, invcadepsilon=1.d14
147 integer :: ixO^L, hxO^L, ix^D, hx^D
148 !-----------------------------------------------------------------------------
149
150 ! Contract indices in idim for output.
151 ixomin^d=ixcmin^d+kr(idims,^d); ixomax^d=ixcmax^d;
152 hxo^l=ixo^l-kr(idims,^d);
153 hx^d=kr(idims,^d);
154
155 ! About the notation: the conventional argument theta (the ratio of slopes)
156 ! would be given by dwC(ixO^S)/dwC(hxO^S). However, in the end one
157 ! multiplies phi(theta) by dwC(hxO^S), which is incorporated in the
158 ! equations below. The minmod limiter can for example be written as:
159 ! A:
160 ! max(0.0d0, min(1.0d0, dwC(ixO^S)/dwC(hxO^S))) * dwC(hxO^S)
161 ! B:
162 ! tmp(ixO^S)*max(0.0d0,min(abs(dwC(ixO^S)),tmp(ixO^S)*dwC(hxO^S)))
163 ! where tmp(ixO^S)=sign(1.0d0,dwC(ixO^S))
164
165 select case (typelim)
166 case (limiter_minmod)
167 ! Minmod limiter eq(3.51e) and (eq.3.38e) with omega=1
168 tmp(ixo^s)=sign(one,dwc(ixo^s))
169 tmp(ixo^s)=tmp(ixo^s)* &
170 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)))
171 if (present(ldw)) ldw(ixo^s) = tmp(ixo^s)
172 if (present(rdw)) rdw(ixo^s) = tmp(ixo^s)
173 case (limiter_woodward)
174 ! Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
175 tmp(ixo^s)=sign(one,dwc(ixo^s))
176 tmp(ixo^s)=2*tmp(ixo^s)* &
177 max(zero,min(abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s),&
178 tmp(ixo^s)*quarter*(dwc(hxo^s)+dwc(ixo^s))))
179 if (present(ldw)) ldw(ixo^s) = tmp(ixo^s)
180 if (present(rdw)) rdw(ixo^s) = tmp(ixo^s)
181 case (limiter_mcbeta)
182 ! Woodward and Collela limiter, with factor beta
183 tmp(ixo^s)=sign(one,dwc(ixo^s))
184 tmp(ixo^s)=tmp(ixo^s)* &
185 max(zero,min(c_mcbeta*abs(dwc(ixo^s)),c_mcbeta*tmp(ixo^s)*dwc(hxo^s),&
186 tmp(ixo^s)*half*(dwc(hxo^s)+dwc(ixo^s))))
187 if (present(ldw)) ldw(ixo^s) = tmp(ixo^s)
188 if (present(rdw)) rdw(ixo^s) = tmp(ixo^s)
189 case (limiter_superbee)
190 ! Roes superbee limiter (eq.3.51i)
191 tmp(ixo^s)=sign(one,dwc(ixo^s))
192 tmp(ixo^s)=tmp(ixo^s)* &
193 max(zero,min(2*abs(dwc(ixo^s)),tmp(ixo^s)*dwc(hxo^s)),&
194 min(abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s)))
195 if (present(ldw)) ldw(ixo^s) = tmp(ixo^s)
196 if (present(rdw)) rdw(ixo^s) = tmp(ixo^s)
197 case (limiter_vanleer)
198 ! van Leer limiter (eq 3.51f), but a missing delta2=1.D-12 is added
199 tmp(ixo^s)=2*max(dwc(hxo^s)*dwc(ixo^s),zero) &
200 /(dwc(ixo^s)+dwc(hxo^s)+qsmall)
201 if (present(ldw)) ldw(ixo^s) = tmp(ixo^s)
202 if (present(rdw)) rdw(ixo^s) = tmp(ixo^s)
203 case (limiter_albada)
204 ! Albada limiter (eq.3.51g) with delta2=1D.-12
205 tmp(ixo^s)=(dwc(hxo^s)*(dwc(ixo^s)**2+qsmall)&
206 +dwc(ixo^s)*(dwc(hxo^s)**2+qsmall))&
207 /(dwc(ixo^s)**2+dwc(hxo^s)**2+qsmall2)
208 if (present(ldw)) ldw(ixo^s) = tmp(ixo^s)
209 if (present(rdw)) rdw(ixo^s) = tmp(ixo^s)
210 case (limiter_koren)
211 tmp(ixo^s)=sign(one,dwc(ixo^s))
212 tmp2(ixo^s)=min(2*abs(dwc(ixo^s)),2*tmp(ixo^s)*dwc(hxo^s))
213 if (present(ldw)) then
214 ldw(ixo^s)=tmp(ixo^s)* &
215 max(zero,min(tmp2(ixo^s),&
216 (dwc(hxo^s)*tmp(ixo^s)+2*abs(dwc(ixo^s)))*third))
217 end if
218 if (present(rdw)) then
219 rdw(ixo^s)=tmp(ixo^s)* &
220 max(zero,min(tmp2(ixo^s),&
221 (2*dwc(hxo^s)*tmp(ixo^s)+abs(dwc(ixo^s)))*third))
222 end if
223 case (limiter_cada)
224 ! This limiter has been rewritten in the usual form, and uses a division
225 ! of the gradients.
226 if (present(ldw)) then
227 ! Cada Left variant
228 ! Compute theta, but avoid division by zero
229 tmp(ixo^s)=dwc(hxo^s)/(dwc(ixo^s) + sign(eps, dwc(ixo^s)))
230 tmp2(ixo^s)=(2+tmp(ixo^s))*third
231 ldw(ixo^s)= max(zero,min(tmp2(ixo^s), &
232 max(-cadalfa*tmp(ixo^s), &
233 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
234 cadgamma)))) * dwc(ixo^s)
235 end if
236
237 if (present(rdw)) then
238 ! Cada Right variant
239 tmp(ixo^s)=dwc(ixo^s)/(dwc(hxo^s) + sign(eps, dwc(hxo^s)))
240 tmp2(ixo^s)=(2+tmp(ixo^s))*third
241 rdw(ixo^s)= max(zero,min(tmp2(ixo^s), &
242 max(-cadalfa*tmp(ixo^s), &
243 min(cadbeta*tmp(ixo^s), tmp2(ixo^s), &
244 cadgamma)))) * dwc(hxo^s)
245 end if
246 case (limiter_cada3)
247 rdelinv=one/(cada3_radius*dxlevel(idims))**2
248 {!DEC$ VECTOR ALWAYS
249 do ix^db=ixomin^db,ixomax^db\}
250 tmpeta(ix^d)=(dwc(ix^d)**2+dwc(ix^d-hx^d)**2)*rdelinv
251 if (present(ldw)) then
252 tmp(ix^d)=dwc(ix^d-hx^d)/(dwc(ix^d)+sign(eps,dwc(ix^d)))
253 ldwa(ix^d)=(two+tmp(ix^d))*third
254 if(tmpeta(ix^d)<=one-cadepsilon) then
255 ldw(ix^d)=ldwa(ix^d)
256 else if(tmpeta(ix^d)>=one+cadepsilon) then
257 ldw(ix^d)=max(zero,min(ldwa(ix^d),max(-cadalfa*tmp(ix^d),&
258 min(cadbeta*tmp(ix^d),ldwa(ix^d),cadgamma))))
259 else
260 tmp2(ix^d)=(tmpeta(ix^d)-one)*invcadepsilon
261 ldw(ix^d)=half*((one-tmp2(ix^d))*ldwa(ix^d)+(one+tmp2(ix^d))*&
262 max(zero,min(ldwa(ix^d),max(-cadalfa*tmp(ix^d),&
263 min(cadbeta*tmp(ix^d),ldwa(ix^d),cadgamma)))))
264 end if
265 ldw(ix^d)=ldw(ix^d)*dwc(ix^d)
266 end if
267 if (present(rdw)) then
268 tmp(ix^d)=dwc(ix^d)/(dwc(ix^d-hx^d)+sign(eps,dwc(ix^d-hx^d)))
269 ldwa(ix^d)=(two+tmp(ix^d))*third
270 if(tmpeta(ix^d)<=one-cadepsilon) then
271 rdw(ix^d)=ldwa(ix^d)
272 else if(tmpeta(ix^d)>=one+cadepsilon) then
273 rdw(ix^d)=max(zero,min(ldwa(ix^d),max(-cadalfa*tmp(ix^d),&
274 min(cadbeta*tmp(ix^d),ldwa(ix^d),cadgamma))))
275 else
276 tmp2(ix^d)=(tmpeta(ix^d)-one)*invcadepsilon
277 rdw(ix^d)=half*((one-tmp2(ix^d))*ldwa(ix^d)+(one+tmp2(ix^d))*&
278 max(zero,min(ldwa(ix^d),max(-cadalfa*tmp(ix^d),&
279 min(cadbeta*tmp(ix^d),ldwa(ix^d),cadgamma)))))
280 end if
281 rdw(ix^d)=rdw(ix^d)*dwc(ix^d-hx^d)
282 end if
283 {end do\}
284 case default
285 call mpistop("Error in dwLimiter: unknown limiter")
286 end select
287
288 end subroutine dwlimiter2
289
290end module mod_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(:), allocatable, parameter d
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
Module with slope/flux limiters.
Definition mod_limiter.t:2
integer, parameter limiter_mcbeta
Definition mod_limiter.t:16
integer, parameter limiter_weno5cu6
Definition mod_limiter.t:35
integer, parameter limiter_mpweno7
Definition mod_limiter.t:38
pure logical function limiter_symmetric(typelim)
integer, parameter limiter_teno5ad
Definition mod_limiter.t:36
integer, parameter limiter_weno3
Definition mod_limiter.t:27
integer, parameter limiter_vanleer
Definition mod_limiter.t:18
integer, parameter limiter_ppm
Definition mod_limiter.t:25
integer function limiter_type(namelim)
Definition mod_limiter.t:43
integer, parameter limiter_cada3
Definition mod_limiter.t:22
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_woodward
Definition mod_limiter.t:15
integer, parameter limiter_superbee
Definition mod_limiter.t:17
integer, parameter limiter_weno5
Definition mod_limiter.t:29
double precision cada3_radius
radius of the asymptotic region [0.001, 10], larger means more accurate in smooth region but more ove...
Definition mod_limiter.t:13
integer, parameter limiter_wenoz5
Definition mod_limiter.t:31
integer, parameter limiter_wenoz5nm
Definition mod_limiter.t:32
integer, parameter limiter_koren
Definition mod_limiter.t:20
integer, parameter limiter_weno7
Definition mod_limiter.t:37
integer, parameter limiter_wenozp5nm
Definition mod_limiter.t:34
integer, parameter limiter_cada
Definition mod_limiter.t:21
integer, parameter limiter_mp5
Definition mod_limiter.t:26
integer, parameter limiter_venk
Definition mod_limiter.t:23
integer, parameter limiter_weno5nm
Definition mod_limiter.t:30
integer, parameter limiter_wenoyc3
Definition mod_limiter.t:28
integer, parameter limiter_albada
Definition mod_limiter.t:19
integer, parameter limiter_minmod
Definition mod_limiter.t:14
Module containing the MP5 (fifth order) flux scheme.
Definition mod_mp5.t:2