preCICE v3.2.0
|
Mapping with radial basis functions using the Petsc library to solve the resulting system. More...
#include <PetRadialBasisFctMapping.hpp>
Public Member Functions | |
PetRadialBasisFctMapping (Mapping::Constraint constraint, int dimensions, const RADIAL_BASIS_FUNCTION_T &function, std::array< bool, 3 > deadAxis, double solverRtol=1e-9, Polynomial polynomial=Polynomial::SEPARATE) | |
Constructor. | |
~PetRadialBasisFctMapping () override | |
Deletes the PETSc objects and the _deadAxis array. | |
void | computeMapping () final override |
Computes the mapping coefficients from the in- and output mesh. | |
void | clear () final override |
Removes a computed mapping. | |
std::string | getName () const final override |
name of the rbf mapping | |
![]() | |
RadialBasisFctBaseMapping (Constraint constraint, int dimensions, const RADIAL_BASIS_FUNCTION_T &function, std::array< bool, 3 > deadAxis, InitialGuessRequirement mappingType) | |
Constructor. | |
~RadialBasisFctBaseMapping () override=default | |
void | tagMeshFirstRound () final override |
Method used by partition. Tags vertices that could be owned by this rank. | |
void | tagMeshSecondRound () final override |
Method used by partition. Tags vertices that can be filtered out. | |
![]() | |
Mapping (Constraint constraint, int dimensions, bool requiresGradientData, InitialGuessRequirement initialGuessRequirement) | |
Constructor, takes mapping constraint. | |
Mapping & | operator= (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::PtrMesh & | getInputMesh () const |
const mesh::PtrMesh & | getOutputMesh () 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. | |
virtual void | mapConservativeAt (const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const Eigen::Ref< const Eigen::MatrixXd > &source, impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > target) |
Just-in-time mapping variant of mapConservative. | |
virtual void | mapConsistentAt (const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > values) |
Just-in-time mapping variant of mapConsistent. | |
virtual void | updateMappingDataCache (impl::MappingDataCache &cache, const Eigen::Ref< const Eigen::VectorXd > &in) |
Allows updating a so-called MappingDataCache for more efficient just-in-time mappings. | |
virtual void | initializeMappingDataCache (impl::MappingDataCache &cache) |
Allocates memory and sets up the data structures inside the MappingDataCache. | |
virtual void | completeJustInTimeMapping (impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > result) |
Completes a just-in-time mapping for conservative constraints. | |
Private Types | |
using | VertexData = std::vector<std::vector<std::pair<int, double>>> |
Stores col -> value for each row. Used to return the already computed values from the preconditioning. | |
Private Member Functions | |
void | mapConservative (const time::Sample &inData, Eigen::VectorXd &outData) final override |
Maps data using a conservative constraint. | |
void | mapConsistent (const time::Sample &inData, Eigen::VectorXd &outData) final override |
Maps data using a consistent constraint. | |
void | printMappingInfo (int dim) const |
Prints an INFO about the current mapping. | |
VertexData | bgPreallocationMatrixC (const mesh::PtrMesh inMesh) |
Preallocate matrix C and saves the coefficients using a boost::geometry spatial tree for neighbor search. | |
VertexData | bgPreallocationMatrixA (const mesh::PtrMesh inMesh, const mesh::PtrMesh outMesh) |
void | loadInitialGuessForDim (int dimension, int allDimensions, petsc::Vector &destination) |
void | storeInitialGuessForDim (int dimension, int allDimensions, petsc::Vector &source) |
Private Attributes | |
logging::Logger | _log {"mapping::PetRadialBasisFctMapping"} |
petsc::Matrix | _matrixC |
Interpolation system matrix. Evaluated basis function on the input mesh. | |
petsc::Matrix | _matrixQ |
Vandermonde Matrix for linear polynomial, constructed from vertices of the input mesh. | |
petsc::Matrix | _matrixA |
Interpolation evaluation matrix. Evaluated basis function on the output mesh. | |
petsc::Matrix | _matrixV |
Coordinates of the output mesh to evaluate the separated polynomial. | |
petsc::Vector | oneInterpolant |
Interpolant of g(x) = 1 evaluated at output sites, used for rescaling. | |
petsc::KSPSolver | _solver |
Used to solve matrixC for the RBF weighting factors. | |
petsc::KSPSolver | _QRsolver |
Used to solve the under-determined system for the separated polynomial. | |
AO | _AOmapping |
Maps globalIndex to local index. | |
const double | _solverRtol |
Polynomial | _polynomial |
Toggles the use of the additional polynomial. | |
bool | useRescaling = true |
Toggles use of rescaled basis functions, only active when Polynomial == SEPARATE. | |
size_t | polyparams |
Number of coefficients for the integrated polynomial. Depends on dimension and number of dead dimensions. | |
size_t | sepPolyparams |
Number of coefficients for the separated polynomial. Depends on dimension and number of dead dimensions. | |
size_t | localPolyparams |
Equals to polyparams if rank == 0. Is 0 everywhere else. | |
utils::Parallel::CommStatePtr | _commState |
The CommState used to communicate. | |
Additional Inherited Members | |
![]() | |
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... | |
![]() | |
int | getPolynomialParameters () const |
Computes the number of polynomial degrees of freedom based on the problem dimension and the dead axis. | |
![]() | |
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 |
![]() | |
RADIAL_BASIS_FUNCTION_T | _basisFunction |
Radial basis function type used in interpolation. | |
std::vector< bool > | _deadAxis |
true if the mapping along some axis should be ignored | |
![]() | |
bool | _hasComputedMapping = false |
Flag to indicate whether computeMapping() has been called. | |
bool | _requiresGradientData |
Flag if gradient data is required for the mapping. | |
Mapping with radial basis functions using the Petsc library to solve the resulting system.
With help of the input data points and values an interpolant is constructed. The interpolant is formed by a weighted sum of conditionally positive radial basis functions and a (low order) polynomial, and evaluated at the output data points.
The radial basis function type has to be given as template parameter.
Definition at line 51 of file PetRadialBasisFctMapping.hpp.
|
private |
Stores col -> value for each row. Used to return the already computed values from the preconditioning.
Definition at line 93 of file PetRadialBasisFctMapping.hpp.
precice::mapping::PetRadialBasisFctMapping< RADIAL_BASIS_FUNCTION_T >::PetRadialBasisFctMapping | ( | Mapping::Constraint | constraint, |
int | dimensions, | ||
const RADIAL_BASIS_FUNCTION_T & | function, | ||
std::array< bool, 3 > | deadAxis, | ||
double | solverRtol = 1e-9, | ||
Polynomial | polynomial = Polynomial::SEPARATE ) |
Constructor.
[in] | constraint | Specifies mapping to be consistent or conservative. |
[in] | dimensions | Dimensionality of the meshes |
[in] | function | Radial basis function used for mapping. |
[in] | xDead,yDead,zDead | Deactivates mapping along an axis |
[in] | solverRtol | Relative tolerance for the linear solver. |
[in] | polynomial | Type of polynomial augmentation |
For description on convergence testing and meaning of solverRtol see http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPConvergedDefault.html#KSPConvergedDefault
Definition at line 169 of file PetRadialBasisFctMapping.hpp.
|
override |
Deletes the PETSc objects and the _deadAxis array.
Definition at line 194 of file PetRadialBasisFctMapping.hpp.
|
private |
Definition at line 922 of file PetRadialBasisFctMapping.hpp.
|
private |
Preallocate matrix C and saves the coefficients using a boost::geometry spatial tree for neighbor search.
Definition at line 842 of file PetRadialBasisFctMapping.hpp.
|
finaloverridevirtual |
Removes a computed mapping.
Implements precice::mapping::RadialBasisFctBaseMapping< RADIAL_BASIS_FUNCTION_T >.
Definition at line 504 of file PetRadialBasisFctMapping.hpp.
|
finaloverridevirtual |
Computes the mapping coefficients from the in- and output mesh.
Implements precice::mapping::RadialBasisFctBaseMapping< RADIAL_BASIS_FUNCTION_T >.
Definition at line 200 of file PetRadialBasisFctMapping.hpp.
|
finaloverridevirtual |
name of the rbf mapping
Implements precice::mapping::Mapping.
Definition at line 520 of file PetRadialBasisFctMapping.hpp.
|
private |
load the initialGuess for a given dimension or allocates the storage for the first iteration
The initialGuess passed to Mapping::map contains one guess for each data dimension as the mapping is evaluated for each dimension. As an examples, a 3D vector thus contains 3 initialGuesses, one for each solver. This extracts it from the overall initialGuess.
Definition at line 526 of file PetRadialBasisFctMapping.hpp.
|
finaloverrideprivatevirtual |
Maps data using a conservative constraint.
[in] | input | Sample to map data from |
[in] | output | Values 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.
Implements precice::mapping::Mapping.
Definition at line 690 of file PetRadialBasisFctMapping.hpp.
|
finaloverrideprivatevirtual |
Maps data using a consistent constraint.
[in] | input | Sample to map data from |
[in] | output | Values 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.
Implements precice::mapping::Mapping.
Definition at line 563 of file PetRadialBasisFctMapping.hpp.
|
private |
Prints an INFO about the current mapping.
Definition at line 820 of file PetRadialBasisFctMapping.hpp.
|
private |
stores the initialGuess for a given dimension
The initialGuess passed to Mapping::map contains one guess for each data dimension as the mapping is evaluated for each dimension. As an examples, a 3D vector thus contains 3 initialGuesses, one for each solver. This stores the initialGuess of a given dimension into the overall initialGuess.
Definition at line 547 of file PetRadialBasisFctMapping.hpp.
|
private |
Maps globalIndex to local index.
Definition at line 119 of file PetRadialBasisFctMapping.hpp.
|
private |
The CommState used to communicate.
Definition at line 142 of file PetRadialBasisFctMapping.hpp.
|
mutableprivate |
Definition at line 95 of file PetRadialBasisFctMapping.hpp.
|
private |
Interpolation evaluation matrix. Evaluated basis function on the output mesh.
Definition at line 104 of file PetRadialBasisFctMapping.hpp.
|
private |
Interpolation system matrix. Evaluated basis function on the input mesh.
Definition at line 98 of file PetRadialBasisFctMapping.hpp.
|
private |
Vandermonde Matrix for linear polynomial, constructed from vertices of the input mesh.
Definition at line 101 of file PetRadialBasisFctMapping.hpp.
|
private |
Coordinates of the output mesh to evaluate the separated polynomial.
Definition at line 107 of file PetRadialBasisFctMapping.hpp.
|
private |
Toggles the use of the additional polynomial.
Definition at line 124 of file PetRadialBasisFctMapping.hpp.
|
private |
Used to solve the under-determined system for the separated polynomial.
Definition at line 116 of file PetRadialBasisFctMapping.hpp.
|
private |
Used to solve matrixC for the RBF weighting factors.
Definition at line 113 of file PetRadialBasisFctMapping.hpp.
|
private |
Definition at line 121 of file PetRadialBasisFctMapping.hpp.
|
private |
Equals to polyparams if rank == 0. Is 0 everywhere else.
Definition at line 136 of file PetRadialBasisFctMapping.hpp.
|
private |
Interpolant of g(x) = 1 evaluated at output sites, used for rescaling.
Definition at line 110 of file PetRadialBasisFctMapping.hpp.
|
private |
Number of coefficients for the integrated polynomial. Depends on dimension and number of dead dimensions.
Definition at line 130 of file PetRadialBasisFctMapping.hpp.
|
private |
Number of coefficients for the separated polynomial. Depends on dimension and number of dead dimensions.
Definition at line 133 of file PetRadialBasisFctMapping.hpp.
|
private |
Toggles use of rescaled basis functions, only active when Polynomial == SEPARATE.
Definition at line 127 of file PetRadialBasisFctMapping.hpp.