DG_implementation Submodule

The implementation of the interfaces defined in DG_structures as well as internal routines that should never be called by a user.


Uses

  • module~~dg_implementation~~UsesGraph module~dg_implementation DG_implementation module~dg_structures DG_structures module~dg_implementation->module~dg_structures module~kinds kinds module~dg_structures->module~kinds

Contents


Functions

function JacobiP(x, alpha, beta, n)

Function to evaluate Jacobi Polynomial of type (with ) at points for order .

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:):: x

The points at which to evaluate the Jacobi Polynomial.

real(kind=wp), intent(in) :: alpha

The value of .

real(kind=wp), intent(in) :: beta

The value of .

integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomial.

Return Value real(kind=wp), dimension(size(x))

Returns .

function Vandermonde1D(n, x)

Initialize the 1D Vandermonde matrix, for the Legendre-Gauss-Lobatto quadrature points.

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomials.

real(kind=wp), intent(in), dimension(:):: x

The points at which to evaluate the Vandermonde matrix.

Return Value real(kind=wp), dimension(size(x),n+1)

The return value is the Vandermonde matrix,

function GradJacobiP(x, alpha, beta, n)

Evaluate the derivative of the Jacobi polynomial of type at points for order .

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:):: x

The points at which to evaluate the derivative of the Jacobi polynomial.

real(kind=wp), intent(in) :: alpha

The value of .

real(kind=wp), intent(in) :: beta

The value of .

integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomials.

Return Value real(kind=wp), dimension(size(x))

The return value is the derivative of the Jacobi polynomial evaluated at ,

function GradVandermonde1D(n, x)

Initialize the gradient of the modal basis at .

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomials.

real(kind=wp), intent(in), dimension(:):: x

The points at which to initialize the gradient of the modal basis.

Return Value real(kind=wp), dimension(size(x),n+1)

The return value is the derivative of the modal basis, .

function Dmatrix1D(n, x, v)

Initialize the differentiation matrix for a reference element of order .

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomials.

real(kind=wp), intent(in), dimension(n+1):: x

The nodal points where the differentiation matrix will provide approximations to the derivative.

real(kind=wp), intent(in), dimension(n+1,n+1):: v

The Vandermonde matrix, .

Return Value real(kind=wp), dimension(n+1,n+1)

The return value is the differentiation matrix, .

function Lift1D(n, v)

Initialize the lift matrix, , used to compute surface integral terms in the Discontinuous Galerkin formulation.

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomials.

real(kind=wp), intent(in), dimension(n+1,n+1):: v

The Vandermonde matrix, .

Return Value real(kind=wp), dimension(n+1,2)

The return value is the lift matrix, , where is a array of the form

function Inverse(a)

Helper function that calculates the inverse of a matrix.

Read more…

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: a

A matrix, .

Return Value real(kind=wp), dimension(size(a,1),size(a,2))

On output contains .

function Filter1D(n, nc, s, v, invv)

Create an exponential filter matrix that can be used to filter out high-frequency noise.

Read more…

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: n

The order of the element.

integer(kind=ip), intent(in) :: nc

The cutoff, below which the low modes are left untouched.

integer(kind=ip), intent(in) :: s

The order of the filter.

real(kind=wp), intent(in), dimension(:,:):: v

The Vandermonde matrix, .

real(kind=wp), intent(in), dimension(:,:):: invv

The inverse of the Vandermonde matric, .

Return Value real(kind=wp), dimension(n+1,n+1)

The return value is the filter matrix, .


Subroutines

subroutine JacobiGQ(alpha, beta, n, x, w)

Compute the n'th order Gauss quadrature points, and weights, , associated with the Jacobi polynomial of type .

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in) :: alpha

The value of .

real(kind=wp), intent(in) :: beta

The value of .

integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomial.

real(kind=wp), intent(out), dimension(n+1):: x

On output the location of the Gauss quadrature points.

real(kind=wp), intent(out), dimension(n+1):: w

On output the weights.

subroutine JacobiGL(alpha, beta, n, r)

Compute the n'th order Gauss Lobatto quadrature points, , associated with the Jacobi polynomial of type .

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in) :: alpha

The value of .

real(kind=wp), intent(in) :: beta

The value of .

integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomial.

real(kind=wp), intent(out), dimension(n+1):: r

On output the location of the Gauss Lobatto quadrature points.

subroutine GaussWeigths(alpha, beta, n, r, w)

