ply_LegPolyProjection_module Module

Module for projection of Q Legendre Polynomials from parent cell to child cells.


Uses

  • module~~ply_legpolyprojection_module~~UsesGraph module~ply_legpolyprojection_module ply_LegPolyProjection_module module~treelmesh_module treelmesh_module module~ply_legpolyprojection_module->module~treelmesh_module module~tem_aux_module tem_aux_module module~ply_legpolyprojection_module->module~tem_aux_module module~env_module env_module module~ply_legpolyprojection_module->module~env_module module~tem_param_module tem_param_module module~ply_legpolyprojection_module->module~tem_param_module

Contents


Variables

TypeVisibilityAttributesNameInitial
integer, private, parameter:: ply_QLegendrePoly_prp =1

Parameter to specify Legendre polynomials as the degrees of freedoms of the elements. The multidimensional polynomias are build as Q-polynomials. The projection is a L2-Projection onto the ansatz function of the finer elements.


Derived Types

type, public :: ply_subsample_type

Components

TypeVisibilityAttributesNameInitial
logical, private :: isActive =.false.

Is subsampling active

integer, private :: sampling_lvl

The current sampling lvl.

integer, private :: projectionType =ply_QLegendrePoly_prp

The type of projection we use to subsample the elemental data.

integer, private :: caplevel =20

Maximal Level down to which subsampling should be done.

integer, private :: minsub =0

Minimal subsampling depth:

integer, private :: maxsub =0

Maximal subsampling depth:

real(kind=rk), private :: eps_osci

Maximum allowed oscillation of the solution. For adaptive subsampling only.

real(kind=rk), private :: dofReducFactor

Factor for the reduction of the degrees of freedom in one subsampling step (per spatial direction).

logical, private :: adaptiveDofReduction

Indicator for limitation of total memory consumption

integer, private :: AbsUpperBoundLevel

Absolute upper bound level to refine to.

type, public :: ply_array_type

Components

TypeVisibilityAttributesNameInitial
real(kind=rk), private, allocatable:: dat(:)

type, private :: ply_ProjCoeff_type

Datatype storing the coefficients arising for the projection of solutions on a parent cell to its children during the subsampling routines.

Components

TypeVisibilityAttributesNameInitial
real(kind=rk), private, allocatable:: projCoeff(:,:,:)

Array holding all the projection coefficients for the projection of a degree of freedom on the parent element (first index) to a degree of freedom on the child element (second index) for a given child element (third index). Therefore the dimension of this array is (nDofs, nDofs, 8).


Functions

private function ply_QLegOneDimCoeff(nDofsOneDim, nChildDofsOneDim) result(projCoeffOneDim)

Routine to create one-dimensional projection coefficient for a coarse element to a fine element.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nDofsOneDim

The number of dofs in one dimension.

integer, intent(in) :: nChildDofsOneDim

The number of dofs in one dimension for the children.

Return Value real(kind=rk),allocatable, (:,:,:)

Projected one-dimensional coefficients.

First index is Legendre polynomial on the parent element, second index is the Legendre polynomial on the child element, third index is left or right projection.

private function ply_QLegSqNorm(polyIndex) result(sqNorm)

Function to calculate the squared L2-Norm of a given Legendre polynomial on the reference element [-1,+1].

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: polyIndex

The Legendre polynomial index to calculate the squared norm for. The first polynomial has index 1.

Return Value real(kind=rk)

The squared L2 Norm of the Legendre polynomial.

private function ply_legVal(points, nPoints, maxPolyDegree) result(val)

Evaluate a given set of Legendre polynomials a given set of 1D points.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: points(:)

A given set of 1D points.

integer, intent(in) :: nPoints

The number of points to evaluate the polynomials at.

integer, intent(in) :: maxPolyDegree

The maximal polynomial degree to evaluate for.

Return Value real(kind=rk),allocatable, (:,:)

Function values for for all Legendre polynomials up to degree maxPolyDegree at all given points. Therefore the dimension of this array is (maxPolyDegree+1, nPoints)


Subroutines

public subroutine ply_QPolyProjection(subsamp, dofReduction, tree, meshData, varDofs, ndims, varcomps, refine_tree, new_refine_tree, newMeshData, newVarDofs)

Subsampling by L2-Projection of the Q-Tensorproduct Legendre polynomials.

Arguments

TypeIntentOptionalAttributesName
type(ply_subsample_type), intent(in) :: subsamp

Parameters for the subsampling.

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

Factor for reducion of degrees of freedom.

type(treelmesh_type), intent(in) :: tree

The tree the data is written for.

type(ply_array_type), intent(in) :: meshData(:)

The data to sub-sample.

