preCICE v3.2.0
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T > Class Template Reference

#include <RadialBasisFctSolver.hpp>

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

Public Types

using DecompositionType = std::conditional_t<RADIAL_BASIS_FUNCTION_T::isStrictlyPositiveDefinite(), Eigen::LLT<Eigen::MatrixXd>, Eigen::ColPivHouseholderQR<Eigen::MatrixXd>>
 
using BASIS_FUNCTION_T = RADIAL_BASIS_FUNCTION_T
 

Public Member Functions

 RadialBasisFctSolver ()=default
 Default constructor.
 
template<typename IndexContainer >
 RadialBasisFctSolver (RADIAL_BASIS_FUNCTION_T basisFunction, const mesh::Mesh &inputMesh, const IndexContainer &inputIDs, const mesh::Mesh &outputMesh, const IndexContainer &outputIDs, std::vector< bool > deadAxis, Polynomial polynomial)
 
Eigen::VectorXd solveConsistent (Eigen::VectorXd &inputData, Polynomial polynomial) const
 Maps the given input data.
 
void computeCacheData (Eigen::MatrixXd &inputData, Polynomial polynomial, Eigen::MatrixXd &polyOut, Eigen::MatrixXd &coeffsOut) const
 
Eigen::VectorXd solveConservative (const Eigen::VectorXd &inputData, Polynomial polynomial) const
 Maps the given input data.
 
void clear ()
 
Eigen::Index getInputSize () const
 
Eigen::Index getOutputSize () const
 
Eigen::Index getNumberOfPolynomials () const
 
template<typename IndexContainer >
Eigen::VectorXd interpolateAt (const mesh::Vertex &v, const Eigen::MatrixXd &poly, const Eigen::MatrixXd &coeffs, const RADIAL_BASIS_FUNCTION_T &function, const IndexContainer &inputIDs, const mesh::Mesh &inMesh) const
 
template<typename IndexContainer >
void addWriteDataToCache (const mesh::Vertex &v, const Eigen::VectorXd &load, Eigen::MatrixXd &epsilon, Eigen::MatrixXd &Au, const RADIAL_BASIS_FUNCTION_T &basisFunction, const IndexContainer &inputIDs, const mesh::Mesh &inMesh) const
 
void evaluateConservativeCache (Eigen::MatrixXd &epsilon, const Eigen::MatrixXd &Au, Eigen::MatrixXd &result) const
 

Private Member Functions

double evaluateRippaLOOCVerror (const Eigen::VectorXd &lambda) const
 

Private Attributes

precice::logging::Logger _log {"mapping::RadialBasisFctSolver"}
 
DecompositionType _decMatrixC
 Decomposition of the interpolation matrix.
 
Eigen::VectorXd _inverseDiagonal
 Diagonal entris of the inverse matrix C, requires for the Rippa scheme.
 
Eigen::ColPivHouseholderQR< Eigen::MatrixXd > _qrMatrixQ
 Decomposition of the polynomial (for separate polynomial)
 
Eigen::MatrixXd _matrixQ
 Polynomial matrix of the input mesh (for separate polynomial)
 
Eigen::MatrixXd _matrixV
 Polynomial matrix of the output mesh (for separate polynomial)
 
Eigen::MatrixXd _matrixA
 Evaluation matrix (output x input)
 
bool computeCrossValidation = false
 
std::array< bool, 3 > _localActiveAxis
 

Detailed Description

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

This class assembles and solves an RBF system, given an input mesh and an output mesh with relevant vertex IDs. The class uses a dense matrix decomposition in order to decompose the resulting system(s) and a backward substitution in order to solve the system at runtime. The functionality uses Eigen and supports only serial execution. In case the polynomial="separate" option is used, the polynomial system is solved using a QR decomposition.

Definition at line 25 of file RadialBasisFctSolver.hpp.

Member Typedef Documentation

◆ BASIS_FUNCTION_T

template<typename RADIAL_BASIS_FUNCTION_T >
using precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::BASIS_FUNCTION_T = RADIAL_BASIS_FUNCTION_T

Definition at line 28 of file RadialBasisFctSolver.hpp.

◆ DecompositionType

template<typename RADIAL_BASIS_FUNCTION_T >
using precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::DecompositionType = std::conditional_t<RADIAL_BASIS_FUNCTION_T::isStrictlyPositiveDefinite(), Eigen::LLT<Eigen::MatrixXd>, Eigen::ColPivHouseholderQR<Eigen::MatrixXd>>

Definition at line 27 of file RadialBasisFctSolver.hpp.

Constructor & Destructor Documentation

◆ RadialBasisFctSolver() [1/2]

template<typename RADIAL_BASIS_FUNCTION_T >
precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::RadialBasisFctSolver ( )
default

Default constructor.

◆ RadialBasisFctSolver() [2/2]

template<typename RADIAL_BASIS_FUNCTION_T >
template<typename IndexContainer >
precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::RadialBasisFctSolver ( RADIAL_BASIS_FUNCTION_T basisFunction,
const mesh::Mesh & inputMesh,
const IndexContainer & inputIDs,
const mesh::Mesh & outputMesh,
const IndexContainer & outputIDs,
std::vector< bool > deadAxis,
Polynomial polynomial )

