ply_poly_transformation_module Module


Uses

  • module~~ply_poly_transformation_module~~UsesGraph module~ply_poly_transformation_module ply_poly_transformation_module module~treelmesh_module treelmesh_module module~ply_poly_transformation_module->module~treelmesh_module module~env_module env_module module~ply_poly_transformation_module->module~env_module

Used by

  • module~~ply_poly_transformation_module~~UsedByGraph module~ply_poly_transformation_module ply_poly_transformation_module module~ply_sampling_module ply_sampling_module module~ply_sampling_module->module~ply_poly_transformation_module module~ply_sampled_tracking_module ply_sampled_tracking_module module~ply_sampled_tracking_module->module~ply_sampling_module program~sdr_harvesting sdr_harvesting program~sdr_harvesting->module~ply_sampled_tracking_module module~sdr_hvs_config_module sdr_hvs_config_module program~sdr_harvesting->module~sdr_hvs_config_module module~sdr_hvs_config_module->module~ply_sampled_tracking_module

Contents


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 :: 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(:)

Functions

private function ply_alpha(mode) result(alpha)

Coefficients from the recursive formulation of legendre polynomials. L_n = alpha * x * L_n-1 + beta * L_n-2

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: mode

The current mode in the polynomial representation.

Return Value real(kind=rk)

Alpha coefficient from the recursive formulation of legendre polynomials.

private function ply_beta(mode) result(beta)

Coefficients from the recursive formulation of legendre polynomials. L_n = alpha * x * L_n-1 + beta * L_n-2

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: mode

The current mode in the polynomial representation.

Return Value real(kind=rk)

Beta coefficient from the recursive formulation of legendre polynomials.

private function ply_alpha_frac(denominator, numerator) result(alpha_frac)

Quotient of two alpha values.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: denominator

Denominator

integer, intent(in) :: numerator

Numerator

Return Value real(kind=rk)

The quotient of two alpha values.

private function ply_alpha_beta(denominator, numerator) result(alpha_beta)

Prodcut of alpha(numerator) * beta(denominator) / alpha(denominator)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: denominator

Denominator

integer, intent(in) :: numerator

Numerator

Return Value real(kind=rk)

The product of alpha(n) * beta(d) / alpha(d)


Subroutines

public subroutine ply_Poly_Transformation(subsamp, dofReduction, mesh, meshData, varDofs, varComps, ndims, refine_tree, new_refine_tree, newMeshData, newVarDofs)

Projection of polynomial data from parent elements to child elements. The projection is done by a direct transformation of the modal coeffiecients to another coordinate system with z=ax+b.

Arguments

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

Parameters for the subsampling

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

Factor for reduction of degrees of freedom.

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

The mesh related to meshData.

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

The data for subsampling.

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

The number of degrees of freedom for every variable.

integer, intent(in) :: varComps(:)

The number of components for every variable.

integer, intent(in) :: ndims

Number of dimensions in the polynomial representation.

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

Logical array that marks elements for refinement of the previous sampling level.

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

Logical array that marks elements for refinement.

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

The subsampled data for new_refine_tree.

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

The number of dofs for the subsampled data.

private subroutine ply_subsampleData(mesh, meshData, nDofs, nChildDofs, nComponents, refine_tree, new_refine_tree, nDims, subsamp, newMeshData)

Arguments

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

The mesh for the data.

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

The data to subsample

integer, intent(in) :: nDofs

The number of degrees of freedom.

integer, intent(in) :: nChildDofs

The number of degrees of freedom for the child elements.

integer, intent(in) :: nComponents

Number of Components.

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

Logical array that marks all elements for refinement for the previous sampling level.

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

Logical array that marks all elements for refinement for the current sampling level.

integer, intent(in) :: nDims

The number of dimensions in the polynomial representation.

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

Parameters for subsampling.

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

The subsampled Data.

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

Subroutine to project element data from a parent cell to its children.

Arguments

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

The polynomial data for a single parent element.

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 for the given variable.

integer, intent(in) :: nDimensions

The number of dimensions.

integer, intent(in) :: nChilds

The number of child elements.

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

The transformation matrix for the linear coordinate transformation.

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

The new data representation for all child cell of the parent cell.

private subroutine ply_transform_matrix(max_modes, v)

Compute the transformation matrix for a projection to the left and right half-interval of Legendre polynomials for the given maximal number of modes.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: max_modes

The maximal number of modes to compute the transformation for.

The resulting matrix v will be max_modes x max_modes large and can be used for the transformation of all polynomials with up to this many modes.

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

The transformation matrix.

Upper triangular matrix is created for shifting and lower triangular for (-1) * shifting. For the right interval we interpret the first index as row index and the second as column. For the left interval this is reverted and we interpret the first index as columns of the matrix.