Routine to calculate the integration weights for the Gauss-Lobatto Quadrature points.

Read more…

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in) :: alpha

The value of .

real(kind=wp), intent(in) :: beta

The value of .

integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomials.

real(kind=wp), intent(in), dimension(:):: r

The points for which the weights should be determined.

real(kind=wp), intent(out), dimension(:):: w

On return the integration weights.

subroutine Jacobian(n, vx, r, dr, x, rx)

Calculate the Jacobian for transforming derivatives from the reference element to the physical element.

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: n

The order of the element.

real(kind=wp), intent(in), dimension(2):: vx

A 1d real array of size 2 that contains the physical coordinates at the boundaries of the element.`

real(kind=wp), intent(in), dimension(n+1):: r

The node locations within the reference element, .

real(kind=wp), intent(in), dimension(:,:):: dr

The derivative matrix for the reference element, .

real(kind=wp), intent(out), dimension(n+1):: x

On output contains the physical coordinates for this physical element, .

real(kind=wp), intent(out), dimension(n+1):: rx

On output contains (\frac{d r_i}{d x_j}.


Module Procedures

module procedure init_ref_element module function init_ref_element(order, sorder, nc)

Arguments

Type IntentOptional AttributesName
integer(kind=ip), intent(in) :: order

The order of the reference element.

integer(kind=ip), optional :: sorder

The order of a smoothing filter. Both this and nc have to be present before a filter is constructed.

integer(kind=ip), optional :: nc

No filtering of orders below this value. Both this and sorder have to be present before a filter is constructed.

Return Value type(ref_element)

The return type has to be a reference element.

module procedure deallocate_ref_element module subroutine deallocate_ref_element(refel)

Arguments

Type IntentOptional AttributesName
type(ref_element) :: refel

The argument has to be a reference element type.

module procedure char_flux_real module function char_flux_real(this, nvar, order, uint, uext, flux, lambda, s, sinv, debug_output)

Arguments

Type IntentOptional AttributesName
class(ref_element), intent(in) :: this

Has to be a reference element type.

integer(kind=ip), intent(in) :: nvar

The number of variables the characteristic flux has to be computed for.

integer(kind=ip), intent(in) :: order

The order of the reference element. Used to declare the size of the return array.

real(kind=wp), intent(in), dimension(2,nvar):: uint

The boundary data internal to this element.

real(kind=wp), intent(in), dimension(2,nvar):: uext

The boundary data external to this element.

real(kind=wp), intent(in), dimension(2,nvar):: flux

The boundary fluxes.

real(kind=wp), intent(in), dimension(2,nvar):: lambda

The boundary characteristic speeds.

real(kind=wp), intent(in), dimension(2,nvar,nvar):: s

The matrix that converts from characteristic to evolved variables.

real(kind=wp), intent(in), dimension(2,nvar,nvar):: sinv

The matrix that converts from evolved to characteristic variables.

logical, intent(in) :: debug_output

Currently ignored for the real version of this routine.

Return Value real(kind=wp), dimension(order+1,nvar)

At return the numerical flux has been lifted to the whole element.

module procedure char_flux_complex module function char_flux_complex(this, nvar, order, uint, uext, flux, lambda, s, sinv, debug_output)

Arguments

Type IntentOptional AttributesName
class(ref_element), intent(in) :: this

Has to be a reference element type.

integer(kind=ip), intent(in) :: nvar

The number of variables the characteristic flux has to be computed for.

integer(kind=ip), intent(in) :: order

The order of the reference element. Used to declare the size of the return array.

complex(kind=wp), intent(in), dimension(2,nvar):: uint

The boundary data internal to this element.

complex(kind=wp), intent(in), dimension(2,nvar):: uext

The boundary data external to this element.

complex(kind=wp), intent(in), dimension(2,nvar):: flux

The boundary fluxes.

real(kind=wp), intent(in), dimension(2,nvar):: lambda

The boundary characteristic speeds.

real(kind=wp), intent(in), dimension(2,nvar,nvar):: s

The matrix that converts from characteristic to evolved variables.

real(kind=wp), intent(in), dimension(2,nvar,nvar):: sinv

The matrix that converts from evolved to characteristic variables.

logical, intent(in) :: debug_output

If .true., produce debug output

Return Value complex(kind=wp), dimension(order+1,nvar)

At return the numerical flux has been lifted to the whole element.