INTERPOLATIONS: Add a couple of limiters.
Use limiter_average() which is not a limiter, just returns average of left and right derivatives, as the default limiter for prolongation. Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
This commit is contained in:
parent
d3c0737c8d
commit
6ffb3f5910
@ -205,8 +205,8 @@ module interpolations
|
||||
!
|
||||
character(len=255) :: sreconstruction = "tvd"
|
||||
character(len=255) :: tlimiter = "mm"
|
||||
character(len=255) :: plimiter = "mm"
|
||||
character(len=255) :: climiter = "mm"
|
||||
character(len=255) :: plimiter = "avg"
|
||||
character(len=255) :: mlp_limiting = "off"
|
||||
character(len=255) :: positivity_fix = "on"
|
||||
character(len=255) :: clip_extrema = "off"
|
||||
@ -476,6 +476,10 @@ module interpolations
|
||||
name_tlim = "van Albada"
|
||||
limiter_tvd => limiter_vanalbada
|
||||
order = max(2, order)
|
||||
case ("ospre")
|
||||
name_tlim = "ospre"
|
||||
limiter_tvd => limiter_ospre
|
||||
order = max(2, order)
|
||||
case default
|
||||
name_tlim = "zero derivative"
|
||||
limiter_tvd => limiter_zero
|
||||
@ -485,6 +489,9 @@ module interpolations
|
||||
! select the prolongation limiter
|
||||
!
|
||||
select case(trim(plimiter))
|
||||
case ("avg", "average", "mean")
|
||||
name_plim = "average"
|
||||
limiter_prol => limiter_average
|
||||
case ("mm", "minmod")
|
||||
name_plim = "minmod"
|
||||
limiter_prol => limiter_minmod
|
||||
@ -497,6 +504,12 @@ module interpolations
|
||||
case ("vl", "vanleer")
|
||||
name_plim = "van Leer"
|
||||
limiter_prol => limiter_vanleer
|
||||
case ("va", "vanalbada")
|
||||
name_plim = "van Albada"
|
||||
limiter_prol => limiter_vanalbada
|
||||
case ("ospre")
|
||||
name_plim = "ospre"
|
||||
limiter_prol => limiter_ospre
|
||||
case default
|
||||
name_plim = "zero derivative"
|
||||
limiter_prol => limiter_zero
|
||||
@ -517,6 +530,9 @@ module interpolations
|
||||
case ("vl", "vanleer")
|
||||
name_clim = "van Leer"
|
||||
limiter_clip => limiter_vanleer
|
||||
case ("ospre")
|
||||
name_clim = "ospre"
|
||||
limiter_clip => limiter_ospre
|
||||
case default
|
||||
name_clim = "zero derivative"
|
||||
limiter_clip => limiter_zero
|
||||
@ -5148,6 +5164,35 @@ module interpolations
|
||||
!
|
||||
!===============================================================================
|
||||
!
|
||||
! function LIMITER_AVERAGE:
|
||||
! ------------------------
|
||||
!
|
||||
! Function returns the average among two arguments.
|
||||
!
|
||||
! Arguments:
|
||||
!
|
||||
! x - scaling factor;
|
||||
! a, b - the input values;
|
||||
!
|
||||
!===============================================================================
|
||||
!
|
||||
function limiter_average(x, a, b) result(c)
|
||||
|
||||
implicit none
|
||||
|
||||
real(kind=8), intent(in) :: x, a, b
|
||||
real(kind=8) :: c
|
||||
!
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
c = 5.0d-01 * (a + b) * x
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
end function limiter_average
|
||||
!
|
||||
!===============================================================================
|
||||
!
|
||||
! function LIMITER_MINMOD:
|
||||
! -----------------------
|
||||
!
|
||||
@ -5290,6 +5335,41 @@ module interpolations
|
||||
!
|
||||
!===============================================================================
|
||||
!
|
||||
! function LIMITER_OSPRE:
|
||||
! ----------------------
|
||||
!
|
||||
! Function returns the minimum module value among two arguments using
|
||||
! OSPRE limiter.
|
||||
!
|
||||
! Arguments:
|
||||
!
|
||||
! x - scaling factor;
|
||||
! a, b - the input values;
|
||||
!
|
||||
!===============================================================================
|
||||
!
|
||||
function limiter_ospre(x, a, b) result(c)
|
||||
|
||||
! local variables are not implicit by default
|
||||
!
|
||||
implicit none
|
||||
|
||||
! input arguments
|
||||
!
|
||||
real(kind=8), intent(in) :: x, a, b
|
||||
real(kind=8) :: c
|
||||
!
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
c = a * b
|
||||
c = 1.5d+00 * x * c * (a + b) / max(eps, a * a + c + b * b)
|
||||
|
||||
!-------------------------------------------------------------------------------
|
||||
!
|
||||
end function limiter_ospre
|
||||
!
|
||||
!===============================================================================
|
||||
!
|
||||
! function LIMITER_VANALBADA:
|
||||
! --------------------------
|
||||
!
|
||||
|
Loading…
x
Reference in New Issue
Block a user