preCICE
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_ON_BOUND_VIOLATION("on-bound-violation"),
52 ATTR_REDUCEDTIMEGRIDQN("reduced-time-grid"),
53 ATTR_IMVJCHUNKSIZE("chunk-size"),
54 ATTR_RSLS_REUSED_TIME_WINDOWS("reused-time-windows-at-restart"),
55 ATTR_RSSVD_TRUNCATIONEPS("truncation-threshold"),
57 ATTR_PRECOND_UPDATE_ON_THRESHOLD("update-on-threshold"),
58 VALUE_IGNORE("ignore"),
59 VALUE_CLAMP("clamp"),
60 VALUE_DISCARD("discard"),
61 VALUE_SCALE_TO_BOUND("scale"),
62 VALUE_CONSTANT("constant"),
63 VALUE_AITKEN("aitken"),
64 VALUE_IQNILS("IQN-ILS"),
65 VALUE_IQNIMVJ("IQN-IMVJ"),
66 VALUE_QR1FILTER("QR1"),
67 VALUE_QR1_ABSFILTER("QR1-absolute"),
68 VALUE_QR2FILTER("QR2"),
69 VALUE_QR3FILTER("QR3"),
74 VALUE_LS_RESTART("RS-LS"),
75 VALUE_ZERO_RESTART("RS-0"),
76 VALUE_SVD_RESTART("RS-SVD"),
77 VALUE_SLIDE_RESTART("RS-SLIDE"),
78 VALUE_NO_RESTART("no-restart"),
79 _meshConfig(meshConfig),
83 _config()
84{
85 PRECICE_ASSERT(meshConfig.get() != nullptr);
86}
87
89{
90 using namespace xml;
91
92 // static int recursionCounter = 0;
93 // recursionCounter++;
94
95 XMLTag::Occurrence occ = XMLTag::OCCUR_NOT_OR_ONCE;
96 std::vector<XMLTag> tags;
97 {
98 XMLTag tag(*this, VALUE_CONSTANT, occ, TAG);
99 tag.setDocumentation("Accelerates coupling data with constant underrelaxation.");
101 tags.push_back(tag);
102 }
103 {
104 XMLTag tag(*this, VALUE_AITKEN, occ, TAG);
105 tag.setDocumentation("Accelerates coupling data with dynamic Aitken under-relaxation.");
107 tags.push_back(tag);
108 }
109 {
110 XMLTag tag(*this, VALUE_IQNILS, occ, TAG);
111 tag.setDocumentation("Accelerates coupling data with the interface quasi-Newton inverse least-squares method.");
112
113 auto reducedTimeGridQN = makeXMLAttribute(ATTR_REDUCEDTIMEGRIDQN, true)
114 .setDocumentation("Whether only the last time step of each time window is used to construct the Jacobian.");
115 tag.addAttribute(reducedTimeGridQN);
116
117 auto onBoundViolation = makeXMLAttribute(ATTR_ON_BOUND_VIOLATION, VALUE_IGNORE)
118 .setOptions({VALUE_IGNORE,
122 .setDocumentation("Defines the strategy to handle updates that violate variable bounds. Use ignore when no special handling is desired. Use clamp to limit the violating components to their bounds. Use discard to skip the QN update when a bound violation occurs. Use scale to scale the QN step with a constant to fit all violating components into the bounds.");
123 tag.addAttribute(onBoundViolation);
124
126 tags.push_back(tag);
127 }
128 {
129 XMLTag tag(*this, VALUE_IQNIMVJ, occ, TAG);
130 tag.setDocumentation("Accelerates coupling data with the interface quasi-Newton inverse multi-vector Jacobian method.");
131
132 auto alwaybuildJacobian = makeXMLAttribute(ATTR_BUILDJACOBIAN, false)
133 .setDocumentation("If set to true, the IMVJ will set up the Jacobian matrix in each coupling iteration, which is inefficient. "
134 "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.");
135 tag.addAttribute(alwaybuildJacobian);
136
137 auto reducedTimeGridQN = makeXMLAttribute(ATTR_REDUCEDTIMEGRIDQN, true)
138 .setDocumentation("Whether only the last time step of each time window is used to construct the Jacobian.");
139 tag.addAttribute(reducedTimeGridQN);
140
141 auto onBoundViolation = makeXMLAttribute(ATTR_ON_BOUND_VIOLATION, VALUE_IGNORE)
142 .setOptions({VALUE_IGNORE,
146 .setDocumentation("Defines the strategy to handle updates that violate variable bounds. Use ignore when no special handling is desired. Use clamp to limit the violating components to their bounds. Use discard to skip the QN update when a bound violation occurs. Use scale to scale the QN step with a constant to fit all violating components into the bounds.");
147 tag.addAttribute(onBoundViolation);
148
150 tags.push_back(tag);
151 }
152
153 for (XMLTag &tag : tags) {
154 parent.addSubtag(tag);
155 }
156}
157
162
164 const xml::ConfigurationContext &context,
165 xml::XMLTag &callingTag)
166{
167 PRECICE_TRACE(callingTag.getFullName());
168
169 if (callingTag.getNamespace() == TAG) {
170 _config.type = callingTag.getName();
171
172 if (_config.type == VALUE_IQNIMVJ) {
173 _config.alwaysBuildJacobian = callingTag.getBooleanAttributeValue(ATTR_BUILDJACOBIAN);
174 }
175
176 if (_config.type == VALUE_IQNIMVJ || _config.type == VALUE_IQNILS) {
177 _config.reducedTimeGridQN = callingTag.getBooleanAttributeValue(ATTR_REDUCEDTIMEGRIDQN);
178 std::string onBoundViolation = callingTag.getStringAttributeValue(ATTR_ON_BOUND_VIOLATION);
179
180 if (onBoundViolation == VALUE_IGNORE) {
182 } else if (onBoundViolation == VALUE_CLAMP) {
184 } else if (onBoundViolation == VALUE_DISCARD) {
186 } else if (onBoundViolation == VALUE_SCALE_TO_BOUND) {
188 } else {
190 }
191 }
192 }
193 if (callingTag.getName() == TAG_RELAX) {
194 _config.relaxationFactor = callingTag.getDoubleAttributeValue(ATTR_VALUE);
195 } else if (callingTag.getName() == TAG_DATA) {
196 std::string dataName = callingTag.getStringAttributeValue(ATTR_NAME);
197 std::string meshName = callingTag.getStringAttributeValue(ATTR_MESH);
198 auto success = _uniqueDataAndMeshNames.emplace(dataName, meshName);
199
200 PRECICE_CHECK(success.second,
201 "You have provided a subtag <data name=\"{}\" mesh=\"{}\"/> more than once in your <acceleration:.../>. "
202 "Please remove the duplicated entry.",
203 dataName, meshName);
204
206 double scaling = 1.0;
207
208 if (_config.type == VALUE_IQNILS || _config.type == VALUE_IQNIMVJ) {
209 scaling = callingTag.getDoubleAttributeValue(ATTR_SCALING);
210 }
211
212 PRECICE_CHECK(_meshConfig->hasMeshName(_meshName) && _meshConfig->getMesh(_meshName)->hasDataName(dataName),
213 "Data with name \"{0}\" associated to mesh \"{1}\" not found on configuration of acceleration. "
214 "Add \"{0}\" to the \"<mesh name={1}>\" tag, or change the data name in the acceleration scheme.",
215 dataName, _meshName);
216
217 const mesh::PtrMesh &mesh = _meshConfig->getMesh(_meshName);
218 const mesh::PtrData &data = mesh->data(dataName);
219 _config.dataIDs.push_back(data->getID());
220 _config.scalings.insert(std::make_pair(data->getID(), scaling));
221
222 _neededMeshes.push_back(_meshName);
223 } else if (callingTag.getName() == TAG_INIT_RELAX) {
224 _userDefinitions.definedRelaxationFactor = true;
225 _config.relaxationFactor = callingTag.getDoubleAttributeValue(ATTR_VALUE);
226 if (callingTag.hasAttribute(ATTR_ENFORCE)) {
227 _config.forceInitialRelaxation = callingTag.getBooleanAttributeValue(ATTR_ENFORCE);
228 } else {
229 _config.forceInitialRelaxation = false;
230 }
231 } else if (callingTag.getName() == TAG_MAX_USED_ITERATIONS) {
232 _userDefinitions.definedMaxIterationsUsed = true;
233 _config.maxIterationsUsed = callingTag.getIntAttributeValue(ATTR_VALUE);
234 } else if (callingTag.getName() == TAG_TIME_WINDOWS_REUSED) {
235 _userDefinitions.definedTimeWindowsReused = true;
236 _config.timeWindowsReused = callingTag.getIntAttributeValue(ATTR_VALUE);
237 } else if (callingTag.getName() == TAG_FILTER) {
238 _userDefinitions.definedFilter = true;
239 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
240 if (f == VALUE_QR1FILTER) {
242 } else if (f == VALUE_QR1_ABSFILTER) {
244 } else if (f == VALUE_QR2FILTER) {
246 } else if (f == VALUE_QR3FILTER) {
248 } else {
249 PRECICE_ASSERT(false);
250 }
251 _config.singularityLimit = callingTag.getDoubleAttributeValue(ATTR_SINGULARITYLIMIT);
252 } else if (callingTag.getName() == TAG_PRECONDITIONER) {
253 _userDefinitions.definedPreconditionerType = true;
254 _config.preconditionerType = callingTag.getStringAttributeValue(ATTR_TYPE);
255 _config.preconditionerUpdateOnThreshold = callingTag.getBooleanAttributeValue(ATTR_PRECOND_UPDATE_ON_THRESHOLD);
256 _config.preconditionerNbNonConstTWindows = callingTag.getIntAttributeValue(ATTR_PRECOND_NONCONST_TIME_WINDOWS);
257 } else if (callingTag.getName() == TAG_IMVJRESTART) {
258 _userDefinitions.defineRestartType = true;
259#ifndef PRECICE_NO_MPI
260 _config.imvjChunkSize = callingTag.getIntAttributeValue(ATTR_IMVJCHUNKSIZE);
261 const auto &f = callingTag.getStringAttributeValue(ATTR_TYPE);
262 PRECICE_CHECK((f == VALUE_NO_RESTART) || (!_config.alwaysBuildJacobian), "IMVJ cannot be in restart mode while parameter always-build-jacobian is set to true. "
263 "Please remove 'always-build-jacobian' from the configuration file or do not run in restart mode.");
264 if (f == VALUE_NO_RESTART) {
266 } else if (f == VALUE_ZERO_RESTART) {
267 _config.imvjRestartType = IQNIMVJAcceleration::RS_ZERO;
268 } else if (f == VALUE_LS_RESTART) {
269 _config.imvjRSLS_reusedTimeWindows = callingTag.getIntAttributeValue(ATTR_RSLS_REUSED_TIME_WINDOWS);
270 _config.imvjRestartType = IQNIMVJAcceleration::RS_LS;
271 } else if (f == VALUE_SVD_RESTART) {
272 _config.imvjRSSVD_truncationEps = callingTag.getDoubleAttributeValue(ATTR_RSSVD_TRUNCATIONEPS);
273 _config.imvjRestartType = IQNIMVJAcceleration::RS_SVD;
274 } else if (f == VALUE_SLIDE_RESTART) {
275 _config.imvjRestartType = IQNIMVJAcceleration::RS_SLIDE;
276 } else {
277 _config.imvjChunkSize = 0;
278 PRECICE_ASSERT(false);
279 }
280#else
281 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
282#endif
283 }
284}
285
287 const xml::ConfigurationContext &context,
288 xml::XMLTag &callingTag)
289{
290 PRECICE_TRACE(callingTag.getName());
291 if (callingTag.getNamespace() == TAG) {
292
293 // create preconditioner
294 if (callingTag.getName() == VALUE_IQNILS || callingTag.getName() == VALUE_AITKEN) {
295 if (_config.preconditionerType == VALUE_CONSTANT_PRECONDITIONER) {
296 _preconditioner = PtrPreconditioner(new ConstantPreconditioner(_config.scalingFactorsInOrder()));
297 } else if (_config.preconditionerType == VALUE_VALUE_PRECONDITIONER) {
298 _preconditioner = PtrPreconditioner(new ValuePreconditioner(_config.preconditionerNbNonConstTWindows));
299 } else if (_config.preconditionerType == VALUE_RESIDUAL_PRECONDITIONER) {
300 _preconditioner = PtrPreconditioner(new ResidualPreconditioner(_config.preconditionerNbNonConstTWindows));
301 } else if (_config.preconditionerType == VALUE_RESIDUAL_SUM_PRECONDITIONER) {
302 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_config.preconditionerNbNonConstTWindows, _config.preconditionerUpdateOnThreshold));
303 } else {
304 // no preconditioner defined
305 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_defaultValuesIQNILS.preconditionerNbNonConstTWindows, _defaultValuesIQNILS.preconditionerUpdateOnThreshold));
306 }
307 }
308
309 PRECICE_CHECK(!_acceleration, "You are trying to define multiple acceleration schemes, which is not allowed. Please remove one of them.");
310 if (callingTag.getName() == VALUE_CONSTANT) {
313 _config.relaxationFactor, _config.dataIDs));
314 } else if (callingTag.getName() == VALUE_AITKEN) {
315 _config.relaxationFactor = (_userDefinitions.definedRelaxationFactor) ? _config.relaxationFactor : _defaultAitkenRelaxationFactor;
318 _config.relaxationFactor, _config.dataIDs, _preconditioner));
319 } else if (callingTag.getName() == VALUE_IQNILS) {
320 _config.relaxationFactor = (_userDefinitions.definedRelaxationFactor) ? _config.relaxationFactor : _defaultValuesIQNILS.relaxationFactor;
321 _config.maxIterationsUsed = (_userDefinitions.definedMaxIterationsUsed) ? _config.maxIterationsUsed : _defaultValuesIQNILS.maxIterationsUsed;
322 _config.timeWindowsReused = (_userDefinitions.definedTimeWindowsReused) ? _config.timeWindowsReused : _defaultValuesIQNILS.timeWindowsReused;
323 _config.filter = (_userDefinitions.definedFilter) ? _config.filter : _defaultValuesIQNILS.filter;
324 _config.singularityLimit = (_userDefinitions.definedFilter) ? _config.singularityLimit : _defaultValuesIQNILS.singularityLimit;
325
328 _config.relaxationFactor,
329 _config.forceInitialRelaxation,
330 _config.maxIterationsUsed,
331 _config.timeWindowsReused,
332 _config.filter, _config.singularityLimit,
333 _config.dataIDs,
334 _config.onBoundViolation,
336 _config.reducedTimeGridQN));
337 } else if (callingTag.getName() == VALUE_IQNIMVJ) {
338#ifndef PRECICE_NO_MPI
339 _config.relaxationFactor = (_userDefinitions.definedRelaxationFactor) ? _config.relaxationFactor : _defaultValuesIQNIMVJ.relaxationFactor;
340 _config.maxIterationsUsed = (_userDefinitions.definedMaxIterationsUsed) ? _config.maxIterationsUsed : _defaultValuesIQNIMVJ.maxIterationsUsed;
341 _config.timeWindowsReused = (_userDefinitions.definedTimeWindowsReused) ? _config.timeWindowsReused : _defaultValuesIQNIMVJ.timeWindowsReused;
342 _config.filter = (_userDefinitions.definedFilter) ? _config.filter : _defaultValuesIQNILS.filter;
343 _config.singularityLimit = (_userDefinitions.definedFilter) ? _config.singularityLimit : _defaultValuesIQNILS.singularityLimit;
344 if (!_config.alwaysBuildJacobian) {
345 _config.imvjRestartType = (_userDefinitions.defineRestartType) ? _config.imvjRestartType : _defaultValuesIQNIMVJ.imvjRestartType;
346 _config.imvjChunkSize = (_userDefinitions.defineRestartType) ? _config.imvjChunkSize : _defaultValuesIQNIMVJ.imvjChunkSize;
347 }
348
349 // create preconditioner
350 // if imvj restart-mode is of type RS-SVD, max number of non-const preconditioned time windows is limited by the chunksize
351 // it is separated from the other acceleration methods, since SVD-restart might be chosen as default here
352 if (_config.imvjRestartType == IQNIMVJAcceleration::RS_SVD)
353 if (_config.preconditionerNbNonConstTWindows > _config.imvjChunkSize)
354 _config.preconditionerNbNonConstTWindows = _config.imvjChunkSize;
355 if (_config.preconditionerType == VALUE_CONSTANT_PRECONDITIONER) {
356 _preconditioner = PtrPreconditioner(new ConstantPreconditioner(_config.scalingFactorsInOrder()));
357 } else if (_config.preconditionerType == VALUE_VALUE_PRECONDITIONER) {
358 _preconditioner = PtrPreconditioner(new ValuePreconditioner(_config.preconditionerNbNonConstTWindows));
359 } else if (_config.preconditionerType == VALUE_RESIDUAL_PRECONDITIONER) {
360 _preconditioner = PtrPreconditioner(new ResidualPreconditioner(_config.preconditionerNbNonConstTWindows));
361 } else if (_config.preconditionerType == VALUE_RESIDUAL_SUM_PRECONDITIONER) {
362 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_config.preconditionerNbNonConstTWindows, _config.preconditionerUpdateOnThreshold));
363 } else {
364 // no preconditioner defined
365 _preconditioner = PtrPreconditioner(new ResidualSumPreconditioner(_defaultValuesIQNILS.preconditionerNbNonConstTWindows, _defaultValuesIQNIMVJ.preconditionerUpdateOnThreshold));
366 }
367
370 _config.relaxationFactor,
371 _config.forceInitialRelaxation,
372 _config.maxIterationsUsed,
373 _config.timeWindowsReused,
374 _config.filter, _config.singularityLimit,
375 _config.dataIDs,
376 _config.onBoundViolation,
378 _config.alwaysBuildJacobian,
379 _config.imvjRestartType,
380 _config.imvjChunkSize,
381 _config.imvjRSLS_reusedTimeWindows,
382 _config.imvjRSSVD_truncationEps,
383 _config.reducedTimeGridQN));
384#else
385 PRECICE_ERROR("Acceleration IQN-IMVJ only works if preCICE is compiled with MPI");
386#endif
387 } else {
388 PRECICE_ASSERT(false);
389 }
390 }
391}
392
400
402{
403 using namespace precice::xml;
404
406 tagData.setDocumentation("The data used to compute the acceleration.");
408 attrName.setDocumentation("The name of the data.");
410 attrMesh.setDocumentation("The name of the mesh which holds the data.");
411 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
412 .setDocumentation(
413 "To improve the performance of a parallel or a multi coupling schemes, "
414 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
415 "We recommend, however, to use an automatic scaling via a preconditioner.");
416 tagData.addAttribute(attrName);
417 tagData.addAttribute(attrMesh);
418 tagData.addAttribute(attrScaling);
419 tag.addSubtag(tagData);
420
422 tagFilter.setDocumentation("Type of filtering technique that is used to "
423 "maintain good conditioning in the least-squares system. Possible filters:\n\n"
424 "- `QR1`: update QR-dec with (relative) test \\\\(R(i,i) < \\epsilon *\\lVert R\\rVert_F\\\\)\n"
425 "- `QR1-absolute`: update QR-dec with (absolute) test \\\\(R(i, i) < \\epsilon\\\\)\n"
426 "- `QR2`: en-block QR-dec with test \\\\(\\lVert v_\\text{orth} \\rVert_2 < \\epsilon * \\lVert v \\rVert_2\\\\)\n\n"
427 "- `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"
428 "Please note that a QR1 is based on Given's rotations whereas QR2 uses modified Gram-Schmidt. "
429 "This can give different results even when no columns are filtered out.");
430 XMLAttribute<double> attrSingularityLimit(ATTR_SINGULARITYLIMIT, 1e-16);
431 attrSingularityLimit.setDocumentation("Limit eps of the filter.");
432 tagFilter.addAttribute(attrSingularityLimit);
433 auto attrFilterName = XMLAttribute<std::string>(ATTR_TYPE)
438 .setDefaultValue(VALUE_QR3FILTER)
439 .setDocumentation("Type of the filter.");
440 tagFilter.addAttribute(attrFilterName);
441 tag.addSubtag(tagFilter);
442}
443
445 xml::XMLTag &tag)
446{
447 using namespace xml;
448 if (tag.getName() == VALUE_CONSTANT) {
449 XMLTag tagRelax(*this, TAG_RELAX, XMLTag::OCCUR_ONCE);
451 attrValue.setDocumentation("Constant relaxation factor.");
452 tagRelax.addAttribute(attrValue);
453 tag.addSubtag(tagRelax);
454 } else if (tag.getName() == VALUE_AITKEN) {
455 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
456 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.5 is used.");
458 attrValue.setDocumentation("Initial relaxation factor.");
459 tagInitRelax.addAttribute(attrValue);
460 tag.addSubtag(tagInitRelax);
461
463 tagData.setDocumentation("The data used to compute the acceleration.");
465 attrName.setDocumentation("The name of the data.");
467 attrMesh.setDocumentation("The name of the mesh which holds the data.");
468 auto attrScaling = makeXMLAttribute(ATTR_SCALING, 1.0)
469 .setDocumentation(
470 "To improve the performance of a parallel or a multi coupling schemes, "
471 "each data set can be manually scaled using this scaling factor with preconditioner type = \"constant\". For all other preconditioner types, the factor is ignored. "
472 "We recommend, however, to use an automatic scaling via a preconditioner.");
473 tagData.addAttribute(attrScaling);
474 tagData.addAttribute(attrName);
475 tagData.addAttribute(attrMesh);
476 tag.addSubtag(tagData);
477
478 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
479 tagPreconditioner.setDocumentation("To improve the numerical stability of multiple data vectors a preconditioner can be applied. "
480 "A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data. "
481 "A value preconditioner scales every acceleration data by the norm of the data in the previous time window. "
482 "A residual preconditioner scales every acceleration data by the current residual. "
483 "A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.");
484 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
489 .setDocumentation("The type of the preconditioner.");
490 tagPreconditioner.addAttribute(attrPreconditionerType);
491 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
492 .setDocumentation(
493 "After the given number of time windows, the preconditioner weights "
494 "are frozen and the preconditioner acts like a constant preconditioner.");
495 tagPreconditioner.addAttribute(nonconstTWindows);
496 tag.addSubtag(tagPreconditioner);
497 } else if (tag.getName() == VALUE_IQNILS) {
498
499 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
500 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
501 XMLAttribute<double> attrDoubleValue(ATTR_VALUE);
502 attrDoubleValue.setDocumentation("Initial relaxation factor.");
503 tagInitRelax.addAttribute(attrDoubleValue);
504 XMLAttribute<bool> attrEnforce(ATTR_ENFORCE, false);
505 attrEnforce.setDocumentation("Enforce initial relaxation in every time window.");
506 tagInitRelax.addAttribute(attrEnforce);
507 tag.addSubtag(tagInitRelax);
508
510 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.");
511 XMLAttribute<int> attrIntValue(ATTR_VALUE);
512 attrIntValue.setDocumentation("The number of columns.");
513 tagMaxUsedIter.addAttribute(attrIntValue);
514 tag.addSubtag(tagMaxUsedIter);
515
516 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
517 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.");
518 XMLAttribute<int> attrNumTimeWindowsReused(ATTR_VALUE);
519 attrNumTimeWindowsReused.setDocumentation("The number of time windows.");
520 tagTimeWindowsReused.addAttribute(attrNumTimeWindowsReused);
521 tag.addSubtag(tagTimeWindowsReused);
522
524
525 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
526 tagPreconditioner.setDocumentation("To improve the performance of a parallel or a multi coupling schemes a preconditioner can be applied.\n\n"
527 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data. \n "
528 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
529 "- A residual preconditioner scales every acceleration data by the current residual.\n"
530 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n\n"
531 "If this tag is not provided, the residual-sum preconditioner is employed.");
532 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
537 .setDocumentation("The type of the preconditioner.");
538 tagPreconditioner.addAttribute(attrPreconditionerType);
539 auto attrpreconditionerUpdateOnThreshold = XMLAttribute<bool>(ATTR_PRECOND_UPDATE_ON_THRESHOLD, true)
540 .setDocumentation("To update the preconditioner weights after the first time window: "
541 "`true`: The preconditioner weights are only updated if the weights will change by more than one order of magnitude. "
542 "`false`: The preconditioner weights are updated after every iteration.");
543 tagPreconditioner.addAttribute(attrpreconditionerUpdateOnThreshold);
544 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
545 .setDocumentation(
546 "After the given number of time windows, the preconditioner weights "
547 "are frozen and the preconditioner acts like a constant preconditioner.");
548 tagPreconditioner.addAttribute(nonconstTWindows);
549 tag.addSubtag(tagPreconditioner);
550
551 } else if (tag.getName() == VALUE_IQNIMVJ) {
552 XMLTag tagInitRelax(*this, TAG_INIT_RELAX, XMLTag::OCCUR_NOT_OR_ONCE);
553 tagInitRelax.setDocumentation("Initial relaxation factor. If this tag is not provided, an initial relaxation of 0.1 is used.");
554 tagInitRelax.addAttribute(
555 XMLAttribute<double>(ATTR_VALUE).setDocumentation("Initial relaxation factor."));
556 tagInitRelax.addAttribute(
557 XMLAttribute<bool>(ATTR_ENFORCE, false).setDocumentation("Enforce initial relaxation in every time window."));
558 tag.addSubtag(tagInitRelax);
559
560 XMLTag tagIMVJRESTART(*this, TAG_IMVJRESTART, XMLTag::OCCUR_NOT_OR_ONCE);
561 auto attrRestartName = XMLAttribute<std::string>(ATTR_TYPE)
567 .setDefaultValue(VALUE_SVD_RESTART)
568 .setDocumentation("Type of the restart mode.");
569 tagIMVJRESTART.addAttribute(attrRestartName);
570 tagIMVJRESTART.setDocumentation("Enable IMVJ Type of IMVJ restart mode that is used: "
571 "`no-restart`: IMVJ runs in normal mode with explicit representation of Jacobian. "
572 "`RS-0`: IMVJ runs in restart mode. After M time windows all Jacobain information is dropped, restart with no information. "
573 "`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. "
574 "`RS-SVD`: IMVJ runs in restart mode. After M time windows a truncated SVD of the Jacobian is updated. "
575 "`RS-SLIDE`: IMVJ runs in sliding window restart mode. "
576 "If this tag is not provided, IMVJ runs in restart mode with SVD-method.");
577 auto attrChunkSize = makeXMLAttribute(ATTR_IMVJCHUNKSIZE, 8)
578 .setDocumentation("Specifies the number of time windows M after which the IMVJ restarts, if run in restart-mode. Default value is M=8.");
579 auto attrReusedTimeWindowsAtRestart = makeXMLAttribute(ATTR_RSLS_REUSED_TIME_WINDOWS, 8)
580 .setDocumentation("If IMVJ restart-mode=RS-LS, the number of reused time windows at restart can be specified.");
581 auto attrRSSVD_truncationEps = makeXMLAttribute(ATTR_RSSVD_TRUNCATIONEPS, 1e-4)
582 .setDocumentation("If IMVJ restart-mode=RS-SVD, the truncation threshold for the updated SVD can be set.");
583 tagIMVJRESTART.addAttribute(attrChunkSize);
584 tagIMVJRESTART.addAttribute(attrReusedTimeWindowsAtRestart);
585 tagIMVJRESTART.addAttribute(attrRSSVD_truncationEps);
586 tag.addSubtag(tagIMVJRESTART);
587
589 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.");
590 XMLAttribute<int> attrIntValue(ATTR_VALUE);
591 attrIntValue.setDocumentation("The number of columns.");
592 tagMaxUsedIter.addAttribute(attrIntValue);
593 tag.addSubtag(tagMaxUsedIter);
594
595 XMLTag tagTimeWindowsReused(*this, TAG_TIME_WINDOWS_REUSED, XMLTag::OCCUR_NOT_OR_ONCE);
596 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.");
597 tagTimeWindowsReused.addAttribute(attrIntValue);
598 tag.addSubtag(tagTimeWindowsReused);
599
601
602 XMLTag tagPreconditioner(*this, TAG_PRECONDITIONER, XMLTag::OCCUR_NOT_OR_ONCE);
603 tagPreconditioner.setDocumentation(
604 "To improve the performance of a parallel or a multi coupling schemes a preconditioner can be applied.\n\n"
605 "- A constant preconditioner scales every acceleration data by a constant value, which you can define as an attribute of data.\n"
606 "- A value preconditioner scales every acceleration data by the norm of the data in the previous time window.\n"
607 "- A residual preconditioner scales every acceleration data by the current residual.\n"
608 "- A residual-sum preconditioner scales every acceleration data by the sum of the residuals from the current time window.\n\n"
609 "If this tag is not provided, the residual-sum preconditioner is employed.");
610 auto attrPreconditionerType = XMLAttribute<std::string>(ATTR_TYPE)
615 .setDocumentation("Type of the preconditioner.");
616 tagPreconditioner.addAttribute(attrPreconditionerType);
617 auto attrpreconditionerUpdateOnThreshold = XMLAttribute<bool>(ATTR_PRECOND_UPDATE_ON_THRESHOLD, true)
618 .setDocumentation("To update the preconditioner weights after the first time window: "
619 "`true`: The preconditioner weights are only updated if the weights will change by more than one order of magnitude. "
620 "`false`: The preconditioner weights are updated after every iteration.");
621 tagPreconditioner.addAttribute(attrpreconditionerUpdateOnThreshold);
622 auto nonconstTWindows = makeXMLAttribute(ATTR_PRECOND_NONCONST_TIME_WINDOWS, -1)
623 .setDocumentation("After the given number of time windows, the preconditioner weights are frozen and the preconditioner acts like a constant preconditioner.");
624 tagPreconditioner.addAttribute(nonconstTWindows);
625 tag.addSubtag(tagPreconditioner);
626
627 } else {
628 PRECICE_ERROR("Acceleration of type \"{}\" is unknown. Please choose a valid acceleration scheme or check the spelling in the configuration file.", tag.getName());
629 }
630}
631
633{
634 std::vector<double> factors;
635 for (int id : dataIDs) {
636 factors.push_back(scalings.at(id));
637 }
638 return factors;
639}
640
641} // 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
#define PRECICE_UNREACHABLE(...)
Definition assertion.hpp:93
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
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)
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21