29template <
typename RADIAL_BASIS_FUNCTION_T>
50 RADIAL_BASIS_FUNCTION_T function,
52 unsigned int verticesPerCluster,
53 double relativeOverlap,
79 void mapConsistentAt(const Eigen::Ref<const Eigen::MatrixXd> &coordinates, const
impl::MappingDataCache &cache, Eigen::Ref<Eigen::MatrixXd> values) final override;
82 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;
135template <
typename RADIAL_BASIS_FUNCTION_T>
139 RADIAL_BASIS_FUNCTION_T function,
141 unsigned int verticesPerCluster,
142 double relativeOverlap,
161template <
typename RADIAL_BASIS_FUNCTION_T>
175 outMesh = this->
input();
177 inMesh = this->
input();
194 meshVertices.clear();
196 _clusters.reserve(centerCandidates.size());
197 for (
const auto &c : centerCandidates) {
200 const VertexID vertexID = meshVertices.size();
205 if (!cluster.
empty()) {
207 meshVertices.emplace_back(std::move(center));
208 _clusters.emplace_back(std::move(cluster));
219 [](
const auto &c) { return c.getNumberOfInputVertices(); }) /
221 PRECICE_DEBUG(
"Maximum number of vertices per cluster {}",
std::max_element(
_clusters.begin(),
_clusters.end(), [](
auto &v1,
auto &v2) { return v1.getNumberOfInputVertices() < v2.getNumberOfInputVertices(); })->getNumberOfInputVertices());
222 PRECICE_DEBUG(
"Minimum number of vertices per cluster {}",
std::min_element(
_clusters.begin(),
_clusters.end(), [](
auto &v1,
auto &v2) { return v1.getNumberOfInputVertices() < v2.getNumberOfInputVertices(); })->getNumberOfInputVertices());
232 for (
const auto &vertex : outMesh->
vertices()) {
236 for (
unsigned int i = 0; i < clusterIDs.size(); ++i) {
238 _clusters[clusterIDs[i]].setNormalizedWeight(normalizedWeights[i], vertex.getID());
254template <
typename RADIAL_BASIS_FUNCTION_T>
271 const auto localNumberOfClusters = clusterIDs.
size();
279 "Output vertex {} of mesh \"{}\" could not be assigned to any cluster in the rbf-pum mapping. This probably means that the meshes do not match well geometry-wise: Visualize the exported preCICE meshes to confirm. "
280 "If the meshes are fine geometry-wise, you can try to increase the number of \"vertices-per-cluster\" (default is 50), the \"relative-overlap\" (default is 0.15), or disable the option \"project-to-input\". "
281 "These options are only valid for the <mapping:rbf-pum-direct/> tag.",
289 std::transform(clusterIDs.cbegin(), clusterIDs.cend(), weights.
begin(), [&](
const auto &ids) { return _clusters[ids].computeWeight(vertex); });
293 if (weightSum <= 0) {
304 return {clusterIDs, weights};
307template <
typename RADIAL_BASIS_FUNCTION_T>
322template <
typename RADIAL_BASIS_FUNCTION_T>
337template <
typename RADIAL_BASIS_FUNCTION_T>
350 for (Eigen::Index v = 0; v < coordinates.cols(); ++v) {
354 for (
std::size_t i = 0; i < clusterIDs.size(); ++i) {
356 auto id = clusterIDs[i];
358 Eigen::VectorXd res = normalizedWeights[i] * source.col(v);
364template <
typename RADIAL_BASIS_FUNCTION_T>
374 if (cache.
p[c].squaredNorm() > 0) {
380template <
typename RADIAL_BASIS_FUNCTION_T>
392template <
typename RADIAL_BASIS_FUNCTION_T>
405template <
typename RADIAL_BASIS_FUNCTION_T>
417 for (Eigen::Index v = 0; v < values.cols(); ++v) {
421 for (
std::size_t i = 0; i < clusterIDs.size(); ++i) {
423 auto id = clusterIDs[i];
426 values.col(v) += localRes;
431template <
typename RADIAL_BASIS_FUNCTION_T>
437 filterMesh = this->
output();
438 outMesh = this->
input();
440 filterMesh = this->
input();
444 if (outMesh->
empty())
483 auto verticesNew = filterMesh->index().getVerticesInsideBox(localBB);
485 std::for_each(verticesNew.begin(), verticesNew.end(), [&filterMesh](
VertexID v) { filterMesh->vertex(v).tag(); });
488template <
typename RADIAL_BASIS_FUNCTION_T>
494template <
typename RADIAL_BASIS_FUNCTION_T>
499 auto dataRadius = centerMesh.
createData(
"radius", 1, -1);
500 auto dataCardinality = centerMesh.
createData(
"number-of-vertices", 1, -1);
503 for (
unsigned int i = 0; i <
_clusters.size(); ++i) {
504 dataCardinality->values()[i] =
static_cast<double>(
_clusters[i].getNumberOfInputVertices());
511 int numberOfVertices = centerMesh.
nVertices();
523 vertexOffsets[0] = centerMesh.
nVertices();
527 int numberOfSecondaryRankVertices = -1;
530 vertexOffsets[secondaryRank] = numberOfSecondaryRankVertices + vertexOffsets[secondaryRank - 1];
539 dataRadius->setSampleAtTime(0,
time::Sample{1, dataRadius->values()});
540 dataCardinality->setSampleAtTime(0,
time::Sample{1, dataCardinality->values()});
545template <
typename RADIAL_BASIS_FUNCTION_T>
555template <
typename RADIAL_BASIS_FUNCTION_T>
558 return "partition-of-unity RBF";
#define PRECICE_DEBUG(...)
#define PRECICE_TRACE(...)
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
VertexContainer & vertices()
Returns modifieable container holding all vertices.
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
std::size_t nVertices() const
Returns the number of vertices.
bool isJustInTime() const
bool empty() const
Does the mesh contain any vertices?
const BoundingBox & getBoundingBox() const
Returns the bounding box of the mesh.
void doExport(int index, double time) final override
Export the mesh and writes files.
mesh::PtrMesh output() const
Returns pointer to output mesh.
Constraint
Specifies additional constraints for a mapping.
int getDimensions() const
Mapping(Constraint constraint, int dimensions, bool requiresGradientData, InitialGuessRequirement initialGuessRequirement)
Constructor, takes mapping constraint.
mesh::PtrMesh input() const
Returns pointer to input mesh.
bool _hasComputedMapping
Flag to indicate whether computeMapping() has been called.
bool isScaledConsistent() const
Returns true if mapping is a form of scaled consistent mapping.
void setInputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the input mesh.
void setOutputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the output mesh.
virtual bool hasConstraint(const Constraint &constraint) const
Checks whether the mapping has the given constraint or not.
InitialGuessRequirement
Specifies whether the mapping requires an initial guess.
void tagMeshSecondRound() final override
void tagMeshFirstRound() final override
std::unique_ptr< mesh::Mesh > _centerMesh
void exportClusterCentersAsVTU(mesh::Mesh ¢ers)
std::pair< std::vector< int >, std::vector< double > > computeNormalizedWeight(const mesh::Vertex &v, std::string_view mesh)
precice::logging::Logger _log
void initializeMappingDataCache(impl::MappingDataCache &cache) final override
const double _relativeOverlap
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
std::vector< SphericalVertexCluster< RBF > > _clusters
const bool _projectToInput
std::string getName() const final override
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a consistent constraint.
PartitionOfUnityMapping(Mapping::Constraint constraint, int dimension, RADIAL_BASIS_FUNCTION_T function, Polynomial polynomial, unsigned int verticesPerCluster, double relativeOverlap, bool projectToInput)
void mapConsistentAt(const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > values) final override
const unsigned int _verticesPerCluster
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a conservative constraint.
void updateMappingDataCache(impl::MappingDataCache &cache, const Eigen::Ref< const Eigen::VectorXd > &in) final override
void computeMapping() final override
void completeJustInTimeMapping(impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > buffer) final override
void clear() final override
Container and creator for meshes.
static constexpr MeshID MESH_ID_UNDEFINED
Use if the id of the mesh is not necessary.
std::size_t nVertices() const
Returns the number of vertices.
PtrData & createData(const std::string &name, int dimension, DataID id, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Create only data for vertex.
const VertexOffsets & getVertexOffsets() const
std::vector< int > VertexOffsets
void setVertexOffsets(VertexOffsets vertexOffsets)
Only used for tests.
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
VertexID getID() const
Returns the unique (among vertices of one mesh on one processor) ID of the vertex.
void setCoords(const VECTOR_T &coordinates)
Sets the coordinates of the vertex.
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
void stop()
Stops a running event.
void addData(std::string_view key, int value)
Adds named integer data, associated to an event.
Class to query the index trees of the mesh.
std::vector< VertexID > getVerticesInsideBox(const mesh::Vertex ¢erVertex, double radius)
Return all the vertices inside the box formed by vertex and radius (boundary exclusive)
static int getSize()
Number of ranks. This includes ranks from both participants, e.g. minimal size is 2.
static Rank getRank()
Current rank.
static bool isPrimary()
True if this process is running the primary rank.
static auto allSecondaryRanks()
Returns an iterable range over salve ranks [1, _size)
static bool isSecondary()
True if this process is running a secondary rank.
static com::PtrCommunication & getCommunication()
Intra-participant communication.
contains the logging framework.
std::tuple< double, Vertices > createClustering(mesh::PtrMesh inMesh, mesh::PtrMesh outMesh, double relativeOverlap, unsigned int verticesPerCluster, bool projectClustersToInput)
Creates a clustering as a collection of Vertices (representing the cluster centers) and a cluster rad...
double estimateClusterRadius(unsigned int verticesPerCluster, mesh::PtrMesh inMesh, const mesh::BoundingBox &bb)
Computes an estimate for the cluster radius, which results in approximately verticesPerCluster vertic...
contains data mapping from points to meshes.
Polynomial
How to handle the polynomial?
provides Mesh, Data and primitives.
std::shared_ptr< Mesh > PtrMesh
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Main namespace of the precice library.
bool syncMode
Enabled further inter- and intra-solver synchronisation.
std::vector< Eigen::MatrixXd > polynomialContributions
std::vector< Eigen::MatrixXd > p
int getDataDimensions() const
Returns the number of data components.