preCICE v3.1.2
Loading...
Searching...
No Matches
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.
219 "The minimal time window size supported by preCICE is {}. "
220 "Please check the <time-window-size value=\"{}\" /> tag "
221 "in the <coupling-scheme:...> of your precice-config.xml and pick an appropriate time window size.",
223 } else {
226 "Time window size value has to be equal to -1 (default), if method=\"first-participant\" is used. "
227 "Please check the <time-window-size value=\"{}\" method=\"{}\" /> "
228 "tag in the <coupling-scheme:...> of your precice-config.xml",
230 }
231 } else if (tag.getName() == TAG_ABS_CONV_MEASURE) {
232 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
233 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
234 double limit = tag.getDoubleAttributeValue(ATTR_LIMIT);
235 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
236 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
238 addAbsoluteConvergenceMeasure(dataName, meshName, limit, suffices, strict);
239 } else if (tag.getName() == TAG_ABS_OR_REL_CONV_MEASURE) {
240 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
241 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
242 double absLimit = tag.getDoubleAttributeValue(ATTR_ABS_LIMIT);
243 double relLimit = tag.getDoubleAttributeValue(ATTR_REL_LIMIT);
244 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
245 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
247 addAbsoluteOrRelativeConvergenceMeasure(dataName, meshName, absLimit, relLimit, suffices, strict);
248 } else if (tag.getName() == TAG_REL_CONV_MEASURE) {
249 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
250 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
251 double limit = tag.getDoubleAttributeValue(ATTR_LIMIT);
252 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
253 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
255 addRelativeConvergenceMeasure(dataName, meshName, limit, suffices, strict);
256 } else if (tag.getName() == TAG_RES_REL_CONV_MEASURE) {
257 const std::string &dataName = tag.getStringAttributeValue(ATTR_DATA);
258 const std::string &meshName = tag.getStringAttributeValue(ATTR_MESH);
259 double limit = tag.getDoubleAttributeValue(ATTR_LIMIT);
260 bool suffices = tag.getBooleanAttributeValue(ATTR_SUFFICES);
261 bool strict = tag.getBooleanAttributeValue(ATTR_STRICT);
263 addResidualRelativeConvergenceMeasure(dataName, meshName, limit, suffices, strict);
264 } else if (tag.getName() == TAG_EXCHANGE) {
267 std::string nameParticipantFrom = tag.getStringAttributeValue(ATTR_FROM);
268 std::string nameParticipantTo = tag.getStringAttributeValue(ATTR_TO);
269 bool initialize = tag.getBooleanAttributeValue(ATTR_INITIALIZE);
270 bool exchangeSubsteps = tag.getBooleanAttributeValue(ATTR_EXCHANGE_SUBSTEPS);
271
272 PRECICE_CHECK(_meshConfig->hasMeshName(nameMesh) && _meshConfig->getMesh(nameMesh)->hasDataName(nameData),
273 "Mesh \"{}\" with data \"{}\" not defined. "
274 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> "
275 "tag in the <coupling-scheme:... /> of your precice-config.xml.",
276 nameMesh, nameData, nameData, nameMesh, nameParticipantFrom, nameParticipantTo);
277
278 mesh::PtrMesh exchangeMesh = _meshConfig->getMesh(nameMesh);
279 PRECICE_ASSERT(exchangeMesh);
280 mesh::PtrData exchangeData = exchangeMesh->data(nameData);
281 PRECICE_ASSERT(exchangeData);
282
283 Config::Exchange newExchange{exchangeData, exchangeMesh, nameParticipantFrom, nameParticipantTo, initialize, exchangeSubsteps};
284 PRECICE_CHECK(!_config.hasExchange(newExchange),
285 R"(Data "{}" of mesh "{}" cannot be exchanged multiple times between participants "{}" and "{}". Please remove one of the exchange tags.)",
286 nameData, nameMesh, nameParticipantFrom, nameParticipantTo);
287
288 _meshConfig->addNeededMesh(nameParticipantFrom, nameMesh);
289 _meshConfig->addNeededMesh(nameParticipantTo, nameMesh);
290 _config.exchanges.emplace_back(std::move(newExchange));
291 } else if (tag.getName() == TAG_MIN_ITERATIONS) {
295 "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.",
297 } else if (tag.getName() == TAG_MAX_ITERATIONS) {
301 "Maximal iteration limit has to be larger than zero. "
302 "Please check the <max-iterations value=\"{0}\" /> subtag in the <coupling-scheme:... /> of your precice-config.xml. "
303 "To disable the iteration limit, remove the <max-iterations value=\"{0}\" /> subtag.",
305 }
306
307 // Additional consistency checks
310 "Maximum iteration limit {1} has to be larger or equal than the minimum iteration limit {0}. "
311 "Please check the <min-iterations value = \"{0}\" /> and <max-iterations value = \"{1}\" /> subtags in the <coupling-scheme:... /> of your precice-config.xml.",
314 }
315}
316
318 const xml::ConfigurationContext &context,
319 xml::XMLTag & tag)
320{
322 if (tag.getNamespace() == TAG) {
324 std::string accessor(_config.participants[0]);
326 addCouplingScheme(scheme, accessor);
327 //_couplingSchemes[accessor] = scheme;
328 accessor = _config.participants[1];
329 scheme = createSerialExplicitCouplingScheme(accessor);
330 addCouplingScheme(scheme, accessor);
331 //_couplingSchemes[accessor] = scheme;
332 _config = Config();
333 } else if (_config.type == VALUE_PARALLEL_EXPLICIT) {
334 std::string accessor(_config.participants[0]);
336 addCouplingScheme(scheme, accessor);
337 //_couplingSchemes[accessor] = scheme;
338 accessor = _config.participants[1];
339 scheme = createParallelExplicitCouplingScheme(accessor);
340 addCouplingScheme(scheme, accessor);
341 //_couplingSchemes[accessor] = scheme;
342 _config = Config();
343 } else if (_config.type == VALUE_SERIAL_IMPLICIT) {
345 std::string accessor(_config.participants[0]);
347 addCouplingScheme(scheme, accessor);
348 //_couplingSchemes[accessor] = scheme;
349 accessor = _config.participants[1];
350 scheme = createSerialImplicitCouplingScheme(accessor);
351 addCouplingScheme(scheme, accessor);
352 //_couplingSchemes[accessor] = scheme;
353 _config = Config();
354 } else if (_config.type == VALUE_PARALLEL_IMPLICIT) {
356 std::string accessor(_config.participants[0]);
358 addCouplingScheme(scheme, accessor);
359 accessor = _config.participants[1];
360 scheme = createParallelImplicitCouplingScheme(accessor);
361 addCouplingScheme(scheme, accessor);
362 _config = Config();
363 } else if (_config.type == VALUE_MULTI) {
366 "One controller per MultiCoupling needs to be defined. "
367 "Please check the <participant name=... /> tags in the <coupling-scheme:... /> of your precice-config.xml. "
368 "Make sure that at least one participant tag provides the attribute <participant name=... control=\"True\"/>.");
369 for (const std::string &accessor : _config.participants) {
371 addCouplingScheme(scheme, accessor);
372 }
373 _config = Config();
374 } else {
376 }
377 }
378}
379
381 const PtrCouplingScheme &cplScheme,
382 const std::string & participantName)
383{
384 PRECICE_TRACE(participantName);
385 if (!utils::contained(participantName, _couplingSchemes)) {
386 PRECICE_DEBUG("No coupling scheme exists for the participant");
387 // Store the new coupling scheme.
388 _couplingSchemes[participantName] = cplScheme;
389 return;
390 }
391 PRECICE_ASSERT(_couplingSchemes.count(participantName) > 0);
392
393 // Create a composition to add the new cplScheme to
394 if (!utils::contained(participantName, _couplingSchemeCompositions)) {
395 PRECICE_DEBUG("Creating a compositional coupling scheme for the participant");
396 auto composition = std::make_shared<CompositionalCouplingScheme>();
397 composition->addCouplingScheme(_couplingSchemes[participantName]);
398 _couplingSchemeCompositions[participantName] = composition.get();
399 _couplingSchemes[participantName] = std::move(composition);
400 }
401
402 PRECICE_ASSERT(_couplingSchemeCompositions.count(participantName) > 0);
403
404 // Add the new scheme to the composition
405 auto composition = _couplingSchemeCompositions.at(participantName);
406 PRECICE_CHECK(!cplScheme->isImplicitCouplingScheme() || !composition->isImplicitCouplingScheme(),
407 "You attempted to define a second implicit coupling-scheme for the participant \"{}\", which is not allowed. "
408 "Please use a multi coupling-scheme for true implicit coupling of multiple participants.",
409 participantName);
410 _couplingSchemeCompositions[participantName]->addCouplingScheme(cplScheme);
411}
412
414 const std::string &type,
415 xml::XMLTag & tag)
416{
417 PRECICE_TRACE(type);
418 addTransientLimitTags(type, tag);
419 _config.type = type;
420 //_config.name = name;
421
422 if (type == VALUE_SERIAL_EXPLICIT) {
424 addTagExchange(tag);
425 } else if (type == VALUE_PARALLEL_EXPLICIT) {
427 addTagExchange(tag);
428 } else if (type == VALUE_PARALLEL_IMPLICIT) {
430 addTagExchange(tag);
438 } else if (type == VALUE_MULTI) {
440 addTagExchange(tag);
448 } else if (type == VALUE_SERIAL_IMPLICIT) {
450 addTagExchange(tag);
458 } else {
459 // 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.
460 PRECICE_ASSERT(false, "Unknown coupling scheme.");
461 }
462}
463
465 const std::string &type,
466 xml::XMLTag & tag)
467{
468 using namespace xml;
469 XMLTag tagMaxTime(*this, TAG_MAX_TIME, XMLTag::OCCUR_NOT_OR_ONCE);
470 tagMaxTime.setDocumentation("Defined the end of the simulation as total time.");
471
472 XMLAttribute<double> attrValueMaxTime(ATTR_VALUE);
473 attrValueMaxTime.setDocumentation("The value of the maximum simulation time.");
474 tagMaxTime.addAttribute(attrValueMaxTime);
475 tag.addSubtag(tagMaxTime);
476
477 XMLTag tagMaxTimeWindows(*this, TAG_MAX_TIME_WINDOWS, XMLTag::OCCUR_NOT_OR_ONCE);
478 tagMaxTimeWindows.setDocumentation("Defined the end of the simulation as a total count of time windows.");
479 XMLAttribute<int> attrValueMaxTimeWindows(ATTR_VALUE);
480 attrValueMaxTimeWindows.setDocumentation("The maximum count of time windows.");
481 tagMaxTimeWindows.addAttribute(attrValueMaxTimeWindows);
482 tag.addSubtag(tagMaxTimeWindows);
483
484 XMLTag tagTimeWindowSize(*this, TAG_TIME_WINDOW_SIZE, XMLTag::OCCUR_ONCE);
485 tagTimeWindowSize.setDocumentation("Defines the size of the time window.");
486 auto attrValueTimeWindowSize = makeXMLAttribute(ATTR_VALUE, CouplingScheme::UNDEFINED_TIME_WINDOW_SIZE)
487 .setDocumentation("The maximum time window size.");
488 tagTimeWindowSize.addAttribute(attrValueTimeWindowSize);
489 if (type == VALUE_SERIAL_EXPLICIT || type == VALUE_SERIAL_IMPLICIT) {
490 // method="first-participant" is only allowed for serial coupling schemes
491 auto attrMethod = makeXMLAttribute(ATTR_METHOD, VALUE_FIXED)
493 .setDocumentation("The method used to determine the time window size. Use `fixed` to fix the time window size for the participants.");
494 tagTimeWindowSize.addAttribute(attrMethod);
495 } else {
496 tagTimeWindowSize.addAttributeHint(ATTR_METHOD, "This feature is only available for serial coupling schemes.");
497 }
498 tag.addSubtag(tagTimeWindowSize);
499}
500
502 xml::XMLTag &tag)
503{
504 using namespace xml;
505 XMLTag tagParticipants(*this, TAG_PARTICIPANTS, XMLTag::OCCUR_ONCE);
506 tagParticipants.setDocumentation("Defines the participants of the coupling scheme.");
507 XMLAttribute<std::string> attrFirst(ATTR_FIRST);
508 attrFirst.setDocumentation("First participant to run the solver.");
509 tagParticipants.addAttribute(attrFirst);
510 XMLAttribute<std::string> attrSecond(ATTR_SECOND);
511 attrSecond.setDocumentation("Second participant to run the solver.");
512 tagParticipants.addAttribute(attrSecond);
513 tag.addSubtag(tagParticipants);
514}
515
517 xml::XMLTag &tag)
518{
519 using namespace xml;
520 XMLTag tagParticipant(*this, TAG_PARTICIPANT, XMLTag::OCCUR_ONCE_OR_MORE);
521 XMLAttribute<std::string> attrName(ATTR_NAME);
522 attrName.setDocumentation("Name of the participant.");
523 tagParticipant.addAttribute(attrName);
524 XMLAttribute<bool> attrControl(ATTR_CONTROL, false);
525 attrControl.setDocumentation("Does this participant control the coupling?");
526 tagParticipant.addAttribute(attrControl);
527 tag.addSubtag(tagParticipant);
528}
529
531 xml::XMLTag &tag)
532{
533 using namespace xml;
534 XMLTag tagExchange(*this, TAG_EXCHANGE, XMLTag::OCCUR_ONCE_OR_MORE);
535 tagExchange.setDocumentation("Defines the flow of data between meshes of participants.");
536
537 auto attrData = XMLAttribute<std::string>(ATTR_DATA).setDocumentation("The data to exchange.");
538 tagExchange.addAttribute(attrData);
539 auto attrMesh = XMLAttribute<std::string>(ATTR_MESH).setDocumentation("The mesh which uses the data.");
540 tagExchange.addAttribute(attrMesh);
541 auto participantFrom = XMLAttribute<std::string>(ATTR_FROM).setDocumentation("The participant sending the data.");
542 tagExchange.addAttribute(participantFrom);
543 auto participantTo = XMLAttribute<std::string>(ATTR_TO).setDocumentation("The participant receiving the data.");
544 tagExchange.addAttribute(participantTo);
545 auto attrInitialize = XMLAttribute<bool>(ATTR_INITIALIZE, false).setDocumentation("Should this data be initialized during initialize?");
546 tagExchange.addAttribute(attrInitialize);
547 auto attrExchangeSubsteps = XMLAttribute<bool>(ATTR_EXCHANGE_SUBSTEPS, false).setDocumentation("Should this data exchange substeps?");
548 tagExchange.addAttribute(attrExchangeSubsteps);
549 tag.addSubtag(tagExchange);
550}
551
553 xml::XMLTag &tag)
554{
555 using namespace xml;
556 XMLTag tagConvergenceMeasure(*this, TAG_ABS_CONV_MEASURE, XMLTag::OCCUR_ARBITRARY);
557 tagConvergenceMeasure.setDocumentation(
558 "Absolute convergence criterion based on the two-norm difference of data values between iterations.\n"
559 "\\$$\\left\\lVert H(x^k) - x^k \\right\\rVert_2 < \\text{limit}\\$$");
560 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
561 XMLAttribute<double> attrLimit(ATTR_LIMIT);
562 attrLimit.setDocumentation("Limit under which the measure is considered to have converged. Must be in \\((0, 1]\\).");
563 tagConvergenceMeasure.addAttribute(attrLimit);
564 tag.addSubtag(tagConvergenceMeasure);
565}
566
568 xml::XMLTag &tag)
569{
570 using namespace xml;
571 XMLTag tagConvergenceMeasure(*this, TAG_ABS_OR_REL_CONV_MEASURE, XMLTag::OCCUR_ARBITRARY);
572 tagConvergenceMeasure.setDocumentation(
573 "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."
574 "\\$$\\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} \\$$ ");
575 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
576 XMLAttribute<double> attrAbsLimit(ATTR_ABS_LIMIT);
577 attrAbsLimit.setDocumentation(R"(Absolute limit under which the measure is considered to have converged.)");
578 tagConvergenceMeasure.addAttribute(attrAbsLimit);
579 XMLAttribute<double> attrRelLimit(ATTR_REL_LIMIT);
580 attrAbsLimit.setDocumentation(R"(Relative limit under which the measure is considered to have converged. Must be in \\‍((0, 1]\\).)");
581 tagConvergenceMeasure.addAttribute(attrRelLimit);
582 tag.addSubtag(tagConvergenceMeasure);
583}
584
586 xml::XMLTag &tag)
587{
588 using namespace xml;
589 XMLTag tagConvergenceMeasure(*this, TAG_RES_REL_CONV_MEASURE,
590 XMLTag::OCCUR_ARBITRARY);
591 tagConvergenceMeasure.setDocumentation(
592 "Relative convergence criterion comparing the currently measured residual to the residual of the first iteration in the time window.\n"
593 "\\$$\\frac{\\left\\lVert H(x^k) - x^k \\right\\rVert_2}{\\left\\lVert H(x^0) - x^0 \\right\\rVert_2} < \\text{limit}\\$$");
594 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
595 XMLAttribute<double> attrLimit(ATTR_LIMIT);
596 attrLimit.setDocumentation("Limit under which the measure is considered to have converged. Must be in \\((0, 1]\\).");
597 tagConvergenceMeasure.addAttribute(attrLimit);
598 tag.addSubtag(tagConvergenceMeasure);
599}
600
602 xml::XMLTag &tag)
603{
604 using namespace xml;
605 XMLTag tagConvergenceMeasure(*this, TAG_REL_CONV_MEASURE, XMLTag::OCCUR_ARBITRARY);
606 tagConvergenceMeasure.setDocumentation(
607 "Relative convergence criterion based on the relative two-norm difference of data values between iterations.\n"
608 "\\$$\\frac{\\left\\lVert H(x^k) - x^k \\right\\rVert_2}{\\left\\lVert H(x^k) \\right\\rVert_2} < \\text{limit} \\$$");
609 addBaseAttributesTagConvergenceMeasure(tagConvergenceMeasure);
610 XMLAttribute<double> attrLimit(ATTR_LIMIT);
611 attrLimit.setDocumentation(R"(Limit under which the measure is considered to have converged. Must be in \\‍((0, 1]\\).)");
612 tagConvergenceMeasure.addAttribute(attrLimit);
613 tag.addSubtag(tagConvergenceMeasure);
614}
615
617 xml::XMLTag &tag)
618{
619 using namespace xml;
620 auto attrData = XMLAttribute<std::string>(ATTR_DATA)
621 .setDocumentation("Data to be measured.");
622 tag.addAttribute(attrData);
623 auto attrMesh = XMLAttribute<std::string>(ATTR_MESH)
624 .setDocumentation("Mesh holding the data.");
625 tag.addAttribute(attrMesh);
626 auto attrSuffices = makeXMLAttribute(ATTR_SUFFICES, false)
627 .setDocumentation("If true, convergence of this measure is sufficient for overall convergence.");
628 tag.addAttribute(attrSuffices);
629 auto attrStrict = makeXMLAttribute(ATTR_STRICT, false)
630 .setDocumentation("If true, non-convergence of this measure ends the simulation. \"strict\" overrules \"suffices\".");
631 tag.addAttribute(attrStrict);
632}
633
635 xml::XMLTag &tag)
636{
637 using namespace xml;
638 XMLTag tagMinIterations(*this, TAG_MIN_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
639 tagMinIterations.setDocumentation("Allows to specify a minimum amount of iterations that must be performed per time window.");
640 XMLAttribute<int> attrValue(ATTR_VALUE);
641 attrValue.setDocumentation("The minimum amount of iterations.");
642 tagMinIterations.addAttribute(attrValue);
643 tag.addSubtag(tagMinIterations);
644}
645
647 xml::XMLTag &tag)
648{
649 using namespace xml;
650 XMLTag tagMaxIterations(*this, TAG_MAX_ITERATIONS, XMLTag::OCCUR_NOT_OR_ONCE);
651 tagMaxIterations.setDocumentation("Allows to specify a maximum amount of iterations per time window.");
652 XMLAttribute<int> attrValue(ATTR_VALUE);
653 attrValue.setDocumentation("The maximum value of iterations.");
654 tagMaxIterations.addAttribute(attrValue);
655 tag.addSubtag(tagMaxIterations);
656}
657
659 xml::XMLTag &tag)
660{
662 if (!_accelerationConfig) {
663 _accelerationConfig = std::make_shared<acceleration::AccelerationConfiguration>(
665 }
666 _accelerationConfig->connectTags(tag);
667}
668
670 const std::string &dataName,
671 const std::string &meshName,
672 double limit,
673 bool suffices,
674 bool strict)
675{
677 PRECICE_CHECK(math::greater(limit, 0.0),
678 "Absolute convergence limit has to be greater than zero. "
679 "Please check the <absolute-convergence-measure limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
680 "in your <coupling-scheme ... /> in the preCICE configuration file.",
681 limit, dataName, meshName);
683 ConvergenceMeasureDefintion convMeasureDef;
684 convMeasureDef.data = getData(dataName, meshName);
685 convMeasureDef.suffices = suffices;
686 convMeasureDef.strict = strict;
687 convMeasureDef.meshName = meshName;
688 convMeasureDef.measure = std::move(measure);
689 convMeasureDef.doesLogging = true;
690 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
691}
692
694 const std::string &dataName,
695 const std::string &meshName,
696 double absLimit,
697 double relLimit,
698 bool suffices,
699 bool strict)
700{
702 PRECICE_CHECK(math::greater(absLimit, 0.0),
703 "Absolute convergence limit has to be greater than zero. "
704 "Please check the <absolute-or-relative-convergence-measure abs-limit=\"{}\" rel-limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
705 "in your <coupling-scheme ... /> in the preCICE configuration file.",
706 absLimit, relLimit, dataName, meshName);
707 PRECICE_CHECK(math::greater(relLimit, 0.0) && math::greaterEquals(1.0, relLimit),
708 "Relative convergence limit has to be in ]0;1]. "
709 "Please check the <absolute-or-relative-convergence-measure abs-limit=\"{}\" rel-limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
710 "in your <coupling-scheme ... /> in the preCICE configuration file.",
711 absLimit, relLimit, dataName, meshName);
713 ConvergenceMeasureDefintion convMeasureDef;
714 convMeasureDef.data = getData(dataName, meshName);
715 convMeasureDef.suffices = suffices;
716 convMeasureDef.strict = strict;
717 convMeasureDef.meshName = meshName;
718 convMeasureDef.measure = std::move(measure);
719 convMeasureDef.doesLogging = true;
720 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
721}
722
724 const std::string &dataName,
725 const std::string &meshName,
726 double limit,
727 bool suffices,
728 bool strict)
729{
731 PRECICE_CHECK(math::greater(limit, 0.0) && math::greaterEquals(1.0, limit),
732 "Relative convergence limit has to be in ]0;1]. "
733 "Please check the <relative-convergence-measure limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
734 "in your <coupling-scheme ... /> in the preCICE configuration file.",
735 limit, dataName, meshName);
738 "The relative convergence limit=\"{}\" is close to the hard-coded numerical resolution=\"{}\" of preCICE. "
739 "This may lead to instabilities. The minimum relative convergence limit should be > \"{}\" ",
741
743 ConvergenceMeasureDefintion convMeasureDef;
744 convMeasureDef.data = getData(dataName, meshName);
745 convMeasureDef.suffices = suffices;
746 convMeasureDef.strict = strict;
747 convMeasureDef.meshName = meshName;
748 convMeasureDef.measure = std::move(measure);
749 convMeasureDef.doesLogging = true;
750 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
751}
752
754 const std::string &dataName,
755 const std::string &meshName,
756 double limit,
757 bool suffices,
758 bool strict)
759{
761 PRECICE_CHECK(math::greater(limit, 0.0) && math::greaterEquals(1.0, limit),
762 "Relative convergence limit has to be in ]0;1]. "
763 "Please check the <residul-relative-convergence-measure limit=\"{}\" data=\"{}\" mesh=\"{}\" /> subtag "
764 "in your <coupling-scheme ... /> in the preCICE configuration file.",
765 limit, dataName, meshName);
768 "The relative convergence limit=\"{}\" is close to the hard-coded numerical resolution=\"{}\" of preCICE. "
769 "This may lead to instabilities. The minimum relative convergence limit should be > \"{}\" ",
771
773 ConvergenceMeasureDefintion convMeasureDef;
774 convMeasureDef.data = getData(dataName, meshName);
775 convMeasureDef.suffices = suffices;
776 convMeasureDef.strict = strict;
777 convMeasureDef.meshName = meshName;
778 convMeasureDef.measure = std::move(measure);
779 convMeasureDef.doesLogging = true;
780 _config.convergenceMeasureDefinitions.push_back(convMeasureDef);
781}
782
784 const std::string &dataName,
785 const std::string &meshName) const
786{
787 PRECICE_CHECK(_meshConfig->hasMeshName(meshName) && _meshConfig->getMesh(meshName)->data(dataName),
788 "Data \"{}\" used by mesh \"{}\" is not configured.", dataName, meshName);
789 const mesh::PtrMesh &mesh = _meshConfig->getMesh(meshName);
790 return mesh->data(dataName);
791}
792
794 int ID) const
795{
796 for (const mesh::PtrMesh &mesh : _meshConfig->meshes()) {
797 if (mesh->hasDataID(ID)) {
798 return mesh->data(ID);
799 }
800 }
801 return nullptr;
802}
803
805 const std::string &accessor) const
806{
807 PRECICE_TRACE(accessor);
808 m2n::PtrM2N m2n = _m2nConfig->getM2N(
811
812 for (const auto &exchange : _config.exchanges) {
813 if ((exchange.from == _config.participants[1]) && exchange.exchangeSubsteps) {
815 "You enabled exchanging substeps in the serial-explicit coupling between the second participant \"{}\" and first participant \"{}\". "
816 "This is inefficient as these substeps will never be used.",
817 exchange.from, exchange.to);
818 break;
819 }
820 }
821
822 addDataToBeExchanged(*scheme, accessor);
823
824 return PtrCouplingScheme(scheme);
825}
826
828 const std::string &accessor) const
829{
830 PRECICE_TRACE(accessor);
831 m2n::PtrM2N m2n = _m2nConfig->getM2N(
834
835 for (const auto &exchange : _config.exchanges) {
836 if (exchange.exchangeSubsteps) {
838 "You enabled exchanging substeps in the parallel-explicit coupling between \"{}\" and \"{}\". "
839 "This is inefficient as these substeps will never be used.",
840 exchange.from, exchange.to);
841 break;
842 }
843 }
844
845 addDataToBeExchanged(*scheme, accessor);
846
847 return PtrCouplingScheme(scheme);
848}
849
851 const std::string &accessor) const
852{
853 PRECICE_TRACE(accessor);
854
855 const auto first = _config.participants[0];
856 const auto second = _config.participants[1];
857
858 m2n::PtrM2N m2n = _m2nConfig->getM2N(
859 first, second);
861
862 addDataToBeExchanged(*scheme, accessor);
864 "No send data configured. "
865 "Use explicit scheme for one-way coupling. "
866 "Please check your <coupling-scheme ... /> and make sure that you provide at least one <exchange .../> subtag, "
867 "where from=\"{}\".",
868 accessor);
869
870 // Add convergence measures
873
874 // Set acceleration
875 setSerialAcceleration(scheme, first, second);
876
877 if (scheme->doesFirstStep() && _accelerationConfig->getAcceleration() && not _accelerationConfig->getAcceleration()->getDataIDs().empty()) {
878 DataID dataID = *(_accelerationConfig->getAcceleration()->getDataIDs().begin());
879 PRECICE_CHECK(not scheme->hasSendData(dataID),
880 "In case of serial coupling, acceleration can be defined for data of second participant only!");
881 }
882
883 return PtrCouplingScheme(scheme);
884}
885
887 const std::string &accessor) const
888{
889 PRECICE_TRACE(accessor);
890 m2n::PtrM2N m2n = _m2nConfig->getM2N(
893
894 addDataToBeExchanged(*scheme, accessor);
896 "No send data configured. Use explicit scheme for one-way coupling. "
897 "Please check your <coupling-scheme ... /> and make sure that you provide at least one <exchange .../> subtag, "
898 "where from=\"{}\".",
899 accessor);
900
901 // Add convergence measures
904
905 // Set acceleration
907
908 return PtrCouplingScheme(scheme);
909}
910
912 const std::string &accessor) const
913{
914 PRECICE_TRACE(accessor);
915
916 BaseCouplingScheme *scheme;
917
919 for (const std::string &participant : _config.participants) {
920 if (_m2nConfig->isM2NConfigured(accessor, participant)) {
921 m2ns[participant] = _m2nConfig->getM2N(accessor, participant);
922 }
923 }
924
925 scheme = new MultiCouplingScheme(
927
928 MultiCouplingScheme *castedScheme = dynamic_cast<MultiCouplingScheme *>(scheme);
929 PRECICE_ASSERT(castedScheme, "The dynamic cast of CouplingScheme failed.");
930 addMultiDataToBeExchanged(*castedScheme, accessor);
931
933 "No send data configured. Use explicit scheme for one-way coupling. "
934 "Please check your <coupling-scheme ... /> and make sure that you provide at least one "
935 "<exchange .../> subtag, where from=\"{}\".",
936 accessor);
937
938 // Add convergence measures
940 if (accessor == _config.controller) {
942 }
943
944 // Set acceleration
946
948 not scheme->doesFirstStep() && _accelerationConfig->getAcceleration() && _accelerationConfig->getAcceleration()->getDataIDs().size() < 3,
949 "Due to numerical reasons, for multi coupling, the number of coupling data vectors should be at least 3, not: {}. "
950 "Please check the <data .../> subtags in your <acceleration:.../> and make sure that you have at least 3.",
951 _accelerationConfig->getAcceleration()->getDataIDs().size());
952 return PtrCouplingScheme(scheme);
953}
954
957 const std::string &method) const
958{
959 PRECICE_TRACE(method);
960 if (method == VALUE_FIXED) {
962 } else if (method == VALUE_FIRST_PARTICIPANT) {
964 } else {
965 // We should never reach this point.
966 PRECICE_UNREACHABLE("Unknown timestepping method '{}'.", method);
967 }
968}
969
976
978{
981 "Not defining convergence measures without providing a maximum iteration limit is forbidden."
982 "Please define a convergence measure or set a maximum iteration limit using <max-iterations value=\"...\" />.");
983
984 PRECICE_INFO("No convergence measures were defined for an implicit coupling scheme. "
985 "It will always iterate the maximum amount iterations, which is {}."
986 "You may want to add a convergence measure in your <coupling-scheme:.../> in your configuration.",
988 }
989}
990
992{
993 const auto &participant = _participantConfig->getParticipant(exchange.to);
994
995 const auto &meshPtr = participant->findMesh(exchange.data->getName()); // related to https://github.com/precice/precice/issues/1694
996
997 if (meshPtr == nullptr) {
998 // Only warn, because might be valid configuration, if summation action is used. See Integration/Serial/SummationActionTwoSources.
999 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.",
1000 exchange.data->getName(), exchange.to, exchange.to, exchange.data->getName());
1001 return; // skip checks below
1002 }
1003
1004 const auto &readDataContext = participant->readDataContext(meshPtr->getName(), exchange.data->getName());
1005 if (readDataContext.getWaveformDegree() == 0) {
1006 PRECICE_CHECK(!exchange.exchangeSubsteps,
1007 "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.",
1008 readDataContext.getDataName(), readDataContext.getWaveformDegree(), exchange.data->getName(), exchange.mesh->getName(), exchange.from, exchange.to);
1009 } else if (readDataContext.getWaveformDegree() >= 2) {
1010 PRECICE_CHECK(exchange.exchangeSubsteps,
1011 "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.",
1012 readDataContext.getDataName(), readDataContext.getWaveformDegree(), exchange.data->getName(), exchange.mesh->getName(), exchange.from, exchange.to);
1013 } else { // For first degree there is no restriction for exchange of substeps
1014 PRECICE_ASSERT(readDataContext.getWaveformDegree() == 1);
1015 }
1016}
1017
1019 BiCouplingScheme & scheme,
1020 const std::string &accessor) const
1021{
1022 PRECICE_TRACE();
1023 for (const Config::Exchange &exchange : _config.exchanges) {
1024 const std::string &from = exchange.from;
1025 const std::string &to = exchange.to;
1026 const std::string &dataName = exchange.data->getName();
1027 const std::string &meshName = exchange.mesh->getName();
1028
1029 PRECICE_CHECK(to != from,
1030 "You cannot define an exchange from and to the same participant. "
1031 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1032 dataName, meshName, from, to);
1033
1035 "Participant \"{}\" is not configured for coupling scheme. "
1036 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1037 from, dataName, meshName, from, to);
1038
1040 "Participant \"{}\" is not configured for coupling scheme. "
1041 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1042 to, dataName, meshName, from, to);
1043
1044 const bool requiresInitialization = exchange.requiresInitialization;
1045
1047 !(requiresInitialization && _participantConfig->getParticipant(from)->isDirectAccessAllowed(exchange.mesh->getName())),
1048 "Participant \"{}\" cannot initialize data of the directly-accessed mesh \"{}\" from the participant\"{}\". "
1049 "Either disable the initialization in the <exchange /> tag or use a locally provided mesh instead.",
1050 from, meshName, to);
1051
1052 const bool exchangeSubsteps = exchange.exchangeSubsteps;
1053
1054 if (from == accessor) {
1055 scheme.addDataToSend(exchange.data, exchange.mesh, requiresInitialization, exchangeSubsteps);
1056 } else if (to == accessor) {
1058 scheme.addDataToReceive(exchange.data, exchange.mesh, requiresInitialization, exchangeSubsteps);
1059 } else {
1061 }
1062 }
1064}
1065
1067 MultiCouplingScheme &scheme,
1068 const std::string & accessor) const
1069{
1070 PRECICE_TRACE();
1071 for (const Config::Exchange &exchange : _config.exchanges) {
1072 const std::string &from = exchange.from;
1073 const std::string &to = exchange.to;
1074 const std::string &dataName = exchange.data->getName();
1075 const std::string &meshName = exchange.mesh->getName();
1076
1077 PRECICE_CHECK(to != from,
1078 "You cannot define an exchange from and to the same participant. "
1079 "Please check the <exchange data=\"{}\" mesh=\"{}\" from=\"{}\" to=\"{}\" /> tag in the <coupling-scheme:... /> of your precice-config.xml.",
1080 dataName, meshName, from, to);
1081
1083 "Participant \"{}\" is not configured for coupling scheme",
1084 from);
1085
1087 "Participant \"{}\" is not configured for coupling scheme", to);
1088
1089 const bool initialize = exchange.requiresInitialization;
1090 const bool exchangeSubsteps = exchange.exchangeSubsteps;
1091
1092 if (from == accessor) {
1093 scheme.addDataToSend(exchange.data, exchange.mesh, initialize, exchangeSubsteps, to);
1094 } else if (to == accessor) {
1095 scheme.addDataToReceive(exchange.data, exchange.mesh, initialize, exchangeSubsteps, from);
1096 }
1097 }
1099}
1100
1102 DataID dataID) const
1103{
1104 const auto match = std::find_if(_config.exchanges.begin(),
1105 _config.exchanges.end(),
1106 [dataID](const Config::Exchange &exchange) { return exchange.data->getID() == dataID; });
1107 if (match != _config.exchanges.end()) {
1108 return;
1109 }
1110
1111 // Data is not being exchanged
1112 std::string dataName = "";
1113 auto dataptr = findDataByID(dataID);
1114 if (dataptr) {
1115 dataName = dataptr->getName();
1116 }
1117
1118 PRECICE_ERROR("You need to exchange every data that you use for convergence measures and/or the iteration acceleration. "
1119 "Data \"{}\" is currently not exchanged over the respective mesh on which it is used for convergence measures and/or iteration acceleration. "
1120 "Please check the <exchange ... /> and <...-convergence-measure ... /> tags in the <coupling-scheme:... /> of your precice-config.xml.",
1121 dataName);
1122}
1123
1125 int dataID,
1126 const std::string &first,
1127 const std::string &second) const
1128{
1129 checkIfDataIsExchanged(dataID);
1130 const auto match = std::find_if(_config.exchanges.begin(),
1131 _config.exchanges.end(),
1132 [dataID](const Config::Exchange &exchange) { return exchange.data->getID() == dataID; });
1133 PRECICE_ASSERT(match != _config.exchanges.end());
1134 const auto &exchange = *match;
1135
1136 // In serial implicit cplschemes, data is only accelerated on the second participant.
1137 if (second == exchange.from) {
1138 return;
1139 }
1140
1141 std::string dataName = "";
1142 auto dataptr = findDataByID(dataID);
1143 if (dataptr) {
1144 dataName = dataptr->getName();
1145 }
1146
1148 "You configured acceleration data \"{}\" in the serial implicit coupling scheme between participants \"{}\" and \"{}\". "
1149 "For serial implicit coupling schemes, only data exchanged from the second to the first participant can be used for acceleration. "
1150 "Here, from \"{}\" to \"{}\". "
1151 "However, you configured data \"{}\" for acceleration, which is exchanged from \"{}\" to \"{}\". "
1152 "Please remove this acceleration data tag or switch to a parallel implicit coupling scheme.",
1153 dataName, first, second, second, first, dataName, first, second);
1154}
1155
1157 BaseCouplingScheme * scheme,
1158 const std::string & participant,
1159 const std::vector<ConvergenceMeasureDefintion> &convergenceMeasureDefinitions) const
1160{
1161 for (auto &elem : convergenceMeasureDefinitions) {
1162 _meshConfig->addNeededMesh(participant, elem.meshName);
1163 checkIfDataIsExchanged(elem.data->getID());
1164 scheme->addConvergenceMeasure(elem.data->getID(), elem.suffices, elem.strict, elem.measure, elem.doesLogging);
1165 }
1166}
1167
1169 BaseCouplingScheme *scheme,
1170 const std::string & first,
1171 const std::string & second) const
1172{
1173 if (_accelerationConfig->getAcceleration().get() != nullptr) {
1174 for (std::string &neededMesh : _accelerationConfig->getNeededMeshes()) {
1175 _meshConfig->addNeededMesh(second, neededMesh);
1176 }
1177 for (const DataID dataID : _accelerationConfig->getAcceleration()->getDataIDs()) {
1178 checkSerialImplicitAccelerationData(dataID, first, second);
1179 }
1180 scheme->setAcceleration(_accelerationConfig->getAcceleration());
1181 }
1182}
1183
1185 BaseCouplingScheme *scheme,
1186 const std::string & participant) const
1187{
1188 if (_accelerationConfig->getAcceleration().get() != nullptr) {
1189 for (std::string &neededMesh : _accelerationConfig->getNeededMeshes()) {
1190 _meshConfig->addNeededMesh(participant, neededMesh);
1191 }
1192 for (const DataID dataID : _accelerationConfig->getAcceleration()->getDataIDs()) {
1193 checkIfDataIsExchanged(dataID);
1194 }
1195 scheme->setAcceleration(_accelerationConfig->getAcceleration());
1196
1198 dynamic_cast<acceleration::AitkenAcceleration *>(_accelerationConfig->getAcceleration().get()) != nullptr,
1199 "You configured participant \"{}\" in a parallel-implicit coupling scheme with \"Aitken\" "
1200 "acceleration, which is known to perform bad in parallel coupling schemes. "
1201 "See https://precice.org/configuration-acceleration.html#dynamic-aitken-under-relaxation for details."
1202 "Consider switching to a serial-implicit coupling scheme or changing the acceleration method.",
1203 participant);
1204 }
1205}
1206
1207} // namespace precice::cplscheme
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:15
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:21
#define PRECICE_WARN(...)
Definition LogMacros.hpp:11
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_INFO(...)
Definition LogMacros.hpp:13
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:35
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
#define PRECICE_UNREACHABLE(...)
Definition assertion.hpp:95
T begin(T... args)
Abstract base class for standard coupling schemes.
void addConvergenceMeasure(int dataID, bool suffices, bool strict, impl::PtrConvergenceMeasure measure, bool doesLogging)
Adds a measure to determine the convergence of coupling iterations.
void setAcceleration(const acceleration::PtrAcceleration &acceleration)
Set an acceleration technique.
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 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.
virtual void xmlEndTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag)
Callback method required when using xml::XMLTag.
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)
virtual void xmlTagCallback(const xml::ConfigurationContext &context, xml::XMLTag &callingTag)
Callback method required when using xml::XMLTag.
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.
Measures the convergence from an old data set to a new one.
Measures the convergence from an old data set to a new one.
Measures the convergence from an old data set to a new one.
Measures the convergence from an old data set to a new one.
XMLAttribute & setOptions(std::vector< ATTRIBUTE_T > options)
XMLAttribute & setDocumentation(std::string documentation)
Sets a documentation string for the attribute.
Represents an XML tag to be configured automatically.
Definition XMLTag.hpp:31
const std::string & getNamespace() const
Returns xml namespace.
Definition XMLTag.hpp:159
std::string getStringAttributeValue(const std::string &name, std::optional< std::string > default_value=std::nullopt) const
Definition XMLTag.cpp:142
bool getBooleanAttributeValue(const std::string &name, std::optional< bool > default_value=std::nullopt) const
Definition XMLTag.cpp:155
XMLTag & setDocumentation(const std::string &documentation)
Adds a description of the purpose of this XML tag.
Definition XMLTag.cpp:29
const std::string & getName() const
Returns name (without namespace).
Definition XMLTag.hpp:153
const std::string & getFullName() const
Returns full name consisting of xml namespace + ":" + name.
Definition XMLTag.hpp:170
int getIntAttributeValue(const std::string &name, std::optional< int > default_value=std::nullopt) const
Definition XMLTag.cpp:129
double getDoubleAttributeValue(const std::string &name, std::optional< double > default_value=std::nullopt) const
Definition XMLTag.cpp:116
XMLTag & addAttribute(const XMLAttribute< double > &attribute)
Removes the XML subtag with given name.
Definition XMLTag.cpp:53
XMLTag & addSubtag(const XMLTag &tag)
Adds an XML tag as subtag by making a copy of the given tag.
Definition XMLTag.cpp:41
T compare(T... args)
T end(T... args)
T find(T... args)
T get(T... args)
vector< double > getData()
Definition mainA.cpp:19
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:39
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:24