MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_limiter.t
Go to the documentation of this file.
1 !> Module with slope/flux limiters
2 module mod_limiter
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  double precision :: schmid_rad^d
15  integer, parameter :: limiter_minmod = 1
16  integer, parameter :: limiter_woodward = 2
17  integer, parameter :: limiter_mcbeta = 3
18  integer, parameter :: limiter_superbee = 4
19  integer, parameter :: limiter_vanleer = 5
20  integer, parameter :: limiter_albada = 6
21  integer, parameter :: limiter_koren = 7
22  integer, parameter :: limiter_cada = 8
23  integer, parameter :: limiter_cada3 = 9
24  integer, parameter :: limiter_schmid = 10
25  integer, parameter :: limiter_venk = 11
26  ! Special cases
27  integer, parameter :: limiter_ppm = 12
28  integer, parameter :: limiter_mp5 = 13
29  integer, parameter :: limiter_weno3 = 14
30  integer, parameter :: limiter_wenoyc3 = 15
31  integer, parameter :: limiter_weno5 = 16
32  integer, parameter :: limiter_weno5nm = 17
33  integer, parameter :: limiter_wenoz5 = 18
34  integer, parameter :: limiter_wenoz5nm = 19
35  integer, parameter :: limiter_wenozp5 = 20
36  integer, parameter :: limiter_wenozp5nm = 21
37  integer, parameter :: limiter_weno5cu6 = 22
38  integer, parameter :: limiter_teno5ad = 23
39  integer, parameter :: limiter_weno7 = 24
40  integer, parameter :: limiter_mpweno7 = 25
41 
42 contains
43 
44  integer function limiter_type(namelim)
45  use mod_comm_lib, only: mpistop
46  character(len=*), intent(in) :: namelim
47 
48  select case (namelim)
49  case ('minmod')
51  case ('woodward')
53  case ('mcbeta')
55  case ('superbee')
57  case ('vanleer')
59  case ('albada')
61  case ('koren')
63  case ('cada')
65  case ('cada3')
67  case ('schmid1')
69  case ('schmid2')
71  case('venk')
73  case ('ppm')
75  case ('mp5')
77  case ('weno3')
79  case ('wenoyc3')
81  case ('weno5')
83  case ('weno5nm')
85  case ('wenoz5')
87  case ('wenoz5nm')
89  case ('wenozp5')
91  case ('wenozp5nm')
93  case ('weno5cu6')
95  case ('teno5ad')
97  case ('weno7')
99  case ('mpweno7')
101 
102  case default
103  limiter_type = -1
104  write(*,*) 'Unknown limiter: ', namelim
105  call mpistop("No such limiter")
106  end select
107  end function limiter_type
108 
109  pure logical function limiter_symmetric(typelim)
110  integer, intent(in) :: typelim
111 
112  select case (typelim)
114  limiter_symmetric = .false.
115  case default
116  limiter_symmetric = .true.
117  end select
118  end function limiter_symmetric
119 
120  !> Limit the centered dwC differences within ixC for iw in direction idim.
121  !> The limiter is chosen according to typelim.
122  !>
123  !> Note that this subroutine is called from upwindLR (hence from methods
124  !> like tvdlf, hancock, hll(c) etc) or directly from tvd.t,
125  !> but also from the gradientS and divvectorS subroutines in geometry.t
126  !> Accordingly, the typelim here corresponds to one of limiter
127  !> or one of gradient_limiter.
128  subroutine dwlimiter2(dwC,ixI^L,ixC^L,idims,typelim,ldw,rdw,a2max)
129 
131  use mod_comm_lib, only: mpistop
132 
133  integer, intent(in) :: ixI^L, ixC^L, idims
134  double precision, intent(in) :: dwC(ixI^S)
135  integer, intent(in) :: typelim
136  !> Result using left-limiter (same as right for symmetric)
137  double precision, intent(out), optional :: ldw(ixI^S)
138  !> Result using right-limiter (same as left for symmetric)
139  double precision, intent(out), optional :: rdw(ixI^S)
140  double precision, intent(in), optional :: a2max
141 
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))
145 
146  ! mcbeta limiter parameter value
147  double precision, parameter :: c_mcbeta=1.4d0
148  ! cada limiter parameter values
149  double precision, parameter :: cadalfa=0.5d0, cadbeta=2.0d0, cadgamma=1.6d0
150  ! full third order cada limiter
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
155  !-----------------------------------------------------------------------------
156 
157  ! Contract indices in idim for output.
158  ixomin^d=ixcmin^d+kr(idims,^d); ixomax^d=ixcmax^d;
159  hxo^l=ixo^l-kr(idims,^d);
160 
161  ! About the notation: the conventional argument theta (the ratio of slopes)
162  ! would be given by dwC(ixO^S)/dwC(hxO^S). However, in the end one
163  ! multiplies phi(theta) by dwC(hxO^S), which is incorporated in the
164  ! equations below. The minmod limiter can for example be written as:
165  ! A:
166  ! max(0.0d0, min(1.0d0, dwC(ixO^S)/dwC(hxO^S))) * dwC(hxO^S)
167  ! B:
168  ! tmp(ixO^S)*max(0.0d0,min(abs(dwC(ixO^S)),tmp(ixO^S)*dwC(hxO^S)))
169  ! where tmp(ixO^S)=sign(1.0d0,dwC(ixO^S))
170 
171  select case (typelim)
172  case (limiter_minmod)
173  ! Minmod limiter eq(3.51e) and (eq.3.38e) with omega=1
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)
179  case (limiter_woodward)
180  ! Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
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)
187  case (limiter_mcbeta)
188  ! Woodward and Collela limiter, with factor beta
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)
195  case (limiter_superbee)
196  ! Roes superbee limiter (eq.3.51i)
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)
203  case (limiter_vanleer)
204  ! van Leer limiter (eq 3.51f), but a missing delta2=1.D-12 is added
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)
209  case (limiter_albada)
210  ! Albada limiter (eq.3.51g) with delta2=1D.-12
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)
216  case (limiter_koren)
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))
223  end if
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))
228  end if
229  case (limiter_cada)
230  ! This limiter has been rewritten in the usual form, and uses a division
231  ! of the gradients.
232  if (present(ldw)) then
233  ! Cada Left variant
234  ! Compute theta, but avoid division by zero
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)
241  end if
242 
243  if (present(rdw)) then
244  ! Cada Right variant
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)
251  end if
252  case (limiter_cada3)
253  rdelinv=one/(cada3_radius*dxlevel(idims))**2
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)
264  elsewhere
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))
270  endwhere
271  ldw(ixo^s)=ldw(ixo^s) * dwc(ixo^s)
272  end if
273 
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)
283  elsewhere
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))
289  endwhere
290  rdw(ixo^s)=rdw(ixo^s) * dwc(hxo^s)
291  end if
292  case(limiter_schmid)
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)
304  else where
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))
310  end where
311  ldw(ixo^s)=ldw(ixo^s)*dwc(ixo^s)
312  end if
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)
322  else where
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))
328  end where
329  rdw(ixo^s)=rdw(ixo^s)*dwc(hxo^s)
330  end if
331  case default
332  call mpistop("Error in dwLimiter: unknown limiter")
333  end select
334 
335  end subroutine dwlimiter2
336 
337 end module mod_limiter
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
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.
Definition: mod_limiter.t:2
integer, parameter limiter_mcbeta
Definition: mod_limiter.t:17
integer, parameter limiter_weno5cu6
Definition: mod_limiter.t:37
integer, parameter limiter_mpweno7
Definition: mod_limiter.t:40
integer, parameter limiter_schmid
Definition: mod_limiter.t:24
pure logical function limiter_symmetric(typelim)
Definition: mod_limiter.t:110
integer, parameter limiter_teno5ad
Definition: mod_limiter.t:38
integer, parameter limiter_weno3
Definition: mod_limiter.t:29
integer, parameter limiter_vanleer
Definition: mod_limiter.t:19
integer, parameter limiter_ppm
Definition: mod_limiter.t:27
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...
Definition: mod_limiter.t:129
integer function limiter_type(namelim)
Definition: mod_limiter.t:45
integer, parameter limiter_cada3
Definition: mod_limiter.t:23
double precision d
Definition: mod_limiter.t:14
integer, parameter limiter_wenozp5
Definition: mod_limiter.t:35
integer, parameter limiter_woodward
Definition: mod_limiter.t:16
integer, parameter limiter_superbee
Definition: mod_limiter.t:18
integer, parameter limiter_weno5
Definition: mod_limiter.t:31
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:33
integer, parameter limiter_wenoz5nm
Definition: mod_limiter.t:34
integer, parameter limiter_koren
Definition: mod_limiter.t:21
integer, parameter limiter_weno7
Definition: mod_limiter.t:39
integer, parameter limiter_wenozp5nm
Definition: mod_limiter.t:36
integer, parameter limiter_cada
Definition: mod_limiter.t:22
integer, parameter limiter_mp5
Definition: mod_limiter.t:28
double precision schmid_rad
Definition: mod_limiter.t:14
integer, parameter limiter_venk
Definition: mod_limiter.t:25
integer, parameter limiter_weno5nm
Definition: mod_limiter.t:32
integer, parameter limiter_wenoyc3
Definition: mod_limiter.t:30
integer, parameter limiter_albada
Definition: mod_limiter.t:20
integer, parameter limiter_minmod
Definition: mod_limiter.t:15
Module containing the MP5 (fifth order) flux scheme.
Definition: mod_mp5.t:2
Definition: mod_ppm.t:1