preCICE
Loading...
Searching...
No Matches
PartitionOfUnityMapping.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <numeric>
5
7#include "io/ExportVTU.hpp"
12#include "mesh/Filter.hpp"
14#include "profiling/Event.hpp"
15#include "query/Index.hpp"
16#include "utils/IntraComm.hpp"
17
18namespace precice {
19extern bool syncMode;
20
21namespace mapping {
22
30template <typename RADIAL_BASIS_FUNCTION_T>
32public:
49 Mapping::Constraint constraint,
50 int dimension,
51 RADIAL_BASIS_FUNCTION_T function,
52 Polynomial polynomial,
53 unsigned int verticesPerCluster,
54 double relativeOverlap,
55 bool projectToInput,
57 bool computeEvaluationOffline = false);
58
68 void computeMapping() final override;
69
71 void clear() final override;
72
74 void tagMeshFirstRound() final override;
75
77 void tagMeshSecondRound() final override;
78
80 std::string getName() const final override;
81
82 void mapConsistentAt(const Eigen::Ref<const Eigen::MatrixXd> &coordinates, const impl::MappingDataCache &cache, Eigen::Ref<Eigen::MatrixXd> values) final override;
83
85 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;
86
87 void updateMappingDataCache(impl::MappingDataCache &cache, const Eigen::Ref<const Eigen::VectorXd> &in) final override;
88
89 void completeJustInTimeMapping(impl::MappingDataCache &cache, Eigen::Ref<Eigen::MatrixXd> buffer) final override;
90
91 void initializeMappingDataCache(impl::MappingDataCache &cache) final override;
92
93private:
95 precice::logging::Logger _log{"mapping::PartitionOfUnityMapping"};
96
97 void _computeCPU(mesh::PtrMesh inMesh, mesh::PtrMesh outMesh, double clusterRadius, const std::vector<mesh::Vertex> &centerCandidates);
98
100 std::vector<SphericalVertexCluster<RADIAL_BASIS_FUNCTION_T>> _clusters;
101
103 RADIAL_BASIS_FUNCTION_T _basisFunction;
104
106
108 const unsigned int _verticesPerCluster;
109
111 const double _relativeOverlap;
112
114 const bool _projectToInput;
115
117 double _clusterRadius = 0;
118
121
123 std::unique_ptr<BatchedRBFSolver<RADIAL_BASIS_FUNCTION_T>> _batchedSolver;
124
126
127 std::unique_ptr<mesh::Mesh> _centerMesh;
128
130 void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) override;
131
133 void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) override;
134
138 std::pair<std::vector<int>, std::vector<double>> computeNormalizedWeight(const mesh::Vertex &v, std::string_view mesh);
139
143
144 // Currently only valid for the Batched RBF solver, maybe move it to dedicated config struct
146};
147
148template <typename RADIAL_BASIS_FUNCTION_T>
150 Mapping::Constraint constraint,
151 int dimension,
152 RADIAL_BASIS_FUNCTION_T function,
153 Polynomial polynomial,
154 unsigned int verticesPerCluster,
155 double relativeOverlap,
156 bool projectToInput,
158 bool computeEvaluationOffline)
159 : Mapping(constraint, dimension, false, Mapping::InitialGuessRequirement::None),
160 _basisFunction(function), _verticesPerCluster(verticesPerCluster), _relativeOverlap(relativeOverlap),
161 _projectToInput(projectToInput), _polynomial(polynomial), _useBatchedSolver(ginkgoParameter.executor != "cpu"),
162 _ginkgoParameter(ginkgoParameter), _computeEvaluationOffline(computeEvaluationOffline)
163{
164 PRECICE_ASSERT(this->getDimensions() <= 3);
165 PRECICE_ASSERT(_polynomial != Polynomial::ON, "Integrated polynomial is not supported for partition of unity data mappings.");
166 PRECICE_ASSERT(_relativeOverlap < 1, "The relative overlap has to be smaller than one.");
167 PRECICE_ASSERT(_verticesPerCluster > 0, "The number of vertices per cluster has to be greater zero.");
168#ifdef PRECICE_NO_KOKKOS_KERNELS
169 PRECICE_ASSERT(_useBatchedSolver == false, "Not implemented");
170#endif
171 if (isScaledConsistent()) {
174 } else {
177 }
178}
179
180template <typename RADIAL_BASIS_FUNCTION_T>
182{
184
185 precice::profiling::Event e("map.pou.computeMapping.From" + this->input()->getName() + "To" + this->output()->getName(), profiling::Synchronize);
186
187 // Recompute the whole clustering
188 PRECICE_ASSERT(!this->_hasComputedMapping, "Please clear the mapping before recomputing.");
189
190 mesh::PtrMesh inMesh;
191 mesh::PtrMesh outMesh;
193 inMesh = this->output();
194 outMesh = this->input();
195 } else { // Consistent or scaled consistent
196 inMesh = this->input();
197 outMesh = this->output();
198 }
199
200 precice::profiling::Event eClusters("map.pou.computeMapping.createClustering.From" + this->input()->getName() + "To" + this->output()->getName());
201 // Step 1: get a tentative clustering consisting of centers and a radius from one of the available algorithms
202 auto [clusterRadius, centerCandidates] = impl::createClustering(inMesh, outMesh, _relativeOverlap, _verticesPerCluster, _projectToInput);
203 eClusters.stop();
204 // Due to the aggressive filtering in the createClustering stage, the canddidates here are already final
205 e.addData("n clusters", centerCandidates.size());
206
207 _clusterRadius = clusterRadius;
208 PRECICE_ASSERT(_clusterRadius > 0 || inMesh->nVertices() == 0 || outMesh->nVertices() == 0);
209
210 if (_useBatchedSolver) {
211 PRECICE_CHECK(!(outMesh->isJustInTime() || inMesh->isJustInTime()), "Just-in-time mappings are not implemented for Kokkos- or Ginkgo-based solvers.");
212 precice::profiling::Event eBatched("map.pou.computeMapping.batchedSolver");
213 _batchedSolver = std::make_unique<BatchedRBFSolver<RADIAL_BASIS_FUNCTION_T>>(_basisFunction, inMesh, outMesh,
214 centerCandidates, _clusterRadius,
216
217 // For the batched solver, we don't register the _centerMesh as such
218 PRECICE_ASSERT(!_centerMesh, "The centerMesh is only utilized for the CPU variant");
219 } else {
220 _computeCPU(inMesh, outMesh, _clusterRadius, centerCandidates);
221 }
222 this->_hasComputedMapping = true;
223}
224
225template <typename RADIAL_BASIS_FUNCTION_T>
226void PartitionOfUnityMapping<RADIAL_BASIS_FUNCTION_T>::_computeCPU(mesh::PtrMesh inMesh, mesh::PtrMesh outMesh, double clusterRadius, const std::vector<mesh::Vertex> &centerCandidates)
227{
228 // Step 1: check, which of the resulting clusters are non-empty and register the cluster centers in a mesh
229 // Here, the VertexCluster computes the matrix decompositions directly in case the cluster is non-empty
230 _centerMesh = std::make_unique<mesh::Mesh>("pou-centers-" + inMesh->getName(), this->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
231 auto &meshVertices = _centerMesh->vertices();
232
233 meshVertices.clear();
234 _clusters.clear();
235 _clusters.reserve(centerCandidates.size());
236 for (const auto &c : centerCandidates) {
237 // We cannot simply copy the vertex from the container in order to fill the vertices of the centerMesh, as the vertexID of each center needs to match the index
238 // of the cluster within the _clusters vector. That's required for the indexing further down and asserted below
239 const VertexID vertexID = meshVertices.size();
240 mesh::Vertex center(c.getCoords(), vertexID);
242
243 // Consider only non-empty clusters (more of a safeguard here)
244 if (!cluster.empty()) {
245 PRECICE_ASSERT(center.getID() == static_cast<int>(_clusters.size()), center.getID(), _clusters.size());
246 meshVertices.emplace_back(std::move(center));
247 _clusters.emplace_back(std::move(cluster));
248 }
249 }
250
251 // Log the average number of resulting clusters
252 PRECICE_DEBUG("Partition of unity data mapping between mesh \"{}\" and mesh \"{}\": mesh \"{}\" on rank {} was decomposed into {} clusters.", this->input()->getName(), this->output()->getName(), inMesh->getName(), utils::IntraComm::getRank(), _clusters.size());
253
254 if (_clusters.size() > 0) {
255 PRECICE_DEBUG("Average number of vertices per cluster {}", std::transform_reduce(
256 _clusters.begin(), _clusters.end(), size_t{0}, std::plus<>(),
257 [](const auto &c) { return c.getNumberOfInputVertices(); }) /
258 _clusters.size());
259 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());
260 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());
261 }
262
263 precice::profiling::Event eWeights("map.pou.computeMapping.computeWeights");
264 // Log a bounding box of the center mesh
265 _centerMesh->computeBoundingBox();
266 PRECICE_DEBUG("Bounding Box of the cluster centers {}", _centerMesh->getBoundingBox());
267
268 // Step 2: Determine PU weights
269 PRECICE_DEBUG("Computing cluster-vertex association");
270 for (const auto &vertex : outMesh->vertices()) {
271 // we use a helper function, as we need the same functionality for just-in-time mapping
272 auto [clusterIDs, normalizedWeights] = computeNormalizedWeight(vertex, outMesh->getName());
273 // Step 3: store the normalized weight in all associated clusters
274 for (unsigned int i = 0; i < clusterIDs.size(); ++i) {
275 PRECICE_ASSERT(clusterIDs[i] < static_cast<int>(_clusters.size()));
276 _clusters[clusterIDs[i]].setNormalizedWeight(normalizedWeights[i], vertex.getID());
277 }
278 }
279 eWeights.stop();
280
281 // Uncomment to add a VTK export of the cluster center distribution for visualization purposes
282 // exportClusterCentersAsVTU(*_centerMesh);
283
284 // we need the center mesh index data structure
285 if (!outMesh->isJustInTime()) {
286 _centerMesh.reset();
287 }
288}
289
290template <typename RADIAL_BASIS_FUNCTION_T>
291std::pair<std::vector<int>, std::vector<double>> PartitionOfUnityMapping<RADIAL_BASIS_FUNCTION_T>::computeNormalizedWeight(const mesh::Vertex &vertex, std::string_view mesh)
292{
293
294 // Step 1: index the clusters / the center mesh in order to define the output vertex -> cluster ownership
295 // the ownership is required to compute the normalized partition of unity weights (Step 2)
296 // query::Index clusterIndex(*_centerMesh.get());
298 query::Index &clusterIndex = _centerMesh->index();
299
300 // Step 2: find all clusters the output vertex lies in, i.e., find all cluster centers which have the distance of a cluster radius from the given output vertex
301 // Here, we do this using the RTree on the centerMesh: VertexID (queried from the centersMesh) == clusterID, by construction above. The loop uses
302 // the vertices to compute the weights required for the partition of unity data mapping.
303 // Note: this could also be done on-the-fly in the map data phase for dynamic queries, which would require to make the mesh as well as the indexTree member variables.
304
305 // Step 2a: get the relevant clusters for the output vertex
306 auto clusterIDs = clusterIndex.getVerticesInsideBox(vertex, _clusterRadius);
307 const auto localNumberOfClusters = clusterIDs.size();
308
309 // Consider the case where we didn't find any cluster (meshes don't match very well)
310 //
311 // In principle, we could assign the vertex to the closest cluster using clusterIDs.emplace_back(clusterIndex.getClosestVertex(vertex.getCoords()).index);
312 // However, this leads to a conflict with weights already set in the corresponding cluster, since we insert the ID and, later on, map the ID to a local weight index
313 // Of course, we could rearrange the weights, but we want to avoid the case here anyway, i.e., prefer to abort.
314 PRECICE_CHECK(localNumberOfClusters > 0,
315 "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. "
316 "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\". "
317 "These options are only valid for the <mapping:rbf-pum-direct/> tag.",
318 vertex.getCoords(), mesh);
319
320 // Next we compute the normalized weights of each output vertex for each partition
321 PRECICE_ASSERT(localNumberOfClusters > 0, "No cluster found for vertex {}", vertex.getCoords());
322
323 // Step 2b: compute the weight in each partition individually and store them in 'weights'
324 std::vector<double> weights(localNumberOfClusters);
325 std::transform(clusterIDs.cbegin(), clusterIDs.cend(), weights.begin(), [&](const auto &ids) { return _clusters[ids].computeWeight(vertex); });
326 double weightSum = std::accumulate(weights.begin(), weights.end(), static_cast<double>(0.));
327 // TODO: This covers the edge case of vertices being at the edge of (several) clusters
328 // In case the sum is equal to zero, we assign equal weights for all clusters
329 if (weightSum <= 0) {
330 PRECICE_ASSERT(weights.size() > 0);
331 std::for_each(weights.begin(), weights.end(), [&weights](auto &w) { w = 1. / weights.size(); });
332 weightSum = 1;
333 }
334 PRECICE_ASSERT(weightSum > 0);
335
336 // Step 2c: Normalize weights
337 std::transform(weights.begin(), weights.end(), weights.begin(), [weightSum](double w) { return w / weightSum; });
338
339 // Return both the cluster IDs and the normalized weights
340 return {clusterIDs, weights};
341}
342
343template <typename RADIAL_BASIS_FUNCTION_T>
345{
347
348 precice::profiling::Event e("map.pou.mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
349
350 // Execute the actual mapping evaluation in all clusters
351 // 1. Assert that all output data values were reset, as we accumulate data in all clusters independently
352 PRECICE_ASSERT(outData.isZero());
353
354 if (_useBatchedSolver) {
355 PRECICE_ASSERT(false, "Not implemented");
356 } else {
357 // 2. Iterate over all clusters and accumulate the result in the output data
358 std::for_each(_clusters.begin(), _clusters.end(), [&](auto &cluster) { cluster.mapConservative(inData, outData); });
359 }
360}
361
362template <typename RADIAL_BASIS_FUNCTION_T>
364{
366
367 precice::profiling::Event e("map.pou.mapData.From" + input()->getName() + "To" + output()->getName(), profiling::Synchronize);
368
369 // Execute the actual mapping evaluation in all clusters
370 // 1. Assert that all output data values were reset, as we accumulate data in all clusters independently
371 PRECICE_ASSERT(outData.isZero());
372
373 if (_useBatchedSolver) {
374 PRECICE_ASSERT(_batchedSolver, "Not initialized");
375 _batchedSolver->solveConsistent(inData, outData);
376 } else {
377 // 2. Execute the actual mapping evaluation in all vertex clusters and accumulate the data
378 std::for_each(_clusters.begin(), _clusters.end(), [&](auto &clusters) { clusters.mapConsistent(inData, outData); });
379 }
380}
381
382template <typename RADIAL_BASIS_FUNCTION_T>
383void PartitionOfUnityMapping<RADIAL_BASIS_FUNCTION_T>::mapConservativeAt(const Eigen::Ref<const Eigen::MatrixXd> &coordinates, const Eigen::Ref<const Eigen::MatrixXd> &source,
384 impl::MappingDataCache &cache, Eigen::Ref<Eigen::MatrixXd>)
385{
386 precice::profiling::Event e("map.pou.mapConservativeAt.From" + input()->getName());
387 // @todo: it would most probably be more efficient to first group the vertices we receive here according to the clusters and then compute the solution
388
391 PRECICE_ASSERT(cache.p.size() == _clusters.size());
392 PRECICE_ASSERT(cache.polynomialContributions.size() == _clusters.size());
393 PRECICE_ASSERT(!_useBatchedSolver, "Not implemented"); // The PRECICE_CHECK is earlier
394
395 mesh::Vertex vertex(coordinates.col(0), -1);
396 for (Eigen::Index v = 0; v < coordinates.cols(); ++v) {
397 vertex.setCoords(coordinates.col(v));
398 auto [clusterIDs, normalizedWeights] = computeNormalizedWeight(vertex, this->input()->getName());
399 // Use the weight to interpolate the solution
400 for (std::size_t i = 0; i < clusterIDs.size(); ++i) {
401 PRECICE_ASSERT(clusterIDs[i] < static_cast<int>(_clusters.size()));
402 auto id = clusterIDs[i];
403 // the input mesh refers here to a consistent constraint
404 Eigen::VectorXd res = normalizedWeights[i] * source.col(v);
405 _clusters[id].addWriteDataToCache(vertex, res, cache.polynomialContributions[id], cache.p[id], *this->output().get());
406 }
407 }
408}
409
410template <typename RADIAL_BASIS_FUNCTION_T>
412{
414 PRECICE_ASSERT(!cache.p.empty());
416 precice::profiling::Event e("map.pou.completeJustInTimeMapping.From" + input()->getName());
417
418 for (std::size_t c = 0; c < _clusters.size(); ++c) {
419 // If there is no contribution, we don't have to evaluate
420 if (cache.p[c].squaredNorm() > 0) {
421 _clusters[c].evaluateConservativeCache(cache.polynomialContributions[c], cache.p[c], buffer);
422 }
423 }
424}
425
426template <typename RADIAL_BASIS_FUNCTION_T>
428{
431 cache.p.resize(_clusters.size());
432 cache.polynomialContributions.resize(_clusters.size());
433 for (std::size_t c = 0; c < _clusters.size(); ++c) {
434 _clusters[c].initializeCacheData(cache.polynomialContributions[c], cache.p[c], cache.getDataDimensions());
435 }
436}
437
438template <typename RADIAL_BASIS_FUNCTION_T>
440{
441 // We cannot synchronize this event, as the call to this function is rank-local only
442 precice::profiling::Event e("map.pou.updateMappingDataCache.From" + input()->getName());
443 PRECICE_ASSERT(cache.p.size() == _clusters.size());
444 PRECICE_ASSERT(cache.polynomialContributions.size() == _clusters.size());
445 Eigen::Map<const Eigen::MatrixXd> inMatrix(in.data(), cache.getDataDimensions(), in.size() / cache.getDataDimensions());
446 for (std::size_t c = 0; c < _clusters.size(); ++c) {
447 _clusters[c].computeCacheData(inMatrix, cache.polynomialContributions[c], cache.p[c]);
448 }
449}
450
451template <typename RADIAL_BASIS_FUNCTION_T>
452void PartitionOfUnityMapping<RADIAL_BASIS_FUNCTION_T>::mapConsistentAt(const Eigen::Ref<const Eigen::MatrixXd> &coordinates, const impl::MappingDataCache &cache, Eigen::Ref<Eigen::MatrixXd> values)
453{
454 precice::profiling::Event e("map.pou.mapConsistentAt.From" + input()->getName());
455 // @todo: it would most probably be more efficient to first group the vertices we receive here according to the clusters and then compute the solution
458 PRECICE_ASSERT(!_useBatchedSolver, "Not implemented"); // The PRECICE_CHECK is earlier
459
460 // First, make sure that everything is reset before we start
461 values.setZero();
462
463 mesh::Vertex vertex(coordinates.col(0), -1);
464 for (Eigen::Index v = 0; v < values.cols(); ++v) {
465 vertex.setCoords(coordinates.col(v));
466 auto [clusterIDs, normalizedWeights] = computeNormalizedWeight(vertex, this->output()->getName());
467 // Use the weight to interpolate the solution
468 for (std::size_t i = 0; i < clusterIDs.size(); ++i) {
469 PRECICE_ASSERT(clusterIDs[i] < static_cast<int>(_clusters.size()));
470 auto id = clusterIDs[i];
471 // the input mesh refers here to a consistent constraint
472 Eigen::VectorXd localRes = normalizedWeights[i] * _clusters[id].interpolateAt(vertex, cache.polynomialContributions[id], cache.p[id], *this->input().get());
473 values.col(v) += localRes;
474 }
475 }
476}
477
478template <typename RADIAL_BASIS_FUNCTION_T>
480{
482 mesh::PtrMesh filterMesh, outMesh;
484 filterMesh = this->output(); // remote
485 outMesh = this->input(); // local
486 } else {
487 filterMesh = this->input(); // remote
488 outMesh = this->output(); // local
489 }
490
491 if (outMesh->empty())
492 return; // Ranks not at the interface should never hold interface vertices
493
494 // The geometric filter of the repartitioning is always disabled for the PU-RBF.
495 // The main rationale: if we use only a fraction of the mesh then we might end up
496 // with too few vertices per rank and we cannot prevent too few vertices from being
497 // tagged, if we have filtered too much vertices beforehand. When using the filtering,
498 // the user could increase the safety-factor or disable the filtering, but that's
499 // a bit hard to understand for users. When no geometric filter is applid,
500 // vertices().size() is here the same as getGlobalNumberOfVertices. Hence, it is much
501 // safer to make use of the unfiltered mesh for the parallel tagging.
502 //
503 // Drawback: the "estimateClusterRadius" below makes use of the mesh R* index tree, and
504 // constructing the tree on the (unfiltered) global mesh is computationally expensive (O( N logN)).
505 // We could pre-filter the global mesh to a local fraction (using our own geometric filtering,
506 // maybe with an increased safety margin or even an iterative increase of the safety margin),
507 // but then there is again the question on how to do this in a safe way, without risking
508 // failures depending on the partitioning. So we stick here to the computationally more
509 // demanding, but safer version.
510 // See also https://github.com/precice/precice/pull/1912#issuecomment-2551143620
511
512 // Get the local bounding boxes
513 auto localBB = outMesh->getBoundingBox();
514 // we cannot check for empty'ness here, as a single output mesh vertex
515 // would lead to a 0D box with zero volume (considered empty). Thus, we
516 // simply check here for default'ness, which is equivalent to outMesh->empty()
517 // further above
518 PRECICE_ASSERT(!localBB.isDefault());
519
520 if (_clusterRadius == 0)
522
523 PRECICE_DEBUG("Cluster radius estimate: {}", _clusterRadius);
525
526 // Now we extend the bounding box by the radius
527 localBB.expandBy(2 * _clusterRadius);
528
529 // ... and tag all affected vertices
530 auto verticesNew = filterMesh->index().getVerticesInsideBox(localBB);
531
532 std::for_each(verticesNew.begin(), verticesNew.end(), [&filterMesh](VertexID v) { filterMesh->vertex(v).tag(); });
533}
534
535template <typename RADIAL_BASIS_FUNCTION_T>
537{
538 // Nothing to be done here. There is no global ownership for matrix entries required and we tag all potentially locally relevant vertices already in the first round.
539}
540
541template <typename RADIAL_BASIS_FUNCTION_T>
543{
545
546 auto dataRadius = centerMesh.createData("radius", 1, -1);
547 auto dataCardinality = centerMesh.createData("number-of-vertices", 1, -1);
548 centerMesh.allocateDataValues();
549 dataRadius->values().fill(_clusterRadius);
550 for (unsigned int i = 0; i < _clusters.size(); ++i) {
551 dataCardinality->values()[i] = static_cast<double>(_clusters[i].getNumberOfInputVertices());
552 }
553
554 // We have to create the global offsets in order to export things in parallel
556 // send number of vertices
557 PRECICE_DEBUG("Send number of vertices: {}", centerMesh.nVertices());
558 int numberOfVertices = centerMesh.nVertices();
559 utils::IntraComm::getCommunication()->send(numberOfVertices, 0);
560
561 // receive vertex offsets
562 mesh::Mesh::VertexOffsets vertexOffsets;
563 utils::IntraComm::getCommunication()->broadcast(vertexOffsets, 0);
564 PRECICE_DEBUG("Vertex offsets: {}", vertexOffsets);
565 PRECICE_ASSERT(centerMesh.getVertexOffsets().empty());
566 centerMesh.setVertexOffsets(std::move(vertexOffsets));
567 } else if (utils::IntraComm::isPrimary()) {
568
570 vertexOffsets[0] = centerMesh.nVertices();
571
572 // receive number of secondary vertices and fill vertex offsets
573 for (int secondaryRank : utils::IntraComm::allSecondaryRanks()) {
574 int numberOfSecondaryRankVertices = -1;
575 utils::IntraComm::getCommunication()->receive(numberOfSecondaryRankVertices, secondaryRank);
576 PRECICE_ASSERT(numberOfSecondaryRankVertices >= 0);
577 vertexOffsets[secondaryRank] = numberOfSecondaryRankVertices + vertexOffsets[secondaryRank - 1];
578 }
579
580 // broadcast vertex offsets
581 PRECICE_DEBUG("Vertex offsets: {}", centerMesh.getVertexOffsets());
582 utils::IntraComm::getCommunication()->broadcast(vertexOffsets);
583 centerMesh.setVertexOffsets(std::move(vertexOffsets));
584 }
585
586 dataRadius->setSampleAtTime(0, time::Sample{1, dataRadius->values()});
587 dataCardinality->setSampleAtTime(0, time::Sample{1, dataCardinality->values()});
589 exporter.doExport(0, 0.0);
590}
591
592template <typename RADIAL_BASIS_FUNCTION_T>
594{
596 _clusters.clear();
597 // TODO: Don't reset this here
598 _clusterRadius = 0;
599 this->_hasComputedMapping = false;
600}
601
602template <typename RADIAL_BASIS_FUNCTION_T>
604{
605 return "partition-of-unity RBF";
606}
607} // namespace mapping
608} // namespace precice
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
VertexContainer & vertices()
Returns modifieable container holding all vertices.
Definition Mesh.cpp:54
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:242
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:64
bool isJustInTime() const
Definition Mesh.hpp:340
bool empty() const
Does the mesh contain any vertices?
Definition Mesh.hpp:88
const BoundingBox & getBoundingBox() const
Returns the bounding box of the mesh.
Definition Mesh.cpp:410
void doExport(int index, double time) final override
Export the mesh and writes files.
Definition ExportXML.cpp:33
mesh::PtrMesh output() const
Returns pointer to output mesh.
Definition Mapping.cpp:92
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:30
Mapping(Constraint constraint, int dimensions, bool requiresGradientData, InitialGuessRequirement initialGuessRequirement)
Constructor, takes mapping constraint.
Definition Mapping.cpp:12
mesh::PtrMesh input() const
Returns pointer to input mesh.
Definition Mapping.cpp:87
bool _hasComputedMapping
Flag to indicate whether computeMapping() has been called.
Definition Mapping.hpp:307
bool isScaledConsistent() const
Returns true if mapping is a form of scaled consistent mapping.
Definition Mapping.cpp:258
void setInputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the input mesh.
Definition Mapping.cpp:97
void setOutputRequirement(MeshRequirement requirement)
Sets the mesh requirement for the output mesh.
Definition Mapping.cpp:103
virtual bool hasConstraint(const Constraint &constraint) const
Checks whether the mapping has the given constraint or not.
Definition Mapping.cpp:248
InitialGuessRequirement
Specifies whether the mapping requires an initial guess.
Definition Mapping.hpp:64
std::unique_ptr< BatchedRBFSolver< RBF > > _batchedSolver
PartitionOfUnityMapping(Mapping::Constraint constraint, int dimension, RADIAL_BASIS_FUNCTION_T function, Polynomial polynomial, unsigned int verticesPerCluster, double relativeOverlap, bool projectToInput, MappingConfiguration::GinkgoParameter ginkgoParameter=MappingConfiguration::GinkgoParameter(), bool computeEvaluationOffline=false)
std::pair< std::vector< int >, std::vector< double > > computeNormalizedWeight(const mesh::Vertex &v, std::string_view mesh)
void initializeMappingDataCache(impl::MappingDataCache &cache) final override
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
MappingConfiguration::GinkgoParameter _ginkgoParameter
std::vector< SphericalVertexCluster< RBF > > _clusters
std::string getName() const final override
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) override
Maps data using a consistent constraint.
void _computeCPU(mesh::PtrMesh inMesh, mesh::PtrMesh outMesh, double clusterRadius, const std::vector< mesh::Vertex > &centerCandidates)
void mapConsistentAt(const Eigen::Ref< const Eigen::MatrixXd > &coordinates, const impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > values) final override
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 completeJustInTimeMapping(impl::MappingDataCache &cache, Eigen::Ref< Eigen::MatrixXd > buffer) final override
Container and creator for meshes.
Definition Mesh.hpp:38
static constexpr MeshID MESH_ID_UNDEFINED
Use if the id of the mesh is not necessary.
Definition Mesh.hpp:57
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:64
PtrData & createData(const std::string &name, int dimension, DataID id, int waveformDegree=time::Time::DEFAULT_WAVEFORM_DEGREE)
Create only data for vertex.
Definition Mesh.cpp:152
const VertexOffsets & getVertexOffsets() const
Definition Mesh.hpp:259
std::vector< int > VertexOffsets
Definition Mesh.hpp:54
void setVertexOffsets(VertexOffsets vertexOffsets)
Only used for tests.
Definition Mesh.hpp:268
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
Definition Mesh.cpp:257
Vertex of a mesh.
Definition Vertex.hpp:16
VertexID getID() const
Returns the unique (among vertices of one mesh on one processor) ID of the vertex.
Definition Vertex.hpp:109
void setCoords(const VECTOR_T &coordinates)
Sets the coordinates of the vertex.
Definition Vertex.hpp:101
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
Definition Vertex.hpp:114
void stop()
Stops a running event.
Definition Event.cpp:51
void addData(std::string_view key, int value)
Adds named integer data, associated to an event.
Definition Event.cpp:63
Class to query the index trees of the mesh.
Definition Index.hpp:64
std::vector< VertexID > getVerticesInsideBox(const mesh::Vertex &centerVertex, double radius)
Return all the vertices inside the box formed by vertex and radius (boundary exclusive)
Definition Index.cpp:225
static int getSize()
Number of ranks. This includes ranks from both participants, e.g. minimal size is 2.
Definition IntraComm.cpp:47
static Rank getRank()
Current rank.
Definition IntraComm.cpp:42
static bool isPrimary()
True if this process is running the primary rank.
Definition IntraComm.cpp:52
static auto allSecondaryRanks()
Returns an iterable range over salve ranks [1, _size)
Definition IntraComm.hpp:37
static bool isSecondary()
True if this process is running a secondary rank.
Definition IntraComm.cpp:57
static com::PtrCommunication & getCommunication()
Intra-participant communication.
Definition IntraComm.hpp:31
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.
Definition Event.hpp:28
Main namespace of the precice library.
int VertexID
Definition Types.hpp:13
bool syncMode
Enabled further inter- and intra-solver synchronisation.
STL namespace.
std::vector< Eigen::MatrixXd > polynomialContributions
int getDataDimensions() const
Returns the number of data components.