preCICE v3.1.2
Loading...
Searching...
No Matches
Namespaces | Classes | Typedefs | Enumerations | Functions | Variables
precice::mapping Namespace Reference

contains data mapping from points to meshes. More...

Namespaces

namespace  impl
 
namespace  tests
 

Classes

class  AxialGeoMultiscaleMapping
 Geometric multiscale mapping in axial direction. More...
 
class  BarycentricBaseMapping
 Base class for interpolation based mappings, where mapping is done using a geometry-based linear combination of input values. Subclasses differ by the way computeMapping() fills the _interpolations and by mesh tagging. Mapping itself is shared. More...
 
class  CompactPolynomialC0
 Wendland radial basis function with compact support. More...
 
class  CompactPolynomialC2
 Wendland radial basis function with compact support. More...
 
class  CompactPolynomialC4
 Wendland radial basis function with compact support. More...
 
class  CompactPolynomialC6
 Wendland radial basis function with compact support. More...
 
class  CompactPolynomialC8
 Wendland radial basis function with compact support. More...
 
struct  CompactSupportBase
 Base class for RBF with compact support. More...
 
class  CompactThinPlateSplinesC2
 Radial basis function with compact support. More...
 
struct  DefiniteFunction
 Base class for RBF functions to distinguish positive definite functions. More...
 
class  Gaussian
 Radial basis function with global and compact support. More...
 
class  GinkgoRadialBasisFctSolver
 
class  InverseMultiquadrics
 Radial basis function with global support. More...
 
class  LinearCellInterpolationMapping
 Mapping using orthogonal projection to nearest triangle/edge/vertex and linear interpolation from projected point. More...
 
class  Mapping
 Abstract base class for mapping of data from one mesh to another. More...
 
class  MappingConfiguration
 Performs XML configuration and holds configured mappings. More...
 
class  Multiquadrics
 Radial basis function with global support. More...
 
class  NearestNeighborBaseMapping
 
class  NearestNeighborGradientMapping
 Mapping using nearest neighboring vertices and their local gradient values. More...
 
class  NearestNeighborMapping
 Mapping using nearest neighboring vertices. More...
 
class  NearestProjectionMapping
 Mapping using orthogonal projection to nearest triangle/edge/vertex and linear interpolation from projected point. More...
 
struct  NoCompactSupportBase
 Base class for RBF without compact support. More...
 
class  PartitionOfUnityMapping
 
class  PetRadialBasisFctMapping
 Mapping with radial basis functions using the Petsc library to solve the resulting system. More...
 
class  Polation
 Calculates the barycentric coordinates of a coordinate on the given vertex/edge/triangle and stores the corresponding weights If all barycentric coordinates are positive, the operation is interpolation. If not, it is an extrapolation. More...
 
class  RadialBasisFctBaseMapping
 Mapping with radial basis functions. More...
 
class  RadialBasisFctMapping
 Mapping with radial basis functions. More...
 
class  RadialBasisFctSolver
 
struct  RadialBasisParameters
 Wrapper struct that is used to transfer RBF-specific parameters to the GPU. More...
 
class  RadialGeoMultiscaleMapping
 Geometric multiscale mapping in radial direction. More...
 
class  SphericalVertexCluster
 
class  ThinPlateSplines
 Radial basis function with global support. More...
 
class  VolumeSplines
 Radial basis function with global support. More...
 
struct  WeightedElement
 Struct that contains weight and index of a vertex. More...
 

Typedefs

using PtrMapping = std::shared_ptr<Mapping>
 
using PtrMappingConfiguration = std::shared_ptr<MappingConfiguration>
 

Enumerations

enum class  Polynomial { ON , OFF , SEPARATE }
 How to handle the polynomial? More...
 
enum class  BasisFunction {
  WendlandC0 , WendlandC2 , WendlandC4 , WendlandC6 ,
  WendlandC8 , ThinPlateSplines , Multiquadrics , InverseMultiquadrics ,
  VolumeSplines , Gaussian , CompactThinPlateSplinesC2
}
 
enum class  GinkgoSolverType { CG , GMRES , QR }
 
enum class  GinkgoPreconditionerType { Jacobi , Cholesky , None }
 

Functions

bool operator< (Mapping::MeshRequirement lhs, Mapping::MeshRequirement rhs)
 
