preCICE v3.1.2
Loading...
Searching...
No Matches
AccelerationConfiguration.cpp
Go to the documentation of this file.
2#include <algorithm>
3#include <memory>
4#include <ostream>
5#include <stdexcept>
6#include <utility>
7#include <vector>
18#include "logging/LogMacros.hpp"
19#include "mesh/Data.hpp"
20#include "mesh/Mesh.hpp"
22#include "utils/assertion.hpp"
23#include "xml/ConfigParser.hpp"
24#include "xml/XMLAttribute.hpp"
25#include "xml/XMLTag.hpp"
26
27namespace precice::acceleration {
28
29using namespace precice::acceleration::impl;
30
32 const mesh::PtrMeshConfiguration &meshConfig)
33 : TAG("acceleration"),
34 TAG_RELAX("relaxation"),
35 TAG_INIT_RELAX("initial-relaxation"),
36 TAG_MAX_USED_ITERATIONS("max-used-iterations"),
37 TAG_TIME_WINDOWS_REUSED("time-windows-reused"),
38 TAG_DATA("data"),
39 TAG_FILTER("filter"),
40 TAG_ESTIMATEJACOBIAN("estimate-jacobian"),
41 TAG_PRECONDITIONER("preconditioner"),
42 TAG_IMVJRESTART("imvj-restart-mode"),
43 ATTR_NAME("name"),
44 ATTR_MESH("mesh"),
45 ATTR_SCALING("scaling"),
46 ATTR_VALUE("value"),
47 ATTR_ENFORCE("enforce"),
48 ATTR_SINGULARITYLIMIT("limit"),
49 ATTR_TYPE("type"),
50 ATTR_BUILDJACOBIAN("always-build-jacobian"),
51 ATTR_IMVJCHUNKSIZE("chunk-size"),
52 ATTR_RSLS_REUSED_TIME_WINDOWS("reused-time-windows-at-restart"),
53 ATTR_RSSVD_TRUNCATIONEPS("truncation-threshold"),
54 ATTR_PRECOND_NONCONST_TIME_WINDOWS("freeze-after"),
55 VALUE_CONSTANT("constant"),
56 VALUE_AITKEN("aitken"),
57 VALUE_IQNILS("IQN-ILS"),
58 VALUE_IQNIMVJ("IQN-IMVJ"),
59 VALUE_QR1FILTER("QR1"),
60 VALUE_QR1_ABSFILTER("QR1-absolute"),
61 VALUE_QR2FILTER("QR2"),
62 VALUE_CONSTANT_PRECONDITIONER("constant"),
63 VALUE_VALUE_PRECONDITIONER("value"),
64 VALUE_RESIDUAL_PRECONDITIONER("residual"),
65 VALUE_RESIDUAL_SUM_PRECONDITIONER("residual-sum"),
66 VALUE_LS_RESTART("RS-LS"),
67 VALUE_ZERO_RESTART("RS-0"),
68 VALUE_SVD_RESTART("RS-SVD"),
69 VALUE_SLIDE_RESTART("RS-SLIDE"),
70 VALUE_NO_RESTART("no-restart"),
71 _meshConfig(meshConfig),
72 _acceleration(),
73 _neededMeshes(),
74 _preconditioner(),
75 _config()
76{
77 PRECICE_ASSERT(meshConfig.get() != nullptr);
78}
79
81{
82 using namespace xml;
83
84 // static int recursionCounter = 0;
85 // recursionCounter++;
86
87 XMLTag::Occurrence occ = XMLTag::OCCUR_NOT_OR_ONCE;
89 {
90 XMLTag tag(*this, VALUE_CONSTANT, occ, TAG);
91 tag.setDocumentation("Accelerates coupling data with constant underrelaxation.");
93 tags.push_back(tag);
94 }
95 {
96 XMLTag tag(*this, VALUE_AITKEN, occ, TAG);
97 tag.setDocumentation("Accelerates coupling data with dynamic Aitken under-relaxation.");
99 tags.push_back(tag);
100 }
101 {
102 XMLTag tag(*this, VALUE_IQNILS, occ, TAG);
103 tag.setDocumentation("Accelerates coupling data with the interface quasi-Newton inverse least-squares method.");
105 tags.push_back(tag);
106 }
107 {
108 XMLTag tag(*this, VALUE_IQNIMVJ, occ, TAG);
109 tag.setDocumentation("Accelerates coupling data with the interface quasi-Newton inverse multi-vector Jacobian method.");
110
111 auto alwaybuildJacobian = makeXMLAttribute(ATTR_BUILDJACOBIAN, false)
112 .setDocumentation("If set to true, the IMVJ will set up the Jacobian matrix"
113 " in each coupling iteration, which is inefficient. If set to false (or not set)"
114 " the Jacobian is only build in the last iteration and the updates are computed using (relatively) cheap MATVEC products.");
115 tag.addAttribute(alwaybuildJacobian);
116
118 tags.push_back(tag);
119 }
120
121 for (XMLTag &tag : tags) {
122 parent.addSubtag(tag);
123 }
124}
125
130
132 const xml::ConfigurationContext &context,
133 xml::XMLTag & callingTag)
134{
135 PRECICE_TRACE(callingTag.getFullName());
136
137 if (callingTag.getNamespace() == TAG) {
138 _config.type = callingTag.getName();
139
142 }
143 if (callingTag.getName() == TAG_RELAX) {
145 } else if (callingTag.getName() == TAG_DATA) {
146 std::string dataName = callingTag.getStringAttributeValue(ATTR_NAME);
147 std::string meshName = callingTag.getStringAttributeValue(ATTR_MESH);
148 auto success = _uniqueDataAndMeshNames.emplace(dataName, meshName);
149 PRECICE_CHECK(success.second,
150 "You have provided a subtag <data name=\"{}\" mesh=\"{}\"/> more than once in your <acceleration:.../>. "
151 "Please remove the duplicated entry.",
152 dataName, meshName);
153
155 double scaling = 1.0;
157 scaling = callingTag.getDoubleAttributeValue(ATTR_SCALING);
158 }
159
160 PRECICE_CHECK(_meshConfig->hasMeshName(_meshName) && _meshConfig->getMesh(_meshName)->hasDataName(dataName),
161 "Data with name \"{0}\" associated to mesh \"{1}\" not found on configuration of acceleration. "
162 "Add \"{0}\" to the \"<mesh name={1}>\" tag, or change the data name in the acceleration scheme.",
163 dataName, _meshName);
164
165 const mesh::PtrMesh &mesh = _meshConfig->getMesh(_meshName);
166 const mesh::PtrData &data = mesh->data(dataName);
167 _config.dataIDs.push_back(data->getID());
168 _config.scalings.insert(std::make_pair(data->getID(), scaling));
169
171 } else if (callingTag.getName() == TAG_INIT_RELAX) {
174 if (callingTag.hasAttribute(ATTR_ENFORCE)) {
176 } else {
178 }
179 } else if (callingTag.getName() == TAG_MAX_USED_ITERATIONS) {
182 } else if (callingTag.getName() == TAG_TIME_WINDOWS_REUSED) {
185 } else if (callingTag.getName() == TAG_FILTER) {
187 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
188 if (f == VALUE_QR1FILTER) {
190 } else if (f == VALUE_QR1_ABSFILTER) {
192 } else if (f == VALUE_QR2FILTER) {
194 } else {
195 PRECICE_ASSERT(false);
196 }
198 } else if (callingTag.getName() == TAG_PRECONDITIONER) {
202 } else if (callingTag.getName() == TAG_IMVJRESTART) {
203
204#ifndef PRECICE_NO_MPI
206 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
207 PRECICE_CHECK((f == VALUE_NO_RESTART) || (!_config.alwaysBuildJacobian), "IMVJ cannot be in restart mode while parameter always-build-jacobian is set to true. "
208 "Please remove 'always-build-jacobian' from the configuration file or do not run in restart mode.");
209 if (f == VALUE_NO_RESTART) {
211 } else if (f == VALUE_ZERO_RESTART) {
213 } else if (f == VALUE_LS_RESTART) {
216 } else if (f == VALUE_SVD_RESTART) {
219 } else if (f == VALUE_SLIDE_RESTART) {
221 } else {
223 PRECICE_ASSERT(false);
224 }
225#else
226 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
227#endif
228 }
229}
230
232 const xml::ConfigurationContext &context,
233 xml::XMLTag & callingTag)
234{
235 PRECICE_TRACE(callingTag.getName());
236 if (callingTag.getNamespace() == TAG) {
237
238 //create preconditioner
239 if (callingTag.getName() == VALUE_IQNILS || callingTag.getName() == VALUE_IQNIMVJ || callingTag.getName() == VALUE_AITKEN) {
240
241 // if imvj restart-mode is of type RS-SVD, max number of non-const preconditioned time windows is limited by the chunksize
242 if (callingTag.getName() == VALUE_IQNIMVJ && _config.imvjRestartType > 0)
253 } else {
254 // no preconditioner defined
256 }
257 }
258
259 if (callingTag.getName() == VALUE_CONSTANT) {
263 } else if (callingTag.getName() == VALUE_AITKEN) {
267 } else if (callingTag.getName() == VALUE_IQNILS) {
282 } else if (callingTag.getName() == VALUE_IQNIMVJ) {
283#ifndef PRECICE_NO_MPI
303#else
304 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
305#endif
306 } else {
307 PRECICE_ASSERT(false);
308 }
309 }
310}
311
318
320{
321 using namespace precice::xml;
322
323 XMLTag tagData(*this, TAG_DATA, XMLTag::OCCUR_ONCE_OR_MORE);
324 tagData.setDocumentation("The data used to compute the acceleration.");
326 attrName.setDocumentation("The name of the data.");
328 attrMesh.setDocumentation("The name of the mesh which holds the data.");
329 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
330 .setDocumentation(
331 "To improve the performance of a parallel or a multi coupling schemes, "
332 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
333 "We recommend, however, to use an automatic scaling via a preconditioner.");
334 tagData.addAttribute(attrScaling);
335 tagData.addAttribute(attrName);
336 tagData.addAttribute(attrMesh);
337 tag.addSubtag(tagData);
338
339 XMLTag tagFilter(*this, TAG_FILTER, XMLTag::OCCUR_NOT_OR_ONCE);
340 tagFilter.setDocumentation("Type of filtering technique that is used to "
341 "maintain good conditioning in the least-squares system. Possible filters:\n"
342 " - `QR1-filter`: updateQR-dec with (relative) test \\\\(R(i,i) < \\epsilon *\\lVert R\\rVert_F\\\\)\n"
343 " - `QR1_absolute-filter`: updateQR-dec with (absolute) test \\\\(R(i, i) < \\epsilon\\\\)\n"
344 " - `QR2-filter`: en-block QR-dec with test \\\\(\\lVert v_\\text{orth} \\rVert_2 < \\epsilon * \\lVert v \\rVert_2\\\\)\n\n"
345 "Please note that a QR1 is based on Given's rotations whereas QR2 uses "
346 "modified Gram-Schmidt. This can give different results even when no columns "
347 "are filtered out.\n"
348 "When this tag is not provided, the QR2-filter with the limit value 1e-2 is used.");
349 XMLAttribute<double> attrSingularityLimit(ATTR_SINGULARITYLIMIT, 1e-16);
350 attrSingularityLimit.setDocumentation("Limit eps of the filter.");
351 tagFilter.addAttribute(attrSingularityLimit);
352 auto attrFilterName = XMLAttribute<std::string>(ATTR_TYPE)
356 .setDocumentation("Type of the filter.");
357 tagFilter.addAttribute(attrFilterName);
358 tag.addSubtag(tagFilter);
359}
360
362 xml::XMLTag &tag)
363{
364 using namespace xml;
365 if (tag.getName() == VALUE_CONSTANT) {
366 XMLTag tagRelax(*this, TAG_RELAX, XMLTag::OCCUR_ONCE);
368 attrValue.setDocumentation("Constant relaxation factor.");
369 tagRelax.addAttribute(attrValue);
370 tag.addSubtag(tagRelax);
371 } else if (tag.getName() == VALUE_AITKEN) {
372 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_ONCE);
373 tagInitRelax.setDocumentation("Initial relaxation factor.");
375 attrValue.setDocumentation("Initial relaxation factor.");
376 tagInitRelax.addAttribute(attrValue);
378
379 XMLTag tagData(*this, TAG_DATA, XMLTag::OCCUR_ONCE_OR_MORE);
380 tagData.setDocumentation("The data used to compute the acceleration.");
382 attrName.setDocumentation("The name of the data.");
384 attrMesh.setDocumentation("The name of the mesh which holds the data.");
385 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
387 "To improve the performance of a parallel or a multi coupling schemes, "
388 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
389 "We recommend, however, to use an automatic scaling via a preconditioner.");
390 tagData.addAttribute(attrScaling);
391 tagData.addAttribute(attrName);
392 tagData.addAttribute(attrMesh);
393 tag.addSubtag(tagData);
394
395 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
396 tagPreconditioner.setDocumentation("To improve the numerical stability of multiple data vectors a preconditioner"
397 " can be applied. A constant preconditioner scales every acceleration data by a constant value, which you can define as"
398 " an attribute of data. "
399 " A value preconditioner scales every acceleration data by the norm of the data in the previous time window."
400 " A residual preconditioner scales every acceleration data by the current residual."
401 " A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.");
407 .setDocumentation("The type of the preconditioner.");
409 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
411 "After the given number of time windows, the preconditioner weights "
412 "are frozen and the preconditioner acts like a constant preconditioner.");
415 } else if (tag.getName() == VALUE_IQNILS) {
416
417 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
418 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
420 attrDoubleValue.setDocumentation("Initial relaxation factor.");
421 tagInitRelax.addAttribute(attrDoubleValue);
423 attrEnforce.setDocumentation("Enforce initial relaxation in every time window.");
424 tagInitRelax.addAttribute(attrEnforce);
426
427 XMLTag tagMaxUsedIter(*this, TAG_MAX_USED_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
428 tagMaxUsedIter.setDocumentation("Maximum number of columns used in low-rank approximation of Jacobian. If this tag is not provided, the attribute value of 100 is used.");
430 attrIntValue.setDocumentation("The number of columns.");
431 tagMaxUsedIter.addAttribute(attrIntValue);
433
434 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
435 tagTimeWindowsReused.setDocumentation("Number of past time windows from which columns are used to approximate Jacobian. If this tag is not provided, the default attribute value of 10 is used.");
437 attrNumTimeWindowsReused.setDocumentation("The number of time windows.");
440
442
443 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
444 tagPreconditioner.setDocumentation("To improve the performance of a parallel or a multi coupling schemes a preconditioner"
445 " can be applied. "
446 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data. \n "
447 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
448 "- A residual preconditioner scales every acceleration data by the current residual.\n"
449 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n"
450 " If this tag is not provided, the residual-sum preconditioner is employed.");
456 .setDocumentation("The type of the preconditioner.");
458 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
460 "After the given number of time windows, the preconditioner weights "
461 "are frozen and the preconditioner acts like a constant preconditioner.");
464
465 } else if (tag.getName() == VALUE_IQNIMVJ) {
466 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
467 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
468 tagInitRelax.addAttribute(
469 XMLAttribute<double>(ATTR_VALUE).setDocumentation("Initial relaxation factor."));
470 tagInitRelax.addAttribute(
471 XMLAttribute<bool>(ATTR_ENFORCE, false).setDocumentation("Enforce initial relaxation in every time window."));
473
474 XMLTag tagIMVJRESTART(*this, TAG_IMVJRESTART, XMLTag::OCCUR_NOT_OR_ONCE);
481 .setDefaultValue(VALUE_SVD_RESTART)
482 .setDocumentation("Type of the restart mode.");
483 tagIMVJRESTART.addAttribute(attrRestartName);
484 tagIMVJRESTART.setDocumentation("Type of IMVJ restart mode that is used:\n"
485 "- `no-restart`: IMVJ runs in normal mode with explicit representation of Jacobian\n"
486 "- `RS-ZERO`: IMVJ runs in restart mode. After M time windows all Jacobain information is dropped, restart with no information\n"
487 "- `RS-LS`: IMVJ runs in restart mode. After M time windows a IQN-LS like approximation for the initial guess of the Jacobian is computed.\n"
488 "- `RS-SVD`: IMVJ runs in restart mode. After M time windows a truncated SVD of the Jacobian is updated.\n"
489 "- `RS-SLIDE`: IMVJ runs in sliding window restart mode.\n"
490 "If this tag is not provided, IMVJ runs in normal mode with explicit representation of Jacobian.");
491 auto attrChunkSize = makeXMLAttribute(ATTR_IMVJCHUNKSIZE, 8)
492 .setDocumentation("Specifies the number of time windows M after which the IMVJ restarts, if run in restart-mode. Default value is M=8.");
494 .setDocumentation("If IMVJ restart-mode=RS-LS, the number of reused time windows at restart can be specified.");
495 auto attrRSSVD_truncationEps = makeXMLAttribute(ATTR_RSSVD_TRUNCATIONEPS, 1e-4)
496 .setDocumentation("If IMVJ restart-mode=RS-SVD, the truncation threshold for the updated SVD can be set.");
497 tagIMVJRESTART.addAttribute(attrChunkSize);
501
502 XMLTag tagMaxUsedIter(*this, TAG_MAX_USED_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
503 tagMaxUsedIter.setDocumentation("Maximum number of columns used in low-rank approximation of Jacobian. If this tag is not provided, the default attribute value of 20 is used.");
505 attrIntValue.setDocumentation("The number of columns.");
506 tagMaxUsedIter.addAttribute(attrIntValue);
508
509 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
510 tagTimeWindowsReused.setDocumentation("Number of past time windows from which columns are used to approximate Jacobian. If this tag is not provided, the attribute value of 0 is used.");
513
515
516 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
518 "To improve the performance of a parallel or a multi coupling schemes a preconditioner can be applied."
519 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data.\n"
520 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
521 "- A residual preconditioner scales every acceleration data by the current residual.\n"
522 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n"
523 " If this tag is not provided, the residual-sum preconditioner is employed.");
529 .setDocumentation("Type of the preconditioner.");
531 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
532 .setDocumentation("After the given number of time windows, the preconditioner weights are frozen and the preconditioner acts like a constant preconditioner.");
535
536 } else {
537 PRECICE_ERROR("Acceleration of type \"{}\" is unknown. Please choose a valid acceleration scheme or check the spelling in the configuration file.", tag.getName());
538 }
539}
540
549
550} // namespace precice::acceleration
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:15
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:35
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
T at(T... args)
virtual void xmlEndTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag)
Callback method required when using xml::XMLTag.
std::set< std::pair< std::string, std::string > > _uniqueDataAndMeshNames
struct precice::acceleration::AccelerationConfiguration::UserDefinitions _userDefinitions
struct precice::acceleration::AccelerationConfiguration::ConfigurationData _config
const struct precice::acceleration::AccelerationConfiguration::DefaultValuesIMVJ _defaultValuesIQNIMVJ
PtrAcceleration getAcceleration()
Returns the configured coupling scheme.
AccelerationConfiguration(const mesh::PtrMeshConfiguration &meshConfig)
const struct precice::acceleration::AccelerationConfiguration::DefaultValuesIQN _defaultValuesIQNILS
virtual void xmlTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag)
Callback method required when using xml::XMLTag.
Interface quasi-Newton with interface least-squares approximation.
Multi vector quasi-Newton update scheme.
Preconditioner that uses the constant user-defined factors to scale the quasi-Newton system.
Preconditioner that uses the recent residual to scale the quasi-Newton system.
Preconditioner that uses the residuals of all iterations of the current time window summed up to scal...
Preconditioner that uses the values from the previous time window to scale the quasi-Newton system.
XMLAttribute & setOptions(std::vector< ATTRIBUTE_T > options)
XMLAttribute & setDocumentation(std::string documentation)
Sets a documentation string for the attribute.
Represents an XML tag to be configured automatically.
Definition XMLTag.hpp:31
const std::string & getNamespace() const
Returns xml namespace.
Definition XMLTag.hpp:159
bool hasAttribute(const std::string &attributeName)
Definition XMLTag.cpp:111
std::string getStringAttributeValue(const std::string &name, std::optional< std::string > default_value=std::nullopt) const
Definition XMLTag.cpp:142
bool getBooleanAttributeValue(const std::string &name, std::optional< bool > default_value=std::nullopt) const
Definition XMLTag.cpp:155
XMLTag & setDocumentation(const std::string &documentation)
Adds a description of the purpose of this XML tag.
Definition XMLTag.cpp:29
const std::string & getName() const
Returns name (without namespace).
Definition XMLTag.hpp:153
const std::string & getFullName() const
Returns full name consisting of xml namespace + ":" + name.
Definition XMLTag.hpp:170
int getIntAttributeValue(const std::string &name, std::optional< int > default_value=std::nullopt) const
Definition XMLTag.cpp:129
double getDoubleAttributeValue(const std::string &name, std::optional< double > default_value=std::nullopt) const
Definition XMLTag.cpp:116
XMLTag & addAttribute(const XMLAttribute< double > &attribute)
Removes the XML subtag with given name.
Definition XMLTag.cpp:53
XMLTag & addSubtag(const XMLTag &tag)
Adds an XML tag as subtag by making a copy of the given tag.
Definition XMLTag.cpp:41
T clear(T... args)
T emplace(T... args)
T get(T... args)
T insert(T... args)
T make_pair(T... args)
std::shared_ptr< Preconditioner > PtrPreconditioner
contains implementations of acceleration schemes.
std::shared_ptr< Acceleration > PtrAcceleration
contains the XML configuration parser.
T push_back(T... args)
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:24