3#include <Eigen/Cholesky>
29template <
typename SOLVER_T,
typename... Args>
53 void clear() final override;
59 precice::logging::Logger
_log{
"mapping::RadialBasisFctMapping"};
79template <
typename SOLVER_T,
typename... Args>
88 _polynomial(polynomial),
89 optionalArgs(
std::make_tuple(
std::forward<Args>(args)...))
91 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.");
94template <
typename SOLVER_T,
typename... Args>
101 PRECICE_ASSERT(this->input()->getDimensions() == this->output()->getDimensions(),
102 this->input()->getDimensions(), this->output()->getDimensions());
103 PRECICE_ASSERT(this->getDimensions() == this->output()->getDimensions(),
104 this->getDimensions(), this->output()->getDimensions());
110 inMesh = this->output();
111 outMesh = this->input();
113 inMesh = this->input();
114 outMesh = this->output();
137 globalInMesh.
addMesh(filteredInMesh);
138 globalOutMesh.
addMesh(*outMesh);
145 globalInMesh.
addMesh(secondaryInMesh);
149 globalOutMesh.
addMesh(secondaryOutMesh);
154 globalOutMesh.
addMesh(*outMesh);
159 _rbfSolver = std::make_unique<SOLVER_T>(this->_basisFunction, globalInMesh, boost::irange<Eigen::Index>(0, globalInMesh.
nVertices()),
160 globalOutMesh, boost::irange<Eigen::Index>(0, globalOutMesh.
nVertices()), this->_deadAxis, _polynomial, std::get<0>(optionalArgs));
162 _rbfSolver = std::make_unique<SOLVER_T>(this->_basisFunction, globalInMesh, boost::irange<Eigen::Index>(0, globalInMesh.
nVertices()),
163 globalOutMesh, boost::irange<Eigen::Index>(0, globalOutMesh.
nVertices()), this->_deadAxis, _polynomial);
166 this->_hasComputedMapping =
true;
170template <
typename SOLVER_T,
typename... Args>
175 this->_hasComputedMapping =
false;
178template <
typename SOLVER_T,
typename... Args>
182 auto param = std::get<0>(optionalArgs);
184 if (param.solver ==
"qr-solver") {
185 return "global-direct RBF (" + exec +
")";
187 return "global-iterative RBF (" + exec +
")";
190 return "global-direct RBF (cpu-executor)";
194template <
typename SOLVER_T,
typename... Args>
206 const auto &localInData = inData.
values;
208 int localOutputSize = 0;
209 for (
const auto &vertex : this->output()->vertices()) {
210 if (vertex.isOwner()) {
225 const auto &localInData = inData.
values;
226 globalInValues.
insert(globalInValues.
begin(), localInData.data(), localInData.data() + localInData.size());
228 int localOutputSize = 0;
229 for (
const auto &vertex : this->output()->vertices()) {
230 if (vertex.isOwner()) {
237 outputValueSizes.
push_back(localOutputSize);
241 int secondaryOutputValueSize;
244 globalInValues.
insert(globalInValues.
end(), secondaryBuffer.
begin(), secondaryBuffer.
end());
247 outputValueSizes.
push_back(secondaryOutputValueSize);
251 const int valueDim = inData.
dataDims;
254 Eigen::Map<Eigen::VectorXd> inputValues(globalInValues.
data(), globalInValues.
size());
255 Eigen::VectorXd outputValues((this->output()->getGlobalNumberOfVertices()) * valueDim);
258 in.resize(_rbfSolver->getOutputSize());
260 outputValues.setZero();
262 for (
int dim = 0; dim < valueDim; dim++) {
263 for (
int i = 0; i < in.size(); i++) {
264 in[i] = inputValues(i * valueDim + dim);
268 out = _rbfSolver->solveConservative(in, _polynomial);
271 for (
int i = 0; i < this->output()->getGlobalNumberOfVertices(); i++) {
272 outputValues[i * valueDim + dim] =
out[i];
280 int outputCounter = 0;
281 for (
int i = 0; i < static_cast<int>(this->output()->nVertices()); ++i) {
282 if (this->output()->vertex(i).isOwner()) {
283 for (
int dim = 0; dim < valueDim; ++dim) {
284 outData[i * valueDim + dim] = outputValues(outputCounter);
291 int beginPoint = outputValueSizes.
at(0);
295 beginPoint += outputValueSizes.
at(rank);
298 outData = outputValues;
304 const int valueDim = inData.
dataDims;
306 int outputCounter = 0;
307 for (
int i = 0; i < static_cast<int>(this->output()->nVertices()); ++i) {
308 if (this->output()->vertex(i).isOwner()) {
309 for (
int dim = 0; dim < valueDim; ++dim) {
310 outData[i * valueDim + dim] = receivedValues.
at(outputCounter);
318template <
typename SOLVER_T,
typename... Args>
330 auto localInDataFiltered = this->input()->getOwnedVertexData(inData.
values);
331 int localOutputSize = outData.size();
339 const int valueDim = inData.
dataDims;
347 const auto &localInData = this->input()->getOwnedVertexData(inData.
values);
348 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.
begin());
351 int inputSizeCounter = localInData.size();
352 int secondaryOutDataSize{0};
357 inputSizeCounter += secondaryBuffer.
size();
360 outValuesSize.
push_back(secondaryOutDataSize);
364 const auto &localInData = inData.
values;
365 std::copy(localInData.data(), localInData.data() + localInData.size(), globalInValues.
begin());
371 in.resize(_rbfSolver->getInputSize());
375 Eigen::Map<Eigen::VectorXd> inputValues(globalInValues.
data(), globalInValues.
size());
377 Eigen::VectorXd outputValues;
378 outputValues.resize((_rbfSolver->getOutputSize()) * valueDim);
381 outputValues.setZero();
384 for (
int dim = 0; dim < valueDim; dim++) {
386 for (
int i = 0; i < this->input()->getGlobalNumberOfVertices(); i++) {
387 in[i] = inputValues[i * valueDim + dim];
390 out = _rbfSolver->solveConsistent(in, _polynomial);
393 for (
int i = 0; i <
out.size(); i++) {
394 outputValues[i * valueDim + dim] =
out[i];
398 outData = Eigen::Map<Eigen::VectorXd>(outputValues.data(), outValuesSize.
at(0));
401 int beginPoint = outValuesSize.
at(0);
407 beginPoint += outValuesSize.
at(rank);
413 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.
Constraint
Specifies additional constraints for a mapping.
InitialGuessRequirement
Specifies whether the mapping requires an initial guess.
Mapping with radial basis functions.
Mapping with radial basis functions.
typename SOLVER_T::BASIS_FUNCTION_T RADIAL_BASIS_FUNCTION_T
Polynomial _polynomial
Treatment of the polynomial.
std::unique_ptr< SOLVER_T > _rbfSolver
precice::logging::Logger _log
std::tuple< Args... > optionalArgs
Optional constructor arguments for the solver class.
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
Removes a computed mapping.
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
name of the rbf mapping
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.
constexpr pointer data() const noexcept
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)
void receiveMesh(Communication &communication, int rankSender, mesh::Mesh &mesh)
Polynomial
How to handle the polynomial?
void filterMesh(Mesh &destination, const Mesh &source, UnaryPredicate p)
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Main namespace of the precice library.
int dataDims
The dimensionality of the data.