ply_sampling_adaptive_module Module

Adaptive sampling of polynomial data.

This module implements the sampling of polynomials with data dependent refinement. Elements, where the polynomials vary above a certain threshold, will be split into their eight children and the polynomial data will be projected onto those. The polynomials in the children can be restricted in their order to limit the memory consumption. In the end only one degree of freedom will be returned for each (refined) element, these are always the mean of the solution in those (refined) elements.

The module provides a data type to describe the configuration of the adaptive sampling: ply_sampling_adaptive_type, one routine to load this configuration ply_sampling_adaptive_load and one routine to actually perform the adaptive sampling ply_sample_adaptive.


Uses

  • module~~ply_sampling_adaptive_module~~UsesGraph module~ply_sampling_adaptive_module ply_sampling_adaptive_module module~tem_refining_module tem_refining_module module~ply_sampling_adaptive_module->module~tem_refining_module module~env_module env_module module~ply_sampling_adaptive_module->module~env_module module~aotus_module aotus_module module~ply_sampling_adaptive_module->module~aotus_module module~ply_split_element_module ply_split_element_module module~ply_sampling_adaptive_module->module~ply_split_element_module iso_c_binding iso_c_binding module~ply_sampling_adaptive_module->iso_c_binding module~ply_sampling_varsys_module ply_sampling_varsys_module module~ply_sampling_adaptive_module->module~ply_sampling_varsys_module mpi mpi module~ply_sampling_adaptive_module->mpi module~tem_topology_module tem_topology_module module~ply_sampling_adaptive_module->module~tem_topology_module module~tem_varsys_module tem_varSys_module module~ply_sampling_adaptive_module->module~tem_varsys_module module~tem_bc_prop_module tem_bc_prop_module module~ply_sampling_adaptive_module->module~tem_bc_prop_module module~tem_subtree_type_module tem_subTree_type_module module~ply_sampling_adaptive_module->module~tem_subtree_type_module module~tem_tracking_module tem_tracking_module module~ply_sampling_adaptive_module->module~tem_tracking_module module~treelmesh_module treelmesh_module module~ply_sampling_adaptive_module->module~treelmesh_module module~tem_aux_module tem_aux_module module~ply_sampling_adaptive_module->module~tem_aux_module module~ply_filter_element_module ply_filter_element_module module~ply_sampling_adaptive_module->module~ply_filter_element_module module~tem_time_module tem_time_module module~ply_sampling_adaptive_module->module~tem_time_module module~tem_subtree_module tem_subTree_module module~ply_sampling_adaptive_module->module~tem_subtree_module module~tem_tools_module tem_tools_module module~ply_sampling_adaptive_module->module~tem_tools_module module~tem_logging_module tem_logging_module module~ply_sampling_adaptive_module->module~tem_logging_module module~ply_split_element_module->module~env_module module~ply_split_legendre_module ply_split_legendre_module module~ply_split_element_module->module~ply_split_legendre_module module~ply_modg_basis_module ply_modg_basis_module module~ply_split_element_module->module~ply_modg_basis_module module~ply_sampling_varsys_module->module~env_module module~ply_sampling_varsys_module->module~tem_topology_module module~ply_sampling_varsys_module->module~tem_varsys_module module~ply_sampling_varsys_module->module~tem_tracking_module module~ply_sampling_varsys_module->module~treelmesh_module module~ply_sampling_varsys_module->module~tem_time_module module~ply_filter_element_module->module~env_module module~ply_filter_element_module->module~aotus_module module~ply_filter_element_module->module~tem_aux_module module~ply_filter_element_module->module~tem_tools_module module~ply_filter_element_module->module~tem_logging_module module~aot_table_module aot_table_module module~ply_filter_element_module->module~aot_table_module module~aot_err_module aot_err_module module~ply_filter_element_module->module~aot_err_module module~ply_split_legendre_module->module~env_module module~ply_modg_basis_module->module~env_module module~ply_dof_module ply_dof_module module~ply_modg_basis_module->module~ply_dof_module module~ply_space_integration_module ply_space_integration_module module~ply_modg_basis_module->module~ply_space_integration_module module~ply_dof_module->module~env_module module~ply_space_integration_module->module~env_module module~tem_param_module tem_param_module module~ply_space_integration_module->module~tem_param_module

Used by

  • module~~ply_sampling_adaptive_module~~UsedByGraph module~ply_sampling_adaptive_module ply_sampling_adaptive_module module~ply_sampling_module ply_sampling_module module~ply_sampling_module->module~ply_sampling_adaptive_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


Variables

TypeVisibilityAttributesNameInitial
integer, private, parameter:: redux_factor =1

