preCICE v3.2.0
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,
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{
37 if (doesFirstStep()) {
40 _participantSetsTimeWindowSize = true; // not allowed to call setTimeWindowSize anymore.
42 } else {
45 }
46 }
47}
48
50 double maxTime,
51 int maxTimeWindows,
52 double timeWindowSize,
53 const std::string &firstParticipant,
54 const std::string &secondParticipant,
58 CouplingMode cplMode)
59 : SerialCouplingScheme(maxTime, maxTimeWindows, timeWindowSize, firstParticipant, secondParticipant, localParticipant, std::move(m2n), dtMethod, cplMode, UNDEFINED_MAX_ITERATIONS, UNDEFINED_MAX_ITERATIONS) {};
60
71
73{
77 getM2N()->receive(dt);
78 PRECICE_DEBUG("Received time window size of {}.", dt);
81 PRECICE_ASSERT(not doesFirstStep(), "Only second participant can receive time window size.");
82
83 if (hasTimeWindowSize() && isImplicitCouplingScheme() && not hasConverged()) { // Restriction necessary as long as extrapolation is not implemented. See https://github.com/precice/precice/issues/1770 for details.
84 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.");
85 }
86
88 }
89}
90
92{
93 // F: send, receive, S: receive, send
95 if (doesFirstStep()) {
99 } else {
101 }
102 if (sendsInitializedData()) { // this send/recv pair is only needed, if no substeps are exchanged.
104 }
105 } else { // second participant
106 if (sendsInitializedData()) {
108 }
109 if (receivesInitializedData()) { // this send/recv pair is only needed, if no substeps are exchanged.
111 }
112 // similar to SerialCouplingScheme::exchangeSecondData()
113 PRECICE_DEBUG("Receiving data...");
115 setTimeWindowSize(getNextTimeWindowSize()); // Needed, because second participant just received _timeWindowSize from first participant, if serial coupling scheme using first participant method.
118 }
119}
120
122{
124 if (doesFirstStep()) { // first participant
125 PRECICE_DEBUG("Sending data...");
129 } else { // second participant
130 PRECICE_DEBUG("Sending data...");
133 }
134 } else {
136
137 if (doesFirstStep()) { // first participant
138 PRECICE_DEBUG("Sending data...");
142 } else { // second participant
143 PRECICE_DEBUG("Perform acceleration (only second participant)...");
145 PRECICE_DEBUG("Sending convergence...");
147 PRECICE_DEBUG("Sending data...");
150 }
151 }
152}
153
155{
157 if (doesFirstStep()) { // first participant
158 PRECICE_DEBUG("Receiving data...");
162 }
163
165
166 if (not doesFirstStep()) { // second participant
167 // the second participant does not want new data in the last iteration of the last time window
168 if (isCouplingOngoing()) {
170 PRECICE_DEBUG("Receiving data...");
174 }
175 }
176 } else {
178
179 if (doesFirstStep()) { // first participant
180 PRECICE_DEBUG("Receiving convergence data...");
182 PRECICE_DEBUG("Receiving data...");
186 }
187
188 if (hasConverged()) {
190 }
191
193
194 if (not doesFirstStep()) { // second participant
195 // the second participant does not want new data in the last iteration of the last time window
196 if (isCouplingOngoing() || not hasConverged()) {
198 PRECICE_DEBUG("Receiving data...");
200 if (hasConverged()) {
202 } else {
203 receiveData(getM2N(), getReceiveData()); // receive data for end of window
204 }
206 }
207 }
208 }
209}
210
212{
213 // SerialCouplingSchemes applies acceleration to send data
214 return getSendData();
215}
216
218{
220 return {};
221 }
222
223 const auto isSecond = !doesFirstStep();
224
225 ImplicitData idata;
226 for (auto cpldata : getReceiveData() | boost::adaptors::map_values) {
227 idata.add(cpldata->getDataID(), isSecond);
228 }
229 return idata;
230}
231
232} // namespace precice::cplscheme
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
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.
bool isCouplingOngoing() const final override
Returns true, when the coupled simulation is still ongoing.
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 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
double getTimeWindowSize() const final override
Returns the time window size, if one is given by the coupling scheme.
bool hasConverged() const override
Checks if the implicit cplscheme has converged.
void moveToNextWindow()
finalizes this window's data and initializes data for next window.
void doImplicitStep()
perform a coupling iteration
double getTime() const final override
getter for _time
bool receivesInitializedData() const
Getter for _receivesInitializedData.
bool sendsInitializedData() const final override
Getter for _sendsInitializedData.
double getNextTimeWindowSize() const
Getter for _nextTimeWindowSize.
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.
void receiveData(const m2n::PtrM2N &m2n, const DataMap &receiveData)
Receives data receiveDataIDs given in mapCouplingData with communication.
bool doesFirstStep() const
Getter for _doesFirstStep.
double getTimeWindowStart() const final override
DataMap & getSendData()
Returns all data to be sent.
DataMap & getReceiveData()
Returns all data to be received.
BiCouplingScheme(double maxTime, int maxTimeWindows, double timeWindowSize, std::string firstParticipant, std::string secondParticipant, const std::string &localParticipant, m2n::PtrM2N m2n, int minIterations, int maxIterations, CouplingMode cplMode, constants::TimesteppingMethod dtMethod)
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.
DataMap & getAccelerationData() final override
interface to provide accelerated data, depending on coupling scheme being used
void receiveAndSetTimeWindowSize()
Receives and sets the time window size, if this participant is the one to receive.
void exchangeInitialData() final override
implements functionality for initialize in base class.
bool _participantSetsTimeWindowSize
Determines, if the time window size is set by the participant.
void exchangeSecondData() final override
Exchanges the second set of data.
void exchangeFirstData() final override
Functions needed for advance()
void sendTimeWindowSize()
Sends time window size, if this participant is the one to send.
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.
ImplicitData implicitDataToReceive() const final override
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.
contains implementations of coupling schemes for coupled simulations.
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.
STL namespace.
void add(DataID did, bool toKeep)