Lagrange Derived Type

type, public, extends(Lagrange_t) :: Lagrange


Contents

Source Code


Components

TypeVisibilityAttributesNameInitial
integer, public :: M

The number of target points.

integer, public :: N

The number of control points.

real(kind=prec), public, pointer, contiguous, dimension(:,:):: bMatrix

The boundary interpolation matrix that is used to map a grid of nodal values at the control points to the element boundaries.

type(c_ptr), public :: bMatrix_gpu
real(kind=prec), public, pointer, contiguous, dimension(:):: bWeights

The barycentric weights that are calculated from the controlPoints and used for interpolation.

character(len=3), public :: backend ='gpu'
type(c_ptr), public :: blas_handle =c_null_ptr

A handle for working with hipblas

integer, public :: controlNodeType
real(kind=prec), public, pointer, contiguous, dimension(:):: controlPoints

The set of nodes in one dimension where data is known. To create higher dimension interpolation and differentiation operators, structured grids in two and three dimensions are created by tensor products of the controlPoints. This design decision implies that all spectral element methods supported by the Lagrange class have the same polynomial degree in each computational/spatial dimension. In practice, the controlPoints are the Legendre-Gauss, Legendre-Gauss-Lobatto, Legendre-Gauss-Radau, Chebyshev-Gauss, Chebyshev-Gauss-Lobatto, or Chebyshev-Gauss-Radau quadrature points over the domain [-1,1] (computational space). The Init routine for this class restricts controlPoints to one of these quadrature types or uniform points on [-1,1].

real(kind=prec), public, pointer, contiguous, dimension(:,:):: dMatrix

The derivative matrix for mapping function nodal values to a nodal values of the derivative estimate. The dMatrix is based on a strong form of the derivative.

type(c_ptr), public :: dMatrix_gpu
real(kind=prec), public, pointer, contiguous, dimension(:,:):: dgMatrix

The derivative matrix for mapping function nodal values to a nodal values of the derivative estimate. The dgMatrix is based on a weak form of the derivative. It must be used with bMatrix to account for boundary contributions in the weak form.

type(c_ptr), public :: dgMatrix_gpu
real(kind=prec), public, pointer, contiguous, dimension(:,:):: iMatrix

The interpolation matrix (transpose) for mapping data from the control grid to the target grid.

type(c_ptr), public :: iMatrix_gpu
real(kind=prec), public, pointer, contiguous, dimension(:):: qWeights

The quadrature weights for discrete integration. The quadradture weights depend on the type of controlPoints provided; one of Legendre-Gauss, Legendre-Gauss-Lobatto, Legendre-Gauss-Radau, Chebyshev-Gauss, Chebyshev-Gauss-Lobatto, Chebyshev-Gauss Radau, or Uniform. If Uniform, the quadrature weights are constant .

type(c_ptr), public :: qWeights_gpu
integer, public :: targetNodeType
real(kind=prec), public, pointer, contiguous, dimension(:):: targetPoints

The set of nodes in one dimension where data is to be interpolated to. To create higher dimension interpolation and differentiation operators, structured grids in two and three dimensions are created by tensor products of the targetPoints. In practice, the targetPoints are set to a uniformly distributed set of points between [-1,1] (computational space) to allow for interpolation from unevenly spaced quadrature points to a plotting grid.


Type-Bound Procedures

procedure, public :: CalculateBarycentricWeights

procedure, public :: CalculateDerivativeMatrix

procedure, public :: CalculateInterpolationMatrix

procedure, public :: CalculateLagrangePolynomials

  • public function CalculateLagrangePolynomials(this, sE) result(lAtS)

    Arguments

    TypeIntentOptionalAttributesName
    class(Lagrange_t) :: this
    real(kind=prec) :: sE

    Return Value real(kind=prec)(0:this%N)

procedure, public :: Free => Free_Lagrange

  • public subroutine Free_Lagrange(this)

    Frees all memory (host and device) associated with an instance of the Lagrange class

    Arguments

    TypeIntentOptionalAttributesName
    class(Lagrange), intent(inout) :: this

    Lagrange class instance

procedure, public :: Init => Init_Lagrange

  • public subroutine Init_Lagrange(this, N, controlNodeType, M, targetNodeType)

    Initialize an instance of the Lagrange class On output, all of the attributes for the Lagrange class are allocated and values are initialized according to the number of control points, number of target points, and the types for the control and target nodes. If a GPU is available, device pointers for the Lagrange attributes are allocated and initialized.

    Arguments

    TypeIntentOptionalAttributesName
    class(Lagrange), intent(out) :: this

    Lagrange class instance

    integer, intent(in) :: N

    The number of control points for interpolant

    integer, intent(in) :: controlNodeType

    The integer code specifying the type of control points. Parameters are defined in SELF_Constants.f90. One of GAUSS(=1), GAUSS_LOBATTO(=2), or UNIFORM(=3)

    integer, intent(in) :: M

    The number of target points for the interpolant

    integer, intent(in) :: targetNodeType

    The integer code specifying the type of target points. Parameters are defined in SELF_Constants.f90. One of GAUSS(=1), GAUSS_LOBATTO(=2), or UNIFORM(=3)

procedure, public :: WriteHDF5 => WriteHDF5_Lagrange_t

  • public subroutine WriteHDF5_Lagrange_t(this, fileId)

    Arguments

    TypeIntentOptionalAttributesName
    class(Lagrange_t), intent(in) :: this
    integer(kind=HID_T), intent(in) :: fileId

Source Code

  type,extends(Lagrange_t),public :: Lagrange
    character(3) :: backend = 'gpu'
    type(c_ptr) :: qWeights_gpu
    type(c_ptr) :: iMatrix_gpu
    type(c_ptr) :: dMatrix_gpu
    type(c_ptr) :: dgMatrix_gpu
    type(c_ptr) :: bMatrix_gpu

  contains

    procedure,public :: Init => Init_Lagrange
    procedure,public :: Free => Free_Lagrange

  endtype Lagrange