weno Derived Type

type, public :: weno

WENO class.


Components

Type Visibility Attributes Name Initial
character(len=:), public, allocatable :: msg

error message

integer, public :: ierr = 0

error status

integer, public :: ncells

number of cells

integer, public :: k = 3

order of reconstruction within the cell (k = 1, 2 or 3)

real(kind=rk), public :: eps = 1e-6_rk

numerical smoothing factor

real(kind=rk), public, allocatable :: cnu(:,:,:)

array of constants for a non-uniform grid


Constructor

public interface weno

  • private pure function weno_init(ncells, k, eps, xedges) result(self)

    Initialize weno object.

    Note

    If the grid is not uniform, the user must supply the edges of the grid.

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: ncells

    number of cells

    integer, intent(in), optional :: k

    (2k - 1) is the order of reconstruction within the cell (1 <= k <= 3). By default, k = 3, i.e. 5th order reconstruction.

    real(kind=rk), intent(in), optional :: eps

    numerical smoothing factor (default=1e-6)

    real(kind=rk), intent(in), optional :: xedges(0:)

    vector(0:ncells) of cell edges; xedges(i) is the value of x at the right boundary of cell i, ; xedges(i-1) is the value of x at the left boundary of cell i, .

    Return Value type(weno)


Type-Bound Procedures

procedure, public, pass(self) :: reconstruct => weno_reconstruct

  • private pure subroutine weno_reconstruct(self, v, vl, vr)

    This subroutine implements the (2k-1)th order WENO method for arbitrary (uniform or non-uniform) finite volume/difference schemes described in ICASE 97-65 (Shu, 1997). The method is applicable to scalar as well as multicomponent problems. In the later case, the reconstruction is applied in a component-by-component fashion. The scheme below depics a generic finite volume discretization and the notation used (see Arguments).

     --0--|--1--| ...  |--(i-1)--|-------(i)-------|--(i+1)--| .... |--nc--|--(ncells+1)---
                                  ^               ^
                                  vl(i)           vr(i)
                                  v_{i-1/2}^+     v_{i+1/2}^-
    

    The procedure can equally be used for finite difference methods. In that case, 'v' is not the average cell value, but rather the flux! See section 2.3.2, page 22.

    @note Note For a scalar 1D problem, this procedure is called once per time step. In contrast, for a scalar 2D problem, it is called (ncells1+ncells2) times per step. So, efficiency is very important. The current implementation is rather general in terms of order and grid type, but at the cost of a number of 'select case' constructs. I wonder if it would be wise to make a specific version for k=2 (3rd order) and uniform grids to get maximum performance for multi-dimensional problems.

    Arguments

    Type IntentOptional Attributes Name
    class(weno), intent(in) :: self

    object

    real(kind=rk), intent(in) :: v(:)

    vector(ncells) of average cell values (if finite volume)

    real(kind=rk), intent(out) :: vl(:)

    vector(ncells) of reconstructed v at left boundary of cell i,

    real(kind=rk), intent(out) :: vr(:)

    vector(ncells) of reconstructed v at right boundary of cell i,

Source Code

   type :: weno
   !! WENO class.
      character(:), allocatable :: msg
         !! error message
      integer :: ierr = 0
         !! error status
      integer :: ncells
         !! number of cells
      integer :: k = 3
         !! order of reconstruction within the cell (k = 1, 2 or 3)
      real(rk) :: eps = 1e-6_rk
         !! numerical smoothing factor
      real(rk), allocatable, private :: d(:)
         !! array of constants
      real(rk), allocatable, private :: c(:, :)
         !! array of constants for uniform grid
      logical, private :: uniform_grid
         !! flag for uniform grid
      real(rk), allocatable :: cnu(:, :, :)
         !! array of constants for a non-uniform grid
   contains
      procedure, pass(self) :: reconstruct => weno_reconstruct
      procedure, pass(self), private :: calc_cnu => weno_calc_cnu
   end type weno