preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BaseCouplingScheme.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <boost/range/adaptor/map.hpp>
4#include <cmath>
5#include <cstddef>
6#include <functional>
7#include <iterator>
8#include <limits>
9#include <sstream>
10#include <utility>
11
21#include "io/TXTTableWriter.hpp"
22#include "logging/LogMacros.hpp"
23#include "math/differences.hpp"
24#include "mesh/Data.hpp"
25#include "mesh/Mesh.hpp"
27#include "profiling/Event.hpp"
29#include "utils/Helpers.hpp"
30#include "utils/IntraComm.hpp"
31#include "utils/assertion.hpp"
32
34
36 double maxTime,
37 int maxTimeWindows,
38 double timeWindowSize,
39 std::string localParticipant,
40 int minIterations,
41 int maxIterations,
42 CouplingMode cplMode,
44 : _couplingMode(cplMode),
45 _maxTime(maxTime),
46 _maxTimeWindows(maxTimeWindows),
47 _timeWindows(1),
48 _timeWindowSize(timeWindowSize),
49 _nextTimeWindowSize(timeWindowSize),
50 _minIterations(minIterations),
51 _maxIterations(maxIterations),
52 _iterations(1),
53 _totalIterations(1),
54 _localParticipant(std::move(localParticipant))
55{
56 PRECICE_ASSERT(not((maxTime != UNDEFINED_MAX_TIME) && (maxTime < 0.0)),
57 "Maximum time has to be larger than zero.");
58 PRECICE_ASSERT(not((maxTimeWindows != UNDEFINED_TIME_WINDOWS) && (maxTimeWindows < 0)),
59 "Maximum number of time windows has to be larger than zero.");
60 PRECICE_ASSERT(not(hasTimeWindowSize() && (timeWindowSize < 0.0)),
61 "Time window size has to be larger than zero.");
62 if (dtMethod == constants::FIXED_TIME_WINDOW_SIZE) {
64 "Time window size has to be given when the fixed time window size method is used.");
65 }
66
70 } else {
74
75 PRECICE_ASSERT(minIterations > 0,
76 minIterations,
77 "Minimal iteration limit has to be larger than zero.");
78 PRECICE_ASSERT((maxIterations == INFINITE_MAX_ITERATIONS) || (maxIterations > 0),
79 maxIterations,
80 "Maximal iteration limit has to be larger than zero or -1 (unlimited).");
81 PRECICE_ASSERT((maxIterations == INFINITE_MAX_ITERATIONS) || (minIterations <= maxIterations),
82 "Minimal iteration limit has to be smaller equal compared to the maximal iteration limit.");
83 }
84
85 if (math::equals(maxTime, UNDEFINED_MAX_TIME)) {
87 } else {
88 _time = impl::TimeHandler(maxTime);
89 }
90}
91
97
99{
100 return _hasConverged;
101}
102
104{
106 PRECICE_DEBUG("Sending number or time steps {}...", times.size());
107 m2n->send(static_cast<int>(times.size()));
108 PRECICE_DEBUG("Sending times...");
109 m2n->send(times);
110}
111
112void BaseCouplingScheme::sendData(const m2n::PtrM2N &m2n, const DataMap &sendData)
113{
115 PRECICE_ASSERT(m2n.get() != nullptr);
116 PRECICE_ASSERT(m2n->isConnected());
117
118 profiling::Event e("waitAndSendData", profiling::Fundamental);
119
120 for (const auto &data : sendData | boost::adaptors::map_values) {
121 if (data->exchangeSubsteps()) {
122 const auto &stamples = data->stamples();
123 PRECICE_ASSERT(!stamples.empty());
124
125 int nTimeSteps = data->timeStepsStorage().nTimes();
126 PRECICE_ASSERT(nTimeSteps > 0);
127 const auto timesAscending = data->timeStepsStorage().getTimes();
128 sendTimes(m2n, timesAscending);
129
130 const auto serialized = com::serialize::SerializedStamples::serialize(*data);
131
132 // Data is actually only send if size>0, which is checked in the derived classes implementation
133 m2n->send(serialized.values(), data->getMeshID(), data->getDimensions() * serialized.nTimeSteps());
134
135 if (data->hasGradient()) {
136 m2n->send(serialized.gradients(), data->getMeshID(), data->getDimensions() * data->meshDimensions() * serialized.nTimeSteps());
137 }
138 } else {
139 const auto &sample = data->timeStepsStorage().getSampleAtEnd();
140 if (data->hasGradient()) {
141 // Data is only received on ranks with size>0, which is checked in the derived class implementation
142 m2n->send(sample.values, data->getMeshID(), data->getDimensions());
143 m2n->send(sample.gradients, data->getMeshID(), data->getDimensions() * data->meshDimensions());
144 } else {
145 // Data is only received on ranks with size>0, which is checked in the derived class implementation
146 m2n->send(sample.values, data->getMeshID(), data->getDimensions());
147 }
148 }
149 }
150}
151
153{
155 PRECICE_DEBUG("Receiving number of time steps...");
156 int numberOfTimeSteps;
157 m2n->receive(numberOfTimeSteps);
158 PRECICE_ASSERT(numberOfTimeSteps > 0);
159
160 std::vector<double> times(numberOfTimeSteps);
161 PRECICE_DEBUG("Receiving times....");
162 m2n->receive(times);
163 PRECICE_DEBUG("Received times {}", times);
164 return times;
165}
166
167void BaseCouplingScheme::receiveData(const m2n::PtrM2N &m2n, const DataMap &receiveData)
168{
170 PRECICE_ASSERT(m2n.get());
171 PRECICE_ASSERT(m2n->isConnected());
172 profiling::Event e("waitAndReceiveData", profiling::Fundamental);
173 for (const auto &data : receiveData | boost::adaptors::map_values) {
174
175 if (data->exchangeSubsteps()) {
176 auto timesAscending = receiveTimes(m2n);
177 const auto nTimeSteps = timesAscending.size();
178
179 auto serialized = com::serialize::SerializedStamples::empty(nTimeSteps, *data);
180
181 // Data is only received on ranks with size>0, which is checked in the derived class implementation
182 m2n->receive(serialized.values(), data->getMeshID(), data->getDimensions() * nTimeSteps);
183
184 if (data->hasGradient()) {
185 m2n->receive(serialized.gradients(), data->getMeshID(), data->getDimensions() * data->meshDimensions() * nTimeSteps);
186 }
187
188 serialized.deserializeInto(timesAscending, *data);
189 } else {
190 if (data->hasGradient()) {
191 // Data is only received on ranks with size>0, which is checked in the derived class implementation
192 time::Sample recvSample(data->getDimensions(), data->nVertices(), data->meshDimensions());
193 m2n->receive(recvSample.values, data->getMeshID(), data->getDimensions());
194 m2n->receive(recvSample.gradients, data->getMeshID(), data->getDimensions() * data->meshDimensions());
195 data->setSampleAtTime(getTime(), recvSample);
196 } else {
197 // Data is only received on ranks with size>0, which is checked in the derived class implementation
198 time::Sample recvSample(data->getDimensions(), data->nVertices());
199 m2n->receive(recvSample.values, data->getMeshID(), data->getDimensions());
200 data->setSampleAtTime(getTime(), recvSample);
201 }
202 }
203 }
204}
205
207{
208 // @TODO This is a hack required until https://github.com/precice/precice/issues/1957
209 // buffer current time
210 const impl::TimeHandler oldTime = _time;
211 // set _time state to point to end of this window such that _time in receiveData is at end of window
213 // receive data for end of window
214 this->receiveData(m2n, receiveData);
215 // reset time state;
216 _time = oldTime;
217}
218
220{
221 for (const auto &data : receiveData | boost::adaptors::map_values) {
222 PRECICE_DEBUG("Initialize {} as zero.", data->getDataName());
223 // just store already initialized zero sample to storage.
224 data->initializeWithZeroAtTime(getTime());
225 }
226}
227
228PtrCouplingData BaseCouplingScheme::addCouplingData(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool communicateSubsteps, CouplingData::Direction direction)
229{
230 const DataID id = data->getID();
231
232 // new data
233 if (_allData.count(id) == 0) {
234 auto ptr = std::make_shared<CouplingData>(data, std::move(mesh), requiresInitialization, communicateSubsteps, direction);
235 _allData.emplace(id, ptr);
236 return ptr;
237 }
238
239 // data is already used by another exchange of this coupling scheme, use existing CouplingData
240 auto &existing = _allData[id];
241 PRECICE_CHECK(existing->getDirection() == direction,
242 "Data \"{0}\" cannot be added for sending and for receiving. Please remove either <exchange data=\"{0}\" ... /> tag", data->getName());
244 "Data \"{0}\" cannot be received from multiple sources. Please remove either <exchange data=\"{0}\" ... /> tag", data->getName());
245 PRECICE_CHECK(existing->exchangeSubsteps() == communicateSubsteps,
246 "Data \"{0}\" is configured with substeps enabled and disabled at the same time. Please make the configuration consistent in the <exchange data=\"{0}\" ... substeps=\"True/False\" ... /> tags", data->getName());
247 PRECICE_CHECK(existing->requiresInitialization == requiresInitialization,
248 "Data \"{0}\" is configured with data initialization enabled and disabled at the same time. Please make the configuration consistent in the <exchange data=\"{0}\" ... /> tags", data->getName());
249
250 return existing;
251}
252
258
259void BaseCouplingScheme::setTimeWindowSize(double timeWindowSize)
260{
261 _timeWindowSize = timeWindowSize;
262}
263
265{
266 _nextTimeWindowSize = timeWindowSize;
267}
268
270{
273 PRECICE_ASSERT(_isInitialized, "Called finalize() before initialize().");
274}
275
277{
278 // initialize with zero data here, might eventually be overwritten in exchangeInitialData
280 // Initialize uses the template method pattern (https://en.wikipedia.org/wiki/Template_method_pattern).
283 _time.resetTo(0);
284 _timeWindows = 1;
285 _hasDataBeenReceived = false;
286
288
291 if (not doesFirstStep()) {
292 // reserve memory and initialize data with zero
293 if (_acceleration) {
294 _acceleration->initialize(getAccelerationData());
296 qnAcceleration) {
297 PRECICE_WARN_IF((qnAcceleration->getMaxUsedTimeWindows() == 0) && (qnAcceleration->getMaxUsedIterations() > _maxIterations), "The maximum number of iterations used in the quasi-Newton acceleration scheme is greater than the maximum number of iterations allowed in one time window. When time-windows-reused is set to 0, and therefore no previous time windows are reused, the actual max-used-ietrations is equal to max-iterations.");
298 }
299 }
300 }
303 }
304
306
307 _isInitialized = true;
308}
309
311{
314
316 // overwrite past iteration with new samples
317 for (const auto &data : _allData | boost::adaptors::map_values) {
318 // TODO: reset CouplingData of changed meshes only #2102
319 data->reinitialize();
320 }
321
322 if (not doesFirstStep()) {
323 if (_acceleration) {
324 _acceleration->initialize(getAccelerationData());
325 }
326 }
328 }
329}
330
335
341
343{
346 PRECICE_ASSERT(_isInitialized, "Before calling advance() coupling scheme has to be initialized via initialize().");
347 _hasDataBeenReceived = false;
348 _isTimeWindowComplete = false;
349
351
353
354 _timeWindows += 1; // increment window counter. If not converged, will be decremented again later.
355
357 }
358}
359
364
366{
369 PRECICE_ASSERT(_isInitialized, "Before calling advance() coupling scheme has to be initialized via initialize().");
371
372 // from first phase
374
376
378
379 if (isImplicitCouplingScheme()) { // check convergence
380 if (not hasConverged()) { // repeat window
381 PRECICE_DEBUG("No convergence achieved");
383 // The computed time window part equals the time window size, since the
384 // time window remainder is zero. Subtract the time window size and do another
385 // coupling iteration.
387 _timeWindows -= 1;
388 _isTimeWindowComplete = false;
389 } else { // write output, prepare for next window
390 PRECICE_DEBUG("Convergence achieved");
392 PRECICE_INFO("Time window completed");
394 if (isCouplingOngoing()) {
395 PRECICE_DEBUG("Setting require create checkpoint");
397 }
398 }
399 // update iterations
401 if (not hasConverged()) {
402 _iterations++;
403 } else {
404 _iterations = 1;
405 }
406 } else {
408 PRECICE_INFO("Time window completed");
410 }
411 if (isCouplingOngoing()) {
413 }
414
415 // Update internal time tracking
418 } else {
420 }
422 }
423}
424
426{
428 for (auto &data : _allData | boost::adaptors::map_values) {
429 data->moveToNextWindow();
430 }
431}
432
437
443
448
450{
451 return _isInitialized;
452}
453
455 double timeToAdd)
456{
457 PRECICE_TRACE(timeToAdd, getTime());
458 PRECICE_ASSERT(isCouplingOngoing(), "Invalid call of addComputedTime() after simulation end.");
459
460 // Check validness
461 bool valid = math::greaterEquals(getNextTimeStepMaxSize(), timeToAdd);
462 PRECICE_CHECK(valid,
463 "The time step size given to preCICE in \"advance\" {} exceeds the maximum allowed time step size {} "
464 "in the remaining of this time window. "
465 "Did you restrict your time step size, \"dt = min(preciceDt, solverDt)\"? "
466 "For more information, consult the adapter example in the preCICE documentation.",
467 timeToAdd, getNextTimeStepMaxSize());
468
469 // add time interval that has been computed in the solver to get the correct time remainder
470 _time.progressBy(timeToAdd);
471
472 return reachedEndOfTimeWindow();
473}
474
476 double lastSolverTimeStepSize) const
477{
478 PRECICE_TRACE(lastSolverTimeStepSize);
479 double remainder = getNextTimeStepMaxSize() - lastSolverTimeStepSize;
480 return not math::greater(remainder, 0.0);
481}
482
487
489{
491}
492
494{
495 PRECICE_ASSERT(not _hasDataBeenReceived, "notifyDataHasBeenReceived() may only be called once within one coupling iteration. If this assertion is triggered this probably means that your coupling scheme has a bug.");
497}
498
503
505{
506 _timeWindows = timeWindows;
507}
508
510{
511 return _time.time();
512}
513
515{
516 return _time.windowStart();
517}
518
520{
521 return _timeWindows;
522}
523
525{
526 if (!isCouplingOngoing()) {
527 return 0.0; // if coupling is not ongoing (i.e. coupling scheme reached end of window) the maximum time step size is zero
528 }
529
530 if (hasTimeWindowSize()) {
532 } else {
533 return _time.untilEnd();
534 }
535}
536
538{
539 bool timestepsLeft = (_maxTimeWindows >= _timeWindows) || (_maxTimeWindows == UNDEFINED_TIME_WINDOWS);
540 return timestepsLeft && !_time.reachedEnd();
541}
542
547
549 Action action) const
550{
551 return _requiredActions.count(action) == 1;
552}
553
555 Action action) const
556{
557 return _fulfilledActions.count(action) == 1;
558}
559
566
568 Action action)
569{
570 _requiredActions.insert(action);
571}
572
574{
576 if (isCouplingOngoing()) {
577 os << "iteration: " << _iterations; //_iterations;
579 os << " of " << _maxIterations;
580 }
582 os << " (min " << _minIterations << ")";
583 }
584 os << ", ";
585 }
587 std::string actionsState = printActionsState();
588 if (!actionsState.empty()) {
589 os << ", " << actionsState;
590 }
591 return os.str();
592}
593
595 int timeWindows,
596 double time) const
597{
599 if (isCouplingOngoing()) {
600 os << "time-window: " << timeWindows;
602 os << " of " << _maxTimeWindows;
603 }
604 os << ", time: " << time;
606 os << " of " << _maxTime;
607 }
608 if (hasTimeWindowSize()) {
609 os << ", time-window-size: " << _timeWindowSize;
610 }
612 os << ", max-time-step-size: " << getNextTimeStepMaxSize();
613 }
614 os << ", ongoing: ";
615 isCouplingOngoing() ? os << "yes" : os << "no";
616 os << ", time-window-complete: ";
617 _isTimeWindowComplete ? os << "yes" : os << "no";
618 } else {
619 os << "Reached end at: final time-window: " << (timeWindows - 1) << ", final time: " << time;
620 }
621 return os.str();
622}
623
625{
627 for (auto action : _requiredActions) {
628 os << toString(action) << ' ';
629 }
630 return os.str();
631}
632
634{
636 std::vector<Action> missing;
639 std::back_inserter(missing));
640 if (not missing.empty()) {
641 std::ostringstream stream;
642 for (auto action : missing) {
643 if (not stream.str().empty()) {
644 stream << ", ";
645 }
646 stream << toString(action);
647 }
648 PRECICE_ERROR("The required actions {} are not fulfilled. "
649 "Did you forget to call \"requiresReadingCheckpoint()\" or \"requiresWritingCheckpoint()\"?",
650 stream.str());
651 }
654}
655
657{
659 for (const auto &data : _allData | boost::adaptors::map_values) {
660 if (data->getDirection() == CouplingData::Direction::Receive) {
661 PRECICE_ASSERT(!data->stamples().empty(), "initializeReceiveDataStorage() didn't initialize data correctly");
662 } else {
663 PRECICE_CHECK(!data->stamples().empty(),
664 "Data {0} on mesh {1} doesn't contain any samples while initializing a coupling scheme of participant {2}. "
665 "There are two common configuration issues that may cause this. "
666 "Either, make sure participant {2} specifies data {0} to be written using tag <write-data mesh=\"{1}\" data=\"{0}\"/>. "
667 "Or ensure participant {2} defines a mapping to mesh {1} from a mesh using data {0}.",
668 data->getDataName(), data->getMeshName(), localParticipant());
669 }
670 }
671}
672
674 const acceleration::PtrAcceleration &acceleration)
675{
676 PRECICE_ASSERT(acceleration.get() != nullptr);
677 _acceleration = acceleration;
678}
679
681{
682 return _doesFirstStep;
683}
684
686{
689 PRECICE_ASSERT(convMeasure.measure.get() != nullptr);
690 convMeasure.measure->newMeasurementSeries();
691 }
692}
693
695 int dataID,
696 bool suffices,
697 bool strict,
699{
700 ConvergenceMeasureContext convMeasure;
701 PRECICE_ASSERT(_allData.count(dataID) == 1, "Data with given data ID must exist!");
702 convMeasure.couplingData = _allData.at(dataID);
703 convMeasure.suffices = suffices;
704 convMeasure.strict = strict;
705 convMeasure.measure = std::move(measure);
706 _convergenceMeasures.push_back(convMeasure);
707}
708
710{
714 _convergenceWriter->writeData("TimeWindow", _timeWindows - 1);
715 _convergenceWriter->writeData("Iteration", _iterations);
716 }
717
718 // If no convergence measures are defined, we never converge
719 if (_convergenceMeasures.empty()) {
720 PRECICE_INFO("No converge measures defined.");
721 return false;
722 }
723
724 // There are convergence measures defined, so we need to check them
725 bool allConverged = true;
726 bool oneSuffices = false; // at least one convergence measure suffices and did converge
727 bool oneStrict = false; // at least one convergence measure is strict and did not converge
728
729 const bool reachedMinIterations = _iterations >= _minIterations;
730 for (const auto &convMeasure : _convergenceMeasures) {
731 PRECICE_ASSERT(convMeasure.couplingData != nullptr);
732 PRECICE_ASSERT(convMeasure.measure.get() != nullptr);
733 PRECICE_ASSERT(convMeasure.couplingData->previousIteration().size() == convMeasure.couplingData->values().size(), convMeasure.couplingData->previousIteration().size(), convMeasure.couplingData->values().size(), convMeasure.couplingData->getDataName());
734 convMeasure.measure->measure(convMeasure.couplingData->previousIteration(), convMeasure.couplingData->values());
735
737 _convergenceWriter->writeData(convMeasure.logHeader(), convMeasure.measure->getNormResidual());
738 }
739
740 if (not convMeasure.measure->isConvergence()) {
741 allConverged = false;
742 if (convMeasure.strict) {
744 oneStrict = true;
746 "The strict convergence measure for data \"" + convMeasure.couplingData->getDataName() +
747 "\" did not converge within the maximum allowed iterations, which terminates the simulation. "
748 "To avoid this forced termination do not mark the convergence measure as strict.");
749 }
750 } else if (convMeasure.suffices == true) {
751 oneSuffices = true;
752 }
753
754 PRECICE_INFO(convMeasure.measure->printState(convMeasure.couplingData->getDataName()));
755 }
756
757 std::string messageSuffix;
758 if (not reachedMinIterations) {
759 messageSuffix = " but hasn't yet reached minimal amount of iterations";
760 }
761 if (allConverged) {
762 PRECICE_INFO("All converged{}", messageSuffix);
763 } else if (oneSuffices && not oneStrict) { // strict overrules suffices
764 PRECICE_INFO("Sufficient measures converged{}", messageSuffix);
765 }
766
767 return reachedMinIterations && (allConverged || (oneSuffices && not oneStrict));
768}
769
771{
773
775 if (not doesFirstStep()) {
777 }
778
779 _iterationsWriter->addData("TimeWindow", io::TXTTableWriter::INT);
780 _iterationsWriter->addData("TotalIterations", io::TXTTableWriter::INT);
781 _iterationsWriter->addData("Iterations", io::TXTTableWriter::INT);
782 _iterationsWriter->addData("Convergence", io::TXTTableWriter::INT);
783
784 if (not doesFirstStep()) {
785 _convergenceWriter->addData("TimeWindow", io::TXTTableWriter::INT);
786 _convergenceWriter->addData("Iteration", io::TXTTableWriter::INT);
787 }
788
789 if (not doesFirstStep()) {
791 _convergenceWriter->addData(convMeasure.logHeader(), io::TXTTableWriter::DOUBLE);
792 }
793 if (_acceleration) {
794 _iterationsWriter->addData("QNColumns", io::TXTTableWriter::INT);
795 _iterationsWriter->addData("DeletedQNColumns", io::TXTTableWriter::INT);
796 _iterationsWriter->addData("DroppedQNColumns", io::TXTTableWriter::INT);
797 }
798 }
799 }
800}
801
803{
805
806 _iterationsWriter->writeData("TimeWindow", _timeWindows - 1);
807 _iterationsWriter->writeData("TotalIterations", _totalIterations);
808 _iterationsWriter->writeData("Iterations", _iterations);
809 bool converged = _iterations >= _minIterations && (_maxIterations < 0 || (_iterations < _maxIterations));
810 _iterationsWriter->writeData("Convergence", converged ? 1 : 0);
811
812 if (not doesFirstStep() && _acceleration) {
814 if (qnAcceleration) {
815 // Only write values for additional columns, if using a QN-based acceleration scheme
816 _iterationsWriter->writeData("QNColumns", qnAcceleration->getLSSystemCols());
817 _iterationsWriter->writeData("DeletedQNColumns", qnAcceleration->getDeletedColumns());
818 _iterationsWriter->writeData("DroppedQNColumns", qnAcceleration->getDroppedColumns());
819 } else {
820 // non-breaking implementation uses "0" for these columns (delete columns in the future?)
821 _iterationsWriter->writeData("QNColumns", 0);
822 _iterationsWriter->writeData("DeletedQNColumns", 0);
823 _iterationsWriter->writeData("DroppedQNColumns", 0);
824 }
825 }
826 }
827}
828
830{
831 if (!hasTimeWindowSize()) {
832 return true; // This participant will always do exactly one step to dictate second's time-window-size
833 }
834
836}
837
839{
841
842 for (const auto &data : _allData | boost::adaptors::map_values) {
843 PRECICE_CHECK(!data->stamples().empty(),
844 "Data {0} on mesh {1} didn't contain any samples while attempting to send it to the coupling partner. "
845 "Make sure participant {2} specifies data {0} to be written using tag <write-data mesh=\"{1}\" data=\"{0}\"/>. "
846 "Alternatively, ensure participant {2} defines a mapping to mesh {1} from a mesh using data {0}.",
847 data->getDataName(), data->getMeshName(), localParticipant());
848 data->storeIteration();
849 }
850}
851
859
866
868{
869 return std::any_of(dataMap.begin(), dataMap.end(),
870 [](auto const &map) { return map.second->requiresInitialization; });
871}
872
874{
875 PRECICE_DEBUG("measure convergence of the coupling iteration");
877 // Stop, when maximal iteration count (given in config) is reached
879 _hasConverged = true;
880
881 // coupling iteration converged for current time window. Advance in time.
882 if (_hasConverged) {
883 if (_acceleration) {
885 _acceleration->iterationsConverged(getAccelerationData(), getTimeWindowStart());
886 }
888 } else {
889 // no convergence achieved for the coupling iteration within the current time window
890 if (_acceleration) {
893 }
894 }
895}
896
898{
900 PRECICE_ASSERT(not doesFirstStep(), "For convergence information the sending participant is never the first one.");
901 profiling::Event e("sendConvergence", profiling::Fundamental);
902 m2n->send(_hasConverged);
903}
904
906{
908 PRECICE_ASSERT(doesFirstStep(), "For convergence information the receiving participant is always the first one.");
909 profiling::Event e("receiveConvergence", profiling::Fundamental);
910 m2n->receive(_hasConverged);
911}
912
914{
915 return _time.windowStart();
916}
917
919{
920 if (hasTimeWindowSize()) {
922 } else {
923 return _time.time() + _time.untilEnd();
924 }
925}
926
928{
929 // Global toggle if a single send data uses substeps
930 for (auto cpldata : _allData | boost::adaptors::map_values) {
931 if (cpldata->getDirection() == CouplingData::Direction::Send && cpldata->exchangeSubsteps()) {
932 return true;
933 }
934 }
935 return false;
936}
937
939{
940 // This implementation covers everything except SerialImplicit
942 return {};
943 }
944
945 ImplicitData idata;
946 for (auto cpldata : _allData | boost::adaptors::map_values) {
947 if (cpldata->getDirection() == CouplingData::Direction::Receive) {
948 idata.add(cpldata->getDataID(), false);
949 }
950 }
951 return idata;
952}
953
958
959} // namespace precice::cplscheme
double const * ptr
Definition ExportCSV.cpp:23
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:16
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:18
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_INFO(...)
Definition LogMacros.hpp:14
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
T any_of(T... args)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T at(T... args)
T back_inserter(T... args)
T begin(T... args)
static SerializedStamples serialize(const cplscheme::CouplingData &data)
Serializes a given CouplingData into SerializedStamples.
static SerializedStamples empty(int nTimeSteps, const cplscheme::CouplingData &data)
Create SerializedStamples with allocated buffers according to size of CouplingData.
bool hasTimeWindowSize() const final override
Function to check whether time window size is defined by coupling scheme.
void initializeWithZeroInitialData(const DataMap &receiveData)
Initializes storage in receiveData as zero.
std::string printCouplingState() const override
Returns coupling state information.
bool isCouplingOngoing() const final override
Returns true, when the coupled simulation is still ongoing.
ChangedMeshes secondSynchronization() final override
void setTimeWindowSize(double timeWindowSize)
Setter for _timeWindowSize.
std::vector< ConvergenceMeasureContext > _convergenceMeasures
All convergence measures of coupling iterations.
void determineInitialReceive(DataMap &receiveData)
Sets _receivesInitializedData, if receiveData requires initialization.
bool addComputedTime(double timeToAdd) final override
Adds newly computed time. Has to be called before every advance.
void sendTimes(const m2n::PtrM2N &m2n, precice::span< double const > times)
std::string _localParticipant
Local participant name.
double _timeWindowSize
size of time window; _timeWindowSize <= _maxTime
double getNextTimeStepMaxSize() const final override
Returns the maximal size of the next time step to be computed.
void determineInitialSend(DataMap &sendData)
Sets _sendsInitializedData, if sendData requires initialization.
bool isTimeWindowComplete() const final override
Returns true, when the accessor can advance to the next time window.
virtual void exchangeSecondData()=0
Exchanges the second set of data.
void reinitialize() final override
Reinitializes the coupling scheme, coupling data, and acceleration schemes.
void receiveDataForWindowEnd(const m2n::PtrM2N &m2n, const DataMap &receiveData)
Like receiveData, but temporarily sets window time to end of window.
bool _hasDataBeenReceived
True, if data has been received from other participant. Flag is used to make sure that coupling schem...
void notifyDataHasBeenReceived()
Used to set flag after data has been received using receiveData().
bool _isInitialized
True, if coupling has been initialized.
bool isExplicitCouplingScheme() const
Function to determine whether coupling scheme is an explicit coupling scheme.
int _maxIterations
Limit of iterations during one time window. Continue to next time window, if _iterations == _maxItera...
void sendData(const m2n::PtrM2N &m2n, const DataMap &sendData)
Sends data sendDataIDs given in mapCouplingData with communication.
virtual void exchangeFirstData()=0
Functions needed for advance()
bool _doesFirstStep
True, if local participant is the one starting the explicit scheme.
bool _receivesInitializedData
True, if this participant has to receive initialized data.
void storeIteration()
used for storing all Data at end of doImplicitStep for later reference.
void initializeTXTWriters()
Initialize txt writers for iterations and convergence tracking.
ChangedMeshes firstSynchronization(const ChangedMeshes &changes) final override
bool isImplicitCouplingScheme() const override
Function to determine whether coupling scheme is an implicit coupling scheme.
bool _hasConverged
True if implicit scheme converged.
void setTimeWindows(int timeWindows)
Setter for _timeWindows.
void sendConvergence(const m2n::PtrM2N &m2n)
sends convergence to other participant via m2n
int _minIterations
Lower limit of iterations during one time window. Prevents convergence if _iterations < _minIteration...
double getTimeWindowSize() const final override
Returns the time window size, if one is given by the coupling scheme.
bool _sendsInitializedData
True, if this participant has to send initialized data.
bool _isTimeWindowComplete
True, if _time == _timeWindowStartTime + _timeWindowSize and (coupling has converged or _iterations =...
int _totalIterations
Number of total iterations performed.
bool anyDataRequiresInitialization(DataMap &dataMap) const
Checks whether any CouplingData in dataMap requires initialization.
bool hasConverged() const override
Checks if the implicit cplscheme has converged.
std::shared_ptr< io::TXTTableWriter > _convergenceWriter
Writes out coupling convergence within all time windows.
double _maxTime
Maximum time being computed. End of simulation is reached, if getTime() == _maxTime.
void finalize() final override
Finalizes the coupling scheme.
std::vector< double > receiveTimes(const m2n::PtrM2N &m2n)
void setAcceleration(const acceleration::PtrAcceleration &acceleration)
Set an acceleration technique.
bool reachedEndOfTimeWindow() const
Function to check whether end of time window is reached. Does not check for convergence.
bool hasDataBeenReceived() const final override
getter for _hasDataBeenReceived
bool requiresSubsteps() const final override
Returns true if any send data of the scheme requires substeps.
void moveToNextWindow()
finalizes this window's data and initializes data for next window.
void addConvergenceMeasure(int dataID, bool suffices, bool strict, impl::PtrConvergenceMeasure measure)
Adds a measure to determine the convergence of coupling iterations.
int _maxTimeWindows
Number of time windows that have to be computed. End of simulation is reached, if _timeWindows == _ma...
void doImplicitStep()
perform a coupling iteration
std::string printBasicState(int timeWindows, double time) const
Prints the coupling state.
double getTime() const final override
getter for _time
std::string printActionsState() const
Prints the action state.
bool receivesInitializedData() const
Getter for _receivesInitializedData.
void requireAction(Action action) final override
Sets an action required to be performed by the accessor.
int getTimeWindows() const final override
getter for _timeWindows
bool isInitialized() const final override
getter for _isInitialized
ImplicitData implicitDataToReceive() const override
Returns a vector of implicit data to receive in the next advance.
bool sendsInitializedData() const final override
Getter for _sendsInitializedData.
PtrCouplingData addCouplingData(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps, CouplingData::Direction direction)
Adds CouplingData with given properties to this BaseCouplingScheme and returns a pointer to the Coupl...
int _timeWindows
number of completed time windows; _timeWindows <= _maxTimeWindows
void initialize() final override
Initializes the coupling scheme and establishes a communication connection to the coupling partner....
void markActionFulfilled(Action action) final override
Tells the coupling scheme that the accessor has performed the given action.
void advanceTXTWriters()
Advance txt writers for iterations and convergence tracking.
acceleration::PtrAcceleration _acceleration
Acceleration method to speedup iteration convergence.
bool willDataBeExchanged(double lastSolverTimeStepSize) const final override
Returns true, if data will be exchanged when calling advance().
virtual DataMap & getAccelerationData()=0
interface to provide accelerated data, depending on coupling scheme being used
double getNextTimeWindowSize() const
Getter for _nextTimeWindowSize.
double _nextTimeWindowSize
time window size of next window (acts as buffer for time windows size provided by first participant,...
DataMap _allData
All send and receive data as a map "data ID -> data".
void checkCouplingDataAvailable()
Issues an error if coupling data does not contain stamples.
CouplingMode _couplingMode
Coupling mode used by coupling scheme.
bool isActionRequired(Action action) const final override
Returns true, if the given action has to be performed by the accessor.
std::string localParticipant() const final override
Returns the name of the local participant.
void receiveConvergence(const m2n::PtrM2N &m2n)
receives convergence from other participant via m2n
void setNextTimeWindowSize(double timeWindowSize)
Setter for _nextTimeWindowSize.
std::shared_ptr< io::TXTTableWriter > _iterationsWriter
Responsible for monitoring iteration count over time window.
bool isActionFulfilled(Action action) const final override
Returns true, if the given action has to be performed by the accessor.
virtual void exchangeInitialData()=0
implements functionality for initialize in base class.
int _iterations
Number of iterations in current time window. _iterations <= _maxIterations.
virtual void initializeReceiveDataStorage()=0
Functions needed for initialize()
void newConvergenceMeasurements()
Reset all convergence measurements after convergence.
BaseCouplingScheme(double maxTime, int maxTimeWindows, double timeWindowSize, std::string localParticipant, int minIterations, int maxIterations, CouplingMode cplMode, constants::TimesteppingMethod dtMethod)
void setDoesFirstStep(bool doesFirstStep)
Setter for _doesFirstStep.
bool measureConvergence()
Measure whether coupling scheme has converged or not.
void receiveData(const m2n::PtrM2N &m2n, const DataMap &receiveData)
Receives data receiveDataIDs given in mapCouplingData with communication.
void checkCompletenessRequiredActions()
If any required actions are open, an error message is issued.
bool doesFirstStep() const
Getter for _doesFirstStep.
double getTimeWindowStart() const final override
static const int INFINITE_MAX_ITERATIONS
To be used, when the number of max iterations is infinite (for implicit coupling).
Action
Actions that are required by CouplingSchemes.
@ WriteCheckpoint
Is the participant required to write a checkpoint?
@ ReadCheckpoint
Is the participant required to read a previously written checkpoint?
@ InitializeData
Is the initialization of coupling data required?
static const double UNDEFINED_MAX_TIME
Does not define a time limit for the coupled simulation.
static const int UNDEFINED_MIN_ITERATIONS
To be used, when the number of min iterations is not defined (for explicit coupling).
static const int UNDEFINED_MAX_ITERATIONS
To be used, when the number of max iterations is not defined (for explicit coupling).
static std::string toString(Action action)
static const double UNDEFINED_TIME_WINDOW_SIZE
To be used, when the time window size is determined dynamically during the coupling.
static const int UNDEFINED_TIME_WINDOWS
Does not define limit on time windows for the coupled simulation.
void resetProgress()
Resets the progress of the time window back to 0.
void resetTo(double timeStart)
Resets the handler to the given time.
bool reachedEndOfWindow(double timeWindowSize) const
double time() const
Returns the current time as a double.
double untilWindowEnd(double timeWindowSize) const
Returns the time distance to the possibly truncated end of the current time window.
void completeTimeWindow(double timeWindowSize)
double windowStart() const
Returns the window start as a double.
void progressBy(double dt)
Progress the time window by the given amount.
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
constexpr size_type size() const noexcept
Definition span.hpp:469
static bool isSecondary()
True if this process is running a secondary rank.
Definition IntraComm.cpp:57
T clear(T... args)
T count(T... args)
T emplace(T... args)
T empty(T... args)
T end(T... args)
T get(T... args)
T insert(T... args)
T make_unique(T... args)
contains implementations of coupling schemes for coupled simulations.
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
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)
static constexpr FundamentalTag Fundamental
Convenience instance of the FundamentalTag.
Definition Event.hpp:19
int DataID
Definition Types.hpp:25
STL namespace.
T set_difference(T... args)
T str(T... args)
Holds meta information to perform a convergence measurement.
void add(DataID did, bool toKeep)
Eigen::MatrixXd gradients
The gradients of the data. Use gradients.col(d*i+k) to get the gradient of vertex i,...
Definition Sample.hpp:67
Eigen::VectorXd values
Definition Sample.hpp:64