linear_weno.f90 Source File


Source Code

!> @file linear_weno.f90
!> @brief Linear (optimal-weight) WENO reconstruction for the smooth hybrid path.
!!
!! Applies the WENO5 and WENO7 substencil reconstructions using fixed optimal
!! weights d_k (constants), bypassing the nonlinear beta-based weighting.  In
!! smooth regions the nonlinear WENO weights converge to d_k anyway; using them
!! directly eliminates the residual dissipation from imperfect weight matching
!! and preserves the full formal order of accuracy:
!!
!!   - linear_weno5_reconstruct : 5th-order (same order as WENO5-JS / WENO5-Z)
!!   - linear_weno7_reconstruct : 7th-order (same order as WENO7-JS)
!!
!! Both subroutines match reconstructor_iface and may be used as drop-in
!! replacements for the nonlinear variants in smooth regions.
!!
!! References:
!!   [1] G.-S. Jiang and C.-W. Shu, J. Comput. Phys., 126:202-228, 1996.
!!   [2] D. S. Balsara and C.-W. Shu, J. Comput. Phys., 160:405-452, 2000.

module linear_weno
  use precision, only: wp
  use solver_state, only: neq
  use weno_family, only: d0, d1, d2, d7_0, d7_1, d7_2, d7_3, weno5_substencils
  implicit none
  private

  public :: linear_weno5_reconstruct, linear_weno7_reconstruct

contains

  ! ---------------------------------------------------------------------------
  !> 5th-order linear WENO reconstruction on a 5-point stencil.
  !!
  !! Equivalent to WENO5-JS in the smooth limit: applies the optimal weights
  !! d0=1/10, d1=3/5, d2=3/10 to the three 3rd-order substencil approximations
  !! without computing beta-based nonlinear corrections.  The resulting
  !! combined coefficients are:
  !!
  !!   f_hat = (1/30)*v1 - (13/60)*v2 + (47/60)*v3 + (27/60)*v4 - (3/60)*v5
  !!
  !! See Jiang & Shu (1996), Eq. (19) combined with substencil Eq. (16a-16c).
  !!
  !! @param[in]  f_stencil  Stencil at 5 points (neq x 5)
  !! @param[out] f_hat      Reconstructed value at the interface
  ! ---------------------------------------------------------------------------
  pure subroutine linear_weno5_reconstruct(f_stencil, f_hat)
    real(wp), intent(in) :: f_stencil(:, :)
    real(wp), intent(out) :: f_hat(:)

    real(wp) :: vp0(neq), vp1(neq), vp2(neq)

    call weno5_substencils(f_stencil, vp0, vp1, vp2)
    ! d0=1/10 pairs with vp0 (leftmost substencil),
    ! d1=3/5  with vp1 (central), d2=3/10 with vp2 (rightmost).
    ! See Jiang & Shu (1996), Eq. (19).
    f_hat = d0 * vp0 + d1 * vp1 + d2 * vp2
  end subroutine linear_weno5_reconstruct

  ! ---------------------------------------------------------------------------
  !> 7th-order linear WENO reconstruction on a 7-point stencil.
  !!
  !! Applies the optimal weights d7_0=1/35, d7_1=12/35, d7_2=18/35, d7_3=4/35
  !! to the four 4th-order substencil approximations from Balsara & Shu (2000)
  !! without nonlinear smoothness weighting.
  !!
  !! Substencil approximations (Balsara & Shu 2000):
  !!   q0 = -(1/4)*v1  + (13/12)*v2 - (23/12)*v3 + (25/12)*v4
  !!   q1 =  (1/12)*v2 -  (5/12)*v3 + (13/12)*v4  + (1/4)*v5
  !!   q2 = -(1/12)*v3 +  (7/12)*v4 +  (7/12)*v5 - (1/12)*v6
  !!   q3 =   (1/4)*v4 + (13/12)*v5 -  (5/12)*v6  + (1/12)*v7
  !!
  !! @param[in]  f_stencil  Stencil at 7 points (neq x 7)
  !! @param[out] f_hat      Reconstructed value at the interface
  ! ---------------------------------------------------------------------------
  pure subroutine linear_weno7_reconstruct(f_stencil, f_hat)
    real(wp), intent(in) :: f_stencil(:, :)
    real(wp), intent(out) :: f_hat(:)

    real(wp) :: q0(neq), q1(neq), q2(neq), q3(neq)

    ! Substencil approximations — same formulas as weno7_js_reconstruct.
    q0 = -1.0_wp / 4.0_wp * f_stencil(:, 1) &
        & + 13.0_wp / 12.0_wp * f_stencil(:, 2) &
        & - 23.0_wp / 12.0_wp * f_stencil(:, 3) &
        & + 25.0_wp / 12.0_wp * f_stencil(:, 4)

    q1 = 1.0_wp / 12.0_wp * f_stencil(:, 2) &
        & - 5.0_wp / 12.0_wp * f_stencil(:, 3) &
        & + 13.0_wp / 12.0_wp * f_stencil(:, 4) &
        & + 1.0_wp / 4.0_wp * f_stencil(:, 5)

    q2 = -1.0_wp / 12.0_wp * f_stencil(:, 3) &
        & + 7.0_wp / 12.0_wp * f_stencil(:, 4) &
        & + 7.0_wp / 12.0_wp * f_stencil(:, 5) &
        & - 1.0_wp / 12.0_wp * f_stencil(:, 6)

    q3 = 1.0_wp / 4.0_wp * f_stencil(:, 4) &
        & + 13.0_wp / 12.0_wp * f_stencil(:, 5) &
        & - 5.0_wp / 12.0_wp * f_stencil(:, 6) &
        & + 1.0_wp / 12.0_wp * f_stencil(:, 7)

    f_hat = d7_0 * q0 + d7_1 * q1 + d7_2 * q2 + d7_3 * q3
  end subroutine linear_weno7_reconstruct

end module linear_weno