preCICE v3.1.1
Loading...
Searching...
No Matches
Public Types | Public 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.
 
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
 

Private Attributes

precice::logging::Logger _log {"mapping::RadialBasisFctSolver"}
 
DecompositionType _decMatrixC
 Decomposition of the interpolation matrix.
 
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)
 

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 23 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 26 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 25 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 222 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

Member Function Documentation

◆ clear()

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

Definition at line 360 of file RadialBasisFctSolver.hpp.

◆ getInputSize()

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

Definition at line 367 of file RadialBasisFctSolver.hpp.

◆ getOutputSize()

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

Definition at line 373 of file RadialBasisFctSolver.hpp.

◆ 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 288 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 336 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 60 of file RadialBasisFctSolver.hpp.

◆ _log

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

Definition at line 57 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 72 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 66 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 69 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 63 of file RadialBasisFctSolver.hpp.


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