!> @file weno5_js.f90 !> @brief WENO5-JS (Jiang & Shu 1996) reconstruction scheme. !! !! 5th-order weighted essentially non-oscillatory reconstruction using the !! Jiang-Shu nonlinear weight formula. Substencil reconstructions and !! smoothness indicators are provided by the weno_family module. !! !! References: !! [1] G.-S. Jiang and C.-W. Shu, "Efficient Implementation of Weighted ENO !! Schemes," J. Comput. Phys., 126:202-228, 1996. module weno5_js use precision, only: wp use solver_state, only: neq use weno_family, only: eps_weno, d0, d1, d2, weno5_substencils, weno5_smoothness implicit none private public :: weno5_js_reconstruct contains ! --------------------------------------------------------------------------- !> Perform WENO5-JS reconstruction on the supplied stencil. !! !! Applies WENO5-JS reconstruction component-wise to a 5-point stencil. !! The stencil may be in physical or characteristic space depending on !! whether spatial_discretization applied characteristic projection first. !! !! Nonlinear weight formula (Jiang & Shu 1996, Eq. 22): !! omega_tilde_k = d_k / (eps + beta_k)^2 !! omega_k = omega_tilde_k / sum(omega_tilde) !! !! @param[in] f_stencil Stencil at 5 points (neq x 5) !! @param[out] f_hat Reconstructed value at the interface ! --------------------------------------------------------------------------- pure subroutine weno5_js_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) real(wp) :: beta0(neq), beta1(neq), beta2(neq) real(wp) :: w0(neq), w1(neq), w2(neq), w_sum(neq) integer :: i call weno5_substencils(f_stencil, vp0, vp1, vp2) call weno5_smoothness(f_stencil, beta0, beta1, beta2) do i = 1, neq w0(i) = d0 / (eps_weno + beta0(i))**2 w1(i) = d1 / (eps_weno + beta1(i))**2 w2(i) = d2 / (eps_weno + beta2(i))**2 end do w_sum = w0 + w1 + w2 do i = 1, neq f_hat(i) = (vp0(i) * w0(i) + vp1(i) * w1(i) + vp2(i) * w2(i)) / w_sum(i) end do end subroutine weno5_js_reconstruct end module weno5_js