preCICE v3.1.2
Loading...
Searching...
No Matches
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
20#include "io/TXTTableWriter.hpp"
21#include "logging/LogMacros.hpp"
22#include "math/differences.hpp"
23#include "mesh/Data.hpp"
24#include "mesh/Mesh.hpp"
27#include "utils/Helpers.hpp"
28#include "utils/IntraComm.hpp"
29#include "utils/assertion.hpp"
30
32
34 double maxTime,
35 int maxTimeWindows,
36 double timeWindowSize,
37 std::string localParticipant,
38 int minIterations,
39 int maxIterations,
40 CouplingMode cplMode,
42 : _couplingMode(cplMode),
43 _maxTime(maxTime),
44 _maxTimeWindows(maxTimeWindows),
45 _timeWindows(1),
46 _timeWindowSize(timeWindowSize),
47 _nextTimeWindowSize(timeWindowSize),
48 _minIterations(minIterations),
49 _maxIterations(maxIterations),
50 _iterations(1),
51 _totalIterations(1),
52 _localParticipant(std::move(localParticipant))
53{
54 PRECICE_ASSERT(not((maxTime != UNDEFINED_MAX_TIME) && (maxTime < 0.0)),
55 "Maximum time has to be larger than zero.");
56 PRECICE_ASSERT(not((maxTimeWindows != UNDEFINED_TIME_WINDOWS) && (maxTimeWindows < 0)),
57 "Maximum number of time windows has to be larger than zero.");
58 PRECICE_ASSERT(not(hasTimeWindowSize() && (timeWindowSize < 0.0)),
59 "Time window size has to be larger than zero.");
60 if (dtMethod == constants::FIXED_TIME_WINDOW_SIZE) {
62 "Time window size has to be given when the fixed time window size method is used.");
63 }
64
68 } else {
72
73 PRECICE_ASSERT(minIterations > 0,
74 minIterations,
75 "Minimal iteration limit has to be larger than zero.");
76 PRECICE_ASSERT((maxIterations == INFINITE_MAX_ITERATIONS) || (maxIterations > 0),
77 maxIterations,
78 "Maximal iteration limit has to be larger than zero or -1 (unlimited).");
79 PRECICE_ASSERT((maxIterations == INFINITE_MAX_ITERATIONS) || (minIterations <= maxIterations),
80 "Minimal iteration limit has to be smaller equal compared to the maximal iteration limit.");
81 }
82
83 if (math::equals(maxTime, UNDEFINED_MAX_TIME)) {
85 } else {
86 _time = impl::TimeHandler(maxTime);
87 }
88}
89
95
97{
98 return _hasConverged;
99}
100
101void BaseCouplingScheme::sendNumberOfTimeSteps(const m2n::PtrM2N &m2n, const int numberOfTimeSteps)
102{
104 PRECICE_DEBUG("Sending number or time steps {}...", numberOfTimeSteps);
105 m2n->send(numberOfTimeSteps);
106}
107
108void BaseCouplingScheme::sendTimes(const m2n::PtrM2N &m2n, const Eigen::VectorXd &times)
109{
111 PRECICE_DEBUG("Sending times...");
112 m2n->send(times);
113}
114
115void BaseCouplingScheme::sendData(const m2n::PtrM2N &m2n, const DataMap &sendData)
116{
118 PRECICE_ASSERT(m2n.get() != nullptr);
119 PRECICE_ASSERT(m2n->isConnected());
120
121 for (const auto &data : sendData | boost::adaptors::map_values) {
122 const auto &stamples = data->stamples();
123 PRECICE_ASSERT(stamples.size() > 0);
124
125 int nTimeSteps = data->timeStepsStorage().nTimes();
126 PRECICE_ASSERT(nTimeSteps > 0);
127
128 if (data->exchangeSubsteps()) {
129 const Eigen::VectorXd timesAscending = data->timeStepsStorage().getTimes();
130 sendNumberOfTimeSteps(m2n, nTimeSteps);
131 sendTimes(m2n, timesAscending);
132
133 const auto serialized = com::serialize::SerializedStamples::serialize(data);
134
135 // Data is actually only send if size>0, which is checked in the derived classes implementation
136 m2n->send(serialized.values(), data->getMeshID(), data->getDimensions() * serialized.nTimeSteps());
137
138 if (data->hasGradient()) {
139 m2n->send(serialized.gradients(), data->getMeshID(), data->getDimensions() * data->meshDimensions() * serialized.nTimeSteps());
140 }
141 } else {
142 data->sample() = stamples.back().sample;
143
144 // Data is only received on ranks with size>0, which is checked in the derived class implementation
145 m2n->send(data->values(), data->getMeshID(), data->getDimensions());
146
147 if (data->hasGradient()) {
148 PRECICE_ASSERT(data->hasGradient());
149 m2n->send(data->gradients(), data->getMeshID(), data->getDimensions() * data->meshDimensions());
150 }
151 }
152 }
153}
154
156{
158 PRECICE_DEBUG("Receiving number of time steps...");
159 int numberOfTimeSteps;
160 m2n->receive(numberOfTimeSteps);
161 return numberOfTimeSteps;
162}
163
164Eigen::VectorXd BaseCouplingScheme::receiveTimes(const m2n::PtrM2N &m2n, int nTimeSteps)
165{
167 PRECICE_DEBUG("Receiving times....");
168 Eigen::VectorXd times(nTimeSteps);
169 m2n->receive(times);
170 PRECICE_DEBUG("Received times {}", times);
171 return times;
172}
173
174void BaseCouplingScheme::receiveData(const m2n::PtrM2N &m2n, const DataMap &receiveData)
175{
177 PRECICE_ASSERT(m2n.get());
178 PRECICE_ASSERT(m2n->isConnected());
179 for (const auto &data : receiveData | boost::adaptors::map_values) {
180
181 if (data->exchangeSubsteps()) {
182 const int nTimeSteps = receiveNumberOfTimeSteps(m2n);
183
184 Eigen::VectorXd serializedValues(nTimeSteps * data->getSize());
185 PRECICE_ASSERT(nTimeSteps > 0);
186 const Eigen::VectorXd timesAscending = receiveTimes(m2n, nTimeSteps);
187
188 auto serialized = com::serialize::SerializedStamples::empty(timesAscending, data);
189
190 // Data is only received on ranks with size>0, which is checked in the derived class implementation
191 m2n->receive(serialized.values(), data->getMeshID(), data->getDimensions() * nTimeSteps);
192
193 if (data->hasGradient()) {
194 m2n->receive(serialized.gradients(), data->getMeshID(), data->getDimensions() * data->meshDimensions() * nTimeSteps);
195 }
196
197 serialized.deserializeInto(timesAscending, data);
198 } else {
199 // Data is only received on ranks with size>0, which is checked in the derived class implementation
200 m2n->receive(data->values(), data->getMeshID(), data->getDimensions());
201
202 if (data->hasGradient()) {
203 PRECICE_ASSERT(data->hasGradient());
204 m2n->receive(data->gradients(), data->getMeshID(), data->getDimensions() * data->meshDimensions());
205 }
206 data->setSampleAtTime(getTime(), data->sample());
207 }
208 }
209}
210
212{
213 // @TODO This is a hack required until https://github.com/precice/precice/issues/1957
214 // buffer current time
215 const impl::TimeHandler oldTime = _time;
216 // set _time state to point to end of this window such that _time in receiveData is at end of window
218 // receive data for end of window
219 this->receiveData(m2n, receiveData);
220 // reset time state;
221 _time = oldTime;
222}
223
225{
226 for (const auto &data : receiveData | boost::adaptors::map_values) {
227 PRECICE_DEBUG("Initialize {} as zero.", data->getDataName());
228 // just store already initialized zero sample to storage.
229 data->setSampleAtTime(getTime(), data->sample());
230 }
231}
232
233PtrCouplingData BaseCouplingScheme::addCouplingData(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool communicateSubsteps, CouplingData::Direction direction)
234{
235 int id = data->getID();
236 PtrCouplingData ptrCplData;
237 if (!utils::contained(id, _allData)) { // data is not used by this coupling scheme yet, create new CouplingData
238 ptrCplData = std::make_shared<CouplingData>(data, std::move(mesh), requiresInitialization, communicateSubsteps, direction);
239 _allData.emplace(id, ptrCplData);
240 } else { // data is already used by another exchange of this coupling scheme, use existing CouplingData
241 ptrCplData = _allData[id];
242 PRECICE_CHECK(ptrCplData->getDirection() == direction, "Data \"{0}\" cannot be added for sending and for receiving. Please remove either <exchange data=\"{0}\" ... /> tag", data->getName());
243 }
244 return ptrCplData;
245}
246
252
253void BaseCouplingScheme::setTimeWindowSize(double timeWindowSize)
254{
255 _timeWindowSize = timeWindowSize;
256}
257
259{
260 _nextTimeWindowSize = timeWindowSize;
261}
262
264{
267 PRECICE_ASSERT(_isInitialized, "Called finalize() before initialize().");
268}
269
270void BaseCouplingScheme::initialize(double startTime, int startTimeWindow)
271{
272 // initialize with zero data here, might eventually be overwritten in exchangeInitialData
274 // Initialize uses the template method pattern (https://en.wikipedia.org/wiki/Template_method_pattern).
275 PRECICE_TRACE(startTime, startTimeWindow);
277 PRECICE_ASSERT(math::greaterEquals(startTime, 0.0), startTime);
278 PRECICE_ASSERT(startTimeWindow >= 0, startTimeWindow);
279 _time.resetTo(startTime);
280 _timeWindows = startTimeWindow;
281 _hasDataBeenReceived = false;
282
285 if (not doesFirstStep()) {
286 // reserve memory and initialize data with zero
287 if (_acceleration) {
288 _acceleration->initialize(getAccelerationData());
289 }
290 }
293 }
294
296
297 _isInitialized = true;
298}
299
304
310
312{
315 PRECICE_ASSERT(_isInitialized, "Before calling advance() coupling scheme has to be initialized via initialize().");
316 _hasDataBeenReceived = false;
317 _isTimeWindowComplete = false;
318
320
322
323 _timeWindows += 1; // increment window counter. If not converged, will be decremented again later.
324
326 }
327}
328
333
335{
338 PRECICE_ASSERT(_isInitialized, "Before calling advance() coupling scheme has to be initialized via initialize().");
340
341 // from first phase
343
345
347
348 if (isImplicitCouplingScheme()) { // check convergence
349 if (not hasConverged()) { // repeat window
350 PRECICE_DEBUG("No convergence achieved");
352 // The computed time window part equals the time window size, since the
353 // time window remainder is zero. Subtract the time window size and do another
354 // coupling iteration.
356 _timeWindows -= 1;
357 _isTimeWindowComplete = false;
358 } else { // write output, prepare for next window
359 PRECICE_DEBUG("Convergence achieved");
361 PRECICE_INFO("Time window completed");
363 if (isCouplingOngoing()) {
364 PRECICE_DEBUG("Setting require create checkpoint");
366 }
367 }
368 //update iterations
370 if (not hasConverged()) {
371 _iterations++;
372 } else {
373 _iterations = 1;
374 }
375 } else {
377 PRECICE_INFO("Time window completed");
379 }
380 if (isCouplingOngoing()) {
382 }
383
384 // Update internal time tracking
387 } else {
389 }
391 }
392}
393
395{
397 for (auto &data : _allData | boost::adaptors::map_values) {
398 data->moveToNextWindow();
399 }
400}
401
406
412
417
419{
420 return _isInitialized;
421}
422
424 double timeToAdd)
425{
426 PRECICE_TRACE(timeToAdd, getTime());
427 PRECICE_ASSERT(isCouplingOngoing(), "Invalid call of addComputedTime() after simulation end.");
428
429 // Check validness
430 bool valid = math::greaterEquals(getNextTimeStepMaxSize(), timeToAdd);
431 PRECICE_CHECK(valid,
432 "The time step size given to preCICE in \"advance\" {} exceeds the maximum allowed time step size {} "
433 "in the remaining of this time window. "
434 "Did you restrict your time step size, \"dt = min(preciceDt, solverDt)\"? "
435 "For more information, consult the adapter example in the preCICE documentation.",
436 timeToAdd, getNextTimeStepMaxSize());
437
438 // add time interval that has been computed in the solver to get the correct time remainder
439 _time.progressBy(timeToAdd);
440
441 return reachedEndOfTimeWindow();
442}
443
445 double lastSolverTimeStepSize) const
446{
447 PRECICE_TRACE(lastSolverTimeStepSize);
448 double remainder = getNextTimeStepMaxSize() - lastSolverTimeStepSize;
449 return not math::greater(remainder, 0.0);
450}
451
456
458{
460}
461
463{
464 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.");
466}
467
472
474{
475 _timeWindows = timeWindows;
476}
477
479{
480 return _time.time();
481}
482
484{
485 return _time.windowStart();
486}
487
489{
490 return _timeWindows;
491}
492
494{
495 if (!isCouplingOngoing()) {
496 return 0.0; // if coupling is not ongoing (i.e. coupling scheme reached end of window) the maximum time step size is zero
497 }
498
499 if (hasTimeWindowSize()) {
501 } else {
502 return _time.untilEnd();
503 }
504}
505
507{
508 bool timestepsLeft = (_maxTimeWindows >= _timeWindows) || (_maxTimeWindows == UNDEFINED_TIME_WINDOWS);
509 return timestepsLeft && !_time.reachedEnd();
510}
511
516
518 Action action) const
519{
520 return _requiredActions.count(action) == 1;
521}
522
524 Action action) const
525{
526 return _fulfilledActions.count(action) == 1;
527}
528
535
537 Action action)
538{
539 _requiredActions.insert(action);
540}
541
543{
545 os << "iteration: " << _iterations; //_iterations;
547 os << " of " << _maxIterations;
548 }
550 os << " (min " << _minIterations << ")";
551 }
552 os << ", " << printBasicState(_timeWindows, getTime()) << ", " << printActionsState();
553 return os.str();
554}
555
557 int timeWindows,
558 double time) const
559{
561 os << "time-window: " << timeWindows;
563 os << " of " << _maxTimeWindows;
564 }
565 os << ", time: " << time;
567 os << " of " << _maxTime;
568 }
569 if (hasTimeWindowSize()) {
570 os << ", time-window-size: " << _timeWindowSize;
571 }
573 os << ", max-time-step-size: " << getNextTimeStepMaxSize();
574 }
575 os << ", ongoing: ";
576 isCouplingOngoing() ? os << "yes" : os << "no";
577 os << ", time-window-complete: ";
578 _isTimeWindowComplete ? os << "yes" : os << "no";
579 return os.str();
580}
581
583{
585 for (auto action : _requiredActions) {
586 os << toString(action) << ' ';
587 }
588 return os.str();
589}
590
592{
594 std::vector<Action> missing;
597 std::back_inserter(missing));
598 if (not missing.empty()) {
599 std::ostringstream stream;
600 for (auto action : missing) {
601 if (not stream.str().empty()) {
602 stream << ", ";
603 }
604 stream << toString(action);
605 }
606 PRECICE_ERROR("The required actions {} are not fulfilled. "
607 "Did you forget to call \"requiresReadingCheckpoint()\" or \"requiresWritingCheckpoint()\"?",
608 stream.str());
609 }
612}
613
615 const acceleration::PtrAcceleration &acceleration)
616{
617 PRECICE_ASSERT(acceleration.get() != nullptr);
618 _acceleration = acceleration;
619}
620
622{
623 return _doesFirstStep;
624}
625
627{
630 PRECICE_ASSERT(convMeasure.measure.get() != nullptr);
631 convMeasure.measure->newMeasurementSeries();
632 }
633}
634
636 int dataID,
637 bool suffices,
638 bool strict,
640 bool doesLogging)
641{
642 ConvergenceMeasureContext convMeasure;
643 PRECICE_ASSERT(_allData.count(dataID) == 1, "Data with given data ID must exist!");
644 convMeasure.couplingData = _allData.at(dataID);
645 convMeasure.suffices = suffices;
646 convMeasure.strict = strict;
647 convMeasure.measure = std::move(measure);
648 convMeasure.doesLogging = doesLogging;
649 _convergenceMeasures.push_back(convMeasure);
650}
651
653{
657 _convergenceWriter->writeData("TimeWindow", _timeWindows - 1);
658 _convergenceWriter->writeData("Iteration", _iterations);
659 }
660
661 // If no convergence measures are defined, we never converge
662 if (_convergenceMeasures.empty()) {
663 PRECICE_INFO("No converge measures defined.");
664 return false;
665 }
666
667 // There are convergence measures defined, so we need to check them
668 bool allConverged = true;
669 bool oneSuffices = false; // at least one convergence measure suffices and did converge
670 bool oneStrict = false; // at least one convergence measure is strict and did not converge
671
672 const bool reachedMinIterations = _iterations >= _minIterations;
673 for (const auto &convMeasure : _convergenceMeasures) {
674 PRECICE_ASSERT(convMeasure.couplingData != nullptr);
675 PRECICE_ASSERT(convMeasure.measure.get() != nullptr);
676 PRECICE_ASSERT(convMeasure.couplingData->previousIteration().size() == convMeasure.couplingData->values().size(), convMeasure.couplingData->previousIteration().size(), convMeasure.couplingData->values().size(), convMeasure.couplingData->getDataName());
677 convMeasure.measure->measure(convMeasure.couplingData->previousIteration(), convMeasure.couplingData->values());
678
679 if (not utils::IntraComm::isSecondary() && convMeasure.doesLogging) {
680 _convergenceWriter->writeData(convMeasure.logHeader(), convMeasure.measure->getNormResidual());
681 }
682
683 if (not convMeasure.measure->isConvergence()) {
684 allConverged = false;
685 if (convMeasure.strict) {
687 oneStrict = true;
689 "The strict convergence measure for data \"" + convMeasure.couplingData->getDataName() +
690 "\" did not converge within the maximum allowed iterations, which terminates the simulation. "
691 "To avoid this forced termination do not mark the convergence measure as strict.");
692 }
693 } else if (convMeasure.suffices == true) {
694 oneSuffices = true;
695 }
696
697 PRECICE_INFO(convMeasure.measure->printState(convMeasure.couplingData->getDataName()));
698 }
699
700 std::string messageSuffix;
701 if (not reachedMinIterations) {
702 messageSuffix = " but hasn't yet reached minimal amount of iterations";
703 }
704 if (allConverged) {
705 PRECICE_INFO("All converged{}", messageSuffix);
706 } else if (oneSuffices && not oneStrict) { // strict overrules suffices
707 PRECICE_INFO("Sufficient measures converged{}", messageSuffix);
708 }
709
710 return reachedMinIterations && (allConverged || (oneSuffices && not oneStrict));
711}
712
714{
716
717 _iterationsWriter = std::make_shared<io::TXTTableWriter>("precice-" + _localParticipant + "-iterations.log");
718 if (not doesFirstStep()) {
719 _convergenceWriter = std::make_shared<io::TXTTableWriter>("precice-" + _localParticipant + "-convergence.log");
720 }
721
722 _iterationsWriter->addData("TimeWindow", io::TXTTableWriter::INT);
723 _iterationsWriter->addData("TotalIterations", io::TXTTableWriter::INT);
724 _iterationsWriter->addData("Iterations", io::TXTTableWriter::INT);
725 _iterationsWriter->addData("Convergence", io::TXTTableWriter::INT);
726
727 if (not doesFirstStep()) {
728 _convergenceWriter->addData("TimeWindow", io::TXTTableWriter::INT);
729 _convergenceWriter->addData("Iteration", io::TXTTableWriter::INT);
730 }
731
732 if (not doesFirstStep()) {
734
735 if (convMeasure.doesLogging) {
736 _convergenceWriter->addData(convMeasure.logHeader(), io::TXTTableWriter::DOUBLE);
737 }
738 }
739 if (_acceleration) {
740 _iterationsWriter->addData("QNColumns", io::TXTTableWriter::INT);
741 _iterationsWriter->addData("DeletedQNColumns", io::TXTTableWriter::INT);
742 _iterationsWriter->addData("DroppedQNColumns", io::TXTTableWriter::INT);
743 }
744 }
745 }
746}
747
749{
751
752 _iterationsWriter->writeData("TimeWindow", _timeWindows - 1);
753 _iterationsWriter->writeData("TotalIterations", _totalIterations);
754 _iterationsWriter->writeData("Iterations", _iterations);
755 bool converged = _iterations >= _minIterations && (_maxIterations < 0 || (_iterations < _maxIterations));
756 _iterationsWriter->writeData("Convergence", converged ? 1 : 0);
757
758 if (not doesFirstStep() && _acceleration) {
759 _iterationsWriter->writeData("QNColumns", _acceleration->getLSSystemCols());
760 _iterationsWriter->writeData("DeletedQNColumns", _acceleration->getDeletedColumns());
761 _iterationsWriter->writeData("DroppedQNColumns", _acceleration->getDroppedColumns());
762 }
763 }
764}
765
767{
768 if (!hasTimeWindowSize()) {
769 return true; // This participant will always do exactly one step to dictate second's time-window-size
770 }
771
773}
774
776{
778 for (const auto &data : _allData | boost::adaptors::map_values) {
779 data->storeIteration();
780 }
781}
782
790
797
799{
801 for (const auto &data : dataMap | boost::adaptors::map_values) {
802 if (data->requiresInitialization) {
803 return true;
804 }
805 }
806 return false;
807}
808
810{
811 PRECICE_DEBUG("measure convergence of the coupling iteration");
813 // Stop, when maximal iteration count (given in config) is reached
815 _hasConverged = true;
816
817 // coupling iteration converged for current time window. Advance in time.
818 if (_hasConverged) {
819 if (_acceleration) {
820 _acceleration->iterationsConverged(getAccelerationData());
821 }
823 } else {
824 // no convergence achieved for the coupling iteration within the current time window
825 if (_acceleration) {
826 // Acceleration works on CouplingData::values(), so we retrieve the data from the storage, perform the acceleration and then put the data back into the storage. See also https://github.com/precice/precice/issues/1645.
827 // @todo For acceleration schemes as described in "Rüth, B, Uekermann, B, Mehl, M, Birken, P, Monge, A, Bungartz, H-J. Quasi-Newton waveform iteration for partitioned surface-coupled multiphysics applications. https://doi.org/10.1002/nme.6443" we need a more elaborate implementation.
828
829 // Load from storage into buffer
830 for (auto &data : getAccelerationData() | boost::adaptors::map_values) {
831 const auto &stamples = data->stamples();
832 PRECICE_ASSERT(stamples.size() > 0);
833 data->sample() = stamples.back().sample;
834 }
835
836 _acceleration->performAcceleration(getAccelerationData());
837
838 // Store from buffer
839 // @todo Currently only data at end of window is accelerated. Remaining data in storage stays as it is.
840 for (auto &data : getAccelerationData() | boost::adaptors::map_values) {
841 data->setSampleAtTime(getTime(), data->sample());
842 }
843 }
844 }
845}
846
848{
850 PRECICE_ASSERT(not doesFirstStep(), "For convergence information the sending participant is never the first one.");
851 m2n->send(_hasConverged);
852}
853
855{
857 PRECICE_ASSERT(doesFirstStep(), "For convergence information the receiving participant is always the first one.");
858 m2n->receive(_hasConverged);
859}
860
862{
863 return _time.windowStart();
864}
865
870
872{
873 // Global toggle if a single send data uses substeps
874 for (auto cpldata : _allData | boost::adaptors::map_values) {
875 if (cpldata->getDirection() == CouplingData::Direction::Send && cpldata->exchangeSubsteps()) {
876 return true;
877 }
878 }
879 return false;
880}
881
883{
884 // This implementation covers everything except SerialImplicit
886 return {};
887 }
888
889 ImplicitData idata;
890 for (auto cpldata : _allData | boost::adaptors::map_values) {
891 if (cpldata->getDirection() == CouplingData::Direction::Receive) {
892 idata.add(cpldata->getDataID(), false);
893 }
894 }
895 return idata;
896}
897
902
903} // namespace precice::cplscheme
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:15
#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
T at(T... args)
T back_inserter(T... args)
T begin(T... args)
static SerializedStamples empty(Eigen::VectorXd timeStamps, const cplscheme::PtrCouplingData data)
Create SerializedStamples with allocated buffers according to size of CouplingData.
static SerializedStamples serialize(const cplscheme::PtrCouplingData data)
Serializes a given CouplingData into SerializedStamples.
void initializeWithZeroInitialData(const DataMap &receiveData)
Initializes storage in receiveData as zero.
std::string printCouplingState() const override
Returns coupling state information.
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 isInitialized() const override final
getter for _isInitialized
ChangedMeshes secondSynchronization() override final
std::string _localParticipant
Local participant name.
double _timeWindowSize
size of time window; _timeWindowSize <= _maxTime
void determineInitialSend(DataMap &sendData)
Sets _sendsInitializedData, if sendData requires initialization.
virtual void exchangeSecondData()=0
Exchanges the second set of data.
ChangedMeshes firstSynchronization(const ChangedMeshes &changes) override final
void receiveDataForWindowEnd(const m2n::PtrM2N &m2n, const DataMap &receiveData)
Like receiveData, but temporarily sets window time to end of window.
bool requiresSubsteps() const override final
Returns true if any send data of the scheme requires substeps.
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().
void addConvergenceMeasure(int dataID, bool suffices, bool strict, impl::PtrConvergenceMeasure measure, bool doesLogging)
Adds a measure to determine the convergence of coupling iterations.
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 initialize(double startTime, int startTimeWindow) override final
Initializes the coupling scheme.
void storeIteration()
used for storing all Data at end of doImplicitStep for later reference.
bool sendsInitializedData() const override final
Getter for _sendsInitializedData.
void initializeTXTWriters()
Initialize txt writers for iterations and convergence tracking.
bool addComputedTime(double timeToAdd) override final
Adds newly computed time. Has to be called before every advance.
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...
void requireAction(Action action) override final
Sets an action required to be performed by the accessor.
bool _sendsInitializedData
True, if this participant has to send initialized data.
bool _isTimeWindowComplete
True, if _time == _timeWindowStartTime + _timeWindowSize and (coupling has converged or _iterations =...
bool isCouplingOngoing() const override final
Returns true, when the coupled simulation is still ongoing.
int _totalIterations
Number of total iterations performed.
bool anyDataRequiresInitialization(DataMap &dataMap) const
Checks whether any CouplingData in dataMap requires initialization.
double getTime() const override final
getter for _time
Eigen::VectorXd receiveTimes(const m2n::PtrM2N &m2n, int nTimeSteps)
bool hasDataBeenReceived() const override final
getter for _hasDataBeenReceived
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.
double getTimeWindowSize() const override final
Returns the time window size, if one is given by the coupling scheme.
double getTimeWindowStart() const override final
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.
int receiveNumberOfTimeSteps(const m2n::PtrM2N &m2n)
void moveToNextWindow()
finalizes this window's data and initializes data for next window.
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 getNextTimeStepMaxSize() const override final
Returns the maximal size of the next time step to be computed.
bool isTimeWindowComplete() const override final
Returns true, when the accessor can advance to the next time window.
std::string printActionsState() const
Prints the action state.
bool hasTimeWindowSize() const override final
Function to check whether time window size is defined by coupling scheme.
bool receivesInitializedData() const
Getter for _receivesInitializedData.
int getTimeWindows() const override final
getter for _timeWindows
ImplicitData implicitDataToReceive() const override
Returns a vector of implicit data to receive in the next advance.
void sendNumberOfTimeSteps(const m2n::PtrM2N &m2n, const int numberOfTimeSteps)
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 advanceTXTWriters()
Advance txt writers for iterations and convergence tracking.
acceleration::PtrAcceleration _acceleration
Acceleration method to speedup iteration convergence.
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 sendTimes(const m2n::PtrM2N &m2n, const Eigen::VectorXd &times)
CouplingMode _couplingMode
Coupling mode used by coupling scheme.
bool willDataBeExchanged(double lastSolverTimeStepSize) const override final
Returns true, if data will be exchanged when calling advance().
bool isActionRequired(Action action) const override final
Returns true, if the given action has to be performed by the accessor.
std::string localParticipant() const override final
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.
virtual void exchangeInitialData()=0
implements functionality for initialize in base class.
int _iterations
Number of iterations in current time window. _iterations <= _maxIterations.
void markActionFulfilled(Action action) override final
Tells the coupling scheme that the accessor has performed the given action.
virtual void initializeReceiveDataStorage()=0
Functions needed for initialize()
bool isActionFulfilled(Action action) const override final
Returns true, if the given action has to be performed by the accessor.
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.
void finalize() override final
Finalizes the coupling scheme.
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.
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.
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)
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)
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
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)