preCICE v3.1.2
Loading...
Searching...
No Matches
SerialCouplingScheme.cpp
Go to the documentation of this file.
1#include <boost/range/adaptor/map.hpp>
2#include <cmath>
3#include <memory>
4#include <ostream>
5#include <utility>
6#include <vector>
7
14#include "logging/LogMacros.hpp"
15#include "m2n/M2N.hpp"
16#include "math/differences.hpp"
17#include "utils/assertion.hpp"
18
19namespace precice::cplscheme {
20
22 double maxTime,
23 int maxTimeWindows,
24 double timeWindowSize,
25 const std::string & firstParticipant,
26 const std::string & secondParticipant,
27 const std::string & localParticipant,
28 m2n::PtrM2N m2n,
30 CouplingMode cplMode,
31 int minIterations,
32 int maxIterations)
33 : BiCouplingScheme(maxTime, maxTimeWindows, timeWindowSize, firstParticipant, secondParticipant, localParticipant, std::move(m2n), minIterations, maxIterations, cplMode, dtMethod)
34{
36 if (doesFirstStep()) {
39 _participantSetsTimeWindowSize = true; // not allowed to call setTimeWindowSize anymore.
41 } else {
44 }
45 }
46}
47
49 double maxTime,
50 int maxTimeWindows,
51 double timeWindowSize,
52 const std::string & firstParticipant,
53 const std::string & secondParticipant,
54 const std::string & localParticipant,
55 m2n::PtrM2N m2n,
57 CouplingMode cplMode)
58 : SerialCouplingScheme(maxTime, maxTimeWindows, timeWindowSize, firstParticipant, secondParticipant, localParticipant, std::move(m2n), dtMethod, cplMode, UNDEFINED_MAX_ITERATIONS, UNDEFINED_MAX_ITERATIONS){};
59
70
72{
76 getM2N()->receive(dt);
77 PRECICE_DEBUG("Received time window size of {}.", dt);
80 PRECICE_ASSERT(not doesFirstStep(), "Only second participant can receive time window size.");
81
82 if (hasTimeWindowSize() && isImplicitCouplingScheme() && not hasConverged()) { // Restriction necessary as long as extrapolation is not implemented. See https://github.com/precice/precice/issues/1770 for details.
83 PRECICE_CHECK(dt == getTimeWindowSize(), "May only use a larger time window size in the first iteration of the window. Otherwise old time window size must equal new time window size.");
84 }
85
87 }
88}
89
91{
92 // F: send, receive, S: receive, send
93 if (doesFirstStep()) {
97 } else {
99 }
100 if (sendsInitializedData()) { // this send/recv pair is only needed, if no substeps are exchanged.
102 }
103 } else { // second participant
104 if (sendsInitializedData()) {
106 }
107 if (receivesInitializedData()) { // this send/recv pair is only needed, if no substeps are exchanged.
109 }
110 // similar to SerialCouplingScheme::exchangeSecondData()
111 PRECICE_DEBUG("Receiving data...");
113 setTimeWindowSize(getNextTimeWindowSize()); // Needed, because second participant just received _timeWindowSize from first participant, if serial coupling scheme using first participant method.
116 }
117}
118
120{
122 if (doesFirstStep()) { // first participant
123 PRECICE_DEBUG("Sending data...");
126 } else { // second participant
127 PRECICE_DEBUG("Sending data...");
129 }
130 } else {
132
133 if (doesFirstStep()) { // first participant
134 PRECICE_DEBUG("Sending data...");
137 } else { // second participant
138 PRECICE_DEBUG("Perform acceleration (only second participant)...");
140 PRECICE_DEBUG("Sending convergence...");
142 PRECICE_DEBUG("Sending data...");
144 }
145 }
146}
147
149{
151 if (doesFirstStep()) { // first participant
152 PRECICE_DEBUG("Receiving data...");
155 }
156
158
159 if (not doesFirstStep()) { // second participant
160 // the second participant does not want new data in the last iteration of the last time window
161 if (isCouplingOngoing()) {
163 PRECICE_DEBUG("Receiving data...");
166 }
167 }
168 } else {
170
171 if (doesFirstStep()) { // first participant
172 PRECICE_DEBUG("Receiving convergence data...");
174 PRECICE_DEBUG("Receiving data...");
177 }
178
179 if (hasConverged()) {
181 }
182
184
185 if (not doesFirstStep()) { // second participant
186 // the second participant does not want new data in the last iteration of the last time window
187 if (isCouplingOngoing() || not hasConverged()) {
189 PRECICE_DEBUG("Receiving data...");
190 if (hasConverged()) {
192 } else {
193 receiveData(getM2N(), getReceiveData()); // receive data for end of window
194 }
196 }
197 }
198 }
199}
200
202{
203 // SerialCouplingSchemes applies acceleration to send data
204 return getSendData();
205}
206
208{
210 return {};
211 }
212
213 const auto isSecond = !doesFirstStep();
214
215 ImplicitData idata;
216 for (auto cpldata : getReceiveData() | boost::adaptors::map_values) {
217 idata.add(cpldata->getDataID(), isSecond);
218 }
219 return idata;
220}
221
222} // namespace precice::cplscheme
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:64
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:95
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:35
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
void initializeWithZeroInitialData(const DataMap &receiveData)
Initializes storage in receiveData as zero.
void setTimeWindowSize(double timeWindowSize)
Setter for _timeWindowSize.
void receiveDataForWindowEnd(const m2n::PtrM2N &m2n, const DataMap &receiveData)
Like receiveData, but temporarily sets window time to end of window.
void notifyDataHasBeenReceived()
Used to set flag after data has been received using receiveData().
bool isExplicitCouplingScheme() const
Function to determine whether coupling scheme is an explicit coupling scheme.
void sendData(const m2n::PtrM2N &m2n, const DataMap &sendData)
Sends data sendDataIDs given in mapCouplingData with communication.
void storeIteration()
used for storing all Data at end of doImplicitStep for later reference.
bool sendsInitializedData() const override final
Getter for _sendsInitializedData.
bool isImplicitCouplingScheme() const override
Function to determine whether coupling scheme is an implicit coupling scheme.
void sendConvergence(const m2n::PtrM2N &m2n)
sends convergence to other participant via m2n
bool isCouplingOngoing() const override final
Returns true, when the coupled simulation is still ongoing.
double getTime() const override final
getter for _time
bool hasConverged() const override
Checks if the implicit cplscheme has converged.
double getTimeWindowSize() const override final
Returns the time window size, if one is given by the coupling scheme.
double getTimeWindowStart() const override final
void moveToNextWindow()
finalizes this window's data and initializes data for next window.
void doImplicitStep()
perform a coupling iteration
bool hasTimeWindowSize() const override final
Function to check whether time window size is defined by coupling scheme.
bool receivesInitializedData() const
Getter for _receivesInitializedData.
double getNextTimeWindowSize() const
Getter for _nextTimeWindowSize.
void receiveConvergence(const m2n::PtrM2N &m2n)
receives convergence from other participant via m2n
void setNextTimeWindowSize(double timeWindowSize)
Setter for _nextTimeWindowSize.
void receiveData(const m2n::PtrM2N &m2n, const DataMap &receiveData)
Receives data receiveDataIDs given in mapCouplingData with communication.
bool doesFirstStep() const
Getter for _doesFirstStep.
Abstract base class for coupling schemes with two participants.
DataMap & getSendData()
Returns all data to be sent.
DataMap & getReceiveData()
Returns all data to be received.
static const double UNDEFINED_TIME_WINDOW_SIZE
To be used, when the time window size is determined dynamically during the coupling.
Coupling scheme for serial coupling, i.e. staggered execution of two coupled participants.
void exchangeSecondData() override final
Exchanges the second set of data.
void receiveAndSetTimeWindowSize()
Receives and sets the time window size, if this participant is the one to receive.
bool _participantSetsTimeWindowSize
Determines, if the time window size is set by the participant.
void sendTimeWindowSize()
Sends time window size, if this participant is the one to send.
void exchangeInitialData() override final
implements functionality for initialize in base class.
SerialCouplingScheme(double maxTime, int maxTimeWindows, double timeWindowSize, const std::string &firstParticipant, const std::string &secondParticipant, const std::string &localParticipant, m2n::PtrM2N m2n, constants::TimesteppingMethod dtMethod, CouplingMode cplMode, int minIterations, int maxIterations)
Constructor.
void exchangeFirstData() override final
Functions needed for advance()
ImplicitData implicitDataToReceive() const override final
Returns a vector of implicit data to receive in the next advance.
bool _participantReceivesTimeWindowSize
Determines, if the time window size is received by the participant.
DataMap & getAccelerationData() override final
interface to provide accelerated data, depending on coupling scheme being used
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.
STL namespace.
void add(DataID did, bool toKeep)