preCICE
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
61 std::unique_ptr<SOLVER_T> _rbfSolver;
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
73 std::tuple<Args...> optionalArgs;
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 (sizeof...(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... Args>
178static std::string getNameWithArgs(const std::tuple<Args...> &optionalArgs)
179{
180 if constexpr (sizeof...(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>
198
199template <typename SOLVER_T, typename... Args>
201{
203 precice::profiling::Event e("map.rbf.mapData.From" + this->input()->getName() + "To" + this->output()->getName(), profiling::Synchronize);
204
205 PRECICE_DEBUG("Map conservative using {}", getName());
206
207 // Gather input data
209
210 const auto &localInData = inData.values;
211
212 int localOutputSize = 0;
213 for (const auto &vertex : this->output()->vertices()) {
214 if (vertex.isOwner()) {
215 ++localOutputSize;
216 }
217 }
218
219 localOutputSize *= inData.dataDims;
220
221 utils::IntraComm::getCommunication()->sendRange(localInData, 0);
222 utils::IntraComm::getCommunication()->send(localOutputSize, 0);
223
224 } else { // Parallel Primary rank or Serial case
225
226 std::vector<double> globalInValues;
227 std::vector<double> outputValueSizes;
228 {
229 const auto &localInData = inData.values;
230 globalInValues.insert(globalInValues.begin(), localInData.data(), localInData.data() + localInData.size());
231
232 int localOutputSize = 0;
233 for (const auto &vertex : this->output()->vertices()) {
234 if (vertex.isOwner()) {
235 ++localOutputSize;
236 }
237 }
238
239 localOutputSize *= inData.dataDims;
240
241 outputValueSizes.push_back(localOutputSize);
242 }
243
244 {
245 int secondaryOutputValueSize;
247 std::vector<double> secondaryBuffer = utils::IntraComm::getCommunication()->receiveRange(rank, com::asVector<double>);
248 globalInValues.insert(globalInValues.end(), secondaryBuffer.begin(), secondaryBuffer.end());
249
250 utils::IntraComm::getCommunication()->receive(secondaryOutputValueSize, rank);
251 outputValueSizes.push_back(secondaryOutputValueSize);
252 }
253 }
254
255 const int valueDim = inData.dataDims;
256
257 using RowMatrixXd = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
258 using MapToMatrix = Eigen::Map<RowMatrixXd, Eigen::Unaligned, Eigen::OuterStride<Eigen::Dynamic>>;
259
260 // copy of input data to Eigen matrix format, rows == outputSize
261 Eigen::MatrixXd in = MapToMatrix(globalInValues.data(), _rbfSolver->getOutputSize(), valueDim, Eigen::OuterStride(valueDim));
262
263 // copy all output values without polynomial entries from col-major to row-major format
264 RowMatrixXd out = _rbfSolver->solveConservative(in, _polynomial).block(0, 0, this->output()->getGlobalNumberOfVertices(), valueDim);
265
266 Eigen::Map<Eigen::VectorXd> outputValues(out.data(), (this->output()->getGlobalNumberOfVertices()) * valueDim);
267
268 // Data scattering to secondary ranks
270
271 // Filter data
272 int outputCounter = 0;
273 for (int i = 0; i < static_cast<int>(this->output()->nVertices()); ++i) {
274 if (this->output()->vertex(i).isOwner()) {
275 for (int dim = 0; dim < valueDim; ++dim) {
276 outData[i * valueDim + dim] = outputValues(outputCounter);
277 ++outputCounter;
278 }
279 }
280 }
281
282 // Data scattering to secondary ranks
283 int beginPoint = outputValueSizes.at(0);
285 precice::span<const double> toSend{outputValues.data() + beginPoint, static_cast<size_t>(outputValueSizes.at(rank))};
286 utils::IntraComm::getCommunication()->sendRange(toSend, rank);
287 beginPoint += outputValueSizes.at(rank);
288 }
289 } else { // Serial
290 outData = outputValues;
291 }
292 }
294 std::vector<double> receivedValues = utils::IntraComm::getCommunication()->receiveRange(0, com::asVector<double>);
295
296 const int valueDim = inData.dataDims;
297
298 int outputCounter = 0;
299 for (int i = 0; i < static_cast<int>(this->output()->nVertices()); ++i) {
300 if (this->output()->vertex(i).isOwner()) {
301 for (int dim = 0; dim < valueDim; ++dim) {
302 outData[i * valueDim + dim] = receivedValues.at(outputCounter);
303 ++outputCounter;
304 }
305 }
306 }
307 }
308}
309
310template <typename SOLVER_T, typename... Args>
312{
314 precice::profiling::Event e("map.rbf.mapData.From" + this->input()->getName() + "To" + this->output()->getName(), profiling::Synchronize);
315
316 PRECICE_DEBUG("Map {} using {}", (this->hasConstraint(Mapping::CONSISTENT) ? "consistent" : "scaled-consistent"), getName());
317
318 // Gather input data
320 // Input data is filtered
321 auto localInDataFiltered = this->input()->getOwnedVertexData(inData.values);
322 int localOutputSize = outData.size();
323
324 // Send data and output size
325 utils::IntraComm::getCommunication()->sendRange(localInDataFiltered, 0);
326 utils::IntraComm::getCommunication()->send(localOutputSize, 0);
327
328 } else { // Primary rank or Serial case
329
330 const int valueDim = inData.dataDims;
331
332 std::vector<double> globalInValues(static_cast<std::size_t>(this->input()->getGlobalNumberOfVertices()) * valueDim, 0.0);
333 std::vector<int> outValuesSize;
334
335 if (utils::IntraComm::isPrimary()) { // Parallel case
336
337 // Filter input data
338 const auto &localInData = this->input()->getOwnedVertexData(inData.values);
339 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.begin());
340 outValuesSize.push_back(outData.size());
341
342 int inputSizeCounter = localInData.size();
343 int secondaryOutDataSize{0};
344
346 std::vector<double> secondaryBuffer = utils::IntraComm::getCommunication()->receiveRange(rank, com::asVector<double>);
347 std::copy(secondaryBuffer.begin(), secondaryBuffer.end(), globalInValues.begin() + inputSizeCounter);
348 inputSizeCounter += secondaryBuffer.size();
349
350 utils::IntraComm::getCommunication()->receive(secondaryOutDataSize, rank);
351 outValuesSize.push_back(secondaryOutDataSize);
352 }
353
354 } else { // Serial case
355 const auto &localInData = inData.values;
356 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.begin());
357 outValuesSize.push_back(outData.size());
358 }
359
360 // copy of input data to Eigen matrix format
361 using RowMatrixXd = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
362 using MapToMatrix = Eigen::Map<RowMatrixXd, Eigen::Unaligned, Eigen::OuterStride<Eigen::Dynamic>>;
363
364 // input matrix with polynomial entries that remain zero
365 Eigen::MatrixXd in = Eigen::MatrixXd::Zero(_rbfSolver->getInputSize(), valueDim);
366
367 in.block(0, 0, this->input()->getGlobalNumberOfVertices(), valueDim) = MapToMatrix(globalInValues.data(), this->input()->getGlobalNumberOfVertices(), valueDim, Eigen::OuterStride(valueDim));
368
369 // copy all output values from col-major to row-major format
370 RowMatrixXd out = _rbfSolver->solveConsistent(in, _polynomial);
371 // copy mapped data at correct position to output data values
372 outData = Eigen::Map<Eigen::VectorXd>(out.data(), outValuesSize.at(0));
373
374 // Data scattering to secondary ranks
375 int beginPoint = outValuesSize.at(0);
376
379 precice::span<const double> toSend{out.data() + beginPoint, static_cast<size_t>(outValuesSize.at(rank))};
380 utils::IntraComm::getCommunication()->sendRange(toSend, rank);
381 beginPoint += outValuesSize.at(rank);
382 }
383 }
384 }
386 std::vector<double> receivedValues = utils::IntraComm::getCommunication()->receiveRange(0, com::asVector<double>);
387 outData = Eigen::Map<Eigen::VectorXd>(receivedValues.data(), receivedValues.size());
388 }
389}
390} // 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
int getDimensions() const
Definition Mesh.cpp:99
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:242
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:354
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
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
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.
static std::string getNameWithArgs(const std::tuple< Args... > &optionalArgs)
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:28
Main namespace of the precice library.
int Rank
Definition Types.hpp:37
STL namespace.
int dataDims
The dimensionality of the data.
Definition Sample.hpp:60
Eigen::VectorXd values
Definition Sample.hpp:64