std::ostreamoperator<< (std::ostream &out, Mapping::MeshRequirement val)
 
std::ostreamoperator<< (std::ostream &os, const WeightedElement &w)
 Make the WeightedElement printable.
 
std::ostreamoperator<< (std::ostream &os, const Polation &p)
 Make the Polation class printable.
 
double computeSquaredDifference (const std::array< double, 3 > &u, std::array< double, 3 > v, const std::array< bool, 3 > &activeAxis={{true, true, true}})
 Deletes all dead directions from fullVector and returns a vector of reduced dimensionality.
 
template<typename IndexContainer >
constexpr void reduceActiveAxis (const mesh::Mesh &mesh, const IndexContainer &IDs, std::array< bool, 3 > &axis)
 given the active axis, computes sets the axis with the lowest spatial expansion to dead
 
template<typename IndexContainer >
void fillPolynomialEntries (Eigen::MatrixXd &matrix, const mesh::Mesh &mesh, const IndexContainer &IDs, Eigen::Index startIndex, std::array< bool, 3 > activeAxis)
 
template<typename RADIAL_BASIS_FUNCTION_T , typename IndexContainer >
Eigen::MatrixXd buildMatrixCLU (RADIAL_BASIS_FUNCTION_T basisFunction, const mesh::Mesh &inputMesh, const IndexContainer &inputIDs, std::array< bool, 3 > activeAxis, Polynomial polynomial)
 
template<typename RADIAL_BASIS_FUNCTION_T , typename IndexContainer >
Eigen::MatrixXd buildMatrixA (RADIAL_BASIS_FUNCTION_T basisFunction, const mesh::Mesh &inputMesh, const IndexContainer &inputIDs, const mesh::Mesh &outputMesh, const IndexContainer outputIDs, std::array< bool, 3 > activeAxis, Polynomial polynomial)
 

Variables

const std::map< std::string, GinkgoSolverTypesolverTypeLookup
 
const std::map< std::string, GinkgoPreconditionerTypepreconditionerTypeLookup
 
const std::map< std::string, std::function< std::shared_ptr< gko::Executor >(const unsigned int, const bool)> > ginkgoExecutorLookup
 

Detailed Description

contains data mapping from points to meshes.

Typedef Documentation

◆ PtrMapping

Definition at line 11 of file SharedPointer.hpp.

◆ PtrMappingConfiguration

Definition at line 12 of file SharedPointer.hpp.

Enumeration Type Documentation

◆ BasisFunction

Enumerator
WendlandC0 
WendlandC2 
WendlandC4 
WendlandC6 
WendlandC8 
ThinPlateSplines 
Multiquadrics 
InverseMultiquadrics 
VolumeSplines 
Gaussian 
CompactThinPlateSplinesC2 

Definition at line 17 of file MappingConfigurationTypes.hpp.

◆ GinkgoPreconditionerType

Enumerator
Jacobi 
Cholesky 
None 

Definition at line 60 of file GinkgoRadialBasisFctSolver.hpp.

◆ GinkgoSolverType

Enumerator
CG 
GMRES 
QR 

Definition at line 54 of file GinkgoRadialBasisFctSolver.hpp.

◆ Polynomial

enum class precice::mapping::Polynomial
strong

How to handle the polynomial?

ON: Include it in the system matrix OFF: Omit it altogether SEPARATE: Compute it separately using least-squares QR.

Enumerator
ON 
OFF 
SEPARATE 

Definition at line 11 of file MappingConfigurationTypes.hpp.

Function Documentation

◆ buildMatrixA()

template<typename RADIAL_BASIS_FUNCTION_T , typename IndexContainer >
Eigen::MatrixXd precice::mapping::buildMatrixA ( RADIAL_BASIS_FUNCTION_T basisFunction,
const mesh::Mesh & inputMesh,
const IndexContainer & inputIDs,
const mesh::Mesh & outputMesh,
const IndexContainer outputIDs,
std::array< bool, 3 > activeAxis,
Polynomial polynomial )

Definition at line 186 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ buildMatrixCLU()

template<typename RADIAL_BASIS_FUNCTION_T , typename IndexContainer >
Eigen::MatrixXd precice::mapping::buildMatrixCLU ( RADIAL_BASIS_FUNCTION_T basisFunction,
const mesh::Mesh & inputMesh,
const IndexContainer & inputIDs,
std::array< bool, 3 > activeAxis,
Polynomial polynomial )