integer, intent(in) :: varDofs(:)

The number of degrees of freedom for each scalar variable.

integer, intent(in) :: ndims

Number of dimensions in the polynomial representation.

integer, intent(in) :: varcomps(:)

Number of components

logical, intent(in) :: refine_tree(:)

Logical array that marks elements for refinement from the last sampling

logical, intent(in) :: new_refine_tree(:)

Logical array that marks elements for refinment.

type(ply_array_type), intent(out), allocatable:: newMeshData(:)

The subsampled data.

integer, intent(out), allocatable:: newVarDofs(:)

The number of dofs for the subsampled dofs.

public subroutine ply_gauleg(x1, x2, x, w, nIntP)

subroutine to create gauss points and weights for one-dimensional integration on the interval [x1,x2].

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: x1

lower limit of integration interval

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

upper limit of integration interval

real(kind=rk), intent(inout), allocatable:: x(:)

The coordinates of the gauss points on the interval [-1,1]. The array has the length nIntP.

real(kind=rk), intent(inout), allocatable:: w(:)

The quadrature weights. The array has the length nIntP.

integer, intent(in) :: nIntP

The number of integration points.

private subroutine ply_initQLegProjCoeff(dofType, nDofs, ndims, nChilds, nChildDofs, projection)

Routine to initialize the projection coefficients for a usage in the subsampling routine to project degrees of freedoms of a parent cell to the degrees of freedoms of a child cell if the degrees of freedoms are Q-Legendre polynomials.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: dofType

The type of degrees of freedom we have in our cells.

integer, intent(in) :: nDofs

The number of degrees of freedom for the parent cells.

integer, intent(in) :: ndims

The number of dimensions in the polynomial representation.

integer, intent(in) :: nChilds

The number of child cells.

integer, intent(in) :: nChildDofs

The number of degrees of freedom for the child cells.

type(ply_ProjCoeff_type), intent(out) :: projection

The subsampling coefficients that will be initialized by this routine.

private subroutine ply_dofToQPoly(dof, nDofs, ndims, xAnsFunc, yAnsFunc, zAnsFunc)

Subroutine to convert linearized dof index to ansatz function number for Q-Polynomials.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: dof

The linearized degree of freedom index.

integer, intent(in) :: nDofs

The number of dofs for all directions.

integer, intent(in) :: ndims

The number of Dimensions in the polynomial representation.

integer, intent(out) :: xAnsFunc

The ansatz function number in x direction.

integer, intent(out) :: yAnsFunc

The ansatz function number in y direction.

integer, intent(out) :: zAnsFunc

The ansatz function number in z direction.

private subroutine ply_subsampleData(tree, meshData, nDofs, nChildDofs, nComponents, projection, projection_oneDof, refine_tree, new_refine_tree, ndims, subsamp, newMeshData)

Routine to subsample mesh information for one refinement level.

Arguments

TypeIntentOptionalAttributesName
type(treelmesh_type), intent(in) :: tree

The tree the data is written for.

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

The data to sub-sample.

integer, intent(in) :: nDofs

The number of degrees of freedom for each scalar variable.

integer, intent(in) :: nChildDofs

The number of degrees of freedom per scalar variable on the child elements.

integer, intent(in) :: nComponents

Number of components

type(ply_ProjCoeff_type), intent(in) :: projection

Projection coefficients for the given data.

type(ply_ProjCoeff_type), intent(in), optional :: projection_oneDof

Projection coeffiecients for the the reduction to polynomial degree of 0.

logical, intent(in), optional :: refine_tree(:)

Logical array that marks all elements for refinement from the last sampling lvl.

logical, intent(in) :: new_refine_tree(:)

Logical array that marks all elements for refinement

integer, intent(in) :: ndims

The number of dimensions in the polynomial representation.

type(ply_subsample_type), intent(in) :: subsamp

Parameters for the subsampling.

real(kind=rk), intent(out), allocatable:: newMeshData(:)

The subsampled data.

private subroutine ply_projDataToChild(parentData, nParentDofs, nChildDofs, nComponents, nChilds, projection, childData)

Subroutine to project elemental data from a parent cell to one of its children.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: parentData(:)

Linearized data for a single variable (can have multiple components) and a single degree of freedom of the parent cell.

integer, intent(in) :: nParentDofs

The number of dofs of the parent element.

integer, intent(in) :: nChildDofs

The total number of dofs for the child cells.

integer, intent(in) :: nComponents

The number of componentns of the given variable.

integer, intent(in) :: nChilds

The number of children.

type(ply_ProjCoeff_type), intent(in) :: projection

The information about the projection coefficients for the parent dofs to the child dofs.

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

The created childData.