preCICE v3.2.0
Loading...
Searching...
No Matches
precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T > Class Template Reference

#include <PartitionOfUnityMapping.hpp>

Inheritance diagram for precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >:
[legend]
Collaboration diagram for precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >:
[legend]

Public Member Functions

 PartitionOfUnityMapping (Mapping::Constraint constraint, int dimension, RADIAL_BASIS_FUNCTION_T function, Polynomial polynomial, unsigned int verticesPerCluster, double relativeOverlap, bool projectToInput)
void computeMapping () final override
void clear () final override
 Clears a computed mapping by deleting the content of the _clusters vector.
void tagMeshFirstRound () final override
 tag the vertices required for the mapping
void tagMeshSecondRound () final override
 nothing to do here
std::string getName () const final override
 name of the pum mapping
void mapConsistentAt (const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > values) final override
 Just-in-time mapping variant of mapConsistent.
void mapConservativeAt (const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const Eigen::Ref< const Eigen::MatrixXd > &source, impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > target) final override
 the target values here remain unused, as we store the (intermediate) result directly in the cache
void updateMappingDataCache (impl::MappingDataCache &cache, const Eigen::Ref< const Eigen::VectorXd > &in) final override
 Allows updating a so-called MappingDataCache for more efficient just-in-time mappings.
void completeJustInTimeMapping (impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > buffer) final override
 Completes a just-in-time mapping for conservative constraints.
void initializeMappingDataCache (impl::MappingDataCache &cache) final override
 Allocates memory and sets up the data structures inside the MappingDataCache.
Public Member Functions inherited from precice::mapping::Mapping
 Mapping (Constraint constraint, int dimensions, bool requiresGradientData, InitialGuessRequirement initialGuessRequirement)
 Constructor, takes mapping constraint.
Mappingoperator= (Mapping &&)=delete
virtual ~Mapping ()=default
 Destructor, empty.
void setMeshes (const mesh::PtrMesh &input, const mesh::PtrMesh &output)
 Sets input and output meshes carrying data to be mapped.
const mesh::PtrMeshgetInputMesh () const
const mesh::PtrMeshgetOutputMesh () const
Constraint getConstraint () const
 Returns the constraint (consistent/conservative) of the mapping.
MeshRequirement getInputRequirement () const
 Returns the requirement on the input mesh.
MeshRequirement getOutputRequirement () const
 Returns the requirement on the output mesh.
bool hasComputedMapping () const
 Returns true, if the mapping has been computed.
virtual bool hasConstraint (const Constraint &constraint) const
 Checks whether the mapping has the given constraint or not.
bool isScaledConsistent () const
 Returns true if mapping is a form of scaled consistent mapping.
bool requiresInitialGuess () const
 Return true if the mapping requires an initial guess.
bool isJustInTimeMapping () const
const Eigen::VectorXd & initialGuess () const
 Return the provided initial guess of a mapping using an initialGuess.
Eigen::VectorXd & initialGuess ()
bool hasInitialGuess () const
 True if initialGuess().size() == 0.
void map (int inputDataID, int outputDataID)
void map (int inputDataID, int outputDataID, Eigen::VectorXd &initialGuess)
void map (const time::Sample &input, Eigen::VectorXd &output)
 Maps an input Sample to output data from input mesh to output mesh.
void map (const time::Sample &input, Eigen::VectorXd &output, Eigen::VectorXd &initialGuess)
 Maps an input Sample to output data from input mesh to output mesh, given an initialGuess.
virtual void scaleConsistentMapping (const Eigen::VectorXd &input, Eigen::VectorXd &output, Constraint type) const
 Scales the consistently mapped output data such that the surface integral of the values on input mesh and output mesh are equal.
bool requiresGradientData () const
 Returns whether the mapping requires gradient data.

Private Member Functions

void mapConservative (const time::Sample &inData, Eigen::VectorXd &outData) override
 Maps data using a conservative constraint.
void mapConsistent (const time::Sample &inData, Eigen::VectorXd &outData) override
 Maps data using a consistent constraint.
std::pair< std::vector< int >, std::vector< double > > computeNormalizedWeight (const mesh::Vertex &v, std::string_view mesh)
void exportClusterCentersAsVTU (mesh::Mesh &centers)

Private Attributes

precice::logging::Logger _log {"mapping::PartitionOfUnityMapping"}
 logger, as usual
std::vector< SphericalVertexCluster< RADIAL_BASIS_FUNCTION_T > > _clusters
 main data container storing all the clusters, which need to be solved individually
RADIAL_BASIS_FUNCTION_T _basisFunction
 Radial basis function type used in interpolation.
