mp5 Module

Applies a 5th-order unlimited reconstruction and clips it to a monotonicity-preserving bound computed from the local curvature. Achieves 5th-order accuracy in smooth regions while suppressing oscillations near discontinuities.

References: [1] A. Suresh and H. T. Huynh, "Accurate Monotonicity-Preserving Schemes with Runge-Kutta Time Stepping," J. Comput. Phys., 136:83-99, 1997.



Variables

Type Visibility Attributes Name Initial
real(kind=wp), private, parameter :: alpha_mp5 = 4.0_wp

MP5 monotonicity parameter α. Used in the upper-bound estimate. α = 4 following Suresh & Huynh (1997), Eq. (2.27).

real(kind=wp), private, parameter :: beta_mp5 = 4.0_wp/3.0_wp

MP5 curvature-limiter coefficient β_L. β_L = 4/3 following Suresh & Huynh (1997), Eq. (2.29).


Functions

private pure elemental function minmod2(a, b) result(m)

Two-argument minmod function. Returns the value with the smaller magnitude if a and b have the same sign, else 0.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: a
real(kind=wp), intent(in) :: b

Return Value real(kind=wp)

private pure elemental function minmod4(a, b, c, d) result(m)

Four-argument minmod function. Returns 0 if the arguments do not all share the same sign; otherwise returns the value with the smallest magnitude.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: a
real(kind=wp), intent(in) :: b
real(kind=wp), intent(in) :: c
real(kind=wp), intent(in) :: d

Return Value real(kind=wp)


Subroutines

public pure subroutine mp5_reconstruct(f_stencil, f_hat)

Perform MP5 reconstruction on the supplied stencil.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: f_stencil(:,:)
real(kind=wp), intent(out) :: f_hat(:)