preCICE v3.2.0
Loading...
Searching...
No Matches
RadialBasisFctMapping.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Cholesky>
4#include <Eigen/Core>
5
7#include "com/Extra.hpp"
10#include "mesh/Filter.hpp"
12#include "profiling/Event.hpp"
13#include "utils/IntraComm.hpp"
14
15namespace precice::mapping {
16
28template <typename SOLVER_T, typename... Args>
29class RadialBasisFctMapping : public RadialBasisFctBaseMapping<typename SOLVER_T::BASIS_FUNCTION_T> {
30public:
31 using RADIAL_BASIS_FUNCTION_T = typename SOLVER_T::BASIS_FUNCTION_T;
41 Mapping::Constraint constraint,
42 int dimensions,
44 std::array<bool, 3> deadAxis,
45 Polynomial polynomial,
46 Args... args);
47
49 void computeMapping() final override;
50
52 void clear() final override;
53
55 std::string getName() const final override;
56
57private:
58 precice::logging::Logger _log{"mapping::RadialBasisFctMapping"};
59
60 // The actual solver
62
64 void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) final override;
65
67 void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) final override;
68
71
74};
75
76// --------------------------------------------------- HEADER IMPLEMENTATIONS
77
78template <typename SOLVER_T, typename... Args>
80 Mapping::Constraint constraint,
81 int dimensions,
83 std::array<bool, 3> deadAxis,
84 Polynomial polynomial,
85 Args... args)
86 : RadialBasisFctBaseMapping<RADIAL_BASIS_FUNCTION_T>(constraint, dimensions, function, deadAxis, Mapping::InitialGuessRequirement::None),
87 _polynomial(polynomial),
88 optionalArgs(std::make_tuple(std::forward<Args>(args)...))
89{
90 PRECICE_CHECK(!(RADIAL_BASIS_FUNCTION_T::isStrictlyPositiveDefinite() && polynomial == Polynomial::ON), "The integrated polynomial (polynomial=\"on\") is not supported for the selected radial-basis function. Please select another radial-basis function or change the polynomial configuration.");
91}
92
93template <typename SOLVER_T, typename... Args>
95{
97
98 precice::profiling::Event e("map.rbf.computeMapping.From" + this->input()->getName() + "To" + this->output()->getName(), profiling::Synchronize);
99
100 PRECICE_ASSERT(this->input()->getDimensions() == this->output()->getDimensions(),
101 this->input()->getDimensions(), this->output()->getDimensions());
102 PRECICE_ASSERT(this->getDimensions() == this->output()->getDimensions(),
103 this->getDimensions(), this->output()->getDimensions());
104
105 mesh::PtrMesh inMesh;
106 mesh::PtrMesh outMesh;
107
109 inMesh = this->output();
110 outMesh = this->input();
111 } else { // Consistent or scaled consistent
112 inMesh = this->input();
113 outMesh = this->output();
114 }
115
117
118 // Input mesh may have overlaps
119 mesh::Mesh filteredInMesh("filteredInMesh", inMesh->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
120 mesh::filterMesh(filteredInMesh, *inMesh, [&](const mesh::Vertex &v) { return v.isOwner(); });
121
122 // Send the mesh
125
126 } else { // Parallel Primary rank or Serial
127
128 mesh::Mesh globalInMesh(inMesh->getName(), inMesh->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
129 mesh::Mesh globalOutMesh(outMesh->getName(), outMesh->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
130
132 {
133 // Input mesh may have overlaps
134 mesh::Mesh filteredInMesh("filteredInMesh", inMesh->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
135 mesh::filterMesh(filteredInMesh, *inMesh, [&](const mesh::Vertex &v) { return v.isOwner(); });
136 globalInMesh.addMesh(filteredInMesh);
137 globalOutMesh.addMesh(*outMesh);
138 }
139
140 // Receive mesh
141 for (Rank secondaryRank : utils::IntraComm::allSecondaryRanks()) {
142 mesh::Mesh secondaryInMesh(inMesh->getName(), inMesh->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
143 com::receiveMesh(*utils::IntraComm::getCommunication(), secondaryRank, secondaryInMesh);
144 globalInMesh.addMesh(secondaryInMesh);
145
146 mesh::Mesh secondaryOutMesh(outMesh->getName(), outMesh->getDimensions(), mesh::Mesh::MESH_ID_UNDEFINED);
147 com::receiveMesh(*utils::IntraComm::getCommunication(), secondaryRank, secondaryOutMesh);
148 globalOutMesh.addMesh(secondaryOutMesh);
149 }
150
151 } else { // Serial
152 globalInMesh.addMesh(*inMesh);
153 globalOutMesh.addMesh(*outMesh);
154 }
155
156 // Forwarding the tuples here requires some template magic I don't want to implement
157 if constexpr (std::tuple_size_v < std::tuple < Args... >>> 0) {
158 _rbfSolver = std::make_unique<SOLVER_T>(this->_basisFunction, globalInMesh, boost::irange<Eigen::Index>(0, globalInMesh.nVertices()),
159 globalOutMesh, boost::irange<Eigen::Index>(0, globalOutMesh.nVertices()), this->_deadAxis, _polynomial, std::get<0>(optionalArgs));
160 } else {
161 _rbfSolver = std::make_unique<SOLVER_T>(this->_basisFunction, globalInMesh, boost::irange<Eigen::Index>(0, globalInMesh.nVertices()),
162 globalOutMesh, boost::irange<Eigen::Index>(0, globalOutMesh.nVertices()), this->_deadAxis, _polynomial);
163 }
164 }
165 this->_hasComputedMapping = true;
166 PRECICE_DEBUG("Compute Mapping is Completed.");
167}
168
169template <typename SOLVER_T, typename... Args>
176
177template <typename SOLVER_T, typename... Args>
179{
180 if constexpr (std::tuple_size_v < std::tuple < Args... >>> 0) {
181 auto param = std::get<0>(optionalArgs);
182 std::string exec = param.executor;
183 if (param.solver == "qr-solver") {
184 return "global-direct RBF (" + exec + ")";
185 } else {
186 return "global-iterative RBF (" + exec + ")";
187 }
188 } else {
189 return "global-direct RBF (cpu-executor)";
190 }
191}
192
193template <typename SOLVER_T, typename... Args>
195{
197 precice::profiling::Event e("map.rbf.mapData.From" + this->input()->getName() + "To" + this->output()->getName(), profiling::Synchronize);
198
199 PRECICE_DEBUG("Map conservative using {}", getName());
200
201 // Gather input data
203
204 const auto &localInData = inData.values;
205
206 int localOutputSize = 0;
207 for (const auto &vertex : this->output()->vertices()) {
208 if (vertex.isOwner()) {
209 ++localOutputSize;
210 }
211 }
212
213 localOutputSize *= inData.dataDims;
214
215 utils::IntraComm::getCommunication()->sendRange(localInData, 0);
216 utils::IntraComm::getCommunication()->send(localOutputSize, 0);
217
218 } else { // Parallel Primary rank or Serial case
219
220 std::vector<double> globalInValues;
221 std::vector<double> outputValueSizes;
222 {
223 const auto &localInData = inData.values;
224 globalInValues.insert(globalInValues.begin(), localInData.data(), localInData.data() + localInData.size());
225
226 int localOutputSize = 0;
227 for (const auto &vertex : this->output()->vertices()) {
228 if (vertex.isOwner()) {
229 ++localOutputSize;
230 }
231 }
232
233 localOutputSize *= inData.dataDims;
234
235 outputValueSizes.push_back(localOutputSize);
236 }
237
238 {
239 int secondaryOutputValueSize;
242 globalInValues.insert(globalInValues.end(), secondaryBuffer.begin(), secondaryBuffer.end());
243
244 utils::IntraComm::getCommunication()->receive(secondaryOutputValueSize, rank);
245 outputValueSizes.push_back(secondaryOutputValueSize);
246 }
247 }
248
249 const int valueDim = inData.dataDims;
250
251 // Construct Eigen vectors
252 Eigen::Map<Eigen::VectorXd> inputValues(globalInValues.data(), globalInValues.size());
253 Eigen::VectorXd outputValues((this->output()->getGlobalNumberOfVertices()) * valueDim);
254
255 Eigen::VectorXd in; // rows == outputSize
256 in.resize(_rbfSolver->getOutputSize()); // rows == outputSize
257
258 outputValues.setZero();
259
260 for (int dim = 0; dim < valueDim; dim++) {
261 for (int i = 0; i < in.size(); i++) { // Fill input data values
262 in[i] = inputValues(i * valueDim + dim);
263 }
264
265 Eigen::VectorXd out;
266 out = _rbfSolver->solveConservative(in, _polynomial);
267
268 // Copy mapped data to output data values
269 for (int i = 0; i < this->output()->getGlobalNumberOfVertices(); i++) {
270 outputValues[i * valueDim + dim] = out[i];
271 }
272 }
273
274 // Data scattering to secondary ranks
276
277 // Filter data
278 int outputCounter = 0;
279 for (int i = 0; i < static_cast<int>(this->output()->nVertices()); ++i) {
280 if (this->output()->vertex(i).isOwner()) {
281 for (int dim = 0; dim < valueDim; ++dim) {
282 outData[i * valueDim + dim] = outputValues(outputCounter);
283 ++outputCounter;
284 }
285 }
286 }
287
288 // Data scattering to secondary ranks
289 int beginPoint = outputValueSizes.at(0);
291 precice::span<const double> toSend{outputValues.data() + beginPoint, static_cast<size_t>(outputValueSizes.at(rank))};
292 utils::IntraComm::getCommunication()->sendRange(toSend, rank);
293 beginPoint += outputValueSizes.at(rank);
294 }
295 } else { // Serial
296 outData = outputValues;
297 }
298 }
301
302 const int valueDim = inData.dataDims;
303
304 int outputCounter = 0;
305 for (int i = 0; i < static_cast<int>(this->output()->nVertices()); ++i) {
306 if (this->output()->vertex(i).isOwner()) {
307 for (int dim = 0; dim < valueDim; ++dim) {
308 outData[i * valueDim + dim] = receivedValues.at(outputCounter);
309 ++outputCounter;
310 }
311 }
312 }
313 }
314}
315
316template <typename SOLVER_T, typename... Args>
318{
320 precice::profiling::Event e("map.rbf.mapData.From" + this->input()->getName() + "To" + this->output()->getName(), profiling::Synchronize);
321
322 PRECICE_DEBUG("Map {} using {}", (this->hasConstraint(Mapping::CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
323
324 // Gather input data
326 // Input data is filtered
327 auto localInDataFiltered = this->input()->getOwnedVertexData(inData.values);
328 int localOutputSize = outData.size();
329
330 // Send data and output size
331 utils::IntraComm::getCommunication()->sendRange(localInDataFiltered, 0);
332 utils::IntraComm::getCommunication()->send(localOutputSize, 0);
333
334 } else { // Primary rank or Serial case
335
336 const int valueDim = inData.dataDims;
337
338 std::vector<double> globalInValues(static_cast<std::size_t>(this->input()->getGlobalNumberOfVertices()) * valueDim, 0.0);
339 std::vector<int> outValuesSize;
340
341 if (utils::IntraComm::isPrimary()) { // Parallel case
342
343 // Filter input data
344 const auto &localInData = this->input()->getOwnedVertexData(inData.values);
345 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.begin());
346 outValuesSize.push_back(outData.size());
347
348 int inputSizeCounter = localInData.size();
349 int secondaryOutDataSize{0};
350
353 std::copy(secondaryBuffer.begin(), secondaryBuffer.end(), globalInValues.begin() + inputSizeCounter);
354 inputSizeCounter += secondaryBuffer.size();
355
356 utils::IntraComm::getCommunication()->receive(secondaryOutDataSize, rank);
357 outValuesSize.push_back(secondaryOutDataSize);
358 }
359
360 } else { // Serial case
361 const auto &localInData = inData.values;
362 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.begin());
363 outValuesSize.push_back(outData.size());
364 }
365
366 Eigen::VectorXd in;
367
368 in.resize(_rbfSolver->getInputSize()); // rows == n
369 in.setZero();
370
371 // Construct Eigen vectors
372 Eigen::Map<Eigen::VectorXd> inputValues(globalInValues.data(), globalInValues.size());
373
374 Eigen::VectorXd outputValues;
375 outputValues.resize((_rbfSolver->getOutputSize()) * valueDim);
376
377 Eigen::VectorXd out;
378 outputValues.setZero();
379
380 // For every data dimension, perform mapping
381 for (int dim = 0; dim < valueDim; dim++) {
382 // Fill input from input data values (last polyparams entries remain zero)
383 for (int i = 0; i < this->input()->getGlobalNumberOfVertices(); i++) {
384 in[i] = inputValues[i * valueDim + dim];
385 }
386
387 out = _rbfSolver->solveConsistent(in, _polynomial);
388
389 // Copy mapped data to output data values
390 for (int i = 0; i < out.size(); i++) {
391 outputValues[i * valueDim + dim] = out[i];
392 }
393 }
394
395 outData = Eigen::Map<Eigen::VectorXd>(outputValues.data(), outValuesSize.at(0));
396
397 // Data scattering to secondary ranks
398 int beginPoint = outValuesSize.at(0);
399
402 precice::span<const double> toSend{outputValues.data() + beginPoint, static_cast<size_t>(outValuesSize.at(rank))};
403 utils::IntraComm::getCommunication()->sendRange(toSend, rank);
404 beginPoint += outValuesSize.at(rank);
405 }
406 }
407 }
410 outData = Eigen::Map<Eigen::VectorXd>(receivedValues.data(), receivedValues.size());
411 }
412}
413} // namespace precice::mapping
#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
T at(T... args)
T begin(T... args)
Abstract base class for mapping of data from one mesh to another.
Definition Mapping.hpp:16
mesh::PtrMesh output() const
Returns pointer to output mesh.
Definition Mapping.cpp:92
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:30
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
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
RadialBasisFctBaseMapping(Constraint constraint, int dimensions, const SOLVER_T::BASIS_FUNCTION_T &function, std::array< bool, 3 > deadAxis, InitialGuessRequirement mappingType)
typename RadialBasisFctSolver< RBF >::BASIS_FUNCTION_T RADIAL_BASIS_FUNCTION_T
void computeMapping() final override
Computes the mapping coefficients from the in- and output mesh.
void mapConsistent(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a consistent constraint.
RadialBasisFctMapping(Mapping::Constraint constraint, int dimensions, RADIAL_BASIS_FUNCTION_T function, std::array< bool, 3 > deadAxis, Polynomial polynomial, Args... args)
Constructor.
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a conservative constraint.
Container and creator for meshes.
Definition Mesh.hpp:38
void addMesh(Mesh &deltaMesh)
Definition Mesh.cpp:332
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:65
Vertex of a mesh.
Definition Vertex.hpp:16
bool isOwner() const
Definition Vertex.cpp:22
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
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
T copy(T... args)
T data(T... args)
T end(T... args)
T insert(T... args)
T make_unique(T... args)
void sendMesh(Communication &communication, int rankReceiver, const mesh::Mesh &mesh)
Definition Extra.cpp:8
constexpr auto asVector
Allows to use Communication::AsVectorTag in a less verbose way.
void receiveMesh(Communication &communication, int rankSender, mesh::Mesh &mesh)
Definition Extra.cpp:13
contains the logging framework.
contains data mapping from points to meshes.
Polynomial
How to handle the polynomial?
void filterMesh(Mesh &destination, const Mesh &source, UnaryPredicate p)
Definition Filter.hpp:16
std::shared_ptr< Mesh > PtrMesh
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:21
Main namespace of the precice library.
int Rank
Definition Types.hpp:37
STL namespace.
T push_back(T... args)
T resize(T... args)
T size(T... args)
int dataDims
The dimensionality of the data.
Definition Sample.hpp:60
Eigen::VectorXd values
Definition Sample.hpp:64
T tuple_size_v