muscl Module

Provides five MUSCL reconstructors, one per TVD limiter. The dispatcher in the reconstruction module selects the appropriate procedure pointer at initialisation time based on the limiter name supplied by the user.

Stencil layout (half_stencil=2, stencil_width=3): col 1 = f_{i-1}, col 2 = f_i (centre), col 3 = f_{i+1}

Formula (component-wise, left-biased stencil): ΔL = s(k,2) - s(k,1) ΔR = s(k,3) - s(k,2) r = ΔR / ΔL (when |ΔL| > tiny, else 0 → phi = 0) f_hat(k) = s(k,2) + 0.5 · phi(r) · ΔL

The right-biased (reversed) stencil for F^- is assembled by spatial_discretization in mirror order; the same formula gives the correct downwind reconstruction automatically.

References: [1] B. van Leer, "Towards the Ultimate Conservative Difference Scheme. V. A Second-Order Sequel to Godunov's Method," J. Comput. Phys., 32:101-136, 1979.


Uses


Abstract Interfaces

abstract interface

Abstract interface for a pure scalar TVD limiter function phi(r).

  • private pure function limiter_iface(r) result(phi)

    Arguments

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

    Return Value real(kind=wp)


Functions

private pure function lim_minmod(r) result(phi)

Minmod TVD limiter. Most dissipative of the standard limiters; guaranteed TVD for any CFL within the stability limit.

Read more…

Arguments

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

Return Value real(kind=wp)

private pure function lim_superbee(r) result(phi)

Superbee TVD limiter. Least dissipative TVD limiter; may over-sharpen contact discontinuities on smooth problems.

Read more…

Arguments

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

Return Value real(kind=wp)

private pure function lim_mc(r) result(phi)

Monotonized Central (MC) TVD limiter. Good balance between dissipation and resolution; recommended for most shock-capturing applications.

Read more…

Arguments

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

Return Value real(kind=wp)

private pure function lim_van_leer(r) result(phi)

van Leer smooth TVD limiter.

Read more…

Arguments

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

Return Value real(kind=wp)

private pure function lim_koren(r) result(phi)

Koren TVD limiter. Third-order accurate for smooth data; preserves monotonicity near discontinuities.

Read more…

Arguments

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

Return Value real(kind=wp)


Subroutines

private pure subroutine muscl_reconstruct(f_stencil, f_hat, lim_fn)

Core MUSCL reconstruction kernel (component-wise, left-biased stencil).

Read more…

Arguments

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

TVD limiter function phi(r)

public pure subroutine muscl_minmod_reconstruct(f_stencil, f_hat)

MUSCL reconstruction with the Minmod limiter.

Arguments

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

public pure subroutine muscl_superbee_reconstruct(f_stencil, f_hat)

MUSCL reconstruction with the Superbee limiter.

Arguments

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

public pure subroutine muscl_mc_reconstruct(f_stencil, f_hat)

MUSCL reconstruction with the Monotonized Central (MC) limiter.

Arguments

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

public pure subroutine muscl_van_leer_reconstruct(f_stencil, f_hat)

MUSCL reconstruction with the van Leer smooth limiter.

Arguments

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

public pure subroutine muscl_koren_reconstruct(f_stencil, f_hat)

MUSCL reconstruction with the Koren limiter.

Arguments

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