40void addSubtagsToParents(std::list<xml::XMLTag> &subtags,
41 std::list<xml::XMLTag> &parents)
43 for (
auto &p : parents) {
44 p.addSubtags(subtags);
50using variant_t = std::variant<xml::XMLAttribute<double>, xml::XMLAttribute<std::string>, xml::XMLAttribute<bool>, xml::XMLAttribute<int>>;
51template <
typename TagStorage>
52void addAttributes(TagStorage &storage,
const std::vector<variant_t> &attributes)
54 for (
auto &s : storage) {
55 for (
auto &a : attributes)
56 std::visit([&s](
auto &&arg) { s.addAttribute(arg); }, a);
61enum struct RBFBackend {
70template <RBFBackend Backend,
typename RBF>
71struct BackendSelector {
76template <
typename RBF>
77struct BackendSelector<RBFBackend::Eigen, RBF> {
78 typedef mapping::RadialBasisFctMapping<RadialBasisFctSolver<RBF>> type;
82#ifndef PRECICE_NO_PETSC
83template <
typename RBF>
84struct BackendSelector<RBFBackend::PETSc, RBF> {
85 typedef mapping::PetRadialBasisFctMapping<RBF> type;
90#ifndef PRECICE_NO_GINKGO
91template <
typename RBF>
92struct BackendSelector<RBFBackend::Ginkgo, RBF> {
93 typedef mapping::RadialBasisFctMapping<GinkgoRadialBasisFctSolver<RBF>, MappingConfiguration::GinkgoParameter> type;
97template <
typename RBF>
98struct BackendSelector<RBFBackend::PUM, RBF> {
99 typedef mapping::PartitionOfUnityMapping<RBF> type;
103using rbf_variant_t = std::variant<CompactPolynomialC0, CompactPolynomialC2, CompactPolynomialC4, CompactPolynomialC6, CompactPolynomialC8, CompactThinPlateSplinesC2, ThinPlateSplines, VolumeSplines, Multiquadrics, InverseMultiquadrics, Gaussian>;
106template <RBFBackend T,
typename RADIAL_BASIS_FUNCTION_T,
typename... Args>
110 return PtrMapping(
new typename BackendSelector<T, RADIAL_BASIS_FUNCTION_T>::type(constraint, dimension, function,
std::forward<Args>(args)...));
114rbf_variant_t constructRBF(
BasisFunction functionType,
double supportRadius,
double shapeParameter)
117 switch (functionType) {
119 return mapping::CompactPolynomialC0(supportRadius);
122 return mapping::CompactPolynomialC2(supportRadius);
125 return mapping::CompactPolynomialC4(supportRadius);
128 return mapping::CompactPolynomialC6(supportRadius);
131 return mapping::CompactPolynomialC8(supportRadius);
154 }
catch (std::invalid_argument &e) {
155 logging::Logger
_log{
"MappingConfiguration"};
163template <RBFBackend T,
typename... Args>
168 auto functionVariant = constructRBF(functionType, supportRadius, shapeParameter);
170 return std::visit([&](
auto &&func) {
return instantiateRBFMapping<T>(constraint, dimension, func,
std::forward<Args>(args)...); }, functionVariant);
183 XMLTag::Occurrence occ = XMLTag::OCCUR_ARBITRARY;
185 XMLTag{*
this,
TYPE_NEAREST_NEIGHBOR, occ,
TAG}.setDocumentation(
"Nearest-neighbour mapping which uses a rstar-spacial index tree to index meshes and run nearest-neighbour queries."),
186 XMLTag{*
this,
TYPE_NEAREST_PROJECTION, occ,
TAG}.setDocumentation(
"Nearest-projection mapping which uses a rstar-spacial index tree to index meshes and locate the nearest projections."),
187 XMLTag{*
this,
TYPE_NEAREST_NEIGHBOR_GRADIENT, occ,
TAG}.setDocumentation(
"Nearest-neighbor-gradient mapping which uses nearest-neighbor mapping with an additional linear approximation using gradient data."),
188 XMLTag{*
this,
TYPE_LINEAR_CELL_INTERPOLATION, occ,
TAG}.setDocumentation(
"Linear cell interpolation mapping which uses a rstar-spacial index tree to index meshes and locate the nearest cell. Only supports 2D meshes.")};
190 XMLTag{*
this,
TYPE_RBF_GLOBAL_DIRECT, occ,
TAG}.setDocumentation(
"Radial-basis-function mapping using a direct solver with a gather-scatter parallelism.")};
192 XMLTag{*
this,
TYPE_RBF_GLOBAL_ITERATIVE, occ,
TAG}.setDocumentation(
"Radial-basis-function mapping using an iterative solver with a distributed parallelism.")};
194 XMLTag{*
this,
TYPE_RBF_PUM_DIRECT, occ,
TAG}.setDocumentation(
"Radial-basis-function mapping using a partition of unity method, which supports a distributed parallelism.")};
196 XMLTag{*
this,
TYPE_RBF_ALIAS, occ,
TAG}.setDocumentation(
"Alias tag, which auto-selects a radial-basis-function mapping depending on the simulation parameter,")};
199 XMLTag{*
this,
TYPE_RADIAL_GEOMETRIC_MULTISCALE, occ,
TAG}.setDocumentation(
"Radial geometric multiscale mapping between multiple 1D and multiple 3D vertices, distributed along a principle axis.")};
204 .setDocumentation(
"Write mappings map written data prior to communication, thus in the same participant who writes the data. "
205 "Read mappings map received data after communication, thus in the same participant who reads the data.");
207 auto attrFromMesh = XMLAttribute<std::string>(
ATTR_FROM,
"")
208 .setDocumentation(
"The mesh to map the data from. The default name is an empty mesh name, which is only valid for a just-in-time mapping (using the API functions \"writeAndMapData\" or \"mapAndReadData\").");
210 auto attrToMesh = XMLAttribute<std::string>(
ATTR_TO,
"")
211 .setDocumentation(
"The mesh to map the data to. The default name is an empty mesh name, which is only valid for a just-in-time mapping (using the API functions \"writeAndMapData\" or \"mapAndReadData\").");
214 .setDocumentation(
"Use conservative to conserve the nodal sum of the data over the interface (needed e.g. for force mapping). Use consistent for normalized quantities such as temperature or pressure. Use scaled-consistent-surface or scaled-consistent-volume for normalized quantities where conservation of integral values (surface or volume) is needed (e.g. velocities when the mass flow rate needs to be conserved). Mesh connectivity is required to use scaled-consistent.")
216 auto attrXDead = makeXMLAttribute(
ATTR_X_DEAD,
false)
217 .setDocumentation(
"If set to true, the x axis will be ignored for the mapping");
218 auto attrYDead = makeXMLAttribute(
ATTR_Y_DEAD,
false)
219 .setDocumentation(
"If set to true, the y axis will be ignored for the mapping");
220 auto attrZDead = makeXMLAttribute(
ATTR_Z_DEAD,
false)
221 .setDocumentation(
"If set to true, the z axis will be ignored for the mapping");
223 .setDocumentation(
"Toggles use of the global polynomial")
226 .setDocumentation(
"Toggles use a local (per cluster) polynomial")
230 .setDocumentation(
"Solver relative tolerance for convergence");
236 .setDocumentation(
"Average number of vertices per cluster (partition) applied in the rbf partition of unity method.");
238 .setDocumentation(
"Value between 0 and 1 indicating the relative overlap between clusters. A value of 0.15 is usually a good trade-off between accuracy and efficiency.");
240 .setDocumentation(
"If enabled, places the cluster centers at the closest vertex of the input mesh. Should be enabled in case of non-uniform point distributions such as for shell structures.");
243 .setDocumentation(
"Type of geometric multiscale mapping. Either 'spread' or 'collect'.")
246 .setDocumentation(
"Principle axis along which geometric multiscale mapping is performed.")
249 .setDocumentation(
"Radius of the circular interface between the 1D and 3D participant.");
252 addAttributes(projectionTags, {attrFromMesh, attrToMesh, attrDirection, attrConstraint});
253 addAttributes(rbfDirectTags, {attrFromMesh, attrToMesh, attrDirection, attrConstraint, attrPolynomial, attrXDead, attrYDead, attrZDead});
254 addAttributes(rbfIterativeTags, {attrFromMesh, attrToMesh, attrDirection, attrConstraint, attrPolynomial, attrXDead, attrYDead, attrZDead, attrSolverRtol});
255 addAttributes(pumDirectTags, {attrFromMesh, attrToMesh, attrDirection, attrConstraint, attrPumPolynomial, verticesPerCluster, relativeOverlap, projectToInput});
256 addAttributes(rbfAliasTag, {attrFromMesh, attrToMesh, attrDirection, attrConstraint, attrXDead, attrYDead, attrZDead});
257 addAttributes(geoMultiscaleTags, {attrFromMesh, attrToMesh, attrDirection, attrConstraint, attrGeoMultiscaleType, attrGeoMultiscaleAxis, attrGeoMultiscaleRadius});
260 XMLTag::Occurrence once = XMLTag::OCCUR_NOT_OR_ONCE;
262 auto attrDeviceId = makeXMLAttribute(
ATTR_DEVICE_ID,
static_cast<int>(0))
263 .setDocumentation(
"Specifies the ID of the GPU that should be used for the Ginkgo GPU backend.");
264 auto attrNThreads = makeXMLAttribute(
ATTR_N_THREADS,
static_cast<int>(0))
265 .setDocumentation(
"Specifies the number of threads for the OpenMP executor that should be used for the Ginkgo OpenMP backend. If a value of \"0\" is set, preCICE doesn't set the number of threads and the default behavior of OpenMP applies.");
270 XMLTag{*
this,
EXECUTOR_CPU, once,
SUBTAG_EXECUTOR}.setDocumentation(
"The default executor, which uses a single-core CPU with a gather-scatter parallelism.")};
272 XMLTag{*
this,
EXECUTOR_CUDA, once,
SUBTAG_EXECUTOR}.setDocumentation(
"Cuda (Nvidia) executor, which uses cuSolver/Ginkgo and a direct QR decomposition with a gather-scatter parallelism."),
273 XMLTag{*
this,
EXECUTOR_HIP, once,
SUBTAG_EXECUTOR}.setDocumentation(
"Hip (AMD/Nvidia) executor, which uses hipSolver/Ginkgo and a direct QR decomposition with a gather-scatter parallelism.")};
275 addAttributes(deviceExecutors, {attrDeviceId});
276 addSubtagsToParents(cpuExecutor, rbfDirectTags);
277 addSubtagsToParents(deviceExecutors, rbfDirectTags);
282 XMLTag{*
this,
EXECUTOR_CPU, once,
SUBTAG_EXECUTOR}.setDocumentation(
"The default executor relying on PETSc, which uses CPUs and distributed memory parallelism via MPI.")};
284 XMLTag{*
this,
EXECUTOR_CUDA, once,
SUBTAG_EXECUTOR}.setDocumentation(
"Cuda (Nvidia) executor, which uses Ginkgo with a gather-scatter parallelism."),
285 XMLTag{*
this,
EXECUTOR_HIP, once,
SUBTAG_EXECUTOR}.setDocumentation(
"Hip (AMD/Nvidia) executor, which uses hipSolver with a gather-scatter parallelism.")};
287 XMLTag{*
this,
EXECUTOR_OMP, once,
SUBTAG_EXECUTOR}.setDocumentation(
"OpenMP executor, which uses Ginkgo with a gather-scatter parallelism.")};
289 addAttributes(deviceExecutors, {attrDeviceId});
290 addAttributes(ompExecutor, {attrNThreads});
291 addSubtagsToParents(cpuExecutor, rbfIterativeTags);
292 addSubtagsToParents(deviceExecutors, rbfIterativeTags);
293 addSubtagsToParents(ompExecutor, rbfIterativeTags);
297 XMLTag{*
this,
EXECUTOR_CPU, once,
SUBTAG_EXECUTOR}.setDocumentation(
"The default (and currently only) executor using a CPU and a distributed memory parallelism via MPI.")};
298 addSubtagsToParents(cpuExecutor, pumDirectTags);
313 .setDocumentation(
"Support radius of each RBF basis function (global choice).");
315 addAttributes(supportRadiusRBF, {attrSupportRadius});
316 addSubtagsToParents(supportRadiusRBF, rbfIterativeTags);
317 addSubtagsToParents(supportRadiusRBF, rbfDirectTags);
318 addSubtagsToParents(supportRadiusRBF, pumDirectTags);
319 addSubtagsToParents(supportRadiusRBF, rbfAliasTag);
327 .setDocumentation(
"Specific shape parameter for RBF basis function.");
329 addAttributes(shapeParameterRBF, {attrShapeParam});
330 addSubtagsToParents(shapeParameterRBF, rbfIterativeTags);
331 addSubtagsToParents(shapeParameterRBF, rbfDirectTags);
332 addSubtagsToParents(shapeParameterRBF, pumDirectTags);
333 addSubtagsToParents(shapeParameterRBF, rbfAliasTag);
340 addAttributes(GaussRBF, {attrShapeParam, attrSupportRadius});
341 addSubtagsToParents(GaussRBF, rbfIterativeTags);
342 addSubtagsToParents(GaussRBF, rbfDirectTags);
343 addSubtagsToParents(GaussRBF, pumDirectTags);
344 addSubtagsToParents(GaussRBF, rbfAliasTag);
351 addSubtagsToParents(attributelessRBFs, rbfIterativeTags);
352 addSubtagsToParents(attributelessRBFs, rbfDirectTags);
353 addSubtagsToParents(attributelessRBFs, pumDirectTags);
354 addSubtagsToParents(attributelessRBFs, rbfAliasTag);
385 PRECICE_CHECK(!toMesh.
empty() || !fromMesh.
empty(),
"Neither a \"to\" nor a \"from\" mesh was defined in a preCICE mapping. Most data mappings require a \"to\" and a \"from\" mesh. "
386 "For just-in-time mapping, at least one of both attributes has to be specified.");
392 "The mapping from mesh \"{0}\" has no \"to\" mesh, which configures a just-in-time mapping for the mesh \"{0}\". For just-in-time mapping, only read-consistent (direction = \"read\" and no \"to\" mesh) and write-conservative (direction = \"write\" and no \"from\" mesh) are implemented.", fromMesh);
394 "The mapping to mesh \"{0}\" has no \"from\" mesh, which configures a just-in-time mapping for the mesh \"{0}\". For just-in-time mapping, only read-consistent (direction = \"read\" and no \"to\" mesh) and write-conservative (direction = \"write\" and no \"from\" mesh) are implemented.", toMesh);
414 PRECICE_CHECK(
_experimental,
"Axial geometric multiscale is experimental and the configuration can change between minor releases. Set experimental=\"on\" in the precice-configuration tag.");
418 throw std::runtime_error{
"Axial geometric multiscale mapping is not available for parallel participants."};
422 throw std::runtime_error{
"Radial geometric multiscale mapping is not available for parallel participants."};
452 PRECICE_CHECK(
_mappings.back().requiresBasisFunction ==
true,
"A basis-function was defined for the mapping "
453 "from mesh \"{}\" to mesh \"{}\", but no basis-function is required for this mapping type. "
454 "Please remove the basis-function tag or configure an rbf mapping, which requires a basis-function.",
457 PRECICE_CHECK(
_rbfConfig.basisFunctionDefined ==
false,
"More than one basis-function was defined for the mapping "
458 "from mesh \"{}\" to mesh \"{}\".",
471 PRECICE_CHECK(exactlyOneSet,
"The specified parameters for the Gaussian RBF mapping are invalid. Please specify either a \"shape-parameter\" or a \"support-radius\".");
500 bool xDead,
bool yDead,
bool zDead,
502 double verticesPerCluster,
503 double relativeOverlap,
504 bool projectToInput)
const
532 rbfConfig.deadAxis = {{xDead, yDead, zDead}};
535 rbfConfig.verticesPerCluster = verticesPerCluster;
536 rbfConfig.relativeOverlap = relativeOverlap;
537 rbfConfig.projectToInput = projectToInput;
549 const double &multiscaleRadius)
const
558 if (!fromMesh && fromMeshName.
empty()) {
560 "Mesh \"{0}\" was not found while creating a mapping. "
561 "Please correct the to=\"{0}\" attribute.",
565 if (!toMesh && toMeshName.
empty()) {
567 "Mesh \"{0}\" was not found while creating a mapping. "
568 "Please correct the from=\"{0}\" attribute.",
574 "Mesh \"{0}\" was not found while creating a mapping. "
575 "Please correct the from=\"{0}\" attribute.",
578 "Mesh \"{0}\" was not found while creating a mapping. "
579 "Please correct the to=\"{0}\" attribute.",
582 PRECICE_CHECK((!toMesh->isJustInTime() && !fromMesh->isJustInTime()) ||
584 "A just-in-time mapping was configured from mesh \"{}\" to mesh \"{}\" using \"mapping:{}\", which is currently not implemented. "
585 "Available mapping types are \"mapping:{}\", \"mapping:{}\" and \"mapping:{}\".",
589 PRECICE_CHECK(fromMesh->getDimensions() == toMesh->getDimensions(),
590 "Mapping between meshes of different dimensions is not allowed yet. "
591 "Please set the same dimensions attribute to meshes \"{}\" and \"{}\", "
592 "or choose different meshes.",
593 fromMeshName, toMeshName);
595 configuredMapping.
fromMesh = fromMesh;
596 configuredMapping.
toMesh = toMesh;
617 "Nearest-neighbor-gradient mapping is not implemented using a \"conservative\" constraint. "
618 "Please select constraint=\" consistent\" or a different mapping method.");
626 "Axial geometric multiscale mapping is not implemented for the \"conservative\" constraint. "
627 "Please select constraint=\" consistent\" or a different mapping method.");
631 if (geoMultiscaleAxis ==
"x") {
633 }
else if (geoMultiscaleAxis ==
"y") {
635 }
else if (geoMultiscaleAxis ==
"z") {
642 if (geoMultiscaleType ==
"spread") {
644 }
else if (geoMultiscaleType ==
"collect") {
656 "Radial geometric multiscale mapping is not implemented for the \"conservative\" constraint. "
657 "Please select constraint=\" consistent\" or a different mapping method.");
661 if (geoMultiscaleAxis ==
"x") {
663 }
else if (geoMultiscaleAxis ==
"y") {
665 }
else if (geoMultiscaleAxis ==
"z") {
672 if (geoMultiscaleType ==
"spread") {
674 }
else if (geoMultiscaleType ==
"collect") {
685 configuredMapping.
mapping =
nullptr;
691 return configuredMapping;
697 bool sameToMesh =
mapping.toMesh->getName() == configuredMapping.toMesh->getName();
698 bool sameFromMesh =
mapping.fromMesh->getName() == configuredMapping.fromMesh->getName();
699 bool sameMapping = sameToMesh && sameFromMesh;
701 "There cannot be two mappings from mesh \"{}\" to mesh \"{}\". "
702 "Please remove one of the duplicated meshes. ",
706 bool sameToMesh =
mapping.toMesh->getName() == configuredMapping.toMesh->getName();
708 bool sameMapping = sameToMesh && isWrite && (
mapping.fromMesh->isJustInTime() || configuredMapping.fromMesh->isJustInTime());
710 "There cannot be two mappings to mesh \"{}\". "
711 "Here, we have a mixture of just-in-time mapping and a conventional mapping. ",
743#ifndef PRECICE_NO_PETSC
749 PRECICE_CHECK(
false,
"The global-iterative RBF solver on a CPU requires a preCICE build with PETSc enabled.");
758#ifndef PRECICE_NO_GINKGO
764#ifndef PRECICE_WITH_CUDA
765 PRECICE_CHECK(
false,
"The cuda-executor (configured for the mapping from mesh {} to mesh {}) requires a Ginkgo and preCICE build with Cuda enabled.",
mapping.fromMesh->getName(),
mapping.toMesh->getName());
769#ifndef PRECICE_WITH_HIP
770 PRECICE_CHECK(
false,
"The hip-executor (configured for the mapping from mesh {} to mesh {}) requires a Ginkgo and preCICE build with HIP enabled.",
mapping.fromMesh->getName(),
mapping.toMesh->getName());
775#ifndef PRECICE_WITH_OPENMP
776 PRECICE_CHECK(
false,
"The omp-executor (configured for the mapping from mesh {} to mesh {}) requires a Ginkgo and preCICE build with OpenMP enabled.",
mapping.fromMesh->getName(),
mapping.toMesh->getName());
790 PRECICE_CHECK(
false,
"The selected executor for the mapping from mesh {} to mesh {} requires a preCICE build with Ginkgo enabled.",
mapping.fromMesh->getName(),
mapping.toMesh->getName());
832 return basisFunction;
#define PRECICE_ERROR(...)
#define PRECICE_TRACE(...)
#define PRECICE_INFO_IF(condition,...)
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
#define PRECICE_UNREACHABLE(...)
Geometric multiscale mapping in axial direction.
MultiscaleType
Geometric multiscale nature of the mapping (spread or collect).
static constexpr double cutoffThreshold
Below that value the function is supposed to be zero. Defines the support radius if not explicitly gi...
Mapping using orthogonal projection to nearest triangle/edge/vertex and linear interpolation from pro...
const std::string ATTR_TO
RBFConfiguration configureRBFMapping(const std::string &type, const std::string &polynomial, bool xDead, bool yDead, bool zDead, double solverRtol, double verticesPerCluster, double relativeOverlap, bool projectToInput) const
const std::string GEOMETRIC_MULTISCALE_AXIS_X
const std::string SUBTAG_BASIS_FUNCTION
const std::string ATTR_VERTICES_PER_CLUSTER
const std::string TYPE_RBF_PUM_DIRECT
void xmlTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
Callback function required for use of automatic configuration.
const std::string RBF_CPOLYNOMIAL_C2
const std::string ATTR_X_DEAD
const std::string SUBTAG_EXECUTOR
const std::vector< ConfiguredMapping > & mappings()
Returns all configured mappings.
GinkgoParameter _ginkgoParameter
const std::string ATTR_DIRECTION
const std::string RBF_VOLUME_SPLINES
Mapping::Constraint constraintValue
const std::string DIRECTION_WRITE
const std::string ATTR_SOLVER_RTOL
mesh::PtrMeshConfiguration _meshConfig
const std::string EXECUTOR_HIP
const std::string DIRECTION_READ
const std::string POLYNOMIAL_OFF
std::unique_ptr< ExecutorConfiguration > _executorConfig
const std::string ATTR_Y_DEAD
const std::string RBF_MULTIQUADRICS
const std::string RBF_CPOLYNOMIAL_C4
const std::string ATTR_SHAPE_PARAM
const std::string CONSTRAINT_SCALED_CONSISTENT_VOLUME
const std::string EXECUTOR_CPU
const std::string TYPE_RBF_ALIAS
const std::string TYPE_RBF_GLOBAL_DIRECT
const std::string ATTR_PROJECT_TO_INPUT
const std::string GEOMETRIC_MULTISCALE_AXIS_Z
bool requiresBasisFunction(const std::string &mappingType) const
const std::string RBF_CPOLYNOMIAL_C6
const std::string TYPE_AXIAL_GEOMETRIC_MULTISCALE
BasisFunction parseBasisFunctions(const std::string &basisFctName) const
Given a basis function name (as a string), transforms the string into an enum of the BasisFunction.
const std::string RBF_TPS
const std::string ATTR_DEVICE_ID
void finishRBFConfiguration()
const std::string RBF_CPOLYNOMIAL_C0
const std::string CONSTRAINT_CONSISTENT
RBFConfiguration _rbfConfig
const std::string TYPE_RBF_GLOBAL_ITERATIVE
const std::string ATTR_RELATIVE_OVERLAP
const std::string ATTR_GEOMETRIC_MULTISCALE_AXIS
const std::string TYPE_LINEAR_CELL_INTERPOLATION
const std::string GEOMETRIC_MULTISCALE_AXIS_Y
const std::string RBF_GAUSSIAN
std::vector< ConfiguredMapping > _mappings
const std::string ATTR_CONSTRAINT
const std::string TYPE_NEAREST_NEIGHBOR
const std::string ATTR_Z_DEAD
const std::string EXECUTOR_CUDA
const std::string ATTR_GEOMETRIC_MULTISCALE_TYPE
const std::string ATTR_SUPPORT_RADIUS
ConfiguredMapping createMapping(const std::string &direction, const std::string &type, const std::string &fromMeshName, const std::string &toMeshName, const std::string &geoMultiscaleType, const std::string &geoMultiscaleAxis, const double &multiscaleRadius) const
void checkDuplicates(const ConfiguredMapping &mapping)
Check whether a mapping to and from the same mesh already exists.
void xmlEndTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
Callback function required for use of automatic configuration.
const std::string TYPE_NEAREST_PROJECTION
const std::string RBF_CPOLYNOMIAL_C8
const std::string POLYNOMIAL_ON
const std::string ATTR_N_THREADS
const std::string ATTR_POLYNOMIAL
const std::string POLYNOMIAL_SEPARATE
const std::string TYPE_NEAREST_NEIGHBOR_GRADIENT
const std::string RBF_INV_MULTIQUADRICS
const RBFConfiguration & rbfConfig() const
const std::string TYPE_RADIAL_GEOMETRIC_MULTISCALE
const std::string RBF_CTPS_C2
MappingConfiguration(xml::XMLTag &parent, mesh::PtrMeshConfiguration meshConfiguration)
const std::string GEOMETRIC_MULTISCALE_TYPE_SPREAD
const std::string CONSTRAINT_CONSERVATIVE
const std::string ATTR_FROM
const std::string EXECUTOR_OMP
void setExperimental(bool experimental)
const std::string ATTR_GEOMETRIC_MULTISCALE_RADIUS
const std::string CONSTRAINT_SCALED_CONSISTENT_SURFACE
const std::string GEOMETRIC_MULTISCALE_TYPE_COLLECT
Constraint
Specifies additional constraints for a mapping.
@ SCALED_CONSISTENT_VOLUME
@ SCALED_CONSISTENT_SURFACE
Mapping using nearest neighboring vertices and their local gradient values.
Mapping using nearest neighboring vertices.
Mapping using orthogonal projection to nearest triangle/edge/vertex and linear interpolation from pro...
Geometric multiscale mapping in radial direction.
MultiscaleType
Geometric multiscale type of the mapping (spread or collect).
static mesh::PtrMesh getJustInTimeMappingMesh(int dimension)
static CommStatePtr current()
Returns an owning pointer to the current CommState.
static void initialize(utils::Parallel::Communicator comm)
Initializes the Petsc environment.
Represents an XML tag to be configured automatically.
const std::string & getNamespace() const
Returns xml namespace.
void addSubtags(const Container &subtags)
std::string getStringAttributeValue(const std::string &name, std::optional< std::string > default_value=std::nullopt) const
bool getBooleanAttributeValue(const std::string &name, std::optional< bool > default_value=std::nullopt) const
const std::string & getName() const
Returns name (without namespace).
int getIntAttributeValue(const std::string &name, std::optional< int > default_value=std::nullopt) const
double getDoubleAttributeValue(const std::string &name, std::optional< double > default_value=std::nullopt) const
contains data mapping from points to meshes.
@ CompactThinPlateSplinesC2
std::shared_ptr< Mapping > PtrMapping
std::shared_ptr< Mesh > PtrMesh
std::shared_ptr< MeshConfiguration > PtrMeshConfiguration
contains the XML configuration parser.
static precice::logging::Logger _log("precicec")
Tightly coupled to the parameters of Participant()