Afivo
0.3
|
Module containing slope limiters. More...
Functions/Subroutines | |
elemental real(dp) function, public | af_limiter_apply (a, b, limiter) |
Apply one of the limiters. More... | |
elemental real(dp) function, public | af_limiter_koren (a, b) |
Modified implementation of Koren limiter, to avoid division and the min/max functions, which can be problematic / expensive. In most literature, you have r = a / b (ratio of gradients). Then the limiter phi(r) is multiplied with b. With this implementation, you get phi(r) * b. More... | |
elemental real(dp) function, public | af_limiter_vanleer (a, b) |
elemental real(dp) function, public | af_limiter_minmod (a, b) |
elemental real(dp) function, public | af_limiter_gminmod43 (a, b) |
elemental real(dp) function, public | af_limiter_mc (a, b) |
elemental real(dp) function | af_limiter_gminmod (a, b, theta) |
Generalized minmod limiter. The parameter theta controls how dissipative the limiter is, with 1 corresponding to the minmod limiter and 2 to the MC limiter. More... | |
Variables | |
integer, parameter | af_num_limiters = 7 |
Number of limiters. More... | |
integer, parameter, public | af_limiter_none_t = 1 |
integer, parameter, public | af_limiter_vanleer_t = 2 |
van Leer limiter More... | |
integer, parameter, public | af_limiter_koren_t = 3 |
Koren limiter. More... | |
integer, parameter, public | af_limiter_minmod_t = 4 |
Minmod limiter. More... | |
integer, parameter, public | af_limiter_mc_t = 5 |
MC limiter. More... | |
integer, parameter, public | af_limiter_gminmod43_t = 6 |
Generalized minmod limiter with theta = 4/3. More... | |
integer, parameter, public | af_limiter_zero_t = 7 |
All slopes are zero. More... | |
logical, dimension(af_num_limiters), parameter, public | af_limiter_symmetric = [.true., .true., .false., .true., .true., .true., .true.] |
Whether limiters are symmetric. More... | |
Module containing slope limiters.
elemental real(dp) function, public m_af_limiters::af_limiter_apply | ( | real(dp), intent(in) | a, |
real(dp), intent(in) | b, | ||
integer, intent(in) | limiter | ||
) |
Apply one of the limiters.
[in] | a | Slopes from one side |
[in] | b | Slopes from other side |
[in] | limiter | Which limiter to use |
Definition at line 41 of file m_af_limiters.f90.
elemental real(dp) function, public m_af_limiters::af_limiter_koren | ( | real(dp), intent(in) | a, |
real(dp), intent(in) | b | ||
) |
Modified implementation of Koren limiter, to avoid division and the min/max functions, which can be problematic / expensive. In most literature, you have r = a / b (ratio of gradients). Then the limiter phi(r) is multiplied with b. With this implementation, you get phi(r) * b.
[in] | a | Density gradient (numerator) |
[in] | b | Density gradient (denominator) |
Definition at line 72 of file m_af_limiters.f90.
elemental real(dp) function, public m_af_limiters::af_limiter_vanleer | ( | real(dp), intent(in) | a, |
real(dp), intent(in) | b | ||
) |
Definition at line 97 of file m_af_limiters.f90.
elemental real(dp) function, public m_af_limiters::af_limiter_minmod | ( | real(dp), intent(in) | a, |
real(dp), intent(in) | b | ||
) |
Definition at line 112 of file m_af_limiters.f90.
elemental real(dp) function, public m_af_limiters::af_limiter_gminmod43 | ( | real(dp), intent(in) | a, |
real(dp), intent(in) | b | ||
) |
Definition at line 120 of file m_af_limiters.f90.
elemental real(dp) function, public m_af_limiters::af_limiter_mc | ( | real(dp), intent(in) | a, |
real(dp), intent(in) | b | ||
) |
Definition at line 129 of file m_af_limiters.f90.
|
private |
Generalized minmod limiter. The parameter theta controls how dissipative the limiter is, with 1 corresponding to the minmod limiter and 2 to the MC limiter.
Definition at line 139 of file m_af_limiters.f90.
|
private |
Number of limiters.
Definition at line 10 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_none_t = 1 |
Definition at line 13 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_vanleer_t = 2 |
van Leer limiter
Definition at line 15 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_koren_t = 3 |
Koren limiter.
Definition at line 17 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_minmod_t = 4 |
Minmod limiter.
Definition at line 19 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_mc_t = 5 |
MC limiter.
Definition at line 21 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_gminmod43_t = 6 |
Generalized minmod limiter with theta = 4/3.
Definition at line 23 of file m_af_limiters.f90.
integer, parameter, public m_af_limiters::af_limiter_zero_t = 7 |
All slopes are zero.
Definition at line 25 of file m_af_limiters.f90.
logical, dimension(af_num_limiters), parameter, public m_af_limiters::af_limiter_symmetric = [.true., .true., .false., .true., .true., .true., .true.] |
Whether limiters are symmetric.
Definition at line 28 of file m_af_limiters.f90.