teno5 Module

Uses the same 5-point stencil and substencil approximations as WENO5-JS (provided by weno_family), but replaces smooth nonlinear weights with a hard threshold that eliminates stencils containing discontinuities before combining the rest using the optimal (linear) weights. This gives sharper shock resolution than WENO5 while maintaining 5th-order accuracy in smooth regions.

Algorithm (Fu, Hu & Adams 2016, Eq. 3.5-3.7): 1. Compute substencil reconstructions vp0..vp2 and smoothness beta0..beta2 (both shared with WENO5-JS via weno_family). 2. Compute global indicator tau5 = |beta0 - beta2| (same as WENO5-Z). 3. Candidate sensitivity: gamma_k = (1 + tau5/(eps+beta_k))^q, q = q_teno. 4. Normalise: chi_k = gamma_k / sum(gamma). 5. Hard cutoff: delta_k = 0 if chi_k < ct_teno, else 1. 6. Weights: alpha_k = delta_k * d_k (d_k are the WENO5 optimal weights). 7. Fallback (sum(alpha) = 0): use the stencil with the smallest beta_k (ENO). 8. f_hat = sum(alpha_k/sum(alpha) * vp_k).

References: [1] L. Fu, X. Y. Hu, and N. A. Adams, "A family of high-order targeted ENO schemes for compressible-fluid simulations," J. Comput. Phys., 305:333-359, 2016.



Variables

Type Visibility Attributes Name Initial
integer, private, parameter :: q_teno = 6

Exponent for the TENO-5 candidate-stencil sensitivity. q = 6 following Fu, Hu & Adams (2016), Eq. (3.7).

real(kind=wp), private, parameter :: ct_teno = 1.0e-5_wp

Hard-threshold parameter for TENO-5 stencil elimination. Stencil k is discarded when chi_k < ct_teno. Fu, Hu & Adams (2016) recommend 1e-5 to 1e-6; smaller values yield sharper shock resolution but less aggressive oscillation suppression.


Subroutines

public pure subroutine teno5_reconstruct(f_stencil, f_hat)

Perform TENO5 reconstruction on the supplied stencil.

Arguments

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