3#include <Eigen/Cholesky>
28template <
typename SOLVER_T,
typename... Args>
78template <
typename SOLVER_T,
typename... Args>
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);
162 globalOutMesh, boost::irange<Eigen::Index>(0, globalOutMesh.
nVertices()), this->_deadAxis,
_polynomial);
169template <
typename SOLVER_T,
typename... Args>
177template <
typename SOLVER_T,
typename... Args>
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>
204 const auto &localInData = inData.
values;
206 int localOutputSize = 0;
207 for (
const auto &vertex : this->
output()->vertices()) {
208 if (vertex.isOwner()) {
223 const auto &localInData = inData.
values;
224 globalInValues.
insert(globalInValues.
begin(), localInData.data(), localInData.data() + localInData.size());
226 int localOutputSize = 0;
227 for (
const auto &vertex : this->
output()->vertices()) {
228 if (vertex.isOwner()) {
235 outputValueSizes.
push_back(localOutputSize);
239 int secondaryOutputValueSize;
242 globalInValues.
insert(globalInValues.
end(), secondaryBuffer.
begin(), secondaryBuffer.
end());
245 outputValueSizes.
push_back(secondaryOutputValueSize);
249 const int valueDim = inData.
dataDims;
252 Eigen::Map<Eigen::VectorXd> inputValues(globalInValues.
data(), globalInValues.
size());
253 Eigen::VectorXd outputValues((this->
output()->getGlobalNumberOfVertices()) * valueDim);
258 outputValues.setZero();
260 for (
int dim = 0; dim < valueDim; dim++) {
261 for (
int i = 0; i < in.size(); i++) {
262 in[i] = inputValues(i * valueDim + dim);
269 for (
int i = 0; i < this->
output()->getGlobalNumberOfVertices(); i++) {
270 outputValues[i * valueDim + dim] = out[i];
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);
289 int beginPoint = outputValueSizes.
at(0);
293 beginPoint += outputValueSizes.
at(rank);
296 outData = outputValues;
302 const int valueDim = inData.
dataDims;
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);
316template <
typename SOLVER_T,
typename... Args>
327 auto localInDataFiltered = this->
input()->getOwnedVertexData(inData.
values);
328 int localOutputSize = outData.size();
336 const int valueDim = inData.
dataDims;
344 const auto &localInData = this->
input()->getOwnedVertexData(inData.
values);
345 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.
begin());
348 int inputSizeCounter = localInData.size();
349 int secondaryOutDataSize{0};
354 inputSizeCounter += secondaryBuffer.
size();
357 outValuesSize.
push_back(secondaryOutDataSize);
361 const auto &localInData = inData.
values;
362 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.
begin());
372 Eigen::Map<Eigen::VectorXd> inputValues(globalInValues.
data(), globalInValues.
size());
374 Eigen::VectorXd outputValues;
375 outputValues.resize((
_rbfSolver->getOutputSize()) * valueDim);
378 outputValues.setZero();
381 for (
int dim = 0; dim < valueDim; dim++) {
383 for (
int i = 0; i < this->
input()->getGlobalNumberOfVertices(); i++) {
384 in[i] = inputValues[i * valueDim + dim];
390 for (
int i = 0; i < out.size(); i++) {
391 outputValues[i * valueDim + dim] = out[i];
395 outData = Eigen::Map<Eigen::VectorXd>(outputValues.data(), outValuesSize.
at(0));
398 int beginPoint = outValuesSize.
at(0);
404 beginPoint += outValuesSize.
at(rank);
410 outData = Eigen::Map<Eigen::VectorXd>(receivedValues.
data(), receivedValues.
size());
#define PRECICE_DEBUG(...)
#define PRECICE_TRACE(...)
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
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.
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.