Definition at line 141 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ computeSquaredDifference()

double precice::mapping::computeSquaredDifference ( const std::array< double, 3 > & u,
std::array< double, 3 > v,
const std::array< bool, 3 > & activeAxis = {{true, true, true}} )
inline

Deletes all dead directions from fullVector and returns a vector of reduced dimensionality.

Definition at line 78 of file RadialBasisFctSolver.hpp.

◆ fillPolynomialEntries()

template<typename IndexContainer >
void precice::mapping::fillPolynomialEntries ( Eigen::MatrixXd & matrix,
const mesh::Mesh & mesh,
const IndexContainer & IDs,
Eigen::Index startIndex,
std::array< bool, 3 > activeAxis )
inline

Definition at line 117 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

◆ operator<()

bool precice::mapping::operator< ( Mapping::MeshRequirement lhs,
Mapping::MeshRequirement rhs )

Defines an ordering for MeshRequirement in terms of specificality

Parameters
[in]lhsthe left-hand side of the binary operator
[in]rhsthe right-hand side of the binary operator

Definition at line 262 of file Mapping.cpp.

◆ operator<<() [1/3]

std::ostream & precice::mapping::operator<< ( std::ostream & os,
const Polation & p )

Make the Polation class printable.

Definition at line 101 of file Polation.cpp.

Here is the call graph for this function:

◆ operator<<() [2/3]

std::ostream & precice::mapping::operator<< ( std::ostream & os,
const WeightedElement & w )

Make the WeightedElement printable.

Definition at line 96 of file Polation.cpp.

◆ operator<<() [3/3]

std::ostream & precice::mapping::operator<< ( std::ostream & out,
Mapping::MeshRequirement val )

Defines the output operation to streams

Parameters
[in,out]outstream to output to.
[in]valthe value to output.

Definition at line 275 of file Mapping.cpp.

◆ reduceActiveAxis()

template<typename IndexContainer >
constexpr void precice::mapping::reduceActiveAxis ( const mesh::Mesh & mesh,
const IndexContainer & IDs,
std::array< bool, 3 > & axis )
constexpr

given the active axis, computes sets the axis with the lowest spatial expansion to dead

Definition at line 93 of file RadialBasisFctSolver.hpp.

Here is the call graph for this function:

Variable Documentation

◆ ginkgoExecutorLookup

const std::map<std::string, std::function<std::shared_ptr<gko::Executor>(const unsigned int, const bool)> > precice::mapping::ginkgoExecutorLookup
Initial value:
{{"reference-executor", [](auto unused, auto unused2) { return gko::ReferenceExecutor::create(); }},
{"omp-executor", [](auto unused, auto unused2) { return gko::OmpExecutor::create(); }},
{"cuda-executor", [](auto deviceId, auto enableUnifiedMemory) { if(enableUnifiedMemory) return gko::CudaExecutor::create(deviceId, gko::OmpExecutor::create(), true, gko::allocation_mode::unified_global); else return gko::CudaExecutor::create(deviceId, gko::OmpExecutor::create(), true, gko::allocation_mode::device); }},
{"hip-executor", [](auto deviceId, auto unused) { return gko::HipExecutor::create(deviceId, gko::OmpExecutor::create(), true); }}}

Definition at line 78 of file GinkgoRadialBasisFctSolver.hpp.

◆ preconditionerTypeLookup

const std::map<std::string, GinkgoPreconditionerType> precice::mapping::preconditionerTypeLookup
Initial value:
{
{"jacobi-preconditioner", GinkgoPreconditionerType::Jacobi},
{"cholesky-preconditioner", GinkgoPreconditionerType::Cholesky},
{"no-preconditioner", GinkgoPreconditionerType::None}}

Definition at line 73 of file GinkgoRadialBasisFctSolver.hpp.

◆ solverTypeLookup

const std::map<std::string, GinkgoSolverType> precice::mapping::solverTypeLookup
Initial value:
{
{"cg-solver", GinkgoSolverType::CG},
{"gmres-solver", GinkgoSolverType::GMRES},
{"qr-solver", GinkgoSolverType::QR}}

Definition at line 68 of file GinkgoRadialBasisFctSolver.hpp.