Afivo  0.3
m_af_limiters.f90
1 #include "cpp_macros.h"
2 
4  use m_af_types
5 
6  implicit none
7  private
8 
10  integer, parameter :: af_num_limiters = 7
11 
12  ! Constants to identify limiters
13  integer, parameter, public :: af_limiter_none_t = 1
15  integer, parameter, public :: af_limiter_vanleer_t = 2
17  integer, parameter, public :: af_limiter_koren_t = 3
19  integer, parameter, public :: af_limiter_minmod_t = 4
21  integer, parameter, public :: af_limiter_mc_t = 5
23  integer, parameter, public :: af_limiter_gminmod43_t = 6
25  integer, parameter, public :: af_limiter_zero_t = 7
26 
28  logical, parameter, public :: af_limiter_symmetric(af_num_limiters) = &
29  [.true., .true., .false., .true., .true., .true., .true.]
30 
31  public :: af_limiter_apply
32  public :: af_limiter_koren
33  public :: af_limiter_vanleer
34  public :: af_limiter_minmod
35  public :: af_limiter_gminmod43
36  public :: af_limiter_mc
37 
38 contains
39 
41  elemental function af_limiter_apply(a, b, limiter) result(slope)
42  real(dp), intent(in) :: a
43  real(dp), intent(in) :: b
44  integer, intent(in) :: limiter
45  real(dp) :: slope
46 
47  select case (limiter)
48  case (af_limiter_none_t)
49  slope = 0.5_dp * (a + b)
50  case (af_limiter_vanleer_t)
51  slope = af_limiter_vanleer(a, b)
52  case (af_limiter_koren_t)
53  slope = af_limiter_koren(a, b)
54  case (af_limiter_minmod_t)
55  slope = af_limiter_minmod(a, b)
56  case (af_limiter_mc_t)
57  slope = af_limiter_mc(a, b)
58  case (af_limiter_gminmod43_t)
59  slope = af_limiter_gminmod43(a, b)
60  case (af_limiter_zero_t)
61  slope = 0.0_dp
62  case default
63  ! Cannot (yet) stop in elemental function, so use default limiter
64  slope = af_limiter_vanleer(a, b)
65  end select
66  end function af_limiter_apply
67 
72  elemental function af_limiter_koren(a, b) result(bphi)
73  real(dp), intent(in) :: a
74  real(dp), intent(in) :: b
75  real(dp), parameter :: third = 1/3.0_dp
76  real(dp) :: bphi, aa, ab
77 
78  aa = a * a
79  ab = a * b
80 
81  if (ab <= 0) then
82  ! a and b have different sign or one of them is zero, so r is either 0,
83  ! inf or negative (special case a == b == 0 is ignored)
84  bphi = 0
85  else if (aa <= 0.25_dp * ab) then
86  ! 0 < a/b <= 1/4, limiter has value 2*a/b
87  bphi = 2*a
88  else if (aa <= 2.5_dp * ab) then
89  ! 1/4 < a/b <= 2.5, limiter has value (1+2*a/b)/3
90  bphi = third * (b + 2*a)
91  else
92  ! (1+2*a/b)/6 >= 1, limiter has value 2
93  bphi = 2*b
94  end if
95  end function af_limiter_koren
96 
97  elemental function af_limiter_vanleer(a, b) result(phi)
98  real(dp), intent(in) :: a
99  real(dp), intent(in) :: b
100  real(dp) :: ab, phi
101 
102  ab = a * b
103  if (ab > 0) then
104  ! Jannis: I think roundoff errors can lead to values slightly outside the
105  ! desired range
106  phi = 2 * ab / (a + b)
107  else
108  phi = 0
109  end if
110  end function af_limiter_vanleer
111 
112  elemental function af_limiter_minmod(a, b) result(phi)
113  real(dp), intent(in) :: a
114  real(dp), intent(in) :: b
115  real(dp) :: phi
116  phi = af_limiter_gminmod(a, b, 1.0_dp)
117  end function af_limiter_minmod
118 
119  ! Generalized minmod limiter with parameter theta = 4/3
120  elemental function af_limiter_gminmod43(a, b) result(phi)
121  real(dp), intent(in) :: a
122  real(dp), intent(in) :: b
123  real(dp) :: phi
124  real(dp), parameter :: theta = 4/3.0_dp
125  phi = af_limiter_gminmod(a, b, theta)
126  end function af_limiter_gminmod43
127 
128  ! MC (monotonized central) limiter
129  elemental function af_limiter_mc(a, b) result(phi)
130  real(dp), intent(in) :: a
131  real(dp), intent(in) :: b
132  real(dp) :: phi
133  phi = af_limiter_gminmod(a, b, 2.0_dp)
134  end function af_limiter_mc
135 
139  elemental function af_limiter_gminmod(a, b, theta) result(phi)
140  real(dp), intent(in) :: a, b, theta
141  real(dp) :: phi
142 
143  if (a * b > 0) then
144  phi = sign(minval(abs([theta * a, theta * b, &
145  0.5_dp * (a + b)])), a)
146  else
147  phi = 0.0_dp
148  end if
149  end function af_limiter_gminmod
150 
151 end module m_af_limiters
Module containing slope limiters.
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Definition: m_af_types.f90:4