assembles the system matrices and computes the decomposition of the interpolation matrix inputMesh refers to the mesh where the interpolants are built on, i.e., the input mesh for consistent mappings and the output mesh for conservative mappings outputMesh refers to the mesh where we evaluate the interpolants, i.e., the output mesh consistent mappings and the input mesh for conservative mappings

Definition at line 342 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

Member Function Documentation

◆ addWriteDataToCache()

template<typename RADIAL_BASIS_FUNCTION_T >
template<typename IndexContainer >
void precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::addWriteDataToCache ( const mesh::Vertex & v,
const Eigen::VectorXd & load,
Eigen::MatrixXd & epsilon,
Eigen::MatrixXd & Au,
const RADIAL_BASIS_FUNCTION_T & basisFunction,
const IndexContainer & inputIDs,
const mesh::Mesh & inMesh ) const

Definition at line 526 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ clear()

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

Definition at line 578 of file RadialBasisFctSolver.hpp.

◆ computeCacheData()

template<typename RADIAL_BASIS_FUNCTION_T >
void precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::computeCacheData ( Eigen::MatrixXd & inputData,
Polynomial polynomial,
Eigen::MatrixXd & polyOut,
Eigen::MatrixXd & coeffsOut ) const

Definition at line 471 of file RadialBasisFctSolver.hpp.

◆ evaluateConservativeCache()

template<typename RADIAL_BASIS_FUNCTION_T >
void precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::evaluateConservativeCache ( Eigen::MatrixXd & epsilon,
const Eigen::MatrixXd & Au,
Eigen::MatrixXd & result ) const

Definition at line 565 of file RadialBasisFctSolver.hpp.

◆ evaluateRippaLOOCVerror()

template<typename RADIAL_BASIS_FUNCTION_T >
double precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::evaluateRippaLOOCVerror ( const Eigen::VectorXd & lambda) const
private

Definition at line 319 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ getInputSize()

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::Index precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::getInputSize ( ) const

Definition at line 585 of file RadialBasisFctSolver.hpp.

◆ getNumberOfPolynomials()

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::Index precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::getNumberOfPolynomials ( ) const

Definition at line 597 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ getOutputSize()

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::Index precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::getOutputSize ( ) const

Definition at line 591 of file RadialBasisFctSolver.hpp.

◆ interpolateAt()

template<typename RADIAL_BASIS_FUNCTION_T >
template<typename IndexContainer >
Eigen::VectorXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::interpolateAt ( const mesh::Vertex & v,
const Eigen::MatrixXd & poly,
const Eigen::MatrixXd & coeffs,
const RADIAL_BASIS_FUNCTION_T & function,
const IndexContainer & inputIDs,
const mesh::Mesh & inMesh ) const

Definition at line 486 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ solveConservative()

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::VectorXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::solveConservative ( const Eigen::VectorXd & inputData,
Polynomial polynomial ) const

Maps the given input data.

Definition at line 414 of file RadialBasisFctSolver.hpp.

◆ solveConsistent()

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::VectorXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::solveConsistent ( Eigen::VectorXd & inputData,
Polynomial polynomial ) const

Maps the given input data.

Definition at line 442 of file RadialBasisFctSolver.hpp.

Member Data Documentation

◆ _decMatrixC

template<typename RADIAL_BASIS_FUNCTION_T >
DecompositionType precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_decMatrixC
private

Decomposition of the interpolation matrix.

Definition at line 78 of file RadialBasisFctSolver.hpp.

◆ _inverseDiagonal

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::VectorXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_inverseDiagonal
private

Diagonal entris of the inverse matrix C, requires for the Rippa scheme.

Definition at line 81 of file RadialBasisFctSolver.hpp.

◆ _localActiveAxis

template<typename RADIAL_BASIS_FUNCTION_T >
std::array<bool, 3> precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_localActiveAxis
private

Definition at line 96 of file RadialBasisFctSolver.hpp.

◆ _log

template<typename RADIAL_BASIS_FUNCTION_T >
precice::logging::Logger precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_log {"mapping::RadialBasisFctSolver"}
mutableprivate

Definition at line 74 of file RadialBasisFctSolver.hpp.

◆ _matrixA

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::MatrixXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_matrixA
private

Evaluation matrix (output x input)

Definition at line 93 of file RadialBasisFctSolver.hpp.

◆ _matrixQ

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::MatrixXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_matrixQ
private

Polynomial matrix of the input mesh (for separate polynomial)

Definition at line 87 of file RadialBasisFctSolver.hpp.

◆ _matrixV

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::MatrixXd precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_matrixV
private

Polynomial matrix of the output mesh (for separate polynomial)

Definition at line 90 of file RadialBasisFctSolver.hpp.

◆ _qrMatrixQ

template<typename RADIAL_BASIS_FUNCTION_T >
Eigen::ColPivHouseholderQR<Eigen::MatrixXd> precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::_qrMatrixQ
private

Decomposition of the polynomial (for separate polynomial)

Definition at line 84 of file RadialBasisFctSolver.hpp.

◆ computeCrossValidation

template<typename RADIAL_BASIS_FUNCTION_T >
bool precice::mapping::RadialBasisFctSolver< RADIAL_BASIS_FUNCTION_T >::computeCrossValidation = false
private

Definition at line 95 of file RadialBasisFctSolver.hpp.


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