Constant to indicate the factor reduction mode.

integer, private, parameter:: redux_decrement =2

Constant to indicate the decrement reduction mode.


Derived Types

Configuration of the adaptive sampling.

Read more…

Components

TypeVisibilityAttributesNameInitial
integer, private :: max_nlevels =0

Maximal number of levels by which any mesh element should be refined.

Read more…
real(kind=rk), private :: eps_osci

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

integer, private :: reduction_mode

Method to use for the reduction.

Read more…
logical, private :: ignore_highmodes =.false.

Indication whether to filter modes during refinement by ignoring all modes in the parent, that exceed the target polynomial degree of the child elements.

Read more…
integer, private :: dof_decrement =1

Number of modes to cut off in each refinement.

Read more…
real(kind=rk), private :: dofReducFactor

Factor to Reduce dofs for every sampling level. Can be used to avoid too drastic increase of memory consumption. For adaptive subsampling only.

logical, private :: adaptiveDofReduction

Indicator for the limitation of memory consumption.

integer, private :: AbsUpperBoundLevel

Absolute upper bound level to refine to.

type(ply_filter_element_type), private :: filter_element

Filtering the poylnomial modes during adaptive refinement.

Read more…

type, private :: realarray_type

Small helping type to allow arrays of arrays for the variable data.

Components

TypeVisibilityAttributesNameInitial
real(kind=rk), private, pointer:: dat(:)=> NULL()

type, private :: sampled_method_data_type

A container for the method data to hold the data in a scalar pointer for the C-pointer conversion.

Components

TypeVisibilityAttributesNameInitial
type(realarray_type), private, allocatable:: component(:)

Subroutines

public subroutine ply_sampling_adaptive_load(me, conf, parent)

Load the configuration for adaptive subsampling.

Arguments

TypeIntentOptionalAttributesName
type(ply_sampling_adaptive_type), intent(out) :: me

Sampling definition to load.

type(flu_State), intent(in) :: conf

Configuration to read the sampling settings from.

integer, intent(in), optional :: parent

Parent table in which to look for the adaptive sampling settings.

public subroutine ply_sample_adaptive(me, ndims, orig_mesh, orig_bcs, varsys, var_degree, lvl_degree, trackInst, trackConfig, time, new_mesh, resvars)

Sample data described by varsys in orig_mesh according to the tracking object trackInst with adaptive refinements.

Read more…

Arguments

TypeIntentOptionalAttributesName
type(ply_sampling_adaptive_type), intent(in) :: me

A ply_sampling_type to describe the sampling method.

integer, intent(in) :: ndims

Number of dimensions in the polynomial representation.

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

The original mesh to be refined.

type(tem_BC_prop_type), intent(in) :: orig_bcs

Boundary conditions for the original mesh.

type(tem_varSys_type), intent(in) :: varsys

Variable system of the original data to do the sampling on.

integer, intent(in) :: var_degree(:)

Maximal polynomial degree for each variable.

Needs to be matching the variable definition in the variable system.

integer, intent(in) :: lvl_degree(:)

Maximal polynomial degree for each level.

type(tem_tracking_instance_type), intent(in) :: trackInst

Tracking object describing what to sample.

type(tem_tracking_config_type), intent(in) :: trackConfig

Tracking configuration with the geometry to obtain from the overall mesh.

type(tem_time_type), intent(in) :: time

Point in time to get the data for.

type(treelmesh_type), intent(out) :: new_mesh

The new mesh with the refined elements.

type(tem_varSys_type), intent(out) :: resvars

Resulting system of variables describing the data in the arrays of subsampled elements.

private subroutine get_sampled_element(fun, varSys, elempos, time, tree, n, nDofs, res)

Get sampled data.

Read more…

Arguments

TypeIntentOptionalAttributesName
class(tem_varSys_op_type), intent(in) :: fun

Description of the method to obtain the variables, here some preset values might be stored, like the space time function to use or the required variables.

type(tem_varSys_type), intent(in) :: varSys

The variable system to obtain the variable from.

integer, intent(in) :: elempos(:)

TreeID of the element to get the variable for.

type(tem_time_type), intent(in) :: time

Point in time at which to evaluate the variable.

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

global treelm mesh info

integer, intent(in) :: n

Number of elements to obtain for this variable (vectorized access).

integer, intent(in) :: nDofs

Number of degrees of freedom within an element.

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

Resulting values for the requested variable.

Linearized array dimension: (nComponents of resulting variable) x (nDegrees of freedom) x (nElems) Access: (iElem-1)fun%nComponentsnDofs + (iDof-1)*fun%nComponents + iComp