preCICE v3.2.0
Loading...
Searching...
No Matches
MappingConfiguration.cpp
Go to the documentation of this file.
2#include <Eigen/Core>
3#include <algorithm>
4#include <cstring>
5#include <list>
6#include <memory>
7#include <ostream>
8#include <utility>
9#include <variant>
10#include "logging/LogMacros.hpp"
14#include "mapping/Mapping.hpp"
25#include "mesh/Mesh.hpp"
28#include "utils/Parallel.hpp"
29#include "utils/Petsc.hpp"
30#include "utils/assertion.hpp"
31#include "xml/ConfigParser.hpp"
32#include "xml/XMLAttribute.hpp"
33#include "xml/XMLTag.hpp"
34namespace precice::mapping {
35
36namespace {
37
38// given a list of subtags and parent tags, this function adds all subtags to all
39// parent tags
40void addSubtagsToParents(std::list<xml::XMLTag> &subtags,
41 std::list<xml::XMLTag> &parents)
42{
43 for (auto &p : parents) {
44 p.addSubtags(subtags);
45 }
46}
47
48// this function uses std::variant in order to add attributes of any type (double, string, bool)
49// to all tags in the list of tags \p storage.
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)
53{
54 for (auto &s : storage) {
55 for (auto &a : attributes)
56 std::visit([&s](auto &&arg) { s.addAttribute(arg); }, a);
57 }
58}
59
60// Enum required for the RBF instantiations
61enum struct RBFBackend {
62 Eigen,
63 PETSc,
64 Ginkgo,
65 PUM
66};
67
68// Helper in order to resolve the template instantiations.
69// Only the template specializations are of interest
70template <RBFBackend Backend, typename RBF>
71struct BackendSelector {
72 typedef RBF type;
73};
74
75// Specialization for the RBF Eigen backend
76template <typename RBF>
77struct BackendSelector<RBFBackend::Eigen, RBF> {
78 typedef mapping::RadialBasisFctMapping<RadialBasisFctSolver<RBF>> type;
79};
80
81// Specialization for the PETSc RBF backend
82#ifndef PRECICE_NO_PETSC
83template <typename RBF>
84struct BackendSelector<RBFBackend::PETSc, RBF> {
85 typedef mapping::PetRadialBasisFctMapping<RBF> type;
86};
87#endif
88
89// Specialization for the Ginkgo RBF backend
90#ifndef PRECICE_NO_GINKGO
91template <typename RBF>
92struct BackendSelector<RBFBackend::Ginkgo, RBF> {
93 typedef mapping::RadialBasisFctMapping<GinkgoRadialBasisFctSolver<RBF>, MappingConfiguration::GinkgoParameter> type;
94};
95#endif
96// Specialization for the RBF PUM backend
97template <typename RBF>
98struct BackendSelector<RBFBackend::PUM, RBF> {
99 typedef mapping::PartitionOfUnityMapping<RBF> type;
100};
101
102// Variant holding all available RBF classes
103using rbf_variant_t = std::variant<CompactPolynomialC0, CompactPolynomialC2, CompactPolynomialC4, CompactPolynomialC6, CompactPolynomialC8, CompactThinPlateSplinesC2, ThinPlateSplines, VolumeSplines, Multiquadrics, InverseMultiquadrics, Gaussian>;
104
105// The actual instantiation of the mapping class, which is called by the visitor \ref getRBFMapping
106template <RBFBackend T, typename RADIAL_BASIS_FUNCTION_T, typename... Args>
107PtrMapping instantiateRBFMapping(mapping::Mapping::Constraint &constraint, int dimension, RADIAL_BASIS_FUNCTION_T function,
108 Args &&...args)
109{
110 return PtrMapping(new typename BackendSelector<T, RADIAL_BASIS_FUNCTION_T>::type(constraint, dimension, function, std::forward<Args>(args)...));
111}
112
113// Constructs the RBF function based on the functionType
114rbf_variant_t constructRBF(BasisFunction functionType, double supportRadius, double shapeParameter)
115{
116 try {
117 switch (functionType) {
119 return mapping::CompactPolynomialC0(supportRadius);
120 }
122 return mapping::CompactPolynomialC2(supportRadius);
123 }
125 return mapping::CompactPolynomialC4(supportRadius);
126 }
128 return mapping::CompactPolynomialC6(supportRadius);
129 }
131 return mapping::CompactPolynomialC8(supportRadius);
132 }
134 return mapping::CompactThinPlateSplinesC2(supportRadius);
135 }
138 }
140 return mapping::VolumeSplines();
141 }
143 return mapping::Multiquadrics(shapeParameter);
144 }
146 return mapping::InverseMultiquadrics(shapeParameter);
147 }
149 return mapping::Gaussian(shapeParameter);
150 }
151 default:
152 PRECICE_UNREACHABLE("No instantiation was found for the selected basis function.");
153 }
154 } catch (std::invalid_argument &e) {
155 logging::Logger _log{"MappingConfiguration"};
156 PRECICE_ERROR(e.what());
157 }
158}
159
160// The actual instantion helper, which avoids enumerating all mapping implementations (more will come) with all RBF kernels
161// The first three arguments of the constructor are prescribed: constraint, dimension and the RBF function object, all other
162// constructor arguments are just forwarded. The first argument (BasisFunction) indicates then the actual instantiation to return.
163template <RBFBackend T, typename... Args>
164PtrMapping getRBFMapping(BasisFunction functionType, mapping::Mapping::Constraint &constraint, int dimension, double supportRadius, double shapeParameter,
165 Args &&...args)
166{
167 // First, construct the RBF function
168 auto functionVariant = constructRBF(functionType, supportRadius, shapeParameter);
169 // ... and instantiate the corresponding RBF mapping class
170 return std::visit([&](auto &&func) { return instantiateRBFMapping<T>(constraint, dimension, func, std::forward<Args>(args)...); }, functionVariant);
171}
172} // namespace
173
175 xml::XMLTag &parent,
176 mesh::PtrMeshConfiguration meshConfiguration)
177 : _meshConfig(std::move(meshConfiguration))
178{
180 using namespace xml;
181
182 // First, we create the available tags
183 XMLTag::Occurrence occ = XMLTag::OCCUR_ARBITRARY;
184 std::list<XMLTag> projectionTags{
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.")};
189 std::list<XMLTag> rbfDirectTags{
190 XMLTag{*this, TYPE_RBF_GLOBAL_DIRECT, occ, TAG}.setDocumentation("Radial-basis-function mapping using a direct solver with a gather-scatter parallelism.")};
191 std::list<XMLTag> rbfIterativeTags{
192 XMLTag{*this, TYPE_RBF_GLOBAL_ITERATIVE, occ, TAG}.setDocumentation("Radial-basis-function mapping using an iterative solver with a distributed parallelism.")};
193 std::list<XMLTag> pumDirectTags{
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.")};
195 std::list<XMLTag> rbfAliasTag{
196 XMLTag{*this, TYPE_RBF_ALIAS, occ, TAG}.setDocumentation("Alias tag, which auto-selects a radial-basis-function mapping depending on the simulation parameter,")};
197 std::list<XMLTag> geoMultiscaleTags{
198 XMLTag{*this, TYPE_AXIAL_GEOMETRIC_MULTISCALE, occ, TAG}.setDocumentation("Axial geometric multiscale mapping between one 1D and multiple 3D vertices."),
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.")};
200
201 // List of all attributes with corresponding documentation
202 auto attrDirection = XMLAttribute<std::string>(ATTR_DIRECTION)
203 .setOptions({DIRECTION_WRITE, DIRECTION_READ})
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.");
206
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\").");
209
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\").");
212
213 auto attrConstraint = XMLAttribute<std::string>(ATTR_CONSTRAINT)
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");
222 auto attrPolynomial = makeXMLAttribute(ATTR_POLYNOMIAL, POLYNOMIAL_SEPARATE)
223 .setDocumentation("Toggles use of the global polynomial")
225 auto attrPumPolynomial = makeXMLAttribute(ATTR_POLYNOMIAL, POLYNOMIAL_SEPARATE)
226 .setDocumentation("Toggles use a local (per cluster) polynomial")
227 .setOptions({POLYNOMIAL_OFF, POLYNOMIAL_SEPARATE});
228
229 auto attrSolverRtol = makeXMLAttribute(ATTR_SOLVER_RTOL, 1e-9)
230 .setDocumentation("Solver relative tolerance for convergence");
231 // TODO: Discuss whether we wanto to introduce this attribute
232 // auto attrMaxIterations = makeXMLAttribute(ATTR_MAX_ITERATIONS, 1e6)
233 // .setDocumentation("Maximum number of iterations of the solver");
234
235 auto verticesPerCluster = XMLAttribute<int>(ATTR_VERTICES_PER_CLUSTER, 50)
236 .setDocumentation("Average number of vertices per cluster (partition) applied in the rbf partition of unity method.");
237 auto relativeOverlap = makeXMLAttribute(ATTR_RELATIVE_OVERLAP, 0.15)
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.");
239 auto projectToInput = XMLAttribute<bool>(ATTR_PROJECT_TO_INPUT, true)
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.");
241
242 auto attrGeoMultiscaleType = XMLAttribute<std::string>(ATTR_GEOMETRIC_MULTISCALE_TYPE)
243 .setDocumentation("Type of geometric multiscale mapping. Either 'spread' or 'collect'.")
245 auto attrGeoMultiscaleAxis = XMLAttribute<std::string>(ATTR_GEOMETRIC_MULTISCALE_AXIS)
246 .setDocumentation("Principle axis along which geometric multiscale mapping is performed.")
248 auto attrGeoMultiscaleRadius = XMLAttribute<double>(ATTR_GEOMETRIC_MULTISCALE_RADIUS)
249 .setDocumentation("Radius of the circular interface between the 1D and 3D participant.");
250
251 // Add the relevant attributes to the relevant tags
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});
258
259 // Now we take care of the subtag executor. We repeat some of the subtags in order to add individual documentation
260 XMLTag::Occurrence once = XMLTag::OCCUR_NOT_OR_ONCE;
261 // TODO, make type an int
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.");
266
267 // First, we have the executors for the direct solvers
268 {
269 std::list<XMLTag> cpuExecutor{
270 XMLTag{*this, EXECUTOR_CPU, once, SUBTAG_EXECUTOR}.setDocumentation("The default executor, which uses a single-core CPU with a gather-scatter parallelism.")};
271 std::list<XMLTag> deviceExecutors{
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.")};
274
275 addAttributes(deviceExecutors, {attrDeviceId});
276 addSubtagsToParents(cpuExecutor, rbfDirectTags);
277 addSubtagsToParents(deviceExecutors, rbfDirectTags);
278 }
279 // Second, the executors for the iterative solver
280 {
281 std::list<XMLTag> cpuExecutor{
282 XMLTag{*this, EXECUTOR_CPU, once, SUBTAG_EXECUTOR}.setDocumentation("The default executor relying on PETSc, which uses CPUs and distributed memory parallelism via MPI.")};
283 std::list<XMLTag> deviceExecutors{
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.")};
286 std::list<XMLTag> ompExecutor{
287 XMLTag{*this, EXECUTOR_OMP, once, SUBTAG_EXECUTOR}.setDocumentation("OpenMP executor, which uses Ginkgo with a gather-scatter parallelism.")};
288
289 addAttributes(deviceExecutors, {attrDeviceId});
290 addAttributes(ompExecutor, {attrNThreads});
291 addSubtagsToParents(cpuExecutor, rbfIterativeTags);
292 addSubtagsToParents(deviceExecutors, rbfIterativeTags);
293 addSubtagsToParents(ompExecutor, rbfIterativeTags);
294 }
295 {
296 std::list<XMLTag> cpuExecutor{
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);
299 }
300 // The alias tag doesn't receive the subtag at all
301
302 // Now we take care of the subtag basis function
303 // First, we have the tags using a support radius
304 std::list<XMLTag> supportRadiusRBF{
305 XMLTag{*this, RBF_CPOLYNOMIAL_C0, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Wendland C0 function"),
306 XMLTag{*this, RBF_CPOLYNOMIAL_C2, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Wendland C2 function"),
307 XMLTag{*this, RBF_CPOLYNOMIAL_C4, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Wendland C4 function"),
308 XMLTag{*this, RBF_CPOLYNOMIAL_C6, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Wendland C6 function"),
309 XMLTag{*this, RBF_CPOLYNOMIAL_C8, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Wendland C8 function"),
310 XMLTag{*this, RBF_CTPS_C2, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Compact thin-plate-spline C2")};
311
312 auto attrSupportRadius = XMLAttribute<double>(ATTR_SUPPORT_RADIUS)
313 .setDocumentation("Support radius of each RBF basis function (global choice).");
314
315 addAttributes(supportRadiusRBF, {attrSupportRadius});
316 addSubtagsToParents(supportRadiusRBF, rbfIterativeTags);
317 addSubtagsToParents(supportRadiusRBF, rbfDirectTags);
318 addSubtagsToParents(supportRadiusRBF, pumDirectTags);
319 addSubtagsToParents(supportRadiusRBF, rbfAliasTag);
320
321 // Now the tags using a shape parameter
322 std::list<XMLTag> shapeParameterRBF{
323 XMLTag{*this, RBF_MULTIQUADRICS, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Multiquadrics"),
324 XMLTag{*this, RBF_INV_MULTIQUADRICS, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Inverse multiquadrics")};
325
326 auto attrShapeParam = XMLAttribute<double>(ATTR_SHAPE_PARAM)
327 .setDocumentation("Specific shape parameter for RBF basis function.");
328
329 addAttributes(shapeParameterRBF, {attrShapeParam});
330 addSubtagsToParents(shapeParameterRBF, rbfIterativeTags);
331 addSubtagsToParents(shapeParameterRBF, rbfDirectTags);
332 addSubtagsToParents(shapeParameterRBF, pumDirectTags);
333 addSubtagsToParents(shapeParameterRBF, rbfAliasTag);
334
335 // For the Gaussian, we need default values as the user can pass a support radius or a shape parameter
336 std::list<XMLTag> GaussRBF{
337 XMLTag{*this, RBF_GAUSSIAN, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Gaussian basis function accepting a support radius or a shape parameter.")};
338 attrShapeParam.setDefaultValue(std::numeric_limits<double>::quiet_NaN());
339 attrSupportRadius.setDefaultValue(std::numeric_limits<double>::quiet_NaN());
340 addAttributes(GaussRBF, {attrShapeParam, attrSupportRadius});
341 addSubtagsToParents(GaussRBF, rbfIterativeTags);
342 addSubtagsToParents(GaussRBF, rbfDirectTags);
343 addSubtagsToParents(GaussRBF, pumDirectTags);
344 addSubtagsToParents(GaussRBF, rbfAliasTag);
345
346 // tags without an attribute
347 std::list<XMLTag> attributelessRBFs{
348 XMLTag{*this, RBF_TPS, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Thin-plate-splines"),
349 XMLTag{*this, RBF_VOLUME_SPLINES, once, SUBTAG_BASIS_FUNCTION}.setDocumentation("Volume splines")};
350
351 addSubtagsToParents(attributelessRBFs, rbfIterativeTags);
352 addSubtagsToParents(attributelessRBFs, rbfDirectTags);
353 addSubtagsToParents(attributelessRBFs, pumDirectTags);
354 addSubtagsToParents(attributelessRBFs, rbfAliasTag);
355
356 // Add all tags to the mapping tag
357 parent.addSubtags(projectionTags);
358 parent.addSubtags(rbfIterativeTags);
359 parent.addSubtags(rbfDirectTags);
360 parent.addSubtags(pumDirectTags);
361 parent.addSubtags(rbfAliasTag);
362 parent.addSubtags(geoMultiscaleTags);
363}
364
366 bool experimental)
367{
368 _experimental = experimental;
369}
370
372 const xml::ConfigurationContext &context,
373 xml::XMLTag &tag)
374{
375 PRECICE_TRACE(tag.getName());
376 if (tag.getNamespace() == TAG) {
377 // Mandatory tags
381 std::string type = tag.getName();
383
384 // Check that either of the two is provided
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.");
387
388 // Restrict to read-consistent and write-conservative for just-in-time mapping
389
390 // The from mesh cannot be empty due to the check above
391 PRECICE_CHECK(!toMesh.empty() || (toMesh.empty() && dir == DIRECTION_READ && constraint == CONSTRAINT_CONSISTENT),
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);
393 PRECICE_CHECK(!fromMesh.empty() || (fromMesh.empty() && dir == DIRECTION_WRITE && constraint == CONSTRAINT_CONSERVATIVE),
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);
395
396 PRECICE_INFO_IF(toMesh.empty(), "Using just-in-time mapping from mesh \"{}\"", fromMesh);
397 PRECICE_INFO_IF(fromMesh.empty(), "Using just-in-time mapping to mesh \"{}\"", toMesh);
398
399 // optional tags
400 // We set here default values, but their actual value doesn't really matter.
401 // It's just for the mapping methods, which do not use these attributes at all.
402 bool xDead = tag.getBooleanAttributeValue(ATTR_X_DEAD, false);
403 bool yDead = tag.getBooleanAttributeValue(ATTR_Y_DEAD, false);
404 bool zDead = tag.getBooleanAttributeValue(ATTR_Z_DEAD, false);
405 double solverRtol = tag.getDoubleAttributeValue(ATTR_SOLVER_RTOL, 1e-9);
407
408 // geometric multiscale related tags
411 double multiscaleRadius = tag.getDoubleAttributeValue(ATTR_GEOMETRIC_MULTISCALE_RADIUS, 1.0);
412
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.");
415 }
416
417 if (type == TYPE_AXIAL_GEOMETRIC_MULTISCALE && context.size > 1) {
418 throw std::runtime_error{"Axial geometric multiscale mapping is not available for parallel participants."};
419 }
420
421 if (type == TYPE_RADIAL_GEOMETRIC_MULTISCALE && context.size > 1) {
422 throw std::runtime_error{"Radial geometric multiscale mapping is not available for parallel participants."};
423 }
424
425 // pum related tags
426 int verticesPerCluster = tag.getIntAttributeValue(ATTR_VERTICES_PER_CLUSTER, 100);
427 double relativeOverlap = tag.getDoubleAttributeValue(ATTR_RELATIVE_OVERLAP, 0.3);
428 bool projectToInput = tag.getBooleanAttributeValue(ATTR_PROJECT_TO_INPUT, true);
429
430 // Convert raw string into enum types as the constructors take enums
431 if (constraint == CONSTRAINT_CONSERVATIVE) {
433 } else if (constraint == CONSTRAINT_CONSISTENT) {
435 } else if (constraint == CONSTRAINT_SCALED_CONSISTENT_SURFACE) {
437 } else if (constraint == CONSTRAINT_SCALED_CONSISTENT_VOLUME) {
439 } else {
440 PRECICE_UNREACHABLE("Unknown mapping constraint \"{}\".", constraint);
441 }
442
443 ConfiguredMapping configuredMapping = createMapping(dir, type, fromMesh, toMesh, geoMultiscaleType, geoMultiscaleAxis, multiscaleRadius);
444
445 _rbfConfig = configureRBFMapping(type, strPolynomial, xDead, yDead, zDead, solverRtol, verticesPerCluster, relativeOverlap, projectToInput);
446
447 checkDuplicates(configuredMapping);
448 _mappings.push_back(configuredMapping);
449 } else if (tag.getNamespace() == SUBTAG_BASIS_FUNCTION) {
450
451 PRECICE_ASSERT(!_mappings.empty());
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.",
455 _mappings.back().fromMesh->getName(), _mappings.back().toMesh->getName());
456 // We can only set one subtag
457 PRECICE_CHECK(_rbfConfig.basisFunctionDefined == false, "More than one basis-function was defined for the mapping "
458 "from mesh \"{}\" to mesh \"{}\".",
459 _mappings.back().fromMesh->getName(), _mappings.back().toMesh->getName());
460
461 std::string basisFctName = tag.getName();
462 double supportRadius = tag.getDoubleAttributeValue(ATTR_SUPPORT_RADIUS, 0.);
463 double shapeParameter = tag.getDoubleAttributeValue(ATTR_SHAPE_PARAM, 0.);
464
465 _rbfConfig.basisFunction = parseBasisFunctions(basisFctName);
466 _rbfConfig.basisFunctionDefined = true;
467 // The Gaussian RBF is always treated as a shape-parameter RBF. Hence, we have to convert the support radius, if necessary
468 if (_rbfConfig.basisFunction == BasisFunction::Gaussian) {
469 const bool exactlyOneSet = (std::isfinite(supportRadius) && !std::isfinite(shapeParameter)) ||
470 (std::isfinite(shapeParameter) && !std::isfinite(supportRadius));
471 PRECICE_CHECK(exactlyOneSet, "The specified parameters for the Gaussian RBF mapping are invalid. Please specify either a \"shape-parameter\" or a \"support-radius\".");
472
473 if (std::isfinite(supportRadius) && !std::isfinite(shapeParameter)) {
474 shapeParameter = std::sqrt(-std::log(Gaussian::cutoffThreshold)) / supportRadius;
475 }
476 }
477
478 _rbfConfig.supportRadius = supportRadius;
479 _rbfConfig.shapeParameter = shapeParameter;
480 } else if (tag.getNamespace() == SUBTAG_EXECUTOR) {
482
483 if (tag.getName() == EXECUTOR_CPU) {
485 } else if (tag.getName() == EXECUTOR_CUDA) {
487 } else if (tag.getName() == EXECUTOR_HIP) {
489 } else if (tag.getName() == EXECUTOR_OMP) {
491 }
492
495 }
496}
497
499 const std::string &polynomial,
500 bool xDead, bool yDead, bool zDead,
501 double solverRtol,
502 double verticesPerCluster,
503 double relativeOverlap,
504 bool projectToInput) const
505{
507
508 if (type == TYPE_RBF_GLOBAL_ITERATIVE)
510 else if (type == TYPE_RBF_GLOBAL_DIRECT)
512 else if (type == TYPE_RBF_PUM_DIRECT)
514 else {
515 // Rather simple auto-selection (for now), consisting of the PUM Eigen backend
517
518 // A more sophisticated criterion here could take the globalNumberOfVertices into account
519 // (the mesh pointer is stored in the configuredMapping anyway), but this quantity is not yet computed
520 // during the configuration time.
521 }
522
523 if (polynomial == POLYNOMIAL_SEPARATE)
524 rbfConfig.polynomial = Polynomial::SEPARATE;
525 else if (polynomial == POLYNOMIAL_ON)
526 rbfConfig.polynomial = Polynomial::ON;
527 else if (polynomial == POLYNOMIAL_OFF)
528 rbfConfig.polynomial = Polynomial::OFF;
529 else
530 PRECICE_UNREACHABLE("Unknown polynomial configuration.");
531
532 rbfConfig.deadAxis = {{xDead, yDead, zDead}};
533 rbfConfig.solverRtol = solverRtol;
534
535 rbfConfig.verticesPerCluster = verticesPerCluster;
536 rbfConfig.relativeOverlap = relativeOverlap;
537 rbfConfig.projectToInput = projectToInput;
538
539 return rbfConfig;
540}
541
543 const std::string &direction,
544 const std::string &type,
545 const std::string &fromMeshName,
546 const std::string &toMeshName,
547 const std::string &geoMultiscaleType,
548 const std::string &geoMultiscaleAxis,
549 const double &multiscaleRadius) const
550{
551 PRECICE_TRACE(direction, type);
552
553 ConfiguredMapping configuredMapping;
554 mesh::PtrMesh fromMesh(_meshConfig->getMesh(fromMeshName));
555 mesh::PtrMesh toMesh(_meshConfig->getMesh(toMeshName));
556
557 // Handle for just-in-time mapping, we copy over the dimension and leave everything else
558 if (!fromMesh && fromMeshName.empty()) {
559 PRECICE_CHECK(toMesh,
560 "Mesh \"{0}\" was not found while creating a mapping. "
561 "Please correct the to=\"{0}\" attribute.",
562 toMeshName);
563 fromMesh = mesh::MeshConfiguration::getJustInTimeMappingMesh(toMesh->getDimensions());
564 }
565 if (!toMesh && toMeshName.empty()) {
566 PRECICE_CHECK(fromMesh,
567 "Mesh \"{0}\" was not found while creating a mapping. "
568 "Please correct the from=\"{0}\" attribute.",
569 fromMeshName);
570 toMesh = mesh::MeshConfiguration::getJustInTimeMappingMesh(fromMesh->getDimensions());
571 }
572
573 PRECICE_CHECK(fromMesh,
574 "Mesh \"{0}\" was not found while creating a mapping. "
575 "Please correct the from=\"{0}\" attribute.",
576 fromMeshName);
577 PRECICE_CHECK(toMesh,
578 "Mesh \"{0}\" was not found while creating a mapping. "
579 "Please correct the to=\"{0}\" attribute.",
580 toMeshName);
581
582 PRECICE_CHECK((!toMesh->isJustInTime() && !fromMesh->isJustInTime()) ||
583 ((toMesh->isJustInTime() || fromMesh->isJustInTime()) && (type == TYPE_NEAREST_NEIGHBOR || type == TYPE_RBF_PUM_DIRECT || type == TYPE_RBF_ALIAS)),
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:{}\".",
586 fromMesh->getName(), toMesh->getName(), type, TYPE_NEAREST_NEIGHBOR, TYPE_RBF_ALIAS, TYPE_RBF_PUM_DIRECT);
587
588 // Check for compatible mesh dimensions
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);
594
595 configuredMapping.fromMesh = fromMesh;
596 configuredMapping.toMesh = toMesh;
597
598 if (direction == DIRECTION_WRITE) {
599 configuredMapping.direction = WRITE;
600 } else if (direction == DIRECTION_READ) {
601 configuredMapping.direction = READ;
602 } else {
603 PRECICE_UNREACHABLE("Unknown mapping direction type \"{}\".", direction);
604 }
605
606 // Create all projection based mappings
607 if (type == TYPE_NEAREST_NEIGHBOR) {
608 configuredMapping.mapping = PtrMapping(new NearestNeighborMapping(constraintValue, fromMesh->getDimensions()));
609 } else if (type == TYPE_NEAREST_PROJECTION) {
610 configuredMapping.mapping = PtrMapping(new NearestProjectionMapping(constraintValue, fromMesh->getDimensions()));
611 } else if (type == TYPE_LINEAR_CELL_INTERPOLATION) {
612 configuredMapping.mapping = PtrMapping(new LinearCellInterpolationMapping(constraintValue, fromMesh->getDimensions()));
613 } else if (type == TYPE_NEAREST_NEIGHBOR_GRADIENT) {
614
615 // NNG is not applicable with the conservative constraint
617 "Nearest-neighbor-gradient mapping is not implemented using a \"conservative\" constraint. "
618 "Please select constraint=\" consistent\" or a different mapping method.");
619
620 configuredMapping.mapping = PtrMapping(new NearestNeighborGradientMapping(constraintValue, fromMesh->getDimensions()));
621
622 } else if (type == TYPE_AXIAL_GEOMETRIC_MULTISCALE) {
623
624 // Axial geometric multiscale is not applicable with the conservative constraint
626 "Axial geometric multiscale mapping is not implemented for the \"conservative\" constraint. "
627 "Please select constraint=\" consistent\" or a different mapping method.");
628
629 // Convert strings into enums
631 if (geoMultiscaleAxis == "x") {
633 } else if (geoMultiscaleAxis == "y") {
635 } else if (geoMultiscaleAxis == "z") {
637 } else {
638 PRECICE_UNREACHABLE("Unknown geometric multiscale axis \"{}\".", geoMultiscaleAxis);
639 }
640
642 if (geoMultiscaleType == "spread") {
644 } else if (geoMultiscaleType == "collect") {
646 } else {
647 PRECICE_UNREACHABLE("Unknown geometric multiscale type \"{}\".", geoMultiscaleType);
648 }
649
650 configuredMapping.mapping = PtrMapping(new AxialGeoMultiscaleMapping(constraintValue, fromMesh->getDimensions(), multiscaleType, multiscaleAxis, multiscaleRadius));
651
652 } else if (type == TYPE_RADIAL_GEOMETRIC_MULTISCALE) {
653
654 // Radial geometric multiscale is not applicable with the conservative constraint
656 "Radial geometric multiscale mapping is not implemented for the \"conservative\" constraint. "
657 "Please select constraint=\" consistent\" or a different mapping method.");
658
659 // Convert strings into enums
661 if (geoMultiscaleAxis == "x") {
663 } else if (geoMultiscaleAxis == "y") {
665 } else if (geoMultiscaleAxis == "z") {
667 } else {
668 PRECICE_UNREACHABLE("Unknown geometric multiscale axis \"{}\".", geoMultiscaleAxis);
669 }
670
672 if (geoMultiscaleType == "spread") {
674 } else if (geoMultiscaleType == "collect") {
676 } else {
677 PRECICE_UNREACHABLE("Unknown geometric multiscale type \"{}\".", geoMultiscaleType);
678 }
679
680 configuredMapping.mapping = PtrMapping(new RadialGeoMultiscaleMapping(constraintValue, fromMesh->getDimensions(), multiscaleType, multiscaleAxis));
681
682 } else {
683 // We need knowledge about the basis function in order to instantiate the rbf related mapping
685 configuredMapping.mapping = nullptr;
686 configuredMapping.configuredWithAliasTag = type == TYPE_RBF_ALIAS;
687 }
688
689 configuredMapping.requiresBasisFunction = requiresBasisFunction(type);
690
691 return configuredMapping;
692}
693
695{
696 for (const ConfiguredMapping &configuredMapping : _mappings) {
697 bool sameToMesh = mapping.toMesh->getName() == configuredMapping.toMesh->getName();
698 bool sameFromMesh = mapping.fromMesh->getName() == configuredMapping.fromMesh->getName();
699 bool sameMapping = sameToMesh && sameFromMesh;
700 PRECICE_CHECK(!sameMapping,
701 "There cannot be two mappings from mesh \"{}\" to mesh \"{}\". "
702 "Please remove one of the duplicated meshes. ",
703 mapping.fromMesh->getName(), mapping.toMesh->getName());
704 }
705 for (const ConfiguredMapping &configuredMapping : _mappings) {
706 bool sameToMesh = mapping.toMesh->getName() == configuredMapping.toMesh->getName();
707 bool isWrite = mapping.direction == mapping::MappingConfiguration::WRITE && configuredMapping.direction == mapping::MappingConfiguration::WRITE;
708 bool sameMapping = sameToMesh && isWrite && (mapping.fromMesh->isJustInTime() || configuredMapping.fromMesh->isJustInTime());
709 PRECICE_CHECK(!sameMapping,
710 "There cannot be two mappings to mesh \"{}\". "
711 "Here, we have a mixture of just-in-time mapping and a conventional mapping. ",
712 mapping.toMesh->getName());
713 }
714}
715
717{
718 if (tag.getNamespace() == TAG) {
719 if (requiresBasisFunction(tag.getName())) {
720 PRECICE_CHECK(_rbfConfig.basisFunctionDefined, "No basis-function was defined for the \"{}\" mapping from mesh \"{}\" to mesh \"{}\".", tag.getName(), _mappings.back().fromMesh->getName(), _mappings.back().toMesh->getName());
721 if (!_executorConfig) {
723 }
725 _executorConfig.reset();
726 }
727 PRECICE_ASSERT(_mappings.back().mapping != nullptr);
728 }
729}
730
732{
735
736 // Instantiate the RBF mapping classes
737 // We first categorize according to the executor
738 // 1. the CPU executor
741 mapping.mapping = getRBFMapping<RBFBackend::Eigen>(_rbfConfig.basisFunction, constraintValue, mapping.fromMesh->getDimensions(), _rbfConfig.supportRadius, _rbfConfig.shapeParameter, _rbfConfig.deadAxis, _rbfConfig.polynomial);
743#ifndef PRECICE_NO_PETSC
744 // for petsc initialization
746
747 mapping.mapping = getRBFMapping<RBFBackend::PETSc>(_rbfConfig.basisFunction, constraintValue, mapping.fromMesh->getDimensions(), _rbfConfig.supportRadius, _rbfConfig.shapeParameter, _rbfConfig.deadAxis, _rbfConfig.solverRtol, _rbfConfig.polynomial);
748#else
749 PRECICE_CHECK(false, "The global-iterative RBF solver on a CPU requires a preCICE build with PETSc enabled.");
750#endif
752 mapping.mapping = getRBFMapping<RBFBackend::PUM>(_rbfConfig.basisFunction, constraintValue, mapping.fromMesh->getDimensions(), _rbfConfig.supportRadius, _rbfConfig.shapeParameter, _rbfConfig.polynomial, _rbfConfig.verticesPerCluster, _rbfConfig.relativeOverlap, _rbfConfig.projectToInput);
753 } else {
754 PRECICE_UNREACHABLE("Unknown RBF solver.");
755 }
756 // 2. any other executor is configured via Ginkgo
757 } else {
758#ifndef PRECICE_NO_GINKGO
760 _ginkgoParameter.usePreconditioner = false;
761 _ginkgoParameter.deviceId = _executorConfig->deviceId;
763 _ginkgoParameter.executor = "cuda-executor";
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());
766#endif
768 _ginkgoParameter.executor = "hip-executor";
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());
771#endif
773 _ginkgoParameter.executor = "omp-executor";
774 _ginkgoParameter.nThreads = _executorConfig->nThreads;
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());
777#endif
778 }
780 _ginkgoParameter.solver = "qr-solver";
782 _ginkgoParameter.solver = "cg-solver";
783 _ginkgoParameter.residualNorm = _rbfConfig.solverRtol;
784 } else {
785 PRECICE_UNREACHABLE("Unknown solver type.");
786 }
787
788 mapping.mapping = getRBFMapping<RBFBackend::Ginkgo>(_rbfConfig.basisFunction, constraintValue, mapping.fromMesh->getDimensions(), _rbfConfig.supportRadius, _rbfConfig.shapeParameter, _rbfConfig.deadAxis, _rbfConfig.polynomial, _ginkgoParameter);
789#else
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());
791#endif
792 }
793}
794
799
801{
802 return mappingType == TYPE_RBF_PUM_DIRECT || mappingType == TYPE_RBF_GLOBAL_DIRECT || mappingType == TYPE_RBF_GLOBAL_ITERATIVE || mappingType == TYPE_RBF_ALIAS;
803}
804
806{
807 BasisFunction basisFunction{};
808 if (basisFctName == RBF_TPS)
809 basisFunction = BasisFunction::ThinPlateSplines;
810 else if (basisFctName == RBF_MULTIQUADRICS)
811 basisFunction = BasisFunction::Multiquadrics;
812 else if (basisFctName == RBF_INV_MULTIQUADRICS)
814 else if (basisFctName == RBF_VOLUME_SPLINES)
815 basisFunction = BasisFunction::VolumeSplines;
816 else if (basisFctName == RBF_GAUSSIAN)
817 basisFunction = BasisFunction::Gaussian;
818 else if (basisFctName == RBF_CTPS_C2)
820 else if (basisFctName == RBF_CPOLYNOMIAL_C0)
821 basisFunction = BasisFunction::WendlandC0;
822 else if (basisFctName == RBF_CPOLYNOMIAL_C2)
823 basisFunction = BasisFunction::WendlandC2;
824 else if (basisFctName == RBF_CPOLYNOMIAL_C4)
825 basisFunction = BasisFunction::WendlandC4;
826 else if (basisFctName == RBF_CPOLYNOMIAL_C6)
827 basisFunction = BasisFunction::WendlandC6;
828 else if (basisFctName == RBF_CPOLYNOMIAL_C8)
829 basisFunction = BasisFunction::WendlandC8;
830 else
831 PRECICE_UNREACHABLE("Unknown basis function \"{}\".", basisFctName);
832 return basisFunction;
833}
834} // namespace precice::mapping
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:16
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_INFO_IF(condition,...)
Definition LogMacros.hpp:25
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
#define PRECICE_UNREACHABLE(...)
Definition assertion.hpp:93
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...
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
void xmlTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
Callback function required for use of automatic configuration.
const std::vector< ConfiguredMapping > & mappings()
Returns all configured mappings.
std::unique_ptr< ExecutorConfiguration > _executorConfig
bool requiresBasisFunction(const std::string &mappingType) const
BasisFunction parseBasisFunctions(const std::string &basisFctName) const
Given a basis function name (as a string), transforms the string into an enum of the BasisFunction.
std::vector< ConfiguredMapping > _mappings
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 RBFConfiguration & rbfConfig() const
MappingConfiguration(xml::XMLTag &parent, mesh::PtrMeshConfiguration meshConfiguration)
Constraint
Specifies additional constraints for a mapping.
Definition Mapping.hpp:30
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.
Definition Parallel.cpp:147
static void initialize(utils::Parallel::Communicator comm)
Initializes the Petsc environment.
Definition Petsc.cpp:38
Represents an XML tag to be configured automatically.
Definition XMLTag.hpp:28
const std::string & getNamespace() const
Returns xml namespace.
Definition XMLTag.hpp:159
void addSubtags(const Container &subtags)
Definition XMLTag.hpp:143
std::string getStringAttributeValue(const std::string &name, std::optional< std::string > default_value=std::nullopt) const
Definition XMLTag.cpp:145
bool getBooleanAttributeValue(const std::string &name, std::optional< bool > default_value=std::nullopt) const
Definition XMLTag.cpp:159
const std::string & getName() const
Returns name (without namespace).
Definition XMLTag.hpp:153
int getIntAttributeValue(const std::string &name, std::optional< int > default_value=std::nullopt) const
Definition XMLTag.cpp:131
double getDoubleAttributeValue(const std::string &name, std::optional< double > default_value=std::nullopt) const
Definition XMLTag.cpp:117
T empty(T... args)
T forward(T... args)
T isfinite(T... args)
T log(T... args)
T make_unique(T... args)
contains data mapping from points to meshes.
std::shared_ptr< Mapping > PtrMapping
std::shared_ptr< Mesh > PtrMesh
std::shared_ptr< MeshConfiguration > PtrMeshConfiguration
contains the XML configuration parser.
STL namespace.
static precice::logging::Logger _log("precicec")
T quiet_NaN(T... args)
T sqrt(T... args)
Direction direction
Direction of mapping (important to set input and output mesh).
bool configuredWithAliasTag
used the automatic rbf alias tag in order to set the mapping
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21
T visit(T... args)