weno9_js.f90 Source File


Source Code

!> @file weno9_js.f90
!> @brief WENO9-JS finite-difference reconstruction scheme.
!!
!! Classical Jiang-Shu WENO on a 9-point stencil using five 5-point
!! substencils.  The coefficients here follow the finite-difference
!! flux-reconstruction form (not plain pointwise interpolation), matching the
!! existing WENO5/WENO7 implementation strategy in this solver.
module weno9_js
  use precision, only: wp
  use weno_family, only: apply_weno_js_generic
  implicit none
  private

  real(wp), parameter :: candidate_coeffs(5, 5) = reshape([ &
                                                          1.0_wp / 5.0_wp, -21.0_wp / 20.0_wp, 137.0_wp / 60.0_wp, -163.0_wp &
                                                          / 60.0_wp, &
                                                          137.0_wp / 60.0_wp, -1.0_wp / 20.0_wp, 17.0_wp / 60.0_wp, -43.0_wp &
                                                          / 60.0_wp, &
                                                          77.0_wp / 60.0_wp, 1.0_wp / 5.0_wp, 1.0_wp / 30.0_wp, -13.0_wp &
                                                          / 60.0_wp, &
                                                          47.0_wp / 60.0_wp, 9.0_wp / 20.0_wp, -1.0_wp / 20.0_wp, -1.0_wp &
                                                          / 20.0_wp, &
                                                          9.0_wp / 20.0_wp, 47.0_wp / 60.0_wp, -13.0_wp / 60.0_wp, 1.0_wp &
                                                          / 30.0_wp, &
                                                          1.0_wp / 5.0_wp, 77.0_wp / 60.0_wp, -43.0_wp / 60.0_wp, 17.0_wp &
                                                          / 60.0_wp, &
                                                          -1.0_wp / 20.0_wp &
                                                          ], [5, 5])

  real(wp), parameter :: linear_weights(5) = [ &
                         1.0_wp / 126.0_wp, 10.0_wp / 63.0_wp, 10.0_wp / 21.0_wp, 20.0_wp / 63.0_wp, &
                         5.0_wp / 126.0_wp]

  real(wp), parameter :: beta_mats(5, 5, 5) = reshape([ &
                                                      11329.0_wp / 2520.0_wp, -208501.0_wp / 5040.0_wp, 121621.0_wp / 1680.0_wp, &
                                                      -288007.0_wp / 5040.0_wp, &
                                                      86329.0_wp / 5040.0_wp, -208501.0_wp / 5040.0_wp, 482963.0_wp / 5040.0_wp, &
                                                      -142033.0_wp / 420.0_wp, &
                                                      679229.0_wp / 2520.0_wp, -411487.0_wp / 5040.0_wp, 121621.0_wp / 1680.0_wp, &
                                                      -142033.0_wp / 420.0_wp, &
                                                      507131.0_wp / 1680.0_wp, -68391.0_wp / 140.0_wp, 252941.0_wp / 1680.0_wp, &
                                                      -288007.0_wp / 5040.0_wp, &
                                                      679229.0_wp / 2520.0_wp, -68391.0_wp / 140.0_wp, 1020563.0_wp / 5040.0_wp, &
                                                      -649501.0_wp / 5040.0_wp, &
                                                      86329.0_wp / 5040.0_wp, -411487.0_wp / 5040.0_wp, 252941.0_wp / 1680.0_wp, &
                                                      -649501.0_wp / 5040.0_wp, &
                                                      53959.0_wp / 2520.0_wp, 1727.0_wp / 1260.0_wp, -60871.0_wp / 5040.0_wp, &
                                                      33071.0_wp / 1680.0_wp, &
                                                      -70237.0_wp / 5040.0_wp, 18079.0_wp / 5040.0_wp, -60871.0_wp / 5040.0_wp, &
                                                      138563.0_wp / 5040.0_wp, &
                                                      -3229.0_wp / 35.0_wp, 168509.0_wp / 2520.0_wp, -88297.0_wp / 5040.0_wp, &
                                                      33071.0_wp / 1680.0_wp, &
                                                      -3229.0_wp / 35.0_wp, 135431.0_wp / 1680.0_wp, -25499.0_wp / 210.0_wp, &
                                                      55051.0_wp / 1680.0_wp, &
                                                      -70237.0_wp / 5040.0_wp, 168509.0_wp / 2520.0_wp, -25499.0_wp / 210.0_wp, &
                                                      242723.0_wp / 5040.0_wp, &
                                                      -140251.0_wp / 5040.0_wp, 18079.0_wp / 5040.0_wp, -88297.0_wp / 5040.0_wp, &
                                                      55051.0_wp / 1680.0_wp, &
                                                      -140251.0_wp / 5040.0_wp, 11329.0_wp / 2520.0_wp, 1727.0_wp / 1260.0_wp, &
                                                      -51001.0_wp / 5040.0_wp, &
                                                      7547.0_wp / 560.0_wp, -38947.0_wp / 5040.0_wp, 8209.0_wp / 5040.0_wp, &
                                                      -51001.0_wp / 5040.0_wp, &
                                                      104963.0_wp / 5040.0_wp, -24923.0_wp / 420.0_wp, 89549.0_wp / 2520.0_wp, &
                                                      -38947.0_wp / 5040.0_wp, &
                                                      7547.0_wp / 560.0_wp, -24923.0_wp / 420.0_wp, 77051.0_wp / 1680.0_wp, &
                                                      -24923.0_wp / 420.0_wp, &
                                                      7547.0_wp / 560.0_wp, -38947.0_wp / 5040.0_wp, 89549.0_wp / 2520.0_wp, &
                                                      -24923.0_wp / 420.0_wp, &
                                                      104963.0_wp / 5040.0_wp, -51001.0_wp / 5040.0_wp, 8209.0_wp / 5040.0_wp, &
                                                      -38947.0_wp / 5040.0_wp, &
                                                      7547.0_wp / 560.0_wp, -51001.0_wp / 5040.0_wp, 1727.0_wp / 1260.0_wp, &
                                                      11329.0_wp / 2520.0_wp, &
                                                      -140251.0_wp / 5040.0_wp, 55051.0_wp / 1680.0_wp, -88297.0_wp / 5040.0_wp, &
                                                      18079.0_wp / 5040.0_wp, &
                                                      -140251.0_wp / 5040.0_wp, 242723.0_wp / 5040.0_wp, -25499.0_wp / 210.0_wp, &
                                                      168509.0_wp / 2520.0_wp, &
                                                      -70237.0_wp / 5040.0_wp, 55051.0_wp / 1680.0_wp, -25499.0_wp / 210.0_wp, &
                                                      135431.0_wp / 1680.0_wp, &
                                                      -3229.0_wp / 35.0_wp, 33071.0_wp / 1680.0_wp, -88297.0_wp / 5040.0_wp, &
                                                      168509.0_wp / 2520.0_wp, &
                                                      -3229.0_wp / 35.0_wp, 138563.0_wp / 5040.0_wp, -60871.0_wp / 5040.0_wp, &
                                                      18079.0_wp / 5040.0_wp, &
                                                      -70237.0_wp / 5040.0_wp, 33071.0_wp / 1680.0_wp, -60871.0_wp / 5040.0_wp, &
                                                      1727.0_wp / 1260.0_wp, &
                                                      53959.0_wp / 2520.0_wp, -649501.0_wp / 5040.0_wp, 252941.0_wp / 1680.0_wp, &
                                                      -411487.0_wp / 5040.0_wp, &
                                                      86329.0_wp / 5040.0_wp, -649501.0_wp / 5040.0_wp, 1020563.0_wp / 5040.0_wp, &
                                                      -68391.0_wp / 140.0_wp, &
                                                      679229.0_wp / 2520.0_wp, -288007.0_wp / 5040.0_wp, 252941.0_wp / 1680.0_wp, &
                                                      -68391.0_wp / 140.0_wp, &
                                                      507131.0_wp / 1680.0_wp, -142033.0_wp / 420.0_wp, 121621.0_wp / 1680.0_wp, &
                                                      -411487.0_wp / 5040.0_wp, &
                                                      679229.0_wp / 2520.0_wp, -142033.0_wp / 420.0_wp, 482963.0_wp / 5040.0_wp, &
                                                      -208501.0_wp / 5040.0_wp, &
                                                      86329.0_wp / 5040.0_wp, -288007.0_wp / 5040.0_wp, 121621.0_wp / 1680.0_wp, &
                                                      -208501.0_wp / 5040.0_wp, &
                                                      11329.0_wp / 2520.0_wp &
                                                      ], [5, 5, 5])

  public :: weno9_js_reconstruct

contains

  pure subroutine weno9_js_reconstruct(f_stencil, f_hat)
    real(wp), intent(in) :: f_stencil(:, :)
    real(wp), intent(out) :: f_hat(:)

    call apply_weno_js_generic(f_stencil, candidate_coeffs, beta_mats, linear_weights, f_hat)
  end subroutine weno9_js_reconstruct

end module weno9_js