preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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_REDUCEDTIMEGRIDQN("reduced-time-grid"),
52 ATTR_IMVJCHUNKSIZE("chunk-size"),
53 ATTR_RSLS_REUSED_TIME_WINDOWS("reused-time-windows-at-restart"),
54 ATTR_RSSVD_TRUNCATIONEPS("truncation-threshold"),
55 ATTR_PRECOND_NONCONST_TIME_WINDOWS("freeze-after"),
56 ATTR_PRECOND_UPDATE_ON_THRESHOLD("update-on-threshold"),
57 VALUE_CONSTANT("constant"),
58 VALUE_AITKEN("aitken"),
59 VALUE_IQNILS("IQN-ILS"),
60 VALUE_IQNIMVJ("IQN-IMVJ"),
61 VALUE_QR1FILTER("QR1"),
62 VALUE_QR1_ABSFILTER("QR1-absolute"),
63 VALUE_QR2FILTER("QR2"),
64 VALUE_QR3FILTER("QR3"),
65 VALUE_CONSTANT_PRECONDITIONER("constant"),
66 VALUE_VALUE_PRECONDITIONER("value"),
67 VALUE_RESIDUAL_PRECONDITIONER("residual"),
68 VALUE_RESIDUAL_SUM_PRECONDITIONER("residual-sum"),
69 VALUE_LS_RESTART("RS-LS"),
70 VALUE_ZERO_RESTART("RS-0"),
71 VALUE_SVD_RESTART("RS-SVD"),
72 VALUE_SLIDE_RESTART("RS-SLIDE"),
73 VALUE_NO_RESTART("no-restart"),
74 _meshConfig(meshConfig),
75 _acceleration(),
76 _neededMeshes(),
77 _preconditioner(),
78 _config()
79{
80 PRECICE_ASSERT(meshConfig.get() != nullptr);
81}
82
84{
85 using namespace xml;
86
87 // static int recursionCounter = 0;
88 // recursionCounter++;
89
90 XMLTag::Occurrence occ = XMLTag::OCCUR_NOT_OR_ONCE;
92 {
93 XMLTag tag(*this, VALUE_CONSTANT, occ, TAG);
94 tag.setDocumentation("Accelerates coupling data with constant underrelaxation.");
96 tags.push_back(tag);
97 }
98 {
99 XMLTag tag(*this, VALUE_AITKEN, occ, TAG);
100 tag.setDocumentation("Accelerates coupling data with dynamic Aitken under-relaxation.");
102 tags.push_back(tag);
103 }
104 {
105 XMLTag tag(*this, VALUE_IQNILS, occ, TAG);
106 tag.setDocumentation("Accelerates coupling data with the interface quasi-Newton inverse least-squares method.");
107
108 auto reducedTimeGridQN = makeXMLAttribute(ATTR_REDUCEDTIMEGRIDQN, true)
109 .setDocumentation("Whether only the last time step of each time window is used to construct the Jacobian.");
110 tag.addAttribute(reducedTimeGridQN);
111
113 tags.push_back(tag);
114 }
115 {
116 XMLTag tag(*this, VALUE_IQNIMVJ, occ, TAG);
117 tag.setDocumentation("Accelerates coupling data with the interface quasi-Newton inverse multi-vector Jacobian method.");
118
119 auto alwaybuildJacobian = makeXMLAttribute(ATTR_BUILDJACOBIAN, false)
120 .setDocumentation("If set to true, the IMVJ will set up the Jacobian matrix"
121 " in each coupling iteration, which is inefficient. If set to false (or not set)"
122 " the Jacobian is only build in the last iteration and the updates are computed using (relatively) cheap MATVEC products.");
123 tag.addAttribute(alwaybuildJacobian);
124
125 auto reducedTimeGridQN = makeXMLAttribute(ATTR_REDUCEDTIMEGRIDQN, true)
126 .setDocumentation("Whether only the last time step of each time window is used to construct the Jacobian.");
127 tag.addAttribute(reducedTimeGridQN);
128
130 tags.push_back(tag);
131 }
132
133 for (XMLTag &tag : tags) {
134 parent.addSubtag(tag);
135 }
136}
137
142
144 const xml::ConfigurationContext &context,
145 xml::XMLTag &callingTag)
146{
147 PRECICE_TRACE(callingTag.getFullName());
148
149 if (callingTag.getNamespace() == TAG) {
150 _config.type = callingTag.getName();
151
154
157 }
158 if (callingTag.getName() == TAG_RELAX) {
160 } else if (callingTag.getName() == TAG_DATA) {
161 std::string dataName = callingTag.getStringAttributeValue(ATTR_NAME);
162 std::string meshName = callingTag.getStringAttributeValue(ATTR_MESH);
163 auto success = _uniqueDataAndMeshNames.emplace(dataName, meshName);
164 PRECICE_CHECK(success.second,
165 "You have provided a subtag <data name=\"{}\" mesh=\"{}\"/> more than once in your <acceleration:.../>. "
166 "Please remove the duplicated entry.",
167 dataName, meshName);
168
170 double scaling = 1.0;
172 scaling = callingTag.getDoubleAttributeValue(ATTR_SCALING);
173 }
174
175 PRECICE_CHECK(_meshConfig->hasMeshName(_meshName) && _meshConfig->getMesh(_meshName)->hasDataName(dataName),
176 "Data with name \"{0}\" associated to mesh \"{1}\" not found on configuration of acceleration. "
177 "Add \"{0}\" to the \"<mesh name={1}>\" tag, or change the data name in the acceleration scheme.",
178 dataName, _meshName);
179
180 const mesh::PtrMesh &mesh = _meshConfig->getMesh(_meshName);
181 const mesh::PtrData &data = mesh->data(dataName);
182 _config.dataIDs.push_back(data->getID());
183 _config.scalings.insert(std::make_pair(data->getID(), scaling));
184
186 } else if (callingTag.getName() == TAG_INIT_RELAX) {
189 if (callingTag.hasAttribute(ATTR_ENFORCE)) {
191 } else {
193 }
194 } else if (callingTag.getName() == TAG_MAX_USED_ITERATIONS) {
197 } else if (callingTag.getName() == TAG_TIME_WINDOWS_REUSED) {
200 } else if (callingTag.getName() == TAG_FILTER) {
202 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
203 if (f == VALUE_QR1FILTER) {
205 } else if (f == VALUE_QR1_ABSFILTER) {
207 } else if (f == VALUE_QR2FILTER) {
209 } else if (f == VALUE_QR3FILTER) {
211 } else {
212 PRECICE_ASSERT(false);
213 }
215 } else if (callingTag.getName() == TAG_PRECONDITIONER) {
220 } else if (callingTag.getName() == TAG_IMVJRESTART) {
222#ifndef PRECICE_NO_MPI
224 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
225 PRECICE_CHECK((f == VALUE_NO_RESTART) || (!_config.alwaysBuildJacobian), "IMVJ cannot be in restart mode while parameter always-build-jacobian is set to true. "
226 "Please remove 'always-build-jacobian' from the configuration file or do not run in restart mode.");
227 if (f == VALUE_NO_RESTART) {
229 } else if (f == VALUE_ZERO_RESTART) {
231 } else if (f == VALUE_LS_RESTART) {
234 } else if (f == VALUE_SVD_RESTART) {
237 } else if (f == VALUE_SLIDE_RESTART) {
239 } else {
241 PRECICE_ASSERT(false);
242 }
243#else
244 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
245#endif
246 }
247}
248
250 const xml::ConfigurationContext &context,
251 xml::XMLTag &callingTag)
252{
253 PRECICE_TRACE(callingTag.getName());
254 if (callingTag.getNamespace() == TAG) {
255
256 // create preconditioner
257 if (callingTag.getName() == VALUE_IQNILS || callingTag.getName() == VALUE_AITKEN) {
266 } else {
267 // no preconditioner defined
269 }
270 }
271
272 PRECICE_CHECK(!_acceleration, "You are trying to define multiple acceleration schemes, which is not allowed. Please remove one of them.");
273 if (callingTag.getName() == VALUE_CONSTANT) {
277 } else if (callingTag.getName() == VALUE_AITKEN) {
282 } else if (callingTag.getName() == VALUE_IQNILS) {
298 } else if (callingTag.getName() == VALUE_IQNIMVJ) {
299#ifndef PRECICE_NO_MPI
308 }
309
310 // create preconditioner
311 // if imvj restart-mode is of type RS-SVD, max number of non-const preconditioned time windows is limited by the chunksize
312 // it is separated from the other acceleration methods, since SVD-restart might be chosen as default here
324 } else {
325 // no preconditioner defined
327 }
328
344#else
345 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
346#endif
347 } else {
348 PRECICE_ASSERT(false);
349 }
350 }
351}
352
360
362{
363 using namespace precice::xml;
364
365 XMLTag tagData(*this, TAG_DATA, XMLTag::OCCUR_ONCE_OR_MORE);
366 tagData.setDocumentation("The data used to compute the acceleration.");
368 attrName.setDocumentation("The name of the data.");
370 attrMesh.setDocumentation("The name of the mesh which holds the data.");
371 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
372 .setDocumentation(
373 "To improve the performance of a parallel or a multi coupling schemes, "
374 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
375 "We recommend, however, to use an automatic scaling via a preconditioner.");
376 tagData.addAttribute(attrScaling);
377 tagData.addAttribute(attrName);
378 tagData.addAttribute(attrMesh);
379 tag.addSubtag(tagData);
380
381 XMLTag tagFilter(*this, TAG_FILTER, XMLTag::OCCUR_NOT_OR_ONCE);
382 tagFilter.setDocumentation("Type of filtering technique that is used to "
383 "maintain good conditioning in the least-squares system. Possible filters:\n"
384 " - `QR1-filter`: updateQR-dec with (relative) test \\\\(R(i,i) < \\epsilon *\\lVert R\\rVert_F\\\\)\n"
385 " - `QR1_absolute-filter`: updateQR-dec with (absolute) test \\\\(R(i, i) < \\epsilon\\\\)\n"
386 " - `QR2-filter`: en-block QR-dec with test \\\\(\\lVert v_\\text{orth} \\rVert_2 < \\epsilon * \\lVert v \\rVert_2\\\\)\n\n"
387 "Please note that a QR1 is based on Given's rotations whereas QR2 uses "
388 "modified Gram-Schmidt. This can give different results even when no columns "
389 "are filtered out.\n"
390 "When this tag is not provided, the QR2-filter with the limit value 1e-2 is used.");
391 XMLAttribute<double> attrSingularityLimit(ATTR_SINGULARITYLIMIT, 1e-16);
392 attrSingularityLimit.setDocumentation("Limit eps of the filter.");
393 tagFilter.addAttribute(attrSingularityLimit);
394 auto attrFilterName = XMLAttribute<std::string>(ATTR_TYPE)
399 .setDocumentation("Type of the filter.");
400 tagFilter.addAttribute(attrFilterName);
401 tag.addSubtag(tagFilter);
402}
403
405 xml::XMLTag &tag)
406{
407 using namespace xml;
408 if (tag.getName() == VALUE_CONSTANT) {
409 XMLTag tagRelax(*this, TAG_RELAX, XMLTag::OCCUR_ONCE);
411 attrValue.setDocumentation("Constant relaxation factor.");
412 tagRelax.addAttribute(attrValue);
413 tag.addSubtag(tagRelax);
414 } else if (tag.getName() == VALUE_AITKEN) {
415 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
416 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.5 is used.");
418 attrValue.setDocumentation("Initial relaxation factor.");
419 tagInitRelax.addAttribute(attrValue);
420 tag.addSubtag(tagInitRelax);
421
422 XMLTag tagData(*this, TAG_DATA, XMLTag::OCCUR_ONCE_OR_MORE);
423 tagData.setDocumentation("The data used to compute the acceleration.");
425 attrName.setDocumentation("The name of the data.");
427 attrMesh.setDocumentation("The name of the mesh which holds the data.");
428 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
429 .setDocumentation(
430 "To improve the performance of a parallel or a multi coupling schemes, "
431 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
432 "We recommend, however, to use an automatic scaling via a preconditioner.");
433 tagData.addAttribute(attrScaling);
434 tagData.addAttribute(attrName);
435 tagData.addAttribute(attrMesh);
436 tag.addSubtag(tagData);
437
438 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
439 tagPreconditioner.setDocumentation("To improve the numerical stability of multiple data vectors a preconditioner"
440 " can be applied. A constant preconditioner scales every acceleration data by a constant value, which you can define as"
441 " an attribute of data. "
442 " A value preconditioner scales every acceleration data by the norm of the data in the previous time window."
443 " A residual preconditioner scales every acceleration data by the current residual."
444 " A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.");
445 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
450 .setDocumentation("The type of the preconditioner.");
451 tagPreconditioner.addAttribute(attrPreconditionerType);
452 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
453 .setDocumentation(
454 "After the given number of time windows, the preconditioner weights "
455 "are frozen and the preconditioner acts like a constant preconditioner.");
456 tagPreconditioner.addAttribute(nonconstTWindows);
457 tag.addSubtag(tagPreconditioner);
458 } else if (tag.getName() == VALUE_IQNILS) {
459
460 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
461 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
462 XMLAttribute<double> attrDoubleValue(ATTR_VALUE);
463 attrDoubleValue.setDocumentation("Initial relaxation factor.");
464 tagInitRelax.addAttribute(attrDoubleValue);
465 XMLAttribute<bool> attrEnforce(ATTR_ENFORCE, false);
466 attrEnforce.setDocumentation("Enforce initial relaxation in every time window.");
467 tagInitRelax.addAttribute(attrEnforce);
468 tag.addSubtag(tagInitRelax);
469
470 XMLTag tagMaxUsedIter(*this, TAG_MAX_USED_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
471 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.");
472 XMLAttribute<int> attrIntValue(ATTR_VALUE);
473 attrIntValue.setDocumentation("The number of columns.");
474 tagMaxUsedIter.addAttribute(attrIntValue);
475 tag.addSubtag(tagMaxUsedIter);
476
477 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
478 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.");
479 XMLAttribute<int> attrNumTimeWindowsReused(ATTR_VALUE);
480 attrNumTimeWindowsReused.setDocumentation("The number of time windows.");
481 tagTimeWindowsReused.addAttribute(attrNumTimeWindowsReused);
482 tag.addSubtag(tagTimeWindowsReused);
483
485
486 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
487 tagPreconditioner.setDocumentation("To improve the performance of a parallel or a multi coupling schemes a preconditioner"
488 " can be applied. "
489 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data. \n "
490 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
491 "- A residual preconditioner scales every acceleration data by the current residual.\n"
492 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n"
493 " If this tag is not provided, the residual-sum preconditioner is employed.");
494 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
499 .setDocumentation("The type of the preconditioner.");
500 tagPreconditioner.addAttribute(attrPreconditionerType);
501 auto attrpreconditionerUpdateOnThreshold = XMLAttribute<bool>(ATTR_PRECOND_UPDATE_ON_THRESHOLD, true)
502 .setDocumentation("To update the preconditioner weights after the first time window:"
503 "- `true`: The preconditioner weights are only updated if the weights will change by more than one order of magnitude.\n"
504 "- `false`: The preconditioner weights are updated after every iteration.");
505 tagPreconditioner.addAttribute(attrpreconditionerUpdateOnThreshold);
506 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
507 .setDocumentation(
508 "After the given number of time windows, the preconditioner weights "
509 "are frozen and the preconditioner acts like a constant preconditioner.");
510 tagPreconditioner.addAttribute(nonconstTWindows);
511 tag.addSubtag(tagPreconditioner);
512
513 } else if (tag.getName() == VALUE_IQNIMVJ) {
514 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
515 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
516 tagInitRelax.addAttribute(
517 XMLAttribute<double>(ATTR_VALUE).setDocumentation("Initial relaxation factor."));
518 tagInitRelax.addAttribute(
519 XMLAttribute<bool>(ATTR_ENFORCE, false).setDocumentation("Enforce initial relaxation in every time window."));
520 tag.addSubtag(tagInitRelax);
521
522 XMLTag tagIMVJRESTART(*this, TAG_IMVJRESTART, XMLTag::OCCUR_NOT_OR_ONCE);
523 auto attrRestartName = XMLAttribute<std::string>(ATTR_TYPE)
529 .setDefaultValue(VALUE_SVD_RESTART)
530 .setDocumentation("Type of the restart mode.");
531 tagIMVJRESTART.addAttribute(attrRestartName);
532 tagIMVJRESTART.setDocumentation("Type of IMVJ restart mode that is used:\n"
533 "- `no-restart`: IMVJ runs in normal mode with explicit representation of Jacobian\n"
534 "- `RS-0`: IMVJ runs in restart mode. After M time windows all Jacobain information is dropped, restart with no information\n"
535 "- `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"
536 "- `RS-SVD`: IMVJ runs in restart mode. After M time windows a truncated SVD of the Jacobian is updated.\n"
537 "- `RS-SLIDE`: IMVJ runs in sliding window restart mode.\n"
538 "If this tag is not provided, IMVJ runs in restart mode with SVD-method.");
539 auto attrChunkSize = makeXMLAttribute(ATTR_IMVJCHUNKSIZE, 8)
540 .setDocumentation("Specifies the number of time windows M after which the IMVJ restarts, if run in restart-mode. Default value is M=8.");
541 auto attrReusedTimeWindowsAtRestart = makeXMLAttribute(ATTR_RSLS_REUSED_TIME_WINDOWS, 8)
542 .setDocumentation("If IMVJ restart-mode=RS-LS, the number of reused time windows at restart can be specified.");
543 auto attrRSSVD_truncationEps = makeXMLAttribute(ATTR_RSSVD_TRUNCATIONEPS, 1e-4)
544 .setDocumentation("If IMVJ restart-mode=RS-SVD, the truncation threshold for the updated SVD can be set.");
545 tagIMVJRESTART.addAttribute(attrChunkSize);
546 tagIMVJRESTART.addAttribute(attrReusedTimeWindowsAtRestart);
547 tagIMVJRESTART.addAttribute(attrRSSVD_truncationEps);
548 tag.addSubtag(tagIMVJRESTART);
549
550 XMLTag tagMaxUsedIter(*this, TAG_MAX_USED_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
551 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.");
552 XMLAttribute<int> attrIntValue(ATTR_VALUE);
553 attrIntValue.setDocumentation("The number of columns.");
554 tagMaxUsedIter.addAttribute(attrIntValue);
555 tag.addSubtag(tagMaxUsedIter);
556
557 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
558 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.");
559 tagTimeWindowsReused.addAttribute(attrIntValue);
560 tag.addSubtag(tagTimeWindowsReused);
561
563
564 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
565 tagPreconditioner.setDocumentation(
566 "To improve the performance of a parallel or a multi coupling schemes a preconditioner can be applied."
567 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data.\n"
568 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
569 "- A residual preconditioner scales every acceleration data by the current residual.\n"
570 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n"
571 " If this tag is not provided, the residual-sum preconditioner is employed.");
572 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
577 .setDocumentation("Type of the preconditioner.");
578 tagPreconditioner.addAttribute(attrPreconditionerType);
579 auto attrpreconditionerUpdateOnThreshold = XMLAttribute<bool>(ATTR_PRECOND_UPDATE_ON_THRESHOLD, true)
580 .setDocumentation("To update the preconditioner weights after the first time window:"
581 "- `true`: The preconditioner weights are only updated if the weights will change by more than one order of magnitude.\n"
582 "- `false`: The preconditioner weights are updated after every iteration.");
583 tagPreconditioner.addAttribute(attrpreconditionerUpdateOnThreshold);
584 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
585 .setDocumentation("After the given number of time windows, the preconditioner weights are frozen and the preconditioner acts like a constant preconditioner.");
586 tagPreconditioner.addAttribute(nonconstTWindows);
587 tag.addSubtag(tagPreconditioner);
588
589 } else {
590 PRECICE_ERROR("Acceleration of type \"{}\" is unknown. Please choose a valid acceleration scheme or check the spelling in the configuration file.", tag.getName());
591 }
592}
593
595{
596 std::vector<double> factors;
597 for (int id : dataIDs) {
598 factors.push_back(scalings.at(id));
599 }
600 return factors;
601}
602
603} // namespace precice::acceleration
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:16
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T at(T... args)
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.
void xmlTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
Callback method required when using xml::XMLTag.
AccelerationConfiguration(const mesh::PtrMeshConfiguration &meshConfig)
const struct precice::acceleration::AccelerationConfiguration::DefaultValuesIQN _defaultValuesIQNILS
void xmlEndTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
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_view documentation)
Sets a documentation string for the attribute.
Represents an XML tag to be configured automatically.
Definition XMLTag.hpp:28
const std::string & getNamespace() const
Returns xml namespace.
Definition XMLTag.hpp:159
bool hasAttribute(const std::string &attributeName) const
Definition XMLTag.cpp:112
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
XMLTag & setDocumentation(std::string_view 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:131
double getDoubleAttributeValue(const std::string &name, std::optional< double > default_value=std::nullopt) const
Definition XMLTag.cpp:117
XMLTag & addAttribute(const XMLAttribute< double > &attribute)
Adds a XML attribute by making a copy of the given attribute.
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:21