preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CouplingSchemeConfiguration.cpp
Go to the documentation of this file.
1#include <algorithm>
2#include <cstddef>
3#include <memory>
4#include <ostream>
5#include <stdexcept>
6#include <utility>
7#include <vector>
8
25#include "logging/LogMacros.hpp"
26#include "m2n/SharedPointer.hpp"
28#include "mesh/Data.hpp"
29#include "mesh/Mesh.hpp"
34#include "utils/Helpers.hpp"
35#include "utils/assertion.hpp"
36#include "xml/ConfigParser.hpp"
37#include "xml/XMLAttribute.hpp"
38#include "xml/XMLTag.hpp"
39
40namespace precice::cplscheme {
41
42const int CouplingSchemeConfiguration::DEFAULT_MIN_ITERATIONS(1); // min 1 iteration
44
46 xml::XMLTag &parent,
50 : TAG("coupling-scheme"),
51 TAG_PARTICIPANTS("participants"),
52 TAG_PARTICIPANT("participant"),
53 TAG_EXCHANGE("exchange"),
54 TAG_MAX_TIME("max-time"),
55 TAG_MAX_TIME_WINDOWS("max-time-windows"),
56 TAG_TIME_WINDOW_SIZE("time-window-size"),
57 TAG_ABS_CONV_MEASURE("absolute-convergence-measure"),
58 TAG_ABS_OR_REL_CONV_MEASURE("absolute-or-relative-convergence-measure"),
59 TAG_REL_CONV_MEASURE("relative-convergence-measure"),
60 TAG_RES_REL_CONV_MEASURE("residual-relative-convergence-measure"),
61 TAG_MIN_ITERATIONS("min-iterations"),
62 TAG_MAX_ITERATIONS("max-iterations"),
63 ATTR_DATA("data"),
64 ATTR_MESH("mesh"),
65 ATTR_PARTICIPANT("participant"),
66 ATTR_INITIALIZE("initialize"),
67 ATTR_EXCHANGE_SUBSTEPS("substeps"),
68 ATTR_TYPE("type"),
69 ATTR_FIRST("first"),
70 ATTR_SECOND("second"),
71 ATTR_VALUE("value"),
72 ATTR_METHOD("method"),
73 ATTR_LIMIT("limit"),
74 ATTR_ABS_LIMIT("abs-limit"),
75 ATTR_REL_LIMIT("rel-limit"),
76 ATTR_NAME("name"),
77 ATTR_FROM("from"),
78 ATTR_TO("to"),
79 ATTR_SUFFICES("suffices"),
80 ATTR_STRICT("strict"),
81 ATTR_CONTROL("control"),
82 VALUE_SERIAL_EXPLICIT("serial-explicit"),
83 VALUE_PARALLEL_EXPLICIT("parallel-explicit"),
84 VALUE_SERIAL_IMPLICIT("serial-implicit"),
85 VALUE_PARALLEL_IMPLICIT("parallel-implicit"),
86 VALUE_MULTI("multi"),
87 VALUE_FIXED("fixed"),
88 VALUE_FIRST_PARTICIPANT("first-participant"),
89 _config(),
90 _meshConfig(std::move(meshConfig)),
91 _m2nConfig(std::move(m2nConfig)),
92 _participantConfig(participantConfig),
93 _couplingSchemes(),
94 _couplingSchemeCompositions()
95{
96 using namespace xml;
97
98 XMLTag::Occurrence occ = XMLTag::OCCUR_ARBITRARY;
100 {
101 XMLTag tag(*this, VALUE_SERIAL_EXPLICIT, occ, TAG);
102 tag.setDocumentation("Explicit coupling scheme according to conventional serial staggered procedure (CSS).");
104 tags.push_back(tag);
105 }
106 {
107 XMLTag tag(*this, VALUE_PARALLEL_EXPLICIT, occ, TAG);
108 tag.setDocumentation("Explicit coupling scheme according to conventional parallel staggered procedure (CPS).");
110 tags.push_back(tag);
111 }
112 {
113 XMLTag tag(*this, VALUE_SERIAL_IMPLICIT, occ, TAG);
114 tag.setDocumentation("Implicit coupling scheme according to block Gauss-Seidel iterations (S-System). "
115 "Improved implicit iterations are achieved by using a acceleration (recommended!).");
117 tags.push_back(tag);
118 }
119 {
120 XMLTag tag(*this, VALUE_PARALLEL_IMPLICIT, occ, TAG);
121 tag.setDocumentation("Parallel Implicit coupling scheme according to block Jacobi iterations (V-System). "
122 "Improved implicit iterations are achieved by using a acceleration (recommended!).");
124 tags.push_back(tag);
125 }
126 {
127 XMLTag tag(*this, VALUE_MULTI, occ, TAG);
128 tag.setDocumentation("Multi coupling scheme according to block Jacobi iterations. "
129 "Improved implicit iterations are achieved by using a acceleration (recommended!).");
131 tags.push_back(tag);
132 }
133
134 for (XMLTag &tag : tags) {
135 parent.addSubtag(tag);
136 }
137}
138
140 const std::string &participantName) const
141{
142 return utils::contained(participantName, _couplingSchemes);
143}
144
146 const std::string &participantName) const
147{
149 "No coupling scheme defined for participant \"{}\". "
150 "Please make sure to provide at least one <coupling-scheme:TYPE> in your "
151 "precice-config.xml that couples this participant using the <participants .../> tag.",
152 participantName);
153 return _couplingSchemes.find(participantName)->second;
154}
155
157 const xml::ConfigurationContext &context,
158 xml::XMLTag &tag)
159{
161 if (tag.getNamespace() == TAG) {
162 _config.type = tag.getName();
163 _accelerationConfig->clear();
164 } else if (tag.getName() == TAG_PARTICIPANTS) {
167
168 PRECICE_CHECK(_participantConfig->hasParticipant(first),
169 "First participant in coupling-scheme <participants first=\"{}\" second=\"{}\" /> is unknown. {}",
170 first, second, _participantConfig->hintFor(first));
171 PRECICE_CHECK(_participantConfig->hasParticipant(second),
172 "Second participant in coupling-scheme <participants first=\"{}\" second=\"{}\" /> is unknown. {}",
173 first, second, _participantConfig->hintFor(second));
174 PRECICE_CHECK(first != second,
175 "First and second participant in coupling scheme are the same. "
176 "Please choose different in the <participants first=\"{}\" second=\"{}\" /> tag in the <coupling-scheme:...> of your precice-config.xml",
177 first, second);
180 } else if (tag.getName() == TAG_PARTICIPANT) {
182 bool control = tag.getBooleanAttributeValue(ATTR_CONTROL);
183 std::string participantName = tag.getStringAttributeValue(ATTR_NAME);
184 PRECICE_CHECK(_participantConfig->hasParticipant(participantName),
185 "Provided participant in multi coupling-scheme <participant name=\"{}\" ... /> is unknown. {}",
186 participantName, _participantConfig->hintFor(participantName));
188 "Participant \"{0}\" is provided multiple times to multi coupling scheme. "
189 "Please make sure that you do not provide the participant multiple times via the <participant name=\"{0}\" /> "
190 "tag in the <coupling-scheme:...> of your precice-config.xml",
191 participantName);
192 if (control) {
194 "Only one controller per MultiCouplingScheme can be defined. "
195 "Please check the <participant name=\"{}\" control=\"{}\" /> tag in the <coupling-scheme:...> of your precice-config.xml",
196 participantName, control);
197 _config.controller = participantName;
198 _config.setController = true;
199 }
200 _config.participants.push_back(participantName);
201 } else if (tag.getName() == TAG_MAX_TIME) {
204 "Maximum time has to be larger than zero. "
205 "Please check the <max-time value=\"{}\" /> tag in the <coupling-scheme:...> of your precice-config.xml",
207 } else if (tag.getName() == TAG_MAX_TIME_WINDOWS) {
210 "Maximum number of time windows has to be larger than zero. "
211 "Please check the <max-time-windows value=\"{}\" /> tag in the <coupling-scheme:...> of your precice-config.xml",
213 } else if (tag.getName() == TAG_TIME_WINDOW_SIZE) {
215 // Attribute does not exist for parallel coupling schemes as it is always fixed.
217
220 "The minimal time window size supported by preCICE is {}. "
221 "Please check the <time-window-size value=\"{}\" /> tag "
222 "in the <coupling-scheme:...> of your precice-config.xml and pick an appropriate time window size.",
224 } else {
227 "You combined a custom time-window-size of {} with method=\"first-participant\". "
228 "The given time-window-size will be ignored as it is prescribed by the participant.",
231 }
232 } else if (tag.getName() == TAG_ABS_CONV_MEASURE) {
233 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
234 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
235 double limit = tag.getDoubleAttributeValue(ATTR_LIMIT);
236 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
237 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
239 addAbsoluteConvergenceMeasure(dataName, meshName, limit, suffices, strict);
240 } else if (tag.getName() == TAG_ABS_OR_REL_CONV_MEASURE) {
241 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
242 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
243 double absLimit = tag.getDoubleAttributeValue(ATTR_ABS_LIMIT);
244 double relLimit = tag.getDoubleAttributeValue(ATTR_REL_LIMIT);
245 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
246 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
248 addAbsoluteOrRelativeConvergenceMeasure(dataName, meshName, absLimit, relLimit, suffices, strict);
249 } else if (tag.getName() == TAG_REL_CONV_MEASURE) {
250 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
251 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
252 double limit = tag.getDoubleAttributeValue(ATTR_LIMIT);
253 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
254 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
256 addRelativeConvergenceMeasure(dataName, meshName, limit, suffices, strict);
257 } else if (tag.getName() == TAG_RES_REL_CONV_MEASURE) {
258 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
259 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
260 double limit = tag.getDoubleAttributeValue(ATTR_LIMIT);
261 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
262 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
264 addResidualRelativeConvergenceMeasure(dataName, meshName, limit, suffices, strict);
265 } else if (tag.getName() == TAG_EXCHANGE) {
268 std::string nameParticipantFrom = tag.getStringAttributeValue(ATTR_FROM);
269 std::string nameParticipantTo = tag.getStringAttributeValue(ATTR_TO);
270 bool initialize = tag.getBooleanAttributeValue(ATTR_INITIALIZE);
271 bool exchangeSubsteps = tag.getBooleanAttributeValue(ATTR_EXCHANGE_SUBSTEPS);
272
273 PRECICE_CHECK(_meshConfig->hasMeshName(nameMesh) && _meshConfig->getMesh(nameMesh)->hasDataName(nameData),
274 "Mesh \"{}\" with data \"{}\" not defined. "
275 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> "
276 "tag in the <coupling-scheme:... /> of your precice-config.xml.",
277 nameMesh, nameData, nameData, nameMesh, nameParticipantFrom, nameParticipantTo);
278
279 mesh::PtrMesh exchangeMesh = _meshConfig->getMesh(nameMesh);
280 PRECICE_ASSERT(exchangeMesh);
281 mesh::PtrData exchangeData = exchangeMesh->data(nameData);
282 PRECICE_ASSERT(exchangeData);
283
284 Config::Exchange newExchange{exchangeData, exchangeMesh, nameParticipantFrom, nameParticipantTo, initialize, exchangeSubsteps};
285 PRECICE_CHECK(!_config.hasExchange(newExchange),
286 R"(Data "{}" of mesh "{}" cannot be exchanged multiple times between participants "{}" and "{}". Please remove one of the exchange tags.)",
287 nameData, nameMesh, nameParticipantFrom, nameParticipantTo);
288
289 _meshConfig->addNeededMesh(nameParticipantFrom, nameMesh);
290 _meshConfig->addNeededMesh(nameParticipantTo, nameMesh);
291 _config.exchanges.emplace_back(std::move(newExchange));
292 } else if (tag.getName() == TAG_MIN_ITERATIONS) {
296 "Minimum iteration limit has to be larger than zero. Please check the <min-iterations value = \"{}\" /> subtag in the <coupling-scheme:... /> of your precice-config.xml.",
298 } else if (tag.getName() == TAG_MAX_ITERATIONS) {
302 "Maximal iteration limit has to be larger than zero. "
303 "Please check the <max-iterations value=\"{0}\" /> subtag in the <coupling-scheme:... /> of your precice-config.xml. "
304 "To disable the iteration limit, remove the <max-iterations value=\"{0}\" /> subtag.",
306 }
307
308 // Additional consistency checks
311 "Maximum iteration limit {1} has to be larger or equal than the minimum iteration limit {0}. "
312 "Please check the <min-iterations value = \"{0}\" /> and <max-iterations value = \"{1}\" /> subtags in the <coupling-scheme:... /> of your precice-config.xml.",
315 }
316}
317
319 const xml::ConfigurationContext &context,
320 xml::XMLTag &tag)
321{
323 if (tag.getNamespace() == TAG) {
325 PRECICE_CHECK(!_allowRemeshing, "Remeshing is currently incompatible with serial coupling schemes. Try using a parallel or a multi coupling scheme instead.");
326 std::string accessor(_config.participants[0]);
328 addCouplingScheme(scheme, accessor);
329 //_couplingSchemes[accessor] = scheme;
330 accessor = _config.participants[1];
331 scheme = createSerialExplicitCouplingScheme(accessor);
332 addCouplingScheme(scheme, accessor);
333 //_couplingSchemes[accessor] = scheme;
334 _config = Config();
335 } else if (_config.type == VALUE_PARALLEL_EXPLICIT) {
336 std::string accessor(_config.participants[0]);
338 addCouplingScheme(scheme, accessor);
339 //_couplingSchemes[accessor] = scheme;
340 accessor = _config.participants[1];
341 scheme = createParallelExplicitCouplingScheme(accessor);
342 addCouplingScheme(scheme, accessor);
343 //_couplingSchemes[accessor] = scheme;
344 _config = Config();
345 } else if (_config.type == VALUE_SERIAL_IMPLICIT) {
346 PRECICE_CHECK(!_allowRemeshing, "Remeshing is currently incompatible with serial coupling schemes. Try using a parallel or a multi coupling scheme instead.");
348 std::string accessor(_config.participants[0]);
350 addCouplingScheme(scheme, accessor);
351 //_couplingSchemes[accessor] = scheme;
352 accessor = _config.participants[1];
353 scheme = createSerialImplicitCouplingScheme(accessor);
354 addCouplingScheme(scheme, accessor);
355 //_couplingSchemes[accessor] = scheme;
356 _config = Config();
357 } else if (_config.type == VALUE_PARALLEL_IMPLICIT) {
358 PRECICE_INFO_IF(_allowRemeshing, "Remeshing for implicit coupling schemes is in development. Currently, the acceleration data is deleted on remeshing.");
360 std::string accessor(_config.participants[0]);
362 addCouplingScheme(scheme, accessor);
363 accessor = _config.participants[1];
364 scheme = createParallelImplicitCouplingScheme(accessor);
365 addCouplingScheme(scheme, accessor);
366 _config = Config();
367 } else if (_config.type == VALUE_MULTI) {
369 PRECICE_CHECK(!_allowRemeshing, "Remeshing is currently incompatible with multi coupling schemes. Try using a parallel coupling scheme instead.");
370 PRECICE_INFO_IF(_allowRemeshing, "Remeshing for implicit coupling schemes is in development. Currently, the acceleration data is deleted on remeshing.");
373 "One controller per MultiCoupling needs to be defined. "
374 "Please check the <participant name=... /> tags in the <coupling-scheme:... /> of your precice-config.xml. "
375 "Make sure that at least one participant tag provides the attribute <participant name=... control=\"True\"/>.");
376 for (const std::string &accessor : _config.participants) {
378 addCouplingScheme(scheme, accessor);
379 }
380 _config = Config();
381 } else {
383 }
384 }
385}
386
388 const PtrCouplingScheme &cplScheme,
389 const std::string &participantName)
390{
391 PRECICE_TRACE(participantName);
392 if (!utils::contained(participantName, _couplingSchemes)) {
393 PRECICE_DEBUG("No coupling scheme exists for the participant");
394 // Store the new coupling scheme.
395 _couplingSchemes[participantName] = cplScheme;
396 return;
397 }
398 PRECICE_ASSERT(_couplingSchemes.count(participantName) > 0);
399
400 PRECICE_CHECK(!_allowRemeshing, "Remeshing is currently incompatible with compositional coupling schemes. If you need remeshing, try using a multi coupling scheme to compose your participants.");
401
402 // Create a composition to add the new cplScheme to
403 if (!utils::contained(participantName, _couplingSchemeCompositions)) {
404 PRECICE_DEBUG("Creating a compositional coupling scheme for the participant");
406 composition->addCouplingScheme(_couplingSchemes[participantName]);
407 _couplingSchemeCompositions[participantName] = composition.get();
408 _couplingSchemes[participantName] = std::move(composition);
409 }
410
411 PRECICE_ASSERT(_couplingSchemeCompositions.count(participantName) > 0);
412
413 // Add the new scheme to the composition
414 auto composition = _couplingSchemeCompositions.at(participantName);
415 PRECICE_CHECK(!cplScheme->isImplicitCouplingScheme() || !composition->isImplicitCouplingScheme(),
416 "You attempted to define a second implicit coupling-scheme for the participant \"{}\", which is not allowed. "
417 "Please use a multi coupling-scheme for true implicit coupling of multiple participants.",
418 participantName);
419 _couplingSchemeCompositions[participantName]->addCouplingScheme(cplScheme);
420}
421
423 const std::string &type,
424 xml::XMLTag &tag)
425{
426 PRECICE_TRACE(type);
427 addTransientLimitTags(type, tag);
428 _config.type = type;
429 //_config.name = name;
430
431 if (type == VALUE_SERIAL_EXPLICIT) {
433 addTagExchange(tag, false);
434 } else if (type == VALUE_PARALLEL_EXPLICIT) {
436 addTagExchange(tag, false);
437 } else if (type == VALUE_PARALLEL_IMPLICIT) {
439 addTagExchange(tag, true);
447 } else if (type == VALUE_MULTI) {
449 addTagExchange(tag, true);
457 } else if (type == VALUE_SERIAL_IMPLICIT) {
459 addTagExchange(tag, true);
467 } else {
468 // If wrong coupling scheme type is provided, this is already caught by the config parser. If the assertion below is triggered, it's a bug in preCICE, not wrong usage.
469 PRECICE_ASSERT(false, "Unknown coupling scheme.");
470 }
471}
472
474 const std::string &type,
475 xml::XMLTag &tag)
476{
477 using namespace xml;
478 XMLTag tagMaxTime(*this, TAG_MAX_TIME, XMLTag::OCCUR_NOT_OR_ONCE);
479 tagMaxTime.setDocumentation("Defined the end of the simulation as total time.");
480
481 XMLAttribute<double> attrValueMaxTime(ATTR_VALUE);
482 attrValueMaxTime.setDocumentation("The value of the maximum simulation time.");
483 tagMaxTime.addAttribute(attrValueMaxTime);
484 tag.addSubtag(tagMaxTime);
485
486 XMLTag tagMaxTimeWindows(*this, TAG_MAX_TIME_WINDOWS, XMLTag::OCCUR_NOT_OR_ONCE);
487 tagMaxTimeWindows.setDocumentation("Defined the end of the simulation as a total count of time windows.");
488 XMLAttribute<int> attrValueMaxTimeWindows(ATTR_VALUE);
489 attrValueMaxTimeWindows.setDocumentation("The maximum count of time windows.");
490 tagMaxTimeWindows.addAttribute(attrValueMaxTimeWindows);
491 tag.addSubtag(tagMaxTimeWindows);
492
493 XMLTag tagTimeWindowSize(*this, TAG_TIME_WINDOW_SIZE, XMLTag::OCCUR_ONCE);
494 tagTimeWindowSize.setDocumentation("Defines the size of the time window.");
495 auto attrValueTimeWindowSize = makeXMLAttribute(ATTR_VALUE, CouplingScheme::UNDEFINED_TIME_WINDOW_SIZE)
496 .setDocumentation("The maximum time window size.");
497 tagTimeWindowSize.addAttribute(attrValueTimeWindowSize);
498 if (type == VALUE_SERIAL_EXPLICIT || type == VALUE_SERIAL_IMPLICIT) {
499 // method="first-participant" is only allowed for serial coupling schemes
500 auto attrMethod = makeXMLAttribute(ATTR_METHOD, VALUE_FIXED)
502 .setDocumentation("The method used to determine the time window size. Use `fixed` to fix the time window size for the participants.");
503 tagTimeWindowSize.addAttribute(attrMethod);
504 } else {
505 tagTimeWindowSize.addAttributeHint(ATTR_METHOD, "This feature is only available for serial coupling schemes.");
506 }
507 tag.addSubtag(tagTimeWindowSize);
508}
509
511 xml::XMLTag &tag)
512{
513 using namespace xml;
514 XMLTag tagParticipants(*this, TAG_PARTICIPANTS, XMLTag::OCCUR_ONCE);
515 tagParticipants.setDocumentation("Defines the participants of the coupling scheme.");
516 XMLAttribute<std::string> attrFirst(ATTR_FIRST);
517 attrFirst.setDocumentation("First participant to run the solver.");
518 tagParticipants.addAttribute(attrFirst);
519 XMLAttribute<std::string> attrSecond(ATTR_SECOND);
520 attrSecond.setDocumentation("Second participant to run the solver.");
521 tagParticipants.addAttribute(attrSecond);
522 tag.addSubtag(tagParticipants);
523}
524
526 xml::XMLTag &tag)
527{
528 using namespace xml;
529 XMLTag tagParticipant(*this, TAG_PARTICIPANT, XMLTag::OCCUR_ONCE_OR_MORE);
530 XMLAttribute<std::string> attrName(ATTR_NAME);
531 attrName.setDocumentation("Name of the participant.");
532 tagParticipant.addAttribute(attrName);
533 XMLAttribute<bool> attrControl(ATTR_CONTROL, false);
534 attrControl.setDocumentation("Does this participant control the coupling?");
535 tagParticipant.addAttribute(attrControl);
536 tag.addSubtag(tagParticipant);
537}
538
540 xml::XMLTag &tag, bool substepsDefault)
541{
542 using namespace xml;
543 XMLTag tagExchange(*this, TAG_EXCHANGE, XMLTag::OCCUR_ONCE_OR_MORE);
544 tagExchange.setDocumentation("Defines the flow of data between meshes of participants.");
545
546 auto attrData = XMLAttribute<std::string>(ATTR_DATA).setDocumentation("The data to exchange.");
547 tagExchange.addAttribute(attrData);
548 auto attrMesh = XMLAttribute<std::string>(ATTR_MESH).setDocumentation("The mesh which uses the data.");
549 tagExchange.addAttribute(attrMesh);
550 auto participantFrom = XMLAttribute<std::string>(ATTR_FROM).setDocumentation("The participant sending the data.");
551 tagExchange.addAttribute(participantFrom);
552 auto participantTo = XMLAttribute<std::string>(ATTR_TO).setDocumentation("The participant receiving the data.");
553 tagExchange.addAttribute(participantTo);
554 auto attrInitialize = XMLAttribute<bool>(ATTR_INITIALIZE, false).setDocumentation("Should this data be initialized during initialize?");
555 tagExchange.addAttribute(attrInitialize);
556 auto attrExchangeSubsteps = XMLAttribute<bool>(ATTR_EXCHANGE_SUBSTEPS, substepsDefault).setDocumentation("Should this data exchange substeps?");
557 tagExchange.addAttribute(attrExchangeSubsteps);
558 tag.addSubtag(tagExchange);
559}
560
562 xml::XMLTag &tag)
563{
564 using namespace xml;
565 XMLTag tagConvergenceMeasure(*this, TAG_ABS_CONV_MEASURE, XMLTag::OCCUR_ARBITRARY);
566 tagConvergenceMeasure.setDocumentation(
567 "Absolute convergence criterion based on the two-norm difference of data values between iterations.\n"
568 "\\$$\\left\\lVert H(x^k) - x^k \\right\\rVert_2 < \\text{limit}\\$$");
569 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
570 XMLAttribute<double> attrLimit(ATTR_LIMIT);
571 attrLimit.setDocumentation("Limit under which the measure is considered to have converged. Must be in \\((0, 1]\\).");
572 tagConvergenceMeasure.addAttribute(attrLimit);
573 tag.addSubtag(tagConvergenceMeasure);
574}
575
577 xml::XMLTag &tag)
578{
579 using namespace xml;
580 XMLTag tagConvergenceMeasure(*this, TAG_ABS_OR_REL_CONV_MEASURE, XMLTag::OCCUR_ARBITRARY);
581 tagConvergenceMeasure.setDocumentation(
582 "Absolute or relative convergence, which is the disjunction of an absolute criterion based on the two-norm difference of data values between iterations and a relative criterion based on the relative two-norm difference of data values between iterations,i.e. convergence is reached as soon as one of the both criteria is fulfilled."
583 "\\$$\\left\\lVert H(x^k) - x^k \\right\\rVert_2 < \\text{abs-limit}\\quad\\text{or}\\quad\\frac{\\left\\lVert H(x^k) - x^k \\right\\rVert_2}{\\left\\lVert H(x^k) \\right\\rVert_2} < \\text{rel-limit} \\$$ ");
584 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
585 XMLAttribute<double> attrAbsLimit(ATTR_ABS_LIMIT);
586 attrAbsLimit.setDocumentation(R"(Absolute limit under which the measure is considered to have converged.)");
587 tagConvergenceMeasure.addAttribute(attrAbsLimit);
588 XMLAttribute<double> attrRelLimit(ATTR_REL_LIMIT);
589 attrAbsLimit.setDocumentation(R"(Relative limit under which the measure is considered to have converged. Must be in \\‍((0, 1]\\).)");
590 tagConvergenceMeasure.addAttribute(attrRelLimit);
591 tag.addSubtag(tagConvergenceMeasure);
592}
593
595 xml::XMLTag &tag)
596{
597 using namespace xml;
598 XMLTag tagConvergenceMeasure(*this, TAG_RES_REL_CONV_MEASURE,
599 XMLTag::OCCUR_ARBITRARY);
600 tagConvergenceMeasure.setDocumentation(
601 "Relative convergence criterion comparing the currently measured residual to the residual of the first iteration in the time window.\n"
602 "\\$$\\frac{\\left\\lVert H(x^k) - x^k \\right\\rVert_2}{\\left\\lVert H(x^0) - x^0 \\right\\rVert_2} < \\text{limit}\\$$");
603 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
604 XMLAttribute<double> attrLimit(ATTR_LIMIT);
605 attrLimit.setDocumentation("Limit under which the measure is considered to have converged. Must be in \\((0, 1]\\).");
606 tagConvergenceMeasure.addAttribute(attrLimit);
607 tag.addSubtag(tagConvergenceMeasure);
608}
609
611 xml::XMLTag &tag)
612{
613 using namespace xml;
614 XMLTag tagConvergenceMeasure(*this, TAG_REL_CONV_MEASURE, XMLTag::OCCUR_ARBITRARY);
615 tagConvergenceMeasure.setDocumentation(
616 "Relative convergence criterion based on the relative two-norm difference of data values between iterations.\n"
617 "\\$$\\frac{\\left\\lVert H(x^k) - x^k \\right\\rVert_2}{\\left\\lVert H(x^k) \\right\\rVert_2} < \\text{limit} \\$$");
618 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
619 XMLAttribute<double> attrLimit(ATTR_LIMIT);
620 attrLimit.setDocumentation(R"(Limit under which the measure is considered to have converged. Must be in \\‍((0, 1]\\).)");
621 tagConvergenceMeasure.addAttribute(attrLimit);
622 tag.addSubtag(tagConvergenceMeasure);
623}
624
626 xml::XMLTag &tag)
627{
628 using namespace xml;
629 auto attrData = XMLAttribute<std::string>(ATTR_DATA)
630 .setDocumentation("Data to be measured.");
631 tag.addAttribute(attrData);
632 auto attrMesh = XMLAttribute<std::string>(ATTR_MESH)
633 .setDocumentation("Mesh holding the data.");
634 tag.addAttribute(attrMesh);
635 auto attrSuffices = makeXMLAttribute(ATTR_SUFFICES, false)
636 .setDocumentation("If true, convergence of this measure is sufficient for overall convergence.");
637 tag.addAttribute(attrSuffices);
638 auto attrStrict = makeXMLAttribute(ATTR_STRICT, false)
639 .setDocumentation("If true, non-convergence of this measure ends the simulation. \"strict\" overrules \"suffices\".");
640 tag.addAttribute(attrStrict);
641}
642
644 xml::XMLTag &tag)
645{
646 using namespace xml;
647 XMLTag tagMinIterations(*this, TAG_MIN_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
648 tagMinIterations.setDocumentation("Allows to specify a minimum amount of iterations that must be performed per time window.");
649 XMLAttribute<int> attrValue(ATTR_VALUE);
650 attrValue.setDocumentation("The minimum amount of iterations.");
651 tagMinIterations.addAttribute(attrValue);
652 tag.addSubtag(tagMinIterations);
653}
654
656 xml::XMLTag &tag)
657{
658 using namespace xml;
659 XMLTag tagMaxIterations(*this, TAG_MAX_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
660 tagMaxIterations.setDocumentation("Allows to specify a maximum amount of iterations per time window.");
661 XMLAttribute<int> attrValue(ATTR_VALUE);
662 attrValue.setDocumentation("The maximum value of iterations.");
663 tagMaxIterations.addAttribute(attrValue);
664 tag.addSubtag(tagMaxIterations);
665}
666
677
679 const std::string &dataName,
680 const std::string &meshName,
681 double limit,
682 bool suffices,
683 bool strict)
684{
686 PRECICE_CHECK(math::greater(limit, 0.0),
687 "Absolute convergence limit has to be greater than zero. "
688 "Please check the <absolute-convergence-measure limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
689 "in your <coupling-scheme ... /> in the preCICE configuration file.",
690 limit, dataName, meshName);
691 ConvergenceMeasureDefintion convMeasureDef;
692 convMeasureDef.data = getData(dataName, meshName);
693 convMeasureDef.suffices = suffices;
694 convMeasureDef.strict = strict;
695 convMeasureDef.meshName = meshName;
697 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
698}
699
701 const std::string &dataName,
702 const std::string &meshName,
703 double absLimit,
704 double relLimit,
705 bool suffices,
706 bool strict)
707{
709 PRECICE_CHECK(math::greater(absLimit, 0.0),
710 "Absolute convergence limit has to be greater than zero. "
711 "Please check the <absolute-or-relative-convergence-measure abs-limit=\"{}\" rel-limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
712 "in your <coupling-scheme ... /> in the preCICE configuration file.",
713 absLimit, relLimit, dataName, meshName);
714 PRECICE_CHECK(math::greater(relLimit, 0.0) && math::greaterEquals(1.0, relLimit),
715 "Relative convergence limit has to be in ]0;1]. "
716 "Please check the <absolute-or-relative-convergence-measure abs-limit=\"{}\" rel-limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
717 "in your <coupling-scheme ... /> in the preCICE configuration file.",
718 absLimit, relLimit, dataName, meshName);
719 ConvergenceMeasureDefintion convMeasureDef;
720 convMeasureDef.data = getData(dataName, meshName);
721 convMeasureDef.suffices = suffices;
722 convMeasureDef.strict = strict;
723 convMeasureDef.meshName = meshName;
725 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
726}
727
729 const std::string &dataName,
730 const std::string &meshName,
731 double limit,
732 bool suffices,
733 bool strict)
734{
736 PRECICE_CHECK(math::greater(limit, 0.0) && math::greaterEquals(1.0, limit),
737 "Relative convergence limit has to be in ]0;1]. "
738 "Please check the <relative-convergence-measure limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
739 "in your <coupling-scheme ... /> in the preCICE configuration file.",
740 limit, dataName, meshName);
743 "The relative convergence limit=\"{}\" is close to the hard-coded numerical resolution=\"{}\" of preCICE. "
744 "This may lead to instabilities. The minimum relative convergence limit should be > \"{}\" ",
746
747 ConvergenceMeasureDefintion convMeasureDef;
748 convMeasureDef.data = getData(dataName, meshName);
749 convMeasureDef.suffices = suffices;
750 convMeasureDef.strict = strict;
751 convMeasureDef.meshName = meshName;
753 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
754}
755
757 const std::string &dataName,
758 const std::string &meshName,
759 double limit,
760 bool suffices,
761 bool strict)
762{
764 PRECICE_CHECK(math::greater(limit, 0.0) && math::greaterEquals(1.0, limit),
765 "Relative convergence limit has to be in ]0;1]. "
766 "Please check the <residul-relative-convergence-measure limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
767 "in your <coupling-scheme ... /> in the preCICE configuration file.",
768 limit, dataName, meshName);
771 "The relative convergence limit=\"{}\" is close to the hard-coded numerical resolution=\"{}\" of preCICE. "
772 "This may lead to instabilities. The minimum relative convergence limit should be > \"{}\" ",
774
775 ConvergenceMeasureDefintion convMeasureDef;
776 convMeasureDef.data = getData(dataName, meshName);
777 convMeasureDef.suffices = suffices;
778 convMeasureDef.strict = strict;
779 convMeasureDef.meshName = meshName;
781 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
782}
783
785 const std::string &dataName,
786 const std::string &meshName) const
787{
788 PRECICE_CHECK(_meshConfig->hasMeshName(meshName) && _meshConfig->getMesh(meshName)->data(dataName),
789 "Data \"{}\" used by mesh \"{}\" is not configured.", dataName, meshName);
790 const mesh::PtrMesh &mesh = _meshConfig->getMesh(meshName);
791 return mesh->data(dataName);
792}
793
795 int ID) const
796{
797 for (const mesh::PtrMesh &mesh : _meshConfig->meshes()) {
798 if (mesh->hasDataID(ID)) {
799 return mesh->data(ID);
800 }
801 }
802 return nullptr;
803}
804
806 const std::string &accessor) const
807{
808 PRECICE_TRACE(accessor);
809 m2n::PtrM2N m2n = _m2nConfig->getM2N(
812
813 for (const auto &exchange : _config.exchanges) {
814 if ((exchange.from == _config.participants[1]) && exchange.exchangeSubsteps) {
816 "Exchange of substeps is activated in the serial-explicit coupling between the second participant \"{}\" and first participant \"{}\". This is inefficient as these substeps will never be used. You can turn this off in your preCICE configuration setting substeps=\"False\" in <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" substeps=\"False\" />", exchange.from, exchange.to, exchange.data->getName(), exchange.mesh->getName(), exchange.from, exchange.to);
817 }
818 }
819
820 addDataToBeExchanged(*scheme, accessor);
821
822 return PtrCouplingScheme(scheme);
823}
824
826 const std::string &accessor) const
827{
828 PRECICE_TRACE(accessor);
830 m2n::PtrM2N m2n = _m2nConfig->getM2N(
833
834 for (const auto &exchange : _config.exchanges) {
835 if (exchange.exchangeSubsteps) {
837 "Exchange of substeps is activated in the parallel-explicit coupling between \"{}\" and \"{}\". This is inefficient as these substeps will never be used. You can turn this off in your preCICE configuration setting substeps=\"False\" in <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" substeps=\"False\" />", exchange.from, exchange.to, exchange.data->getName(), exchange.mesh->getName(), exchange.from, exchange.to);
838 }
839 }
840
841 addDataToBeExchanged(*scheme, accessor);
842
843 return PtrCouplingScheme(scheme);
844}
845
847 const std::string &accessor) const
848{
849 PRECICE_TRACE(accessor);
850
851 const auto first = _config.participants[0];
852 const auto second = _config.participants[1];
853
854 m2n::PtrM2N m2n = _m2nConfig->getM2N(
855 first, second);
857
858 addDataToBeExchanged(*scheme, accessor);
860 "No send data configured. "
861 "Use explicit scheme for one-way coupling. "
862 "Please check your <coupling-scheme ... /> and make sure that you provide at least one <exchange .../> subtag, "
863 "where from=\"{}\".",
864 accessor);
865
866 // Add convergence measures
869
870 // Set acceleration
871 setSerialAcceleration(scheme, first, second);
872
873 if (scheme->doesFirstStep() && _accelerationConfig->getAcceleration() && not _accelerationConfig->getAcceleration()->getPrimaryDataIDs().empty()) {
874 DataID dataID = *(_accelerationConfig->getAcceleration()->getPrimaryDataIDs().begin());
875 PRECICE_CHECK(not scheme->hasSendData(dataID),
876 "In case of serial coupling, acceleration can be defined for data of second participant only!");
877 }
878
879 return PtrCouplingScheme(scheme);
880}
881
883 const std::string &accessor) const
884{
885 PRECICE_TRACE(accessor);
887 m2n::PtrM2N m2n = _m2nConfig->getM2N(
890
891 addDataToBeExchanged(*scheme, accessor);
893 "No send data configured. Use explicit scheme for one-way coupling. "
894 "Please check your <coupling-scheme ... /> and make sure that you provide at least one <exchange .../> subtag, "
895 "where from=\"{}\".",
896 accessor);
897
898 // Add convergence measures
901
902 // Set acceleration
904
905 return PtrCouplingScheme(scheme);
906}
907
909 const std::string &accessor) const
910{
911 PRECICE_TRACE(accessor);
913
914 BaseCouplingScheme *scheme;
915
917 for (const std::string &participant : _config.participants) {
918 if (_m2nConfig->isM2NConfigured(accessor, participant)) {
919 m2ns[participant] = _m2nConfig->getM2N(accessor, participant);
920 }
921 }
922
923 scheme = new MultiCouplingScheme(
925
926 MultiCouplingScheme *castedScheme = dynamic_cast<MultiCouplingScheme *>(scheme);
927 PRECICE_ASSERT(castedScheme, "The dynamic cast of CouplingScheme failed.");
928 addMultiDataToBeExchanged(*castedScheme, accessor);
929
931 "No send data configured. Use explicit scheme for one-way coupling. "
932 "Please check your <coupling-scheme ... /> and make sure that you provide at least one "
933 "<exchange .../> subtag, where from=\"{}\".",
934 accessor);
935
936 // Add convergence measures
938 if (accessor == _config.controller) {
940 }
941
942 // Set acceleration
944
946 not scheme->doesFirstStep() && _accelerationConfig->getAcceleration() && _accelerationConfig->getAcceleration()->getPrimaryDataIDs().size() < 3,
947 "Due to numerical reasons, for multi coupling, the number of coupling data vectors should be at least 3, not: {}. "
948 "Please check the <data .../> subtags in your <acceleration:.../> and make sure that you have at least 3.",
949 _accelerationConfig->getAcceleration()->getPrimaryDataIDs().size());
950 return PtrCouplingScheme(scheme);
951}
952
955 const std::string &method) const
956{
957 PRECICE_TRACE(method);
958 if (method == VALUE_FIXED) {
960 } else if (method == VALUE_FIRST_PARTICIPANT) {
962 } else {
963 // We should never reach this point.
964 PRECICE_UNREACHABLE("Unknown timestepping method '{}'.", method);
965 }
966}
967
974
976{
979 "Not defining convergence measures without providing a maximum iteration limit is forbidden."
980 "Please define a convergence measure or set a maximum iteration limit using <max-iterations value=\"...\" />.");
981
982 PRECICE_INFO("No convergence measures were defined for an implicit coupling scheme. "
983 "It will always iterate the maximum amount iterations, which is {}."
984 "You may want to add a convergence measure in your <coupling-scheme:.../> in your configuration.",
986 }
987}
988
990{
991 const auto &participant = _participantConfig->getParticipant(exchange.to);
992
993 const auto &meshPtr = participant->findMesh(exchange.data->getName()); // related to https://github.com/precice/precice/issues/1694
994
995 if (meshPtr == nullptr) {
996 // Only warn, because might be valid configuration, if summation action is used. See Integration/Serial/SummationActionTwoSources.
997 PRECICE_WARN("You defined <exchange data=\"{}\" ... to=\"{}\" /> in the <coupling-scheme:... />, but <participant name=\"{}\"> has no corresponding <read-data name=\"{}\" ... />. Usually this means that there is an error in your configuration.",
998 exchange.data->getName(), exchange.to, exchange.to, exchange.data->getName());
999 return; // skip checks below
1000 }
1001
1002 const auto &readDataContext = participant->readDataContext(meshPtr->getName(), exchange.data->getName());
1003 if (readDataContext.getWaveformDegree() == 0) {
1004 PRECICE_CHECK(!exchange.exchangeSubsteps,
1005 "You configured <data:scalar/vector name=\"{}\" waveform-degree=\"{}\" />. Please deactivate exchange of substeps by setting substeps=\"false\" in the following exchange tag of your coupling scheme: <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" />. Reason: For constant interpolation no exchange of data for substeps is needed. Please consider using waveform-degree=\"1\" or higher, if you want to use subcycling.",
1006 readDataContext.getDataName(), readDataContext.getWaveformDegree(), exchange.data->getName(), exchange.mesh->getName(), exchange.from, exchange.to);
1007 } else if (readDataContext.getWaveformDegree() >= 2) {
1008 PRECICE_CHECK(exchange.exchangeSubsteps,
1009 "You configured <data:scalar/vector name=\"{}\" waveform-degree=\"{}\" />. Please activate exchange of substeps by setting substeps=\"true\" in the following exchange tag of your coupling scheme: <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" />. Reason: For higher-order interpolation exchange of data for substeps is required. If you don't want to activate exchange of additional data, please consider using waveform-degree=\"1\". Note that deactivating exchange of substep data might lead to worse results, if you use subcycling.",
1010 readDataContext.getDataName(), readDataContext.getWaveformDegree(), exchange.data->getName(), exchange.mesh->getName(), exchange.from, exchange.to);
1011 } else { // For first degree there is no restriction for exchange of substeps
1012 PRECICE_ASSERT(readDataContext.getWaveformDegree() == 1);
1013 }
1014}
1015
1017 BiCouplingScheme &scheme,
1018 const std::string &accessor) const
1019{
1020 PRECICE_TRACE();
1021 for (const Config::Exchange &exchange : _config.exchanges) {
1022 const std::string &from = exchange.from;
1023 const std::string &to = exchange.to;
1024 const std::string &dataName = exchange.data->getName();
1025 const std::string &meshName = exchange.mesh->getName();
1026
1027 PRECICE_CHECK(to != from,
1028 "You cannot define an exchange from and to the same participant. "
1029 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1030 dataName, meshName, from, to);
1031
1033 "Participant \"{}\" is not configured for coupling scheme. "
1034 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1035 from, dataName, meshName, from, to);
1036
1038 "Participant \"{}\" is not configured for coupling scheme. "
1039 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1040 to, dataName, meshName, from, to);
1041
1042 const bool requiresInitialization = exchange.requiresInitialization;
1043
1045 !(requiresInitialization && _participantConfig->getParticipant(from)->isDirectAccessAllowed(exchange.mesh->getName())),
1046 "Participant \"{}\" cannot initialize data of the directly-accessed mesh \"{}\" from the participant\"{}\". "
1047 "Either disable the initialization in the <exchange /> tag or use a locally provided mesh instead.",
1048 from, meshName, to);
1049
1050 const bool exchangeSubsteps = exchange.exchangeSubsteps;
1051
1052 if (from == accessor) {
1053 scheme.addDataToSend(exchange.data, exchange.mesh, requiresInitialization, exchangeSubsteps);
1054 } else if (to == accessor) {
1056 scheme.addDataToReceive(exchange.data, exchange.mesh, requiresInitialization, exchangeSubsteps);
1057 } else {
1059 }
1060 }
1062}
1063
1065 MultiCouplingScheme &scheme,
1066 const std::string &accessor) const
1067{
1068 PRECICE_TRACE();
1069 for (const Config::Exchange &exchange : _config.exchanges) {
1070 const std::string &from = exchange.from;
1071 const std::string &to = exchange.to;
1072 const std::string &dataName = exchange.data->getName();
1073 const std::string &meshName = exchange.mesh->getName();
1074
1075 PRECICE_CHECK(to != from,
1076 "You cannot define an exchange from and to the same participant. "
1077 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1078 dataName, meshName, from, to);
1079
1081 "Participant \"{}\" is not configured for coupling scheme",
1082 from);
1083
1085 "Participant \"{}\" is not configured for coupling scheme", to);
1086
1087 const bool initialize = exchange.requiresInitialization;
1088 const bool exchangeSubsteps = exchange.exchangeSubsteps;
1089
1090 if (from == accessor) {
1091 scheme.addDataToSend(exchange.data, exchange.mesh, initialize, exchangeSubsteps, to);
1092 } else if (to == accessor) {
1093 scheme.addDataToReceive(exchange.data, exchange.mesh, initialize, exchangeSubsteps, from);
1094 }
1095 }
1097}
1098
1100 DataID dataID) const
1101{
1102 const auto match = std::find_if(_config.exchanges.begin(),
1103 _config.exchanges.end(),
1104 [dataID](const Config::Exchange &exchange) { return exchange.data->getID() == dataID; });
1105 if (match != _config.exchanges.end()) {
1106 return;
1107 }
1108
1109 // Data is not being exchanged
1110 std::string dataName = "";
1111 auto dataptr = findDataByID(dataID);
1112 if (dataptr) {
1113 dataName = dataptr->getName();
1114 }
1115
1116 PRECICE_ERROR("You need to exchange every data that you use for convergence measures and/or the iteration acceleration. "
1117 "Data \"{}\" is currently not exchanged over the respective mesh on which it is used for convergence measures and/or iteration acceleration. "
1118 "Please check the <exchange ... /> and <...-convergence-measure ... /> tags in the <coupling-scheme:... /> of your precice-config.xml.",
1119 dataName);
1120}
1121
1123 int dataID,
1124 const std::string &first,
1125 const std::string &second) const
1126{
1127 checkIfDataIsExchanged(dataID);
1128 const auto match = std::find_if(_config.exchanges.begin(),
1129 _config.exchanges.end(),
1130 [dataID](const Config::Exchange &exchange) { return exchange.data->getID() == dataID; });
1131 PRECICE_ASSERT(match != _config.exchanges.end());
1132 const auto &exchange = *match;
1133
1134 // In serial implicit cplschemes, data is only accelerated on the second participant.
1135 if (second == exchange.from) {
1136 return;
1137 }
1138
1139 std::string dataName = "";
1140 auto dataptr = findDataByID(dataID);
1141 if (dataptr) {
1142 dataName = dataptr->getName();
1143 }
1144
1146 "You configured acceleration data \"{}\" in the serial implicit coupling scheme between participants \"{}\" and \"{}\". "
1147 "For serial implicit coupling schemes, only data exchanged from the second to the first participant can be used for acceleration. "
1148 "Here, from \"{}\" to \"{}\". "
1149 "However, you configured data \"{}\" for acceleration, which is exchanged from \"{}\" to \"{}\". "
1150 "Please remove this acceleration data tag or switch to a parallel implicit coupling scheme.",
1151 dataName, first, second, second, first, dataName, first, second);
1152}
1153
1155 BaseCouplingScheme *scheme,
1156 const std::string &participant,
1157 const std::vector<ConvergenceMeasureDefintion> &convergenceMeasureDefinitions) const
1158{
1159 for (auto &elem : convergenceMeasureDefinitions) {
1160 _meshConfig->addNeededMesh(participant, elem.meshName);
1161 checkIfDataIsExchanged(elem.data->getID());
1162 scheme->addConvergenceMeasure(elem.data->getID(), elem.suffices, elem.strict, elem.measure);
1163 }
1164}
1165
1167 BaseCouplingScheme *scheme,
1168 const std::string &first,
1169 const std::string &second) const
1170{
1171 if (_accelerationConfig->getAcceleration().get() != nullptr) {
1172 for (std::string &neededMesh : _accelerationConfig->getNeededMeshes()) {
1173 _meshConfig->addNeededMesh(second, neededMesh);
1174 }
1175 for (const DataID dataID : _accelerationConfig->getAcceleration()->getPrimaryDataIDs()) {
1176 checkSerialImplicitAccelerationData(dataID, first, second);
1177 }
1178 scheme->setAcceleration(_accelerationConfig->getAcceleration());
1179 }
1180}
1181
1183 BaseCouplingScheme *scheme,
1184 const std::string &participant) const
1185{
1186 if (_accelerationConfig->getAcceleration().get() != nullptr) {
1187 for (std::string &neededMesh : _accelerationConfig->getNeededMeshes()) {
1188 _meshConfig->addNeededMesh(participant, neededMesh);
1189 }
1190 for (const DataID dataID : _accelerationConfig->getAcceleration()->getPrimaryDataIDs()) {
1191 checkIfDataIsExchanged(dataID);
1192 }
1193 scheme->setAcceleration(_accelerationConfig->getAcceleration());
1194
1196 dynamic_cast<acceleration::AitkenAcceleration *>(_accelerationConfig->getAcceleration().get()) != nullptr,
1197 "You configured participant \"{}\" in a parallel-implicit coupling scheme with \"Aitken\" "
1198 "acceleration, which is known to perform bad in parallel coupling schemes. "
1199 "See https://precice.org/configuration-acceleration.html#dynamic-aitken-under-relaxation for details."
1200 "Consider switching to a serial-implicit coupling scheme or changing the acceleration method.",
1201 participant);
1202 }
1203}
1204
1206 bool allowed)
1207{
1208 _allowRemeshing = allowed;
1209}
1210
1211} // namespace precice::cplscheme
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:16
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:18
#define PRECICE_WARN(...)
Definition LogMacros.hpp:12
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_INFO_IF(condition,...)
Definition LogMacros.hpp:25
#define PRECICE_INFO(...)
Definition LogMacros.hpp:14
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
#define PRECICE_UNREACHABLE(...)
Definition assertion.hpp:93
T begin(T... args)
Abstract base class for standard coupling schemes.
void setAcceleration(const acceleration::PtrAcceleration &acceleration)
Set an acceleration technique.
void addConvergenceMeasure(int dataID, bool suffices, bool strict, impl::PtrConvergenceMeasure measure)
Adds a measure to determine the convergence of coupling iterations.
bool doesFirstStep() const
Getter for _doesFirstStep.
Abstract base class for coupling schemes with two participants.
void addDataToReceive(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps)
Adds data to be received on data exchange.
void addDataToSend(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps)
Adds data to be sent on data exchange and possibly be modified during coupling iterations.
void determineInitialDataExchange() override
Determines which data is initialized and therefore has to be exchanged during initialize.
PtrCouplingScheme createSerialExplicitCouplingScheme(const std::string &accessor) const
void checkIterationLimits() const
Helper function to check iteration limits in conjunction with convergence measures.
const PtrCouplingScheme & getCouplingScheme(const std::string &participantName) const
Returns the configured coupling scheme.
precice::config::PtrParticipantConfiguration _participantConfig
void addMultiDataToBeExchanged(MultiCouplingScheme &scheme, const std::string &accessor) const
Adds configured exchange data to be sent or received to scheme. Only used specifically for MultiCoupl...
PtrCouplingScheme createSerialImplicitCouplingScheme(const std::string &accessor) const
void xmlEndTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
Callback method required when using xml::XMLTag.
void xmlTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag) override
Callback method required when using xml::XMLTag.
void setSerialAcceleration(BaseCouplingScheme *scheme, const std::string &first, const std::string &second) const
void addResidualRelativeConvergenceMeasure(const std::string &dataName, const std::string &meshName, double limit, bool suffices, bool strict)
mesh::PtrData getData(const std::string &dataName, const std::string &meshName) const
void addAbsoluteConvergenceMeasure(const std::string &dataName, const std::string &meshName, double limit, bool suffices, bool strict)
void addDataToBeExchanged(BiCouplingScheme &scheme, const std::string &accessor) const
Adds configured exchange data to be sent or received to scheme.
struct precice::cplscheme::CouplingSchemeConfiguration::Config _config
void checkSubstepExchangeWaveformDegree(const Config::Exchange &exchange) const
Helper function to check that waveform-degree and substep exchange are compatible.
constants::TimesteppingMethod getTimesteppingMethod(const std::string &method) const
void checkSerialImplicitAccelerationData(DataID dataID, const std::string &first, const std::string &second) const
acceleration::PtrAccelerationConfiguration _accelerationConfig
std::map< std::string, PtrCouplingScheme > _couplingSchemes
Map from participant name to coupling scheme (composition).
CouplingSchemeConfiguration(xml::XMLTag &parent, mesh::PtrMeshConfiguration meshConfig, m2n::M2NConfiguration::SharedPointer m2nConfig, config::PtrParticipantConfiguration participantConfig)
Constructor.
void addTagExchange(xml::XMLTag &tag, bool substepsDefault)
void addTypespecifcSubtags(const std::string &type, xml::XMLTag &tag)
PtrCouplingScheme createMultiCouplingScheme(const std::string &accessor) const
void addCouplingScheme(const PtrCouplingScheme &cplScheme, const std::string &participantName)
Adds a manually configured coupling scheme for a participant.
bool hasCouplingScheme(const std::string &participantName) const
Check, if a coupling scheme is configured for a participant.
void updateConfigForImplicitCoupling()
Helper to update some configs which may have a different meaning in implicit coupling.
void addAbsoluteOrRelativeConvergenceMeasure(const std::string &dataName, const std::string &meshName, double absLimit, double relLimit, bool suffices, bool strict)
PtrCouplingScheme createParallelImplicitCouplingScheme(const std::string &accessor) const
PtrCouplingScheme createParallelExplicitCouplingScheme(const std::string &accessor) const
void addTransientLimitTags(const std::string &type, xml::XMLTag &tag)
void addRelativeConvergenceMeasure(const std::string &dataName, const std::string &meshName, double limit, bool suffices, bool strict)
std::map< std::string, CompositionalCouplingScheme * > _couplingSchemeCompositions
If a participant has more than one coupling scheme, a composition is created.
void addConvergenceMeasures(BaseCouplingScheme *scheme, const std::string &participant, const std::vector< ConvergenceMeasureDefintion > &convergenceMeasureDefinitions) const
void setParallelAcceleration(BaseCouplingScheme *scheme, const std::string &participant) const
static const int INFINITE_MAX_ITERATIONS
To be used, when the number of max iterations is infinite (for implicit coupling).
static const int UNDEFINED_MAX_ITERATIONS
To be used, when the number of max iterations is not defined (for explicit coupling).
static const double UNDEFINED_TIME_WINDOW_SIZE
To be used, when the time window size is determined dynamically during the coupling.
A coupling scheme with multiple participants.
void addDataToSend(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps, const std::string &to)
Adds data to be sent on data exchange and possibly be modified during coupling iterations.
void determineInitialDataExchange() override
Determines which data is initialized and therefore has to be exchanged during initialize.
void addDataToReceive(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps, const std::string &from)
Adds data to be received on data exchange.
Coupling scheme for parallel coupling, i.e. simultaneous execution of two coupled participants.
Coupling scheme for serial coupling, i.e. staggered execution of two coupled participants.
Represents an XML tag to be configured automatically.
Definition XMLTag.hpp:28
const std::string & getNamespace() const
Returns xml namespace.
Definition XMLTag.hpp:159
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 compare(T... args)
T end(T... args)
T find(T... args)
T get(T... args)
vector< double > getData()
Definition mainA.cpp:19
T make_unique(T... args)
contains implementations of coupling schemes for coupled simulations.
std::shared_ptr< CouplingScheme > PtrCouplingScheme
constexpr double NUMERICAL_ZERO_DIFFERENCE
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greaterEquals(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greater(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
bool contained(const ELEMENT_T &element, const std::vector< ELEMENT_T > &vec)
Returns true, if given element is in vector, otherwise false.
Definition Helpers.hpp:38
int DataID
Definition Types.hpp:25
STL namespace.
T push_back(T... args)
std::vector< ConvergenceMeasureDefintion > convergenceMeasureDefinitions
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21