3#include <Eigen/Cholesky>
28template <
typename SOLVER_T,
typename... Args>
44 std::array<bool, 3> deadAxis,
78template <
typename SOLVER_T,
typename... Args>
83 std::array<bool, 3> deadAxis,
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.");
93template <
typename SOLVER_T,
typename... Args>
110 outMesh = this->
input();
112 inMesh = this->
input();
136 globalInMesh.
addMesh(filteredInMesh);
137 globalOutMesh.
addMesh(*outMesh);
144 globalInMesh.
addMesh(secondaryInMesh);
148 globalOutMesh.
addMesh(secondaryOutMesh);
153 globalOutMesh.
addMesh(*outMesh);
157 if constexpr (
sizeof...(Args) > 0) {
162 globalOutMesh, boost::irange<Eigen::Index>(0, globalOutMesh.
nVertices()), this->_deadAxis,
_polynomial);
169template <
typename SOLVER_T,
typename... Args>
177template <
typename... Args>
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 +
")";
186 return "global-iterative RBF (" + exec +
")";
189 return "global-direct RBF (cpu-executor)";
193template <
typename SOLVER_T,
typename... Args>
199template <
typename SOLVER_T,
typename... Args>
210 const auto &localInData = inData.
values;
212 int localOutputSize = 0;
213 for (
const auto &vertex : this->
output()->vertices()) {
214 if (vertex.isOwner()) {
226 std::vector<double> globalInValues;
227 std::vector<double> outputValueSizes;
229 const auto &localInData = inData.
values;
230 globalInValues.insert(globalInValues.begin(), localInData.data(), localInData.data() + localInData.size());
232 int localOutputSize = 0;
233 for (
const auto &vertex : this->
output()->vertices()) {
234 if (vertex.isOwner()) {
241 outputValueSizes.push_back(localOutputSize);
245 int secondaryOutputValueSize;
248 globalInValues.insert(globalInValues.end(), secondaryBuffer.begin(), secondaryBuffer.end());
251 outputValueSizes.push_back(secondaryOutputValueSize);
255 const int valueDim = inData.
dataDims;
257 using RowMatrixXd = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
258 using MapToMatrix = Eigen::Map<RowMatrixXd, Eigen::Unaligned, Eigen::OuterStride<Eigen::Dynamic>>;
261 Eigen::MatrixXd in = MapToMatrix(globalInValues.data(),
_rbfSolver->getOutputSize(), valueDim, Eigen::OuterStride(valueDim));
264 RowMatrixXd out =
_rbfSolver->solveConservative(in,
_polynomial).block(0, 0, this->
output()->getGlobalNumberOfVertices(), valueDim);
266 Eigen::Map<Eigen::VectorXd> outputValues(out.data(), (this->output()->getGlobalNumberOfVertices()) * valueDim);
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);
283 int beginPoint = outputValueSizes.at(0);
287 beginPoint += outputValueSizes.at(rank);
290 outData = outputValues;
296 const int valueDim = inData.
dataDims;
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);
310template <
typename SOLVER_T,
typename... Args>
321 auto localInDataFiltered = this->
input()->getOwnedVertexData(inData.
values);
322 int localOutputSize = outData.size();
330 const int valueDim = inData.
dataDims;
332 std::vector<double> globalInValues(
static_cast<std::size_t
>(this->
input()->getGlobalNumberOfVertices()) * valueDim, 0.0);
333 std::vector<int> outValuesSize;
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());
342 int inputSizeCounter = localInData.size();
343 int secondaryOutDataSize{0};
347 std::copy(secondaryBuffer.begin(), secondaryBuffer.end(), globalInValues.begin() + inputSizeCounter);
348 inputSizeCounter += secondaryBuffer.size();
351 outValuesSize.push_back(secondaryOutDataSize);
355 const auto &localInData = inData.
values;
356 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.begin());
357 outValuesSize.push_back(outData.size());
361 using RowMatrixXd = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
362 using MapToMatrix = Eigen::Map<RowMatrixXd, Eigen::Unaligned, Eigen::OuterStride<Eigen::Dynamic>>;
365 Eigen::MatrixXd in = Eigen::MatrixXd::Zero(
_rbfSolver->getInputSize(), valueDim);
367 in.block(0, 0, this->
input()->getGlobalNumberOfVertices(), valueDim) = MapToMatrix(globalInValues.data(), this->input()->getGlobalNumberOfVertices(), valueDim, Eigen::OuterStride(valueDim));
372 outData = Eigen::Map<Eigen::VectorXd>(out.data(), outValuesSize.at(0));
375 int beginPoint = outValuesSize.at(0);
381 beginPoint += outValuesSize.at(rank);
387 outData = Eigen::Map<Eigen::VectorXd>(receivedValues.data(), receivedValues.size());
#define PRECICE_DEBUG(...)
#define PRECICE_TRACE(...)
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
int getDimensions() const
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Abstract base class for mapping of data from one mesh to another.
mesh::PtrMesh output() const
Returns pointer to output mesh.
Constraint
Specifies additional constraints for a mapping.
int getDimensions() const
mesh::PtrMesh input() const
Returns pointer to input mesh.
bool _hasComputedMapping
Flag to indicate whether computeMapping() has been called.
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.
RadialBasisFctBaseMapping(Constraint constraint, int dimensions, const SOLVER_T::BASIS_FUNCTION_T &function, std::array< bool, 3 > deadAxis, InitialGuessRequirement mappingType)
SOLVER_T::BASIS_FUNCTION_T _basisFunction
typename RadialBasisFctSolver< RBF >::BASIS_FUNCTION_T RADIAL_BASIS_FUNCTION_T
std::unique_ptr< RadialBasisFctSolver< RBF > > _rbfSolver
precice::logging::Logger _log
std::tuple< Args... > optionalArgs
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.
void clear() final override
RadialBasisFctMapping(Mapping::Constraint constraint, int dimensions, RADIAL_BASIS_FUNCTION_T function, std::array< bool, 3 > deadAxis, Polynomial polynomial, Args... args)
Constructor.
std::string getName() const final override
void mapConservative(const time::Sample &inData, Eigen::VectorXd &outData) final override
Maps data using a conservative constraint.
Container and creator for meshes.
void addMesh(Mesh &deltaMesh)
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.
A C++ 11 implementation of the non-owning C++20 std::span type.
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.
void sendMesh(Communication &communication, int rankReceiver, const mesh::Mesh &mesh)
constexpr auto asVector
Allows to use Communication::AsVectorTag in a less verbose way.
void receiveMesh(Communication &communication, int rankSender, mesh::Mesh &mesh)
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)
std::shared_ptr< Mesh > PtrMesh
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Main namespace of the precice library.
int dataDims
The dimensionality of the data.