weno_family Module

Provides the parameters and helper subroutines that are common across the solver's WENO-family schemes. The original WENO5 helpers are retained for WENO5-JS, WENO5-Z, and TENO5. More general candidate-polynomial and smoothness-evaluation helpers support WENO-CU6, WENO9-JS, and WENO11-JS without duplicating large blocks of algebra.

Exported parameters: eps_weno — small stabiliser to avoid division by zero in smoothness weights. d0, d1, d2 — linear (optimal) WENO5 weights; sum = 1.

Exported subroutines: weno5_substencils — computes 3rd-order substencil reconstructions vp0..vp2. weno5_smoothness — computes Jiang-Shu smoothness indicators beta0..beta2. compute_candidate_reconstructions — evaluates candidate face polynomials from tabulated coefficients. compute_candidate_smoothness — evaluates Jiang-Shu smoothness forms from tabulated quadratic matrices. compute_full_stencil_smoothness — evaluates a quadratic smoothness form on the full stencil (used by CU6). apply_weno_js_generic — generic Jiang-Shu nonlinear weighting for odd-width WENO schemes.

References: [1] G.-S. Jiang and C.-W. Shu, "Efficient Implementation of Weighted ENO Schemes," J. Comput. Phys., 126:202-228, 1996.


Uses


Variables

Type Visibility Attributes Name Initial
real(kind=wp), public, parameter :: eps_weno = 1.0e-6_wp

Small parameter to avoid division by zero in WENO weights. See Jiang & Shu (1996), Eq. (22). Typical literature range: 1e-6 to 1e-4; order-of-magnitude changes have negligible effect in smooth regions but may affect shock sharpness.

real(kind=wp), public, parameter :: d0 = 0.1_wp

Linear (optimal) weights for WENO5 5th-order accuracy. gamma_0 = 1/10, gamma_1 = 3/5, gamma_2 = 3/10. See Jiang & Shu (1996), Eq. (19).

real(kind=wp), public, parameter :: d1 = 0.6_wp
real(kind=wp), public, parameter :: d2 = 0.3_wp
real(kind=wp), public, parameter :: d7_0 = 1.0_wp/35.0_wp

Linear (optimal) weights for WENO7 7th-order accuracy. d7_0=1/35, d7_1=12/35, d7_2=18/35, d7_3=4/35; sum = 1. See Balsara & Shu (2000), Eq. (2.16).

real(kind=wp), public, parameter :: d7_1 = 12.0_wp/35.0_wp
real(kind=wp), public, parameter :: d7_2 = 18.0_wp/35.0_wp
real(kind=wp), public, parameter :: d7_3 = 4.0_wp/35.0_wp

Subroutines

public pure subroutine weno5_substencils(f_stencil, vp0, vp1, vp2)

Compute the three 3rd-order substencil reconstructions for WENO5.

Read more…

Arguments

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

public pure subroutine weno5_smoothness(f_stencil, beta0, beta1, beta2)

Compute the Jiang-Shu smoothness indicators for a 5-point WENO5 stencil.

Read more…

Arguments

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

public pure subroutine compute_candidate_reconstructions(f_stencil, coeffs, q_candidate)

Evaluate candidate reconstruction polynomials from tabulated coefficients.

Read more…

Arguments

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

public pure subroutine compute_candidate_smoothness(f_stencil, beta_mats, beta)

Evaluate candidate Jiang-Shu smoothness indicators from quadratic forms.

Read more…

Arguments

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

public pure subroutine compute_full_stencil_smoothness(f_stencil, beta_mat, beta)

Evaluate a quadratic smoothness form on the full stencil.

Read more…

Arguments

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

public pure subroutine apply_weno_js_generic(f_stencil, coeffs, beta_mats, linear_weights, f_hat)

Generic Jiang-Shu WENO reconstruction for odd-width schemes.

Read more…

Arguments

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