preCICE v3.3.0
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_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"),
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"),
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),
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 in each coupling iteration, which is inefficient. "
121 "If set to false (or not set) the Jacobian is only build in the last iteration and the updates are computed using (relatively) cheap MATVEC products.");
122 tag.addAttribute(alwaybuildJacobian);
123
124 auto reducedTimeGridQN = makeXMLAttribute(ATTR_REDUCEDTIMEGRIDQN, true)
125 .setDocumentation("Whether only the last time step of each time window is used to construct the Jacobian.");
126 tag.addAttribute(reducedTimeGridQN);
127
129 tags.push_back(tag);
130 }
131
132 for (XMLTag &tag : tags) {
133 parent.addSubtag(tag);
134 }
135}
136
141
143 const xml::ConfigurationContext &context,
144 xml::XMLTag &callingTag)
145{
146 PRECICE_TRACE(callingTag.getFullName());
147
148 if (callingTag.getNamespace() == TAG) {
149 _config.type = callingTag.getName();
150
151 if (_config.type == VALUE_IQNIMVJ)
152 _config.alwaysBuildJacobian = callingTag.getBooleanAttributeValue(ATTR_BUILDJACOBIAN);
153
154 if (_config.type == VALUE_IQNIMVJ || _config.type == VALUE_IQNILS)
155 _config.reducedTimeGridQN = callingTag.getBooleanAttributeValue(ATTR_REDUCEDTIMEGRIDQN);
156 }
157 if (callingTag.getName() == TAG_RELAX) {
158 _config.relaxationFactor = callingTag.getDoubleAttributeValue(ATTR_VALUE);
159 } else if (callingTag.getName() == TAG_DATA) {
160 std::string dataName = callingTag.getStringAttributeValue(ATTR_NAME);
161 std::string meshName = callingTag.getStringAttributeValue(ATTR_MESH);
162 auto success = _uniqueDataAndMeshNames.emplace(dataName, meshName);
163 PRECICE_CHECK(success.second,
164 "You have provided a subtag <data name=\"{}\" mesh=\"{}\"/> more than once in your <acceleration:.../>. "
165 "Please remove the duplicated entry.",
166 dataName, meshName);
167
169 double scaling = 1.0;
170 if (_config.type == VALUE_IQNILS || _config.type == VALUE_IQNIMVJ) {
171 scaling = callingTag.getDoubleAttributeValue(ATTR_SCALING);
172 }
173
174 PRECICE_CHECK(_meshConfig->hasMeshName(_meshName) && _meshConfig->getMesh(_meshName)->hasDataName(dataName),
175 "Data with name \"{0}\" associated to mesh \"{1}\" not found on configuration of acceleration. "
176 "Add \"{0}\" to the \"<mesh name={1}>\" tag, or change the data name in the acceleration scheme.",
177 dataName, _meshName);
178
179 const mesh::PtrMesh &mesh = _meshConfig->getMesh(_meshName);
180 const mesh::PtrData &data = mesh->data(dataName);
181 _config.dataIDs.push_back(data->getID());
182 _config.scalings.insert(std::make_pair(data->getID(), scaling));
183
184 _neededMeshes.push_back(_meshName);
185 } else if (callingTag.getName() == TAG_INIT_RELAX) {
186 _userDefinitions.definedRelaxationFactor = true;
187 _config.relaxationFactor = callingTag.getDoubleAttributeValue(ATTR_VALUE);
188 if (callingTag.hasAttribute(ATTR_ENFORCE)) {
189 _config.forceInitialRelaxation = callingTag.getBooleanAttributeValue(ATTR_ENFORCE);
190 } else {
191 _config.forceInitialRelaxation = false;
192 }
193 } else if (callingTag.getName() == TAG_MAX_USED_ITERATIONS) {
194 _userDefinitions.definedMaxIterationsUsed = true;
195 _config.maxIterationsUsed = callingTag.getIntAttributeValue(ATTR_VALUE);
196 } else if (callingTag.getName() == TAG_TIME_WINDOWS_REUSED) {
197 _userDefinitions.definedTimeWindowsReused = true;
198 _config.timeWindowsReused = callingTag.getIntAttributeValue(ATTR_VALUE);
199 } else if (callingTag.getName() == TAG_FILTER) {
200 _userDefinitions.definedFilter = true;
201 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
202 if (f == VALUE_QR1FILTER) {
204 } else if (f == VALUE_QR1_ABSFILTER) {
206 } else if (f == VALUE_QR2FILTER) {
208 } else if (f == VALUE_QR3FILTER) {
210 } else {
211 PRECICE_ASSERT(false);
212 }
213 _config.singularityLimit = callingTag.getDoubleAttributeValue(ATTR_SINGULARITYLIMIT);
214 } else if (callingTag.getName() == TAG_PRECONDITIONER) {
215 _userDefinitions.definedPreconditionerType = true;
216 _config.preconditionerType = callingTag.getStringAttributeValue(ATTR_TYPE);
217 _config.preconditionerUpdateOnThreshold = callingTag.getBooleanAttributeValue(ATTR_PRECOND_UPDATE_ON_THRESHOLD);
218 _config.preconditionerNbNonConstTWindows = callingTag.getIntAttributeValue(ATTR_PRECOND_NONCONST_TIME_WINDOWS);
219 } else if (callingTag.getName() == TAG_IMVJRESTART) {
220 _userDefinitions.defineRestartType = true;
221#ifndef PRECICE_NO_MPI
222 _config.imvjChunkSize = callingTag.getIntAttributeValue(ATTR_IMVJCHUNKSIZE);
223 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
224 PRECICE_CHECK((f == VALUE_NO_RESTART) || (!_config.alwaysBuildJacobian), "IMVJ cannot be in restart mode while parameter always-build-jacobian is set to true. "
225 "Please remove 'always-build-jacobian' from the configuration file or do not run in restart mode.");
226 if (f == VALUE_NO_RESTART) {
228 } else if (f == VALUE_ZERO_RESTART) {
229 _config.imvjRestartType = IQNIMVJAcceleration::RS_ZERO;
230 } else if (f == VALUE_LS_RESTART) {
231 _config.imvjRSLS_reusedTimeWindows = callingTag.getIntAttributeValue(ATTR_RSLS_REUSED_TIME_WINDOWS);
232 _config.imvjRestartType = IQNIMVJAcceleration::RS_LS;
233 } else if (f == VALUE_SVD_RESTART) {
234 _config.imvjRSSVD_truncationEps = callingTag.getDoubleAttributeValue(ATTR_RSSVD_TRUNCATIONEPS);
235 _config.imvjRestartType = IQNIMVJAcceleration::RS_SVD;
236 } else if (f == VALUE_SLIDE_RESTART) {
237 _config.imvjRestartType = IQNIMVJAcceleration::RS_SLIDE;
238 } else {
239 _config.imvjChunkSize = 0;
240 PRECICE_ASSERT(false);
241 }
242#else
243 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
244#endif
245 }
246}
247
249 const xml::ConfigurationContext &context,
250 xml::XMLTag &callingTag)
251{
252 PRECICE_TRACE(callingTag.getName());
253 if (callingTag.getNamespace() == TAG) {
254
255 // create preconditioner
256 if (callingTag.getName() == VALUE_IQNILS || callingTag.getName() == VALUE_AITKEN) {
257 if (_config.preconditionerType == VALUE_CONSTANT_PRECONDITIONER) {
258 _preconditioner = PtrPreconditioner(new ConstantPreconditioner(_config.scalingFactorsInOrder()));
259 } else if (_config.preconditionerType == VALUE_VALUE_PRECONDITIONER) {
260 _preconditioner = PtrPreconditioner(new ValuePreconditioner(_config.preconditionerNbNonConstTWindows));
261 } else if (_config.preconditionerType == VALUE_RESIDUAL_PRECONDITIONER) {
262 _preconditioner = PtrPreconditioner(new ResidualPreconditioner(_config.preconditionerNbNonConstTWindows));
263 } else if (_config.preconditionerType == VALUE_RESIDUAL_SUM_PRECONDITIONER) {
264 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_config.preconditionerNbNonConstTWindows, _config.preconditionerUpdateOnThreshold));
265 } else {
266 // no preconditioner defined
267 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_defaultValuesIQNILS.preconditionerNbNonConstTWindows, _defaultValuesIQNILS.preconditionerUpdateOnThreshold));
268 }
269 }
270
271 PRECICE_CHECK(!_acceleration, "You are trying to define multiple acceleration schemes, which is not allowed. Please remove one of them.");
272 if (callingTag.getName() == VALUE_CONSTANT) {
275 _config.relaxationFactor, _config.dataIDs));
276 } else if (callingTag.getName() == VALUE_AITKEN) {
277 _config.relaxationFactor = (_userDefinitions.definedRelaxationFactor) ? _config.relaxationFactor : _defaultAitkenRelaxationFactor;
280 _config.relaxationFactor, _config.dataIDs, _preconditioner));
281 } else if (callingTag.getName() == VALUE_IQNILS) {
282 _config.relaxationFactor = (_userDefinitions.definedRelaxationFactor) ? _config.relaxationFactor : _defaultValuesIQNILS.relaxationFactor;
283 _config.maxIterationsUsed = (_userDefinitions.definedMaxIterationsUsed) ? _config.maxIterationsUsed : _defaultValuesIQNILS.maxIterationsUsed;
284 _config.timeWindowsReused = (_userDefinitions.definedTimeWindowsReused) ? _config.timeWindowsReused : _defaultValuesIQNILS.timeWindowsReused;
285 _config.filter = (_userDefinitions.definedFilter) ? _config.filter : _defaultValuesIQNILS.filter;
286 _config.singularityLimit = (_userDefinitions.definedFilter) ? _config.singularityLimit : _defaultValuesIQNILS.singularityLimit;
289 _config.relaxationFactor,
290 _config.forceInitialRelaxation,
291 _config.maxIterationsUsed,
292 _config.timeWindowsReused,
293 _config.filter, _config.singularityLimit,
294 _config.dataIDs,
296 _config.reducedTimeGridQN));
297 } else if (callingTag.getName() == VALUE_IQNIMVJ) {
298#ifndef PRECICE_NO_MPI
299 _config.relaxationFactor = (_userDefinitions.definedRelaxationFactor) ? _config.relaxationFactor : _defaultValuesIQNIMVJ.relaxationFactor;
300 _config.maxIterationsUsed = (_userDefinitions.definedMaxIterationsUsed) ? _config.maxIterationsUsed : _defaultValuesIQNIMVJ.maxIterationsUsed;
301 _config.timeWindowsReused = (_userDefinitions.definedTimeWindowsReused) ? _config.timeWindowsReused : _defaultValuesIQNIMVJ.timeWindowsReused;
302 _config.filter = (_userDefinitions.definedFilter) ? _config.filter : _defaultValuesIQNILS.filter;
303 _config.singularityLimit = (_userDefinitions.definedFilter) ? _config.singularityLimit : _defaultValuesIQNILS.singularityLimit;
304 if (!_config.alwaysBuildJacobian) {
305 _config.imvjRestartType = (_userDefinitions.defineRestartType) ? _config.imvjRestartType : _defaultValuesIQNIMVJ.imvjRestartType;
306 _config.imvjChunkSize = (_userDefinitions.defineRestartType) ? _config.imvjChunkSize : _defaultValuesIQNIMVJ.imvjChunkSize;
307 }
308
309 // create preconditioner
310 // if imvj restart-mode is of type RS-SVD, max number of non-const preconditioned time windows is limited by the chunksize
311 // it is separated from the other acceleration methods, since SVD-restart might be chosen as default here
312 if (_config.imvjRestartType == IQNIMVJAcceleration::RS_SVD)
313 if (_config.preconditionerNbNonConstTWindows > _config.imvjChunkSize)
314 _config.preconditionerNbNonConstTWindows = _config.imvjChunkSize;
315 if (_config.preconditionerType == VALUE_CONSTANT_PRECONDITIONER) {
316 _preconditioner = PtrPreconditioner(new ConstantPreconditioner(_config.scalingFactorsInOrder()));
317 } else if (_config.preconditionerType == VALUE_VALUE_PRECONDITIONER) {
318 _preconditioner = PtrPreconditioner(new ValuePreconditioner(_config.preconditionerNbNonConstTWindows));
319 } else if (_config.preconditionerType == VALUE_RESIDUAL_PRECONDITIONER) {
320 _preconditioner = PtrPreconditioner(new ResidualPreconditioner(_config.preconditionerNbNonConstTWindows));
321 } else if (_config.preconditionerType == VALUE_RESIDUAL_SUM_PRECONDITIONER) {
322 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_config.preconditionerNbNonConstTWindows, _config.preconditionerUpdateOnThreshold));
323 } else {
324 // no preconditioner defined
325 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_defaultValuesIQNILS.preconditionerNbNonConstTWindows, _defaultValuesIQNIMVJ.preconditionerUpdateOnThreshold));
326 }
327
330 _config.relaxationFactor,
331 _config.forceInitialRelaxation,
332 _config.maxIterationsUsed,
333 _config.timeWindowsReused,
334 _config.filter, _config.singularityLimit,
335 _config.dataIDs,
337 _config.alwaysBuildJacobian,
338 _config.imvjRestartType,
339 _config.imvjChunkSize,
340 _config.imvjRSLS_reusedTimeWindows,
341 _config.imvjRSSVD_truncationEps,
342 _config.reducedTimeGridQN));
343#else
344 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
345#endif
346 } else {
347 PRECICE_ASSERT(false);
348 }
349 }
350}
351
359
361{
362 using namespace precice::xml;
363
365 tagData.setDocumentation("The data used to compute the acceleration.");
367 attrName.setDocumentation("The name of the data.");
369 attrMesh.setDocumentation("The name of the mesh which holds the data.");
370 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
371 .setDocumentation(
372 "To improve the performance of a parallel or a multi coupling schemes, "
373 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
374 "We recommend, however, to use an automatic scaling via a preconditioner.");
375 tagData.addAttribute(attrScaling);
376 tagData.addAttribute(attrName);
377 tagData.addAttribute(attrMesh);
378 tag.addSubtag(tagData);
379
381 tagFilter.setDocumentation("Type of filtering technique that is used to "
382 "maintain good conditioning in the least-squares system. Possible filters:\n\n"
383 "- `QR1`: update QR-dec with (relative) test \\\\(R(i,i) < \\epsilon *\\lVert R\\rVert_F\\\\)\n"
384 "- `QR1-absolute`: update QR-dec with (absolute) test \\\\(R(i, i) < \\epsilon\\\\)\n"
385 "- `QR2`: en-block QR-dec with test \\\\(\\lVert v_\\text{orth} \\rVert_2 < \\epsilon * \\lVert v \\rVert_2\\\\)\n\n"
386 "- `QR3`: update QR-dec only when the pre-scaling weights have changed or there is one or more columns are to be removed 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 modified Gram-Schmidt. "
388 "This can give different results even when no columns are filtered out.");
389 XMLAttribute<double> attrSingularityLimit(ATTR_SINGULARITYLIMIT, 1e-16);
390 attrSingularityLimit.setDocumentation("Limit eps of the filter.");
391 tagFilter.addAttribute(attrSingularityLimit);
392 auto attrFilterName = XMLAttribute<std::string>(ATTR_TYPE)
397 .setDefaultValue(VALUE_QR3FILTER)
398 .setDocumentation("Type of the filter.");
399 tagFilter.addAttribute(attrFilterName);
400 tag.addSubtag(tagFilter);
401}
402
404 xml::XMLTag &tag)
405{
406 using namespace xml;
407 if (tag.getName() == VALUE_CONSTANT) {
408 XMLTag tagRelax(*this, TAG_RELAX, XMLTag::OCCUR_ONCE);
410 attrValue.setDocumentation("Constant relaxation factor.");
411 tagRelax.addAttribute(attrValue);
412 tag.addSubtag(tagRelax);
413 } else if (tag.getName() == VALUE_AITKEN) {
414 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
415 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.5 is used.");
417 attrValue.setDocumentation("Initial relaxation factor.");
418 tagInitRelax.addAttribute(attrValue);
419 tag.addSubtag(tagInitRelax);
420
422 tagData.setDocumentation("The data used to compute the acceleration.");
424 attrName.setDocumentation("The name of the data.");
426 attrMesh.setDocumentation("The name of the mesh which holds the data.");
427 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
428 .setDocumentation(
429 "To improve the performance of a parallel or a multi coupling schemes, "
430 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
431 "We recommend, however, to use an automatic scaling via a preconditioner.");
432 tagData.addAttribute(attrScaling);
433 tagData.addAttribute(attrName);
434 tagData.addAttribute(attrMesh);
435 tag.addSubtag(tagData);
436
437 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
438 tagPreconditioner.setDocumentation("To improve the numerical stability of multiple data vectors a preconditioner can be applied. "
439 "A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data. "
440 "A value preconditioner scales every acceleration data by the norm of the data in the previous time window. "
441 "A residual preconditioner scales every acceleration data by the current residual. "
442 "A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.");
443 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
448 .setDocumentation("The type of the preconditioner.");
449 tagPreconditioner.addAttribute(attrPreconditionerType);
450 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
451 .setDocumentation(
452 "After the given number of time windows, the preconditioner weights "
453 "are frozen and the preconditioner acts like a constant preconditioner.");
454 tagPreconditioner.addAttribute(nonconstTWindows);
455 tag.addSubtag(tagPreconditioner);
456 } else if (tag.getName() == VALUE_IQNILS) {
457
458 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
459 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
460 XMLAttribute<double> attrDoubleValue(ATTR_VALUE);
461 attrDoubleValue.setDocumentation("Initial relaxation factor.");
462 tagInitRelax.addAttribute(attrDoubleValue);
463 XMLAttribute<bool> attrEnforce(ATTR_ENFORCE, false);
464 attrEnforce.setDocumentation("Enforce initial relaxation in every time window.");
465 tagInitRelax.addAttribute(attrEnforce);
466 tag.addSubtag(tagInitRelax);
467
469 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.");
470 XMLAttribute<int> attrIntValue(ATTR_VALUE);
471 attrIntValue.setDocumentation("The number of columns.");
472 tagMaxUsedIter.addAttribute(attrIntValue);
473 tag.addSubtag(tagMaxUsedIter);
474
475 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
476 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.");
477 XMLAttribute<int> attrNumTimeWindowsReused(ATTR_VALUE);
478 attrNumTimeWindowsReused.setDocumentation("The number of time windows.");
479 tagTimeWindowsReused.addAttribute(attrNumTimeWindowsReused);
480 tag.addSubtag(tagTimeWindowsReused);
481
483
484 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
485 tagPreconditioner.setDocumentation("To improve the performance of a parallel or a multi coupling schemes a preconditioner can be applied.\n\n"
486 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data. \n "
487 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
488 "- A residual preconditioner scales every acceleration data by the current residual.\n"
489 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n\n"
490 "If this tag is not provided, the residual-sum preconditioner is employed.");
491 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
496 .setDocumentation("The type of the preconditioner.");
497 tagPreconditioner.addAttribute(attrPreconditionerType);
498 auto attrpreconditionerUpdateOnThreshold = XMLAttribute<bool>(ATTR_PRECOND_UPDATE_ON_THRESHOLD, true)
499 .setDocumentation("To update the preconditioner weights after the first time window: "
500 "`true`: The preconditioner weights are only updated if the weights will change by more than one order of magnitude. "
501 "`false`: The preconditioner weights are updated after every iteration.");
502 tagPreconditioner.addAttribute(attrpreconditionerUpdateOnThreshold);
503 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
504 .setDocumentation(
505 "After the given number of time windows, the preconditioner weights "
506 "are frozen and the preconditioner acts like a constant preconditioner.");
507 tagPreconditioner.addAttribute(nonconstTWindows);
508 tag.addSubtag(tagPreconditioner);
509
510 } else if (tag.getName() == VALUE_IQNIMVJ) {
511 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
512 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
513 tagInitRelax.addAttribute(
514 XMLAttribute<double>(ATTR_VALUE).setDocumentation("Initial relaxation factor."));
515 tagInitRelax.addAttribute(
516 XMLAttribute<bool>(ATTR_ENFORCE, false).setDocumentation("Enforce initial relaxation in every time window."));
517 tag.addSubtag(tagInitRelax);
518
519 XMLTag tagIMVJRESTART(*this, TAG_IMVJRESTART, XMLTag::OCCUR_NOT_OR_ONCE);
520 auto attrRestartName = XMLAttribute<std::string>(ATTR_TYPE)
526 .setDefaultValue(VALUE_SVD_RESTART)
527 .setDocumentation("Type of the restart mode.");
528 tagIMVJRESTART.addAttribute(attrRestartName);
529 tagIMVJRESTART.setDocumentation("Enable IMVJ Type of IMVJ restart mode that is used: "
530 "`no-restart`: IMVJ runs in normal mode with explicit representation of Jacobian. "
531 "`RS-0`: IMVJ runs in restart mode. After M time windows all Jacobain information is dropped, restart with no information. "
532 "`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. "
533 "`RS-SVD`: IMVJ runs in restart mode. After M time windows a truncated SVD of the Jacobian is updated. "
534 "`RS-SLIDE`: IMVJ runs in sliding window restart mode. "
535 "If this tag is not provided, IMVJ runs in restart mode with SVD-method.");
536 auto attrChunkSize = makeXMLAttribute(ATTR_IMVJCHUNKSIZE, 8)
537 .setDocumentation("Specifies the number of time windows M after which the IMVJ restarts, if run in restart-mode. Default value is M=8.");
538 auto attrReusedTimeWindowsAtRestart = makeXMLAttribute(ATTR_RSLS_REUSED_TIME_WINDOWS, 8)
539 .setDocumentation("If IMVJ restart-mode=RS-LS, the number of reused time windows at restart can be specified.");
540 auto attrRSSVD_truncationEps = makeXMLAttribute(ATTR_RSSVD_TRUNCATIONEPS, 1e-4)
541 .setDocumentation("If IMVJ restart-mode=RS-SVD, the truncation threshold for the updated SVD can be set.");
542 tagIMVJRESTART.addAttribute(attrChunkSize);
543 tagIMVJRESTART.addAttribute(attrReusedTimeWindowsAtRestart);
544 tagIMVJRESTART.addAttribute(attrRSSVD_truncationEps);
545 tag.addSubtag(tagIMVJRESTART);
546
548 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.");
549 XMLAttribute<int> attrIntValue(ATTR_VALUE);
550 attrIntValue.setDocumentation("The number of columns.");
551 tagMaxUsedIter.addAttribute(attrIntValue);
552 tag.addSubtag(tagMaxUsedIter);
553
554 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
555 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.");
556 tagTimeWindowsReused.addAttribute(attrIntValue);
557 tag.addSubtag(tagTimeWindowsReused);
558
560
561 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
562 tagPreconditioner.setDocumentation(
563 "To improve the performance of a parallel or a multi coupling schemes a preconditioner can be applied.\n\n"
564 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data.\n"
565 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
566 "- A residual preconditioner scales every acceleration data by the current residual.\n"
567 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n\n"
568 "If this tag is not provided, the residual-sum preconditioner is employed.");
569 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
574 .setDocumentation("Type of the preconditioner.");
575 tagPreconditioner.addAttribute(attrPreconditionerType);
576 auto attrpreconditionerUpdateOnThreshold = XMLAttribute<bool>(ATTR_PRECOND_UPDATE_ON_THRESHOLD, true)
577 .setDocumentation("To update the preconditioner weights after the first time window: "
578 "`true`: The preconditioner weights are only updated if the weights will change by more than one order of magnitude. "
579 "`false`: The preconditioner weights are updated after every iteration.");
580 tagPreconditioner.addAttribute(attrpreconditionerUpdateOnThreshold);
581 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
582 .setDocumentation("After the given number of time windows, the preconditioner weights are frozen and the preconditioner acts like a constant preconditioner.");
583 tagPreconditioner.addAttribute(nonconstTWindows);
584 tag.addSubtag(tagPreconditioner);
585
586 } else {
587 PRECICE_ERROR("Acceleration of type \"{}\" is unknown. Please choose a valid acceleration scheme or check the spelling in the configuration file.", tag.getName());
588 }
589}
590
592{
593 std::vector<double> factors;
594 for (int id : dataIDs) {
595 factors.push_back(scalings.at(id));
596 }
597 return factors;
598}
599
600} // 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
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 get(T... args)
T make_pair(T... args)
std::shared_ptr< Preconditioner > PtrPreconditioner
contains implementations of acceleration schemes.
std::shared_ptr< Acceleration > PtrAcceleration
provides Mesh, Data and primitives.
std::shared_ptr< Data > PtrData
std::shared_ptr< Mesh > PtrMesh
std::shared_ptr< MeshConfiguration > PtrMeshConfiguration
contains the XML configuration parser.
XMLAttribute< std::string > makeXMLAttribute(std::string name, const char *defaultValue)
T push_back(T... args)
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21