weno5_z Module

5th-order WENO reconstruction with the Z-weight formula. Uses the same substencil reconstructions and smoothness indicators as WENO5-JS but replaces the Jiang-Shu weight formula with the Z-weight formula, which achieves full 5th-order accuracy at smooth extrema where WENO5-JS degrades to 3rd order.

Weight formula (Borges et al. 2008, Eq. 3.4-3.5): tau5 = |beta0 - beta2| (global indicator) alpha_k = d_k * (1 + (tau5/(eps+beta_k))^p), p = p_wenoz = 2 omega_k = alpha_k / sum(alpha)

References: [1] R. Borges, M. Carmona, B. Costa, and W. S. Don, "An improved weighted essentially non-oscillatory scheme for hyperbolic conservation laws," J. Comput. Phys., 227:3191-3211, 2008. [2] G.-S. Jiang and C.-W. Shu, "Efficient Implementation of Weighted ENO Schemes," J. Comput. Phys., 126:202-228, 1996.



Variables

Type Visibility Attributes Name Initial
integer, private, parameter :: p_wenoz = 2

Exponent for the WENO-Z enhancement factor. p = 2 gives better shock sharpness than p = 1. See Borges et al. (2008), Eq. (3.5).


Subroutines

public pure subroutine weno5_z_reconstruct(f_stencil, f_hat)

Perform WENO5-Z reconstruction on the supplied stencil.

Arguments

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