const unsigned int _verticesPerCluster
 Input parameters provided by the user for the clustering algorithm:
const double _relativeOverlap
 overlap of vertex clusters
const bool _projectToInput
 toggles whether we project the cluster centers to the input mesh
double _clusterRadius = 0
 derived parameter based on the input above: the radius of each cluster
Polynomial _polynomial
 polynomial treatment of the RBF system
std::unique_ptr< mesh::Mesh_centerMesh

Additional Inherited Members

Public Types inherited from precice::mapping::Mapping
enum  Constraint { CONSISTENT , CONSERVATIVE , SCALED_CONSISTENT_SURFACE , SCALED_CONSISTENT_VOLUME }
 Specifies additional constraints for a mapping. More...
enum class  MeshRequirement { UNDEFINED = 0 , VERTEX = 1 , FULL = 2 }
 Specifies requirements for the input and output meshes of a mapping. More...
enum class  InitialGuessRequirement : bool { Required = true , None = false }
 Specifies whether the mapping requires an initial guess. More...
Protected Member Functions inherited from precice::mapping::Mapping
mesh::PtrMesh input () const
 Returns pointer to input mesh.
mesh::PtrMesh output () const
 Returns pointer to output mesh.
void setInputRequirement (MeshRequirement requirement)
 Sets the mesh requirement for the input mesh.
void setOutputRequirement (MeshRequirement requirement)
 Sets the mesh requirement for the output mesh.
int getDimensions () const
Protected Attributes inherited from precice::mapping::Mapping
bool _hasComputedMapping = false
 Flag to indicate whether computeMapping() has been called.
bool _requiresGradientData
 Flag if gradient data is required for the mapping.

Detailed Description

template<typename RADIAL_BASIS_FUNCTION_T>
class precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >

Mapping using partition of unity decomposition strategies: The class here inherits from the Mapping class and orchestrates the partitions (called vertex clusters) in order to represent a partition of unity. This means in particular that the class computes the weights for the evaluation vertices and the necessary association between evaluation vertices and the clusters during initialization and traverses through all vertex clusters when evaluating the mapping.

Definition at line 30 of file PartitionOfUnityMapping.hpp.

Constructor & Destructor Documentation

◆ PartitionOfUnityMapping()

template<typename RADIAL_BASIS_FUNCTION_T>
precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::PartitionOfUnityMapping ( Mapping::Constraint constraint,
int dimension,
RADIAL_BASIS_FUNCTION_T function,
Polynomial polynomial,
unsigned int verticesPerCluster,
double relativeOverlap,
bool projectToInput )

Constructor, which mostly sets the mesh connectivity requirements and initializes member variables.

Parameters
[in]constraintSpecifies mapping to be consistent or conservative.
[in]dimensionDimensionality of the meshes
[in]functionRadial basis function type used in interpolation
[in]polynomialThe handling of the polynomial in the RBF system. Valid choices are 'off' and 'separate'
[in]verticesPerClusterTarget number of vertices to be clustered together
[in]relativeOverlapOverlap between clusters: The parameter here determines the distance between two cluster centers, given the cluster radius (already determined through verticesPerCluster ). A value of 1 would correspond to no distance between cluster centers (i.e. completely overlapping clusters), 0 to distance of 2 x radius between clusters centers.
[in]projectToInputif enabled, places the cluster centers at the closest vertex of the input mesh. See also mapping::impl::createClustering()

Definition at line 136 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

Member Function Documentation

◆ clear()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::clear ( )
finaloverridevirtual

Clears a computed mapping by deleting the content of the _clusters vector.

Implements precice::mapping::Mapping.

Definition at line 544 of file PartitionOfUnityMapping.hpp.

◆ completeJustInTimeMapping()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::completeJustInTimeMapping ( impl::MappingDataCache & cache,
Eigen::Ref< Eigen::MatrixXd > result )
finaloverridevirtual

Completes a just-in-time mapping for conservative constraints.

For conservative constraints (only implemented in write direction at the moment), we potentially buffer (partially mapped) data written by the user into the MappingDataCache. Once the user has written all data, we complete the mapping to generate the full output data. In many cases, the initial 'partial processing' or partial mapping is a computationally cheaper operation, whereas the finalization in this function call performs the computationally more demanding operations. To prevent the expensive operations from entering the API functions, we use the MappingDataCache and the additional completeJustInTimeMapping function. The function needs to be called in the ParticipantImpl before calling storeBufferedData (which adds the (mapped and written) output data to the time step storage). The function is called through the (Write)DataContext owning the just-in-time mapping.

This is not relevant for consistent constraints (only implemented in read direction at the moment), as we map before we give the data to the user, i.e., there is no completion to defer.

See also the documentation in impl::MappingDataCache for more information.

Parameters
[in]cacheThe MappingDataCache holding the buffered data we want to evaluate
[out]resultThe data vector we store the output data in.

Reimplemented from precice::mapping::Mapping.

Definition at line 363 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ computeMapping()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::computeMapping ( )
finaloverridevirtual

Computes the clustering for the partition of unity method and fills the _clusters vector, which allows to travers through all vertex cluster computed. Each vertex cluster in the vector directly computes local mapping matrices and matrix decompositions. In addition, the method computes the normalized weights (Shepard's method) for the partition of unity method and stores them directly in each relevant vertex cluster. In debug mode, the function also exports the partition centers as a separate mesh for visualization purpose.

Implements precice::mapping::Mapping.

Definition at line 162 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ computeNormalizedWeight()

template<typename RADIAL_BASIS_FUNCTION_T>
std::pair< std::vector< int >, std::vector< double > > precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::computeNormalizedWeight ( const mesh::Vertex & v,
std::string_view mesh )
private

Given an output vertex, computes the normalized PU weight at the given location The mesh name is only required for logging purposes returns the clusterIDs and all weights for these clusters

Definition at line 252 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ exportClusterCentersAsVTU()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::exportClusterCentersAsVTU ( mesh::Mesh & centers)
private

export the center vertices of all clusters as a mesh with some additional data on it such as vertex count only enabled in debug builds and mainly for debugging purpose

Definition at line 493 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ getName()

template<typename RADIAL_BASIS_FUNCTION_T>
std::string precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::getName ( ) const
finaloverridevirtual

name of the pum mapping

Implements precice::mapping::Mapping.

Definition at line 554 of file PartitionOfUnityMapping.hpp.

◆ initializeMappingDataCache()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::initializeMappingDataCache ( impl::MappingDataCache & cache)
finaloverridevirtual

Allocates memory and sets up the data structures inside the MappingDataCache.

This is only relevant for just-in-time mappings. The MappingDataCache is specific to a Data-Mapping pair and stores intermediate computatons for efficiency reasons. The initialization happens immediately after ParticipantImpl::computeMappings through the DataContext::initializeMappingDataCache, which owns the just-in-time mapping, since only then the size of relevant data structures is known.

See the documentation in impl::MappingDataCache for more information.

Parameters
cachethe cache in use, specific to the mapping-data pair

Reimplemented from precice::mapping::Mapping.

Definition at line 379 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ mapConservative()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::mapConservative ( const time::Sample & inData,
Eigen::VectorXd & outData )
overrideprivatevirtual

Maps data using a conservative constraint.

Parameters
[in]inputSample to map data from
[in]outputValues to map to

If requiresInitialGuess(), then the initial guess is available via initialGuess(). Provide a new initial guess by overwriting it. The mapping has full control over its size.

See also
For mappings requiring an initialGuess: initialGuess() hasInitialGuess()

Implements precice::mapping::Mapping.

Definition at line 306 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ mapConservativeAt()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::mapConservativeAt ( const Eigen::Ref< const Eigen::MatrixXd > & coordinates,
const Eigen::Ref< const Eigen::MatrixXd > & source,
impl::MappingDataCache & cache,
Eigen::Ref< Eigen::MatrixXd > target )
finaloverridevirtual

the target values here remain unused, as we store the (intermediate) result directly in the cache

Reimplemented from precice::mapping::Mapping.

Definition at line 336 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ mapConsistent()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::mapConsistent ( const time::Sample & inData,
Eigen::VectorXd & outData )
overrideprivatevirtual

Maps data using a consistent constraint.

Parameters
[in]inputSample to map data from
[in]outputValues to map to

If requiresInitialGuess(), then the initial guess is available via initialGuess(). Provide a new initial guess by overwriting it. The mapping has full control over its size.

See also
For mappings requiring an initialGuess: initialGuess() hasInitialGuess()

Implements precice::mapping::Mapping.

Definition at line 321 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ mapConsistentAt()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::mapConsistentAt ( const Eigen::Ref< const Eigen::MatrixXd > & coordinates,
const impl::MappingDataCache & cache,
Eigen::Ref< Eigen::MatrixXd > values )
finaloverridevirtual

Just-in-time mapping variant of mapConsistent.

Parameters
coordinates[in]where to compute the mapping
cache[in]the mapping data cache previously computed with updateMappingDataCache
values[out]data buffer passed from the user (needs to have the correct shape)
Note
the default implementation in this class simply aborts the code and actual implementations are in derived classes

Reimplemented from precice::mapping::Mapping.

Definition at line 404 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ tagMeshFirstRound()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::tagMeshFirstRound ( )
finaloverridevirtual

tag the vertices required for the mapping

Implements precice::mapping::Mapping.

Definition at line 430 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

◆ tagMeshSecondRound()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::tagMeshSecondRound ( )
finaloverridevirtual

nothing to do here

Implements precice::mapping::Mapping.

Definition at line 487 of file PartitionOfUnityMapping.hpp.

◆ updateMappingDataCache()

template<typename RADIAL_BASIS_FUNCTION_T>
void precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::updateMappingDataCache ( impl::MappingDataCache & cache,
const Eigen::Ref< const Eigen::VectorXd > & in )
finaloverridevirtual

Allows updating a so-called MappingDataCache for more efficient just-in-time mappings.

To compute a mapping just-in-time (either with mapConservativeAt or mapConsistentAt ), we first sample in time and then (using this time sample) compute an interpolant in space. Since we typically need to evaluate at a specific point in time for many different locations (potentially corresponding to many different calls of mapConservativeAt or mapConsistentAt ) we would compute the same time interpolant many times. Depending on the mapping method selected, we need to compute further derived data structures and values from the time interpolant (e.g. RBF coefficients) which is computationally demanding. The cache enables to store the latest queried time snapshot and (potentially, depending on the mapping) further derived data structures which remain constant for a specific point in time.

In the default implementation, the cache stores the latest waveform sample in time.

This function updates the cache in use to the provided time-interpolated data ( in ). It is called in the ReadDataContext::mapAndReadValues and called whenever the timestamp of the MappingDataCache is different from the actual 'readTime'. See also the documentation of the impl::MappingDataCache for more information.

Parameters
cachethe cache in use, specific to a particular mapping-data pair
inthe time-interpolated data values

Reimplemented from precice::mapping::Mapping.

Definition at line 391 of file PartitionOfUnityMapping.hpp.

Here is the call graph for this function:

Member Data Documentation

◆ _basisFunction

template<typename RADIAL_BASIS_FUNCTION_T>
RADIAL_BASIS_FUNCTION_T precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_basisFunction
private

Radial basis function type used in interpolation.

Definition at line 98 of file PartitionOfUnityMapping.hpp.

◆ _centerMesh

template<typename RADIAL_BASIS_FUNCTION_T>
std::unique_ptr<mesh::Mesh> precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_centerMesh
private

Definition at line 117 of file PartitionOfUnityMapping.hpp.

◆ _clusterRadius

template<typename RADIAL_BASIS_FUNCTION_T>
double precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_clusterRadius = 0
private

derived parameter based on the input above: the radius of each cluster

Definition at line 112 of file PartitionOfUnityMapping.hpp.

◆ _clusters

template<typename RADIAL_BASIS_FUNCTION_T>
std::vector<SphericalVertexCluster<RADIAL_BASIS_FUNCTION_T> > precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_clusters
private

main data container storing all the clusters, which need to be solved individually

Definition at line 95 of file PartitionOfUnityMapping.hpp.

◆ _log

template<typename RADIAL_BASIS_FUNCTION_T>
precice::logging::Logger precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_log {"mapping::PartitionOfUnityMapping"}
private

logger, as usual

Definition at line 92 of file PartitionOfUnityMapping.hpp.

◆ _polynomial

template<typename RADIAL_BASIS_FUNCTION_T>
Polynomial precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_polynomial
private

polynomial treatment of the RBF system

Definition at line 115 of file PartitionOfUnityMapping.hpp.

◆ _projectToInput

template<typename RADIAL_BASIS_FUNCTION_T>
const bool precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_projectToInput
private

toggles whether we project the cluster centers to the input mesh

Definition at line 109 of file PartitionOfUnityMapping.hpp.

◆ _relativeOverlap

template<typename RADIAL_BASIS_FUNCTION_T>
const double precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_relativeOverlap
private

overlap of vertex clusters

Definition at line 106 of file PartitionOfUnityMapping.hpp.

◆ _verticesPerCluster

template<typename RADIAL_BASIS_FUNCTION_T>
const unsigned int precice::mapping::PartitionOfUnityMapping< RADIAL_BASIS_FUNCTION_T >::_verticesPerCluster
private

Input parameters provided by the user for the clustering algorithm:

target number of input vertices for each cluster

Definition at line 103 of file PartitionOfUnityMapping.hpp.


The documentation for this class was generated from the following file: