preCICE v3.3.0
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
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,
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),
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
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
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
212 _time.progressBy(_nextTimeWindowSize);
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{
272 PRECICE_ASSERT(_isInitialized, "Called finalize() before initialize().");
273}
274
276{
277 // initialize with zero data here, might eventually be overwritten in exchangeInitialData
279 // Initialize uses the template method pattern (https://en.wikipedia.org/wiki/Template_method_pattern).
282 _time.resetTo(0);
283 _timeWindows = 1;
284 _hasDataBeenReceived = false;
285
287
290 if (not doesFirstStep()) {
291 // reserve memory and initialize data with zero
292 if (_acceleration) {
293 _acceleration->initialize(getAccelerationData());
295 qnAcceleration) {
296 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.");
297 }
298 }
299 }
302 }
303
305
306 _isInitialized = true;
307}
308
310{
313
315 // overwrite past iteration with new samples
316 for (const auto &data : _allData | boost::adaptors::map_values) {
317 // TODO: reset CouplingData of changed meshes only #2102
318 data->reinitialize();
319 }
320
321 if (not doesFirstStep()) {
322 if (_acceleration) {
323 _acceleration->initialize(getAccelerationData());
324 }
325 }
327 }
328}
329
334
340
342{
345 PRECICE_ASSERT(_isInitialized, "Before calling advance() coupling scheme has to be initialized via initialize().");
346 _hasDataBeenReceived = false;
347 _isTimeWindowComplete = false;
348
350
352
353 _timeWindows += 1; // increment window counter. If not converged, will be decremented again later.
354
356 }
357}
358
363
365{
368 PRECICE_ASSERT(_isInitialized, "Before calling advance() coupling scheme has to be initialized via initialize().");
370
371 // from first phase
373
375
377
378 if (isImplicitCouplingScheme()) { // check convergence
379 if (not hasConverged()) { // repeat window
380 PRECICE_DEBUG("No convergence achieved");
382 // The computed time window part equals the time window size, since the
383 // time window remainder is zero. Subtract the time window size and do another
384 // coupling iteration.
386 _timeWindows -= 1;
387 _isTimeWindowComplete = false;
388 } else { // write output, prepare for next window
389 PRECICE_DEBUG("Convergence achieved");
391 PRECICE_INFO("Time window completed");
393 if (isCouplingOngoing()) {
394 PRECICE_DEBUG("Setting require create checkpoint");
396 }
397 }
398 // update iterations
400 if (not hasConverged()) {
401 _iterations++;
402 } else {
403 _iterations = 1;
404 }
405 } else {
407 PRECICE_INFO("Time window completed");
409 }
410 if (isCouplingOngoing()) {
412 }
413
414 // Update internal time tracking
416 _time.completeTimeWindow(_timeWindowSize);
417 } else {
418 _time.resetProgress();
419 }
421 }
422}
423
425{
427 for (auto &data : _allData | boost::adaptors::map_values) {
428 data->moveToNextWindow();
429 }
430}
431
436
442
447
449{
450 return _isInitialized;
451}
452
454 double timeToAdd)
455{
456 PRECICE_TRACE(timeToAdd, getTime());
457 PRECICE_ASSERT(isCouplingOngoing(), "Invalid call of addComputedTime() after simulation end.");
458
459 if (!hasTimeWindowSize()) {
461 "advance() was called with the max value of double which is not allowed. "
462 "As this participant prescribes the time-window size using <time-window-size method=\"first-participant\" />, directly using getMaxTimeStepSize() is not permitted. "
463 "Make sure to pass your own desired time-step size or use the recommended limiting \"dt = min(solver_dt, getMaxTimeStepSize())\".");
464 }
465
466 // Check validness
467 bool valid = math::greaterEquals(getNextTimeStepMaxSize(), timeToAdd);
468 PRECICE_CHECK(valid,
469 "The time step size given to preCICE in \"advance\" {} exceeds the maximum allowed time step size {} "
470 "in the remaining of this time window. "
471 "Did you restrict your time step size, \"dt = min(preciceDt, solverDt)\"? "
472 "For more information, consult the adapter example in the preCICE documentation.",
473 timeToAdd, getNextTimeStepMaxSize());
474
475 // add time interval that has been computed in the solver to get the correct time remainder
476 _time.progressBy(timeToAdd);
477
478 return reachedEndOfTimeWindow();
479}
480
482 double lastSolverTimeStepSize) const
483{
484 PRECICE_TRACE(lastSolverTimeStepSize);
485 double remainder = getNextTimeStepMaxSize() - lastSolverTimeStepSize;
486 return not math::greater(remainder, 0.0);
487}
488
493
498
500{
501 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.");
503}
504
509
511{
512 _timeWindows = timeWindows;
513}
514
516{
517 return _time.time();
518}
519
521{
522 return _time.windowStart();
523}
524
526{
527 return _time.windowProgress();
528}
529
531{
532 return _timeWindows;
533}
534
536{
537 if (!isCouplingOngoing()) {
538 return 0.0; // if coupling is not ongoing (i.e. coupling scheme reached end of window) the maximum time step size is zero
539 }
540
541 if (hasTimeWindowSize()) {
542 return _time.untilWindowEnd(_timeWindowSize);
543 } else {
544 return _time.untilEnd();
545 }
546}
547
549{
550 bool timestepsLeft = (_maxTimeWindows >= _timeWindows) || (_maxTimeWindows == UNDEFINED_TIME_WINDOWS);
551 return timestepsLeft && !_time.reachedEnd();
552}
553
558
560 Action action) const
561{
562 return _requiredActions.count(action) == 1;
563}
564
566 Action action) const
567{
568 return _fulfilledActions.count(action) == 1;
569}
570
577
583
585{
586 if (!isCouplingOngoing()) {
587 return fmt::format("Reached end at: final time-window: {}, final time: {}", (_timeWindows - 1), getTime());
588 }
589 std::string str;
590 auto out = std::back_inserter(str);
591
595
596 if (hasMax && hasMin) {
597 fmt::format_to(out, "it {} (min: {}, max: {}), ", _iterations, _minIterations, _maxIterations);
598 } else if (hasMax) {
599 fmt::format_to(out, "it {} (max: {}), ", _iterations, _maxIterations);
600 } else if (hasMin) {
601 fmt::format_to(out, "it {} (min: {}), ", _iterations, _minIterations);
602 } else {
603 fmt::format_to(out, "it {}, ", _iterations);
604 }
605 }
606
607 fmt::format_to(out, "time-window {}", _timeWindows);
609 fmt::format_to(out, " (max: {})", _maxTimeWindows);
610 }
611 fmt::format_to(out, ", t {}", getTime());
613 fmt::format_to(out, " (max: {})", _maxTime);
614 }
615 if (hasTimeWindowSize()) {
616 fmt::format_to(out, ", Dt {}, max-dt {}", _timeWindowSize, getNextTimeStepMaxSize());
617 }
618 return str;
619}
620
622{
625 "Data was not initialized. Did you forget to call \"requiresWritingData()\"?");
626
628 "The iteration checkpoint wasn't read. Did you forget to call \"requiresReadingCheckpoint()\"?");
629
631 "The iteration checkpoint wasn't written. Did you forget to call \"requiresWritingCheckpoint()\"?");
632
633 _requiredActions.clear();
634 _fulfilledActions.clear();
635}
636
638{
640 for (const auto &data : _allData | boost::adaptors::map_values) {
641 if (data->getDirection() == CouplingData::Direction::Receive) {
642 PRECICE_ASSERT(!data->stamples().empty(), "initializeReceiveDataStorage() didn't initialize data correctly");
643 } else {
644 PRECICE_CHECK(!data->stamples().empty(),
645 "Data {0} on mesh {1} doesn't contain any samples while initializing a coupling scheme of participant {2}. "
646 "There are two common configuration issues that may cause this. "
647 "Either, make sure participant {2} specifies data {0} to be written using tag <write-data mesh=\"{1}\" data=\"{0}\"/>. "
648 "Or ensure participant {2} defines a mapping to mesh {1} from a mesh using data {0}.",
649 data->getDataName(), data->getMeshName(), localParticipant());
650 }
651 }
652}
653
660
662{
663 return _doesFirstStep;
664}
665
667{
670 PRECICE_ASSERT(convMeasure.measure.get() != nullptr);
671 convMeasure.measure->newMeasurementSeries();
672 }
673}
674
676 int dataID,
677 bool suffices,
678 bool strict,
680{
681 ConvergenceMeasureContext convMeasure;
682 PRECICE_ASSERT(_allData.count(dataID) == 1, "Data with given data ID must exist!");
683 convMeasure.couplingData = _allData.at(dataID);
684 convMeasure.suffices = suffices;
685 convMeasure.strict = strict;
686 convMeasure.measure = std::move(measure);
687 _convergenceMeasures.push_back(convMeasure);
688}
689
691{
695 _convergenceWriter->writeData("TimeWindow", _timeWindows - 1);
696 _convergenceWriter->writeData("Iteration", _iterations);
697 }
698
699 // If no convergence measures are defined, we never converge
700 if (_convergenceMeasures.empty()) {
701 PRECICE_INFO("No converge measures defined.");
702 return false;
703 }
704
705 // There are convergence measures defined, so we need to check them
706 bool allConverged = true;
707 bool oneSuffices = false; // at least one convergence measure suffices and did converge
708 bool oneStrict = false; // at least one convergence measure is strict and did not converge
709
710 const bool reachedMinIterations = _iterations >= _minIterations;
711 for (const auto &convMeasure : _convergenceMeasures) {
712 PRECICE_ASSERT(convMeasure.couplingData != nullptr);
713 PRECICE_ASSERT(convMeasure.measure.get() != nullptr);
714 PRECICE_ASSERT(convMeasure.couplingData->previousIteration().size() == convMeasure.couplingData->values().size(), convMeasure.couplingData->previousIteration().size(), convMeasure.couplingData->values().size(), convMeasure.couplingData->getDataName());
715 convMeasure.measure->measure(convMeasure.couplingData->previousIteration(), convMeasure.couplingData->values());
716
718 _convergenceWriter->writeData(convMeasure.logHeader(), convMeasure.measure->getNormResidual());
719 }
720
721 if (not convMeasure.measure->isConvergence()) {
722 allConverged = false;
723 if (convMeasure.strict) {
725 oneStrict = true;
727 "The strict convergence measure for data \"" + convMeasure.couplingData->getDataName() +
728 "\" did not converge within the maximum allowed iterations, which terminates the simulation. "
729 "To avoid this forced termination do not mark the convergence measure as strict.");
730 }
731 } else if (convMeasure.suffices == true) {
732 oneSuffices = true;
733 }
734
735 PRECICE_INFO(convMeasure.measure->printState(convMeasure.couplingData->getDataName()));
736 }
737
738 std::string messageSuffix;
739 if (not reachedMinIterations) {
740 messageSuffix = " but hasn't yet reached minimal amount of iterations";
741 }
742 if (allConverged) {
743 PRECICE_INFO("All converged{}", messageSuffix);
744 } else if (oneSuffices && not oneStrict) { // strict overrules suffices
745 PRECICE_INFO("Sufficient measures converged{}", messageSuffix);
746 }
747
748 return reachedMinIterations && (allConverged || (oneSuffices && not oneStrict));
749}
750
752{
754
756 if (not doesFirstStep()) {
758 }
759
760 _iterationsWriter->addData("TimeWindow", io::TXTTableWriter::INT);
761 _iterationsWriter->addData("TotalIterations", io::TXTTableWriter::INT);
762 _iterationsWriter->addData("Iterations", io::TXTTableWriter::INT);
763 _iterationsWriter->addData("Convergence", io::TXTTableWriter::INT);
764
765 if (not doesFirstStep()) {
766 _convergenceWriter->addData("TimeWindow", io::TXTTableWriter::INT);
767 _convergenceWriter->addData("Iteration", io::TXTTableWriter::INT);
768 }
769
770 if (not doesFirstStep()) {
772 _convergenceWriter->addData(convMeasure.logHeader(), io::TXTTableWriter::DOUBLE);
773 }
774 if (_acceleration) {
775 _iterationsWriter->addData("QNColumns", io::TXTTableWriter::INT);
776 _iterationsWriter->addData("DeletedQNColumns", io::TXTTableWriter::INT);
777 _iterationsWriter->addData("DroppedQNColumns", io::TXTTableWriter::INT);
778 }
779 }
780 }
781}
782
784{
786
787 _iterationsWriter->writeData("TimeWindow", _timeWindows - 1);
788 _iterationsWriter->writeData("TotalIterations", _totalIterations);
789 _iterationsWriter->writeData("Iterations", _iterations);
790 bool converged = _iterations >= _minIterations && (_maxIterations < 0 || (_iterations < _maxIterations));
791 _iterationsWriter->writeData("Convergence", converged ? 1 : 0);
792
793 if (not doesFirstStep() && _acceleration) {
795 if (qnAcceleration) {
796 // Only write values for additional columns, if using a QN-based acceleration scheme
797 _iterationsWriter->writeData("QNColumns", qnAcceleration->getLSSystemCols());
798 _iterationsWriter->writeData("DeletedQNColumns", qnAcceleration->getDeletedColumns());
799 _iterationsWriter->writeData("DroppedQNColumns", qnAcceleration->getDroppedColumns());
800 } else {
801 // non-breaking implementation uses "0" for these columns (delete columns in the future?)
802 _iterationsWriter->writeData("QNColumns", 0);
803 _iterationsWriter->writeData("DeletedQNColumns", 0);
804 _iterationsWriter->writeData("DroppedQNColumns", 0);
805 }
806 }
807 }
808}
809
811{
812 if (!hasTimeWindowSize()) {
813 return true; // This participant will always do exactly one step to dictate second's time-window-size
814 }
815
816 return _time.reachedEndOfWindow(_timeWindowSize);
817}
818
820{
822
823 for (const auto &data : _allData | boost::adaptors::map_values) {
824 PRECICE_CHECK(!data->stamples().empty(),
825 "Data {0} on mesh {1} didn't contain any samples while attempting to send it to the coupling partner. "
826 "Make sure participant {2} specifies data {0} to be written using tag <write-data mesh=\"{1}\" data=\"{0}\"/>. "
827 "Alternatively, ensure participant {2} defines a mapping to mesh {1} from a mesh using data {0}.",
828 data->getDataName(), data->getMeshName(), localParticipant());
829 data->storeIteration();
830 }
831}
832
840
847
849{
850 return std::any_of(dataMap.begin(), dataMap.end(),
851 [](auto const &map) { return map.second->requiresInitialization; });
852}
853
855{
856 PRECICE_DEBUG("measure convergence of the coupling iteration");
858 // Stop, when maximal iteration count (given in config) is reached
860 _hasConverged = true;
861
862 // coupling iteration converged for current time window. Advance in time.
863 if (_hasConverged) {
864 if (_acceleration) {
866 _acceleration->iterationsConverged(getAccelerationData(), getTimeWindowStart());
867 }
869 } else {
870 // no convergence achieved for the coupling iteration within the current time window
871 if (_acceleration) {
874 }
875 }
876}
877
879{
881 PRECICE_ASSERT(not doesFirstStep(), "For convergence information the sending participant is never the first one.");
882 profiling::Event e("sendConvergence", profiling::Fundamental);
883 m2n->send(_hasConverged);
884}
885
887{
889 PRECICE_ASSERT(doesFirstStep(), "For convergence information the receiving participant is always the first one.");
890 profiling::Event e("receiveConvergence", profiling::Fundamental);
891 m2n->receive(_hasConverged);
892}
893
895{
896 return _time.windowStart();
897}
898
900{
901 if (hasTimeWindowSize()) {
902 return _time.time() + _time.untilWindowEnd(_timeWindowSize);
903 } else {
904 return _time.time() + _time.untilEnd();
905 }
906}
907
909{
910 // Global toggle if a single send data uses substeps
911 for (auto cpldata : _allData | boost::adaptors::map_values) {
912 if (cpldata->getDirection() == CouplingData::Direction::Send && cpldata->exchangeSubsteps()) {
913 return true;
914 }
915 }
916 return false;
917}
918
920{
921 // This implementation covers everything except SerialImplicit
923 return {};
924 }
925
926 ImplicitData idata;
927 for (auto cpldata : _allData | boost::adaptors::map_values) {
928 if (cpldata->getDirection() == CouplingData::Direction::Receive) {
929 idata.add(cpldata->getDataID(), false);
930 }
931 }
932 return idata;
933}
934
939
940} // namespace precice::cplscheme
#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 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
double getTime() const final override
getter for _time
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 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.
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 empty(T... args)
T end(T... args)
T make_shared(T... args)
T max(T... args)
contains implementations of acceleration schemes.
std::shared_ptr< Acceleration > PtrAcceleration
contains actions to modify exchanged data.
Definition Action.hpp:6
std::shared_ptr< ConvergenceMeasure > PtrConvergenceMeasure
contains implementations of coupling schemes for coupled simulations.
std::shared_ptr< CouplingData > PtrCouplingData
std::map< int, PtrCouplingData > DataMap
contains the logic of the parallel communication between participants.
Definition BoundM2N.cpp:12
std::shared_ptr< M2N > PtrM2N
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)
provides Mesh, Data and primitives.
std::shared_ptr< Data > PtrData
std::shared_ptr< Mesh > PtrMesh
static constexpr Group Fundamental
Convenience instance of the Cat::Fundamental.
Definition Event.hpp:22
int DataID
Definition Types.hpp:25
STL namespace.
T dynamic_pointer_cast(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