preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SerialImplicitCouplingSchemeTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <boost/test/tools/old/interface.hpp>
4#include <iterator>
5#include <memory>
6#include <string>
7#include <vector>
12#include "com/SharedPointer.hpp"
22#include "logging/LogMacros.hpp"
24#include "m2n/M2N.hpp"
25#include "m2n/SharedPointer.hpp"
27#include "math/differences.hpp"
28#include "mesh/Data.hpp"
29#include "mesh/Mesh.hpp"
31#include "mesh/Vertex.hpp"
37#include "testing/Testing.hpp"
38#include "xml/XMLTag.hpp"
39
40using namespace precice;
41using namespace precice::cplscheme;
42
43#ifndef PRECICE_NO_MPI
44
45BOOST_AUTO_TEST_SUITE(CplSchemeTests)
46
48 CouplingScheme &cplScheme,
49 const std::string &nameParticipant,
50 const mesh::MeshConfiguration &meshConfig,
51 const std::vector<int> &validIterations)
52{
53 BOOST_REQUIRE(meshConfig.meshes().size() == 1);
54 mesh::PtrMesh mesh = meshConfig.meshes().at(0);
55 BOOST_REQUIRE(mesh->data().size() == 2);
56 BOOST_REQUIRE(!mesh->empty());
57 BOOST_REQUIRE(!validIterations.empty());
58
59 mesh::Vertex &vertex = mesh->vertex(0);
60 int index = vertex.getID();
61 auto &dataValues0 = mesh->data(0)->values();
62 auto &dataValues1 = mesh->data(1)->values();
63 double initialStepsizeData0 = 5.0;
64 double stepsizeData0 = 5.0;
65 Eigen::VectorXd initialStepsizeData1 = Eigen::VectorXd::Constant(3, 5.0);
66 Eigen::VectorXd stepsizeData1 = Eigen::VectorXd::Constant(3, 5.0);
67 double computedTime = 0.0;
68 int computedTimesteps = 0;
69 std::string nameParticipant0("Participant0");
70 std::string nameParticipant1("Participant1");
71 BOOST_TEST(((nameParticipant == nameParticipant0) || (nameParticipant == nameParticipant1)));
72 int iterationCount = 0;
73 std::vector<int>::const_iterator iterValidIterations = validIterations.begin();
74
75 if (nameParticipant == nameParticipant0) {
76 mesh->data(0)->setSampleAtTime(0, time::Sample{mesh->data(0)->getDimensions(), mesh->data(0)->values()});
77 cplScheme.initialize();
78 BOOST_TEST(not cplScheme.isTimeWindowComplete());
79 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
80 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
81 BOOST_TEST(not cplScheme.hasDataBeenReceived());
82
83 // Tells coupling scheme, that a checkpoint has been created.
84 // All required actions have to be performed before calling advance().
85 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
86 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
87
88 while (cplScheme.isCouplingOngoing()) {
89 dataValues0(index) += stepsizeData0;
90 // The max time step size is required to be obeyed.
91 double maxTimeStepSize = cplScheme.getNextTimeStepMaxSize();
92 cplScheme.addComputedTime(maxTimeStepSize);
93 mesh->data(0)->setSampleAtTime(cplScheme.getTime(), time::Sample{mesh->data(0)->getDimensions(), mesh->data(0)->values()});
94 cplScheme.firstSynchronization({});
95 cplScheme.firstExchange();
96 cplScheme.secondSynchronization();
97 cplScheme.secondExchange();
98 iterationCount++;
99 // A coupling timestep is complete, when the coupling iterations are
100 // globally converged and if subcycling steps have filled one global
101 // timestep.
102 if (cplScheme.isTimeWindowComplete()) {
103 // Advance participant time and timestep
104 computedTime += maxTimeStepSize;
105 computedTimesteps++;
106 BOOST_TEST(testing::equals(computedTime, cplScheme.getTime()));
107 BOOST_TEST(computedTimesteps == cplScheme.getTimeWindows() - 1);
108 // The iteration number is enforced by the controlled decrease of the
109 // change of data written
110 BOOST_TEST(iterationCount == *iterValidIterations);
111 if (cplScheme.isCouplingOngoing()) {
112 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
113 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
114 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
115 } else {
116 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
117 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
118 }
119 iterationCount = 0;
120 iterValidIterations++;
121 if (iterValidIterations == validIterations.end()) {
122 BOOST_REQUIRE(not cplScheme.isCouplingOngoing());
123 }
124 // Reset data values, to simulate same convergence behavior of
125 // interface values in next timestep.
126 stepsizeData0 = initialStepsizeData0;
127 } else { // coupling timestep is not yet complete
128 BOOST_TEST(cplScheme.isCouplingOngoing());
129 BOOST_TEST(iterationCount < *iterValidIterations);
130 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
131 cplScheme.markActionFulfilled(CouplingScheme::Action::ReadCheckpoint);
132 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::ReadCheckpoint));
133 // The written data value is decreased in a regular manner, in order
134 // to achieve a predictable convergence.
135 stepsizeData0 -= 1.0;
136 }
137 // the first participant always receives new data
138 // if(cplScheme.isCouplingOngoing())
139 BOOST_TEST(cplScheme.hasDataBeenReceived());
140 }
141 cplScheme.finalize(); // Ends the coupling scheme
142 BOOST_TEST(testing::equals(computedTime, 0.3));
143 BOOST_TEST(computedTimesteps == 3);
144 } else if (nameParticipant == nameParticipant1) {
145 mesh->data(1)->setSampleAtTime(0, time::Sample{mesh->data(1)->getDimensions(), mesh->data(1)->values()});
146 cplScheme.initialize();
147 BOOST_TEST(not cplScheme.isTimeWindowComplete());
148 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
149 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
150 BOOST_TEST(cplScheme.hasDataBeenReceived());
151
152 // Tells coupling scheme, that a checkpoint has been created.
153 // All required actions have to be performed before calling advance().
154 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
155 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
156
157 while (cplScheme.isCouplingOngoing()) {
158 Eigen::VectorXd currentData(3);
159 currentData = dataValues1.segment(index * 3, 3);
160 currentData += stepsizeData1;
161 dataValues1.segment(index * 3, 3) = currentData;
162 // The max time step size is required to be obeyed.
163 double maxTimeStepSize = cplScheme.getNextTimeStepMaxSize();
164 cplScheme.addComputedTime(maxTimeStepSize);
165 mesh->data(1)->setSampleAtTime(cplScheme.getTime(), time::Sample{mesh->data(1)->getDimensions(), mesh->data(1)->values()});
166 cplScheme.firstSynchronization({});
167 cplScheme.firstExchange();
168 cplScheme.secondSynchronization();
169 cplScheme.secondExchange();
170 iterationCount++;
171 // A coupling timestep is complete, when the coupling iterations are
172 // globally converged and if subcycling steps have filled one global
173 // timestep.
174 if (cplScheme.isTimeWindowComplete()) {
175 // Advance participant time and timestep
176 computedTime += maxTimeStepSize;
177 computedTimesteps++;
178 BOOST_TEST(testing::equals(computedTime, cplScheme.getTime()));
179 BOOST_TEST(computedTimesteps == cplScheme.getTimeWindows() - 1);
180 // The iterations are enforced by the controlled decrease of the
181 // change of data written
182 BOOST_TEST(iterationCount == *iterValidIterations);
183 if (cplScheme.isCouplingOngoing()) {
184 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
185 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
186 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
187 } else {
188 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
189 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
190 }
191 iterationCount = 0;
192 iterValidIterations++;
193 if (iterValidIterations == validIterations.end()) {
194 BOOST_REQUIRE(not cplScheme.isCouplingOngoing());
195 }
196 // Reset data values, to simulate same convergence behavior of
197 // interface values in next timestep.
198 stepsizeData1 = initialStepsizeData1;
199 } else { // coupling timestep is not yet complete
200 BOOST_TEST(cplScheme.isCouplingOngoing());
201 BOOST_TEST(iterationCount < *iterValidIterations);
202 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
203 // The load checkpoint action requires to fallback to the cplScheme of the
204 // first implicit iteration of the current timestep/time.
205 cplScheme.markActionFulfilled(CouplingScheme::Action::ReadCheckpoint);
206 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::ReadCheckpoint));
207 // The written data value is decreased in a regular manner, in order
208 // to achieve a predictable convergence.
209 // stepsizeData1 -= 1.0;
210 stepsizeData1 -= Eigen::Vector3d::Constant(1.0);
211 }
212 // only check if data is received
213 if (cplScheme.isCouplingOngoing())
214 BOOST_TEST(cplScheme.hasDataBeenReceived());
215 }
216 cplScheme.finalize(); // Ends the coupling scheme
217 BOOST_TEST(testing::equals(computedTime, 0.3));
218 BOOST_TEST(computedTimesteps == 3);
219 }
220}
221
223 CouplingScheme &cplScheme,
224 const std::string &nameParticipant,
225 const mesh::MeshConfiguration &meshConfig,
226 const std::vector<int> &validIterations)
227{
228 BOOST_REQUIRE(meshConfig.meshes().size() == 1);
229 mesh::PtrMesh mesh = meshConfig.meshes().at(0);
230 BOOST_REQUIRE(mesh->data().size() == 2);
231 BOOST_REQUIRE(!mesh->empty());
232 BOOST_REQUIRE(!validIterations.empty());
233
234 double initialStepsizeData0 = 5.0;
235 double stepsizeData0 = 5.0;
236 Eigen::Vector3d initialStepsizeData1 = Eigen::Vector3d::Constant(5.0);
237 Eigen::Vector3d stepsizeData1 = Eigen::Vector3d::Constant(5.0);
238 double computedTime = 0.0;
239 int computedTimesteps = 0;
240 std::string nameParticipant0("Participant0");
241 std::string nameParticipant1("Participant1");
242 BOOST_TEST(((nameParticipant == nameParticipant0) || (nameParticipant == nameParticipant1)));
243 int iterationCount = 0;
244 std::vector<int>::const_iterator iterValidIterations =
245 validIterations.begin();
246
247 if (nameParticipant == nameParticipant0) {
248 iterationCount++; // different handling due to subcycling
249 mesh->data(0)->setSampleAtTime(0, time::Sample{mesh->data(0)->getDimensions(), mesh->data(0)->values()});
250 cplScheme.initialize();
251 BOOST_TEST(not cplScheme.isTimeWindowComplete());
252 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
253 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
254 BOOST_TEST(not cplScheme.hasDataBeenReceived());
255
256 // Tells coupling scheme, that a checkpoint has been created.
257 // All required actions have to be performed before calling advance().
258 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
259 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
260
261 double maxTimeStepSize = cplScheme.getNextTimeStepMaxSize();
262 double computedTimeStepSize = maxTimeStepSize / 2.0;
263 int subcyclingStep = 0;
264
265 // Clear data for iteration.
266 mesh->data(0)->timeStepsStorage().trim();
267
268 // Main coupling loop
269 while (cplScheme.isCouplingOngoing()) {
270 cplScheme.addComputedTime(computedTimeStepSize);
271 mesh->data(0)->setSampleAtTime(cplScheme.getTime(), time::Sample{mesh->data(0)->getDimensions(), mesh->data(0)->values()});
272 cplScheme.firstSynchronization({});
273 cplScheme.firstExchange();
274 cplScheme.secondSynchronization();
275 cplScheme.secondExchange();
276 // A coupling timestep is complete, when the coupling iterations are
277 // globally converged and if subcycling steps have filled one global
278 // timestep.
279 if (cplScheme.isTimeWindowComplete()) {
280 // Advance participant time and timestep
281 mesh->data(0)->timeStepsStorage().trim();
282 computedTime += maxTimeStepSize;
283 computedTimesteps++;
284 BOOST_TEST(testing::equals(computedTime, cplScheme.getTime()));
285 BOOST_TEST(computedTimesteps == cplScheme.getTimeWindows() - 1);
286 // The iteration number is enforced by the controlled decrease of the
287 // change of data written
288 BOOST_TEST(iterationCount == *iterValidIterations);
289 if (cplScheme.isCouplingOngoing()) {
290 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
291 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
292 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
293 } else {
294 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
295 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
296 }
297 iterationCount = 1;
298 iterValidIterations++;
299 if (iterValidIterations == validIterations.end()) {
300 BOOST_REQUIRE(not cplScheme.isCouplingOngoing());
301 }
302 // Reset data values, to simulate same convergence behavior of
303 // interface values in next time step.
304 stepsizeData0 = initialStepsizeData0;
305 BOOST_TEST(subcyclingStep == 1);
306 subcyclingStep = 0;
307 } else { // coupling timestep is not yet complete
308 BOOST_TEST(cplScheme.isCouplingOngoing());
309 // If end of time window is reached
310 if (cplScheme.hasDataBeenReceived()) {
311 BOOST_TEST(iterationCount <= *iterValidIterations);
312 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
313 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
314 cplScheme.markActionFulfilled(CouplingScheme::Action::ReadCheckpoint);
315 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::ReadCheckpoint));
316 // The written data value is decreased in a regular manner, in order
317 // to achieve a predictable convergence.
318 stepsizeData0 -= 1.0;
319 subcyclingStep = 0; // Subcycling steps
320 iterationCount++; // Implicit coupling iterations
321 } else { // If subcycling
322 BOOST_TEST(iterationCount <= *iterValidIterations);
323 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
324 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
325 BOOST_TEST(subcyclingStep < 2);
326 subcyclingStep++;
327 }
328 }
329 }
330 cplScheme.finalize(); // Ends the coupling scheme
331 BOOST_TEST(testing::equals(computedTime, 0.3));
332 BOOST_TEST(computedTimesteps == 3);
333 BOOST_TEST(stepsizeData0 == 5.0);
334 }
335
336 else if (nameParticipant == nameParticipant1) {
337 iterationCount++; // different handling due to subcycling
338 mesh->data(1)->setSampleAtTime(0, time::Sample{mesh->data(1)->getDimensions(), mesh->data(1)->values()});
339 cplScheme.initialize();
340 BOOST_TEST(not cplScheme.isTimeWindowComplete());
341 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
342 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
343 BOOST_TEST(cplScheme.hasDataBeenReceived());
344
345 // Tells coupling scheme, that a checkpoint has been created.
346 // All required actions have to be performed before calling advance().
347 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
348 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
349
350 double maxTimeStepSize = cplScheme.getNextTimeStepMaxSize();
351 double preferredTimeStepSize = maxTimeStepSize / 2.5;
352 double computedTimeStepSize = preferredTimeStepSize;
353 int subcyclingStep = 0;
354
355 // Clear data for iteration.
356 mesh->data(1)->timeStepsStorage().trim();
357
358 // Main coupling loop
359 while (cplScheme.isCouplingOngoing()) {
360 cplScheme.addComputedTime(computedTimeStepSize);
361 mesh->data(1)->setSampleAtTime(cplScheme.getTime(), time::Sample{mesh->data(1)->getDimensions(), mesh->data(1)->values()});
362 cplScheme.firstSynchronization({});
363 cplScheme.firstExchange();
364 cplScheme.secondSynchronization();
365 cplScheme.secondExchange();
366 computedTimeStepSize =
367 cplScheme.getNextTimeStepMaxSize() < preferredTimeStepSize
368 ? cplScheme.getNextTimeStepMaxSize()
369 : preferredTimeStepSize;
370 // A coupling time step is complete, when the coupling iterations are
371 // globally converged and if subcycling steps have filled one global
372 // time step.
373 if (cplScheme.isTimeWindowComplete()) {
374 mesh->data(1)->timeStepsStorage().trim();
375 // Advance participant time and time step
376 computedTime += maxTimeStepSize;
377 computedTimesteps++;
378 BOOST_TEST(testing::equals(computedTime, cplScheme.getTime()));
379 BOOST_TEST(computedTimesteps == cplScheme.getTimeWindows() - 1);
380 // The iteration number is enforced by the controlled decrease of the
381 // change of data written
382 BOOST_TEST(iterationCount == *iterValidIterations);
383 if (cplScheme.isCouplingOngoing()) {
384 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
385 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
386 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::WriteCheckpoint));
387 } else {
388 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
389 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
390 }
391 iterationCount = 1;
392 iterValidIterations++;
393 if (iterValidIterations == validIterations.end()) {
394 BOOST_REQUIRE(not cplScheme.isCouplingOngoing());
395 }
396 // Reset data values, to simulate same convergence behavior of
397 // interface values in next time step.
398 stepsizeData1 = initialStepsizeData1;
399 BOOST_TEST(subcyclingStep == 2);
400 subcyclingStep = 0;
401 } else { // coupling timestep is not yet complete
402 BOOST_TEST(cplScheme.isCouplingOngoing());
403 // If end of time window is reached
404 if (cplScheme.hasDataBeenReceived()) {
405 BOOST_TEST(iterationCount <= *iterValidIterations);
406 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
407 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
408 cplScheme.markActionFulfilled(CouplingScheme::Action::ReadCheckpoint);
409 BOOST_TEST(cplScheme.isActionFulfilled(CouplingScheme::Action::ReadCheckpoint));
410 // The written data value is decreased in a regular manner, in order
411 // to achieve a predictable convergence.
412 stepsizeData1.array() -= 1.0;
413 subcyclingStep = 0; // Subcycling steps
414 iterationCount++; // Implicit coupling iterations
415 } else { // If subcycling
416 BOOST_TEST(iterationCount <= *iterValidIterations);
417 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint));
418 BOOST_TEST(not cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint));
419 BOOST_TEST(subcyclingStep < 3);
420 subcyclingStep++;
421 }
422 }
423 }
424 cplScheme.finalize(); // Ends the coupling scheme
425 BOOST_TEST(testing::equals(computedTime, 0.3));
426 BOOST_TEST(computedTimesteps == 3);
427 }
428}
429
438
439BOOST_FIXTURE_TEST_SUITE(SerialImplicitCouplingSchemeTests, SerialImplicitCouplingSchemeFixture)
440
441PRECICE_TEST_SETUP(1_rank)
442BOOST_AUTO_TEST_CASE(testParseConfigurationWithRelaxation)
443{
444 PRECICE_TEST();
445 using namespace mesh;
446
447 std::string path(_pathToTests + "serial-implicit-cplscheme-relax-const-config.xml");
448
450 PtrDataConfiguration dataConfig(new DataConfiguration(root));
451 PtrMeshConfiguration meshConfig(new MeshConfiguration(root, dataConfig));
453 new m2n::M2NConfiguration(root));
455 CouplingSchemeConfiguration cplSchemeConfig(root, meshConfig, m2nConfig, participantConfig);
456
458 BOOST_CHECK(cplSchemeConfig.getData("Data0", "Mesh") != cplSchemeConfig.getData("Data1", "Mesh"));
459 BOOST_CHECK(cplSchemeConfig.findDataByID(cplSchemeConfig.getData("Data0", "Mesh")->getID()) != cplSchemeConfig.findDataByID(cplSchemeConfig.getData("Data1", "Mesh")->getID()));
460 BOOST_CHECK(cplSchemeConfig.getData("Data0", "Mesh") == cplSchemeConfig.findDataByID(cplSchemeConfig.getData("Data0", "Mesh")->getID()));
461 BOOST_CHECK(cplSchemeConfig.getData("Data1", "Mesh") == cplSchemeConfig.findDataByID(cplSchemeConfig.getData("Data1", "Mesh")->getID()));
462 BOOST_CHECK(cplSchemeConfig.findDataByID(2) == nullptr); // nullptr, there are only two pieces of data.
463 BOOST_CHECK(cplSchemeConfig._accelerationConfig->getAcceleration().get()); // no nullptr
464}
465
467PRECICE_TEST_SETUP("Participant0"_on(1_rank), "Participant1"_on(1_rank), Require::Events)
468BOOST_AUTO_TEST_CASE(testAbsConvergenceMeasureSynchronized)
469{
470 PRECICE_TEST();
472 options.useOnlyPrimaryCom = true;
473 auto m2n = context.connectPrimaryRanks("Participant0", "Participant1", options);
474
475 using namespace mesh;
476
478 // Create a data configuration, to simplify configuration of data
479 PtrDataConfiguration dataConfig(new DataConfiguration(root));
480 dataConfig->addData("data0", mesh::Data::typeName::SCALAR);
481 dataConfig->addData("data1", mesh::Data::typeName::VECTOR);
482
483 MeshConfiguration meshConfig(root, dataConfig);
484 mesh::PtrMesh mesh(new Mesh("Mesh", 3, testing::nextMeshID()));
485 mesh->createData("data0", 1, 0_dataID);
486 mesh->createData("data1", 3, 1_dataID);
487 mesh->createVertex(Eigen::Vector3d::Zero());
488 mesh->allocateDataValues();
489 meshConfig.insertMeshToMeshDimensionsMap(mesh->getName(), mesh->getDimensions());
490 meshConfig.addMesh(mesh);
491
492 // Create all parameters necessary to create an ImplicitCouplingScheme object
493 const double maxTime = 1.0;
494 const int maxTimeWindows = 3;
495 const double timeWindowSize = 0.1;
496 std::string nameParticipant0("Participant0");
497 std::string nameParticipant1("Participant1");
498 int sendDataIndex = -1;
499 int receiveDataIndex = -1;
500 int convergenceDataIndex = -1;
501 if (context.isNamed(nameParticipant0)) {
502 sendDataIndex = 0;
503 receiveDataIndex = 1;
504 convergenceDataIndex = receiveDataIndex;
505 } else {
506 sendDataIndex = 1;
507 receiveDataIndex = 0;
508 convergenceDataIndex = sendDataIndex;
509 }
510
511 // Create the coupling scheme object
512 const int minIterations = 1;
513 const int maxIterations = 100;
514 cplscheme::SerialCouplingScheme cplScheme(maxTime, maxTimeWindows, timeWindowSize, nameParticipant0, nameParticipant1, context.name, m2n, constants::FIXED_TIME_WINDOW_SIZE, BaseCouplingScheme::Implicit, minIterations, maxIterations);
515 cplScheme.addDataToSend(mesh->data(sendDataIndex), mesh, false, true);
516 cplScheme.addDataToReceive(mesh->data(receiveDataIndex), mesh, false, true);
518
519 double convergenceLimit1 = sqrt(3.0); // when diff_vector = (1.0, 1.0, 1.0)
520 cplscheme::impl::PtrConvergenceMeasure absoluteConvMeasure1(
521 new cplscheme::impl::AbsoluteConvergenceMeasure(convergenceLimit1));
522 cplScheme.addConvergenceMeasure(convergenceDataIndex, false, false, absoluteConvMeasure1);
523
524 // Expected iterations per implicit timesptep
525 std::vector<int> validIterations = {5, 5, 5};
526 runCoupling(cplScheme, context.name, meshConfig, validIterations);
527}
528
529PRECICE_TEST_SETUP("Participant0"_on(1_rank), "Participant1"_on(1_rank), Require::Events)
530BOOST_AUTO_TEST_CASE(testConfiguredAbsConvergenceMeasureSynchronized)
531{
532 PRECICE_TEST();
533
534 using namespace mesh;
535
536 std::string configurationPath(
537 _pathToTests + "serial-implicit-cplscheme-absolute-config.xml");
538
540 PtrDataConfiguration dataConfig(new DataConfiguration(root));
541 PtrMeshConfiguration meshConfig(new MeshConfiguration(root, dataConfig));
544 CouplingSchemeConfiguration cplSchemeConfig(root, meshConfig, m2nConfig, participantConfig);
545
546 xml::configure(root, xml::ConfigurationContext{}, configurationPath);
547 m2n::PtrM2N m2n = m2nConfig->getM2N("Participant0", "Participant1");
548 useOnlyPrimaryCom(m2n) = true;
549
550 // some dummy mesh
551 meshConfig->meshes().at(0)->createVertex(Eigen::Vector3d(1.0, 1.0, 1.0));
552 meshConfig->meshes().at(0)->createVertex(Eigen::Vector3d(2.0, 1.0, -1.0));
553 meshConfig->meshes().at(0)->createVertex(Eigen::Vector3d(3.0, 1.0, 1.0));
554 meshConfig->meshes().at(0)->createVertex(Eigen::Vector3d(4.0, 1.0, -1.0));
555 meshConfig->meshes().at(0)->allocateDataValues();
556
557 std::vector<int> validIterations = {5, 5, 5};
558
559 if (context.isNamed("Participant0")) {
560 m2n->requestPrimaryRankConnection("Participant1", "Participant0", "");
561 } else {
562 m2n->acceptPrimaryRankConnection("Participant1", "Participant0", "");
563 }
564
565 runCoupling(*cplSchemeConfig.getCouplingScheme(context.name),
566 context.name, *meshConfig, validIterations);
567}
568
569PRECICE_TEST_SETUP("Participant0"_on(1_rank), "Participant1"_on(1_rank), Require::Events)
570BOOST_AUTO_TEST_CASE(testMinIterConvergenceMeasureSynchronized)
571{
572 PRECICE_TEST();
574 options.useOnlyPrimaryCom = true;
575 auto m2n = context.connectPrimaryRanks("Participant0", "Participant1", options);
576
578 // Create a data configuration, to simplify configuration of data
580 dataConfig->addData("data0", mesh::Data::typeName::SCALAR);
581 dataConfig->addData("data1", mesh::Data::typeName::VECTOR);
582
583 mesh::MeshConfiguration meshConfig(root, dataConfig);
584 mesh::PtrMesh mesh(new mesh::Mesh("Mesh", 3, testing::nextMeshID()));
585 mesh->createData("data0", 1, 0_dataID);
586 mesh->createData("data1", 3, 1_dataID);
587 mesh->createVertex(Eigen::Vector3d::Zero());
588 mesh->allocateDataValues();
589 meshConfig.insertMeshToMeshDimensionsMap(mesh->getName(), mesh->getDimensions());
590 meshConfig.addMesh(mesh);
591
592 // Create all parameters necessary to create an ImplicitCouplingScheme object
593 const double maxTime = 1.0;
594 const int maxTimeWindows = 3;
595 const double timeWindowSize = 0.1;
596 std::string nameParticipant0("Participant0");
597 std::string nameParticipant1("Participant1");
598 int sendDataIndex = -1;
599 int receiveDataIndex = -1;
600 if (context.isNamed(nameParticipant0)) {
601 sendDataIndex = 0;
602 receiveDataIndex = 1;
603 } else {
604 sendDataIndex = 1;
605 receiveDataIndex = 0;
606 }
607
608 // Create the coupling scheme object
609 const int minIterations = 1;
610 const int maxIterations = 3;
611 cplscheme::SerialCouplingScheme cplScheme(maxTime, maxTimeWindows, timeWindowSize, nameParticipant0, nameParticipant1, context.name, m2n, constants::FIXED_TIME_WINDOW_SIZE, BaseCouplingScheme::Implicit, minIterations, maxIterations);
612 cplScheme.addDataToSend(mesh->data(sendDataIndex), mesh, false, true);
613 cplScheme.addDataToReceive(mesh->data(receiveDataIndex), mesh, false, true);
614 cplScheme.determineInitialDataExchange();
615
616 // Expected iterations per implicit timesptep
617 std::vector<int> validIterations = {3, 3, 3};
618 runCoupling(cplScheme, context.name, meshConfig, validIterations);
619}
620
621PRECICE_TEST_SETUP("Participant0"_on(1_rank), "Participant1"_on(1_rank), Require::Events)
622BOOST_AUTO_TEST_CASE(testMinIterConvergenceMeasureSynchronizedWithSubcycling)
623{
624 PRECICE_TEST();
626 options.useOnlyPrimaryCom = true;
627 auto m2n = context.connectPrimaryRanks("Participant0", "Participant1", options);
628
630 // Create a data configuration, to simplify configuration of data
632 dataConfig->addData("data0", mesh::Data::typeName::SCALAR);
633 dataConfig->addData("data1", mesh::Data::typeName::VECTOR);
634
635 mesh::MeshConfiguration meshConfig(root, dataConfig);
636 mesh::PtrMesh mesh(new mesh::Mesh("Mesh", 3, testing::nextMeshID()));
637 mesh->createData("data0", 1, 0_dataID);
638 mesh->createData("data1", 3, 1_dataID);
639 mesh->createVertex(Eigen::Vector3d::Zero());
640 mesh->allocateDataValues();
641 meshConfig.insertMeshToMeshDimensionsMap(mesh->getName(), mesh->getDimensions());
642 meshConfig.addMesh(mesh);
643
644 // Create all parameters necessary to create an ImplicitCouplingScheme object
645 double maxTime = 1.0;
646 int maxTimeWindows = 3;
647 double timeWindowSize = 0.1;
648 std::string nameParticipant0("Participant0");
649 std::string nameParticipant1("Participant1");
650 int sendDataIndex = -1;
651 int receiveDataIndex = -1;
652 std::vector<int> validIterations;
653 if (context.isNamed(nameParticipant0)) {
654 sendDataIndex = 0;
655 receiveDataIndex = 1;
656 validIterations = {3, 3, 3};
657 } else {
658 sendDataIndex = 1;
659 receiveDataIndex = 0;
660 validIterations = {3, 3, 3};
661 }
662
663 // Create the coupling scheme object
664 const int minIterations = 1;
665 const int maxIterations = 3;
666 cplscheme::SerialCouplingScheme cplScheme(maxTime, maxTimeWindows, timeWindowSize, nameParticipant0, nameParticipant1, context.name, m2n, constants::FIXED_TIME_WINDOW_SIZE, BaseCouplingScheme::Implicit, minIterations, maxIterations);
667 cplScheme.addDataToSend(mesh->data(sendDataIndex), mesh, false, true);
668 cplScheme.addDataToReceive(mesh->data(receiveDataIndex), mesh, false, true);
669 cplScheme.determineInitialDataExchange();
670
672 cplScheme, context.name, meshConfig, validIterations);
673}
674
675PRECICE_TEST_SETUP("Participant0"_on(1_rank), "Participant1"_on(1_rank), Require::Events)
676BOOST_AUTO_TEST_CASE(testInitializeData)
677{
678 PRECICE_TEST();
680 options.useOnlyPrimaryCom = true;
681 auto m2n = context.connectPrimaryRanks("Participant0", "Participant1", options);
682
684
685 // Create a data configuration, to simplify configuration of data
686
688 dataConfig->addData("Data0", mesh::Data::typeName::SCALAR);
689 dataConfig->addData("Data1", mesh::Data::typeName::VECTOR);
690
691 mesh::MeshConfiguration meshConfig(root, dataConfig);
692 mesh::PtrMesh mesh(new mesh::Mesh("Mesh", 3, testing::nextMeshID()));
693 const auto dataID0 = mesh->createData("Data0", 1, 0_dataID)->getID();
694 const auto dataID1 = mesh->createData("Data1", 3, 1_dataID)->getID();
695 mesh->createVertex(Eigen::Vector3d::Zero());
696 mesh->allocateDataValues();
697 meshConfig.insertMeshToMeshDimensionsMap(mesh->getName(), mesh->getDimensions());
698 meshConfig.addMesh(mesh);
699
700 // Create all parameters necessary to create an ImplicitCouplingScheme object
701 const double maxTime = 1.0;
702 const int maxTimeWindows = 3;
703 const double timeWindowSize = 0.1;
704 const double timeStepSize = timeWindowSize; // solver is not subcycling
705 std::string nameParticipant0("Participant0");
706 std::string nameParticipant1("Participant1");
707 int sendDataIndex = -1;
708 int receiveDataIndex = -1;
709 bool dataRequiresInitialization = false;
710 if (context.isNamed(nameParticipant0)) {
711 sendDataIndex = dataID0;
712 receiveDataIndex = dataID1;
713 } else {
714 sendDataIndex = dataID1;
715 receiveDataIndex = dataID0;
716 dataRequiresInitialization = true;
717 }
718
719 // Create the coupling scheme object
720 const int minIterations = 1;
721 const int maxIterations = 3;
722 cplscheme::SerialCouplingScheme cplScheme(maxTime, maxTimeWindows, timeWindowSize, nameParticipant0, nameParticipant1, context.name, m2n, constants::FIXED_TIME_WINDOW_SIZE, BaseCouplingScheme::Implicit, minIterations, maxIterations);
724
725 cplScheme.addDataToSend(mesh->data(sendDataIndex), mesh, dataRequiresInitialization, true);
726 CouplingData *sendCouplingData = Fixture::getSendData(cplScheme, sendDataIndex);
727 cplScheme.addDataToReceive(mesh->data(receiveDataIndex), mesh, not dataRequiresInitialization, true);
728 CouplingData *receiveCouplingData = Fixture::getReceiveData(cplScheme, receiveDataIndex);
729 cplScheme.determineInitialDataExchange();
730
731 if (context.isNamed(nameParticipant0)) {
732 // ensure that read data is uninitialized
733 BOOST_TEST(receiveCouplingData->getSize() == 3);
734 BOOST_TEST(testing::equals(receiveCouplingData->values(), Eigen::Vector3d(0.0, 0.0, 0.0)));
735 // ensure that write data is uninitialized
736 BOOST_TEST(sendCouplingData->getSize() == 1);
737 BOOST_TEST(testing::equals(sendCouplingData->values()(0), 0.0));
738
739 BOOST_TEST(Fixture::isImplicitCouplingScheme(cplScheme));
740 sendCouplingData->setSampleAtTime(0, time::Sample{sendCouplingData->getDimensions(), sendCouplingData->values()});
741 cplScheme.initialize();
742 BOOST_TEST(cplScheme.hasDataBeenReceived());
743 // ensure that initial data was read
744 BOOST_TEST(receiveCouplingData->getSize() == 3);
745 BOOST_TEST(testing::equals(receiveCouplingData->values(), Eigen::Vector3d(1.0, 2.0, 3.0)));
746 BOOST_TEST(receiveCouplingData->getPreviousIterationSize() == 3);
747 BOOST_TEST(testing::equals(receiveCouplingData->previousIteration(), Eigen::Vector3d(0.0, 0.0, 0.0)));
748 // ensure that write data is still uninitialized
749 BOOST_TEST(sendCouplingData->getSize() == 1);
750 BOOST_TEST(testing::equals(sendCouplingData->values()(0), 0.0));
751 BOOST_TEST(sendCouplingData->getPreviousIterationSize() == 1);
752 BOOST_TEST(testing::equals(sendCouplingData->previousIteration()(0), 0.0));
753 BOOST_TEST(sendCouplingData->getPreviousIterationSize() == 1);
754 // set write data
755 while (cplScheme.isCouplingOngoing()) {
756 if (cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint)) {
757 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
758 }
759 BOOST_TEST(cplScheme.getNextTimeStepMaxSize() == timeStepSize);
760 sendCouplingData->setSampleAtTime(cplScheme.getTime() + timeStepSize, time::Sample{sendCouplingData->getDimensions(), Eigen::VectorXd::Constant(sendCouplingData->getSize(), 4.0)});
761 cplScheme.addComputedTime(timeStepSize);
762 sendCouplingData->setSampleAtTime(cplScheme.getTime(), time::Sample{sendCouplingData->getDimensions(), Eigen::VectorXd::Constant(sendCouplingData->getSize(), 4.0)});
763 cplScheme.firstSynchronization({});
764 cplScheme.firstExchange();
765 cplScheme.secondSynchronization();
766 cplScheme.secondExchange();
767 BOOST_TEST(cplScheme.hasDataBeenReceived());
768 if (cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint)) {
769 cplScheme.markActionFulfilled(CouplingScheme::Action::ReadCheckpoint);
770 }
771 }
772 } else {
773 BOOST_TEST(context.isNamed(nameParticipant1));
774 BOOST_TEST(cplScheme.isActionRequired(CouplingScheme::Action::InitializeData));
775 Eigen::VectorXd v(3);
776 v << 1.0, 2.0, 3.0;
777 sendCouplingData->setSampleAtTime(0, time::Sample{sendCouplingData->getDimensions(), v});
778 cplScheme.markActionFulfilled(CouplingScheme::Action::InitializeData);
779 BOOST_TEST(receiveCouplingData->getSize() == 1);
780 BOOST_TEST(testing::equals(receiveCouplingData->values()(0), 0.0));
781 BOOST_TEST(sendCouplingData->getSize() == 3);
782 BOOST_TEST(testing::equals(sendCouplingData->values(), Eigen::Vector3d(1.0, 2.0, 3.0)));
783 cplScheme.initialize();
784 BOOST_TEST(cplScheme.hasDataBeenReceived());
785 BOOST_TEST(receiveCouplingData->getSize() == 1);
786 BOOST_TEST(testing::equals(receiveCouplingData->values()(0), 4.0));
787 BOOST_TEST(receiveCouplingData->getPreviousIterationSize() == 1);
788 BOOST_TEST(testing::equals(receiveCouplingData->previousIteration()(0), 0.0));
789 BOOST_TEST(sendCouplingData->getSize() == 3);
790 BOOST_TEST(testing::equals(sendCouplingData->values(), Eigen::Vector3d(1.0, 2.0, 3.0)));
791 BOOST_TEST(sendCouplingData->getPreviousIterationSize() == 3);
792 BOOST_TEST(testing::equals(sendCouplingData->previousIteration(), Eigen::Vector3d(1.0, 2.0, 3.0)));
793 while (cplScheme.isCouplingOngoing()) {
794 if (cplScheme.isActionRequired(CouplingScheme::Action::WriteCheckpoint)) {
795 cplScheme.markActionFulfilled(CouplingScheme::Action::WriteCheckpoint);
796 }
797 BOOST_TEST(cplScheme.getNextTimeStepMaxSize() == timeStepSize);
798 sendCouplingData->setSampleAtTime(cplScheme.getTime() + timeStepSize, time::Sample{sendCouplingData->getDimensions(), v});
799 cplScheme.addComputedTime(timeStepSize);
800 cplScheme.firstSynchronization({});
801 cplScheme.firstExchange();
802 cplScheme.secondSynchronization();
803 cplScheme.secondExchange();
804 if (cplScheme.isCouplingOngoing()) {
805 BOOST_TEST(cplScheme.hasDataBeenReceived());
806 }
807 if (cplScheme.isActionRequired(CouplingScheme::Action::ReadCheckpoint)) {
808 cplScheme.markActionFulfilled(CouplingScheme::Action::ReadCheckpoint);
809 }
810 }
811 }
812 cplScheme.finalize();
813}
814
817
818#endif // not PRECICE_NO_MPI
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
void runCouplingWithSubcycling(CouplingScheme &cplScheme, const std::string &nameParticipant, const mesh::MeshConfiguration &meshConfig, const std::vector< int > &validIterations)
BOOST_AUTO_TEST_CASE(testParseConfigurationWithRelaxation)
void runCoupling(CouplingScheme &cplScheme, const std::string &nameParticipant, const mesh::MeshConfiguration &meshConfig, const std::vector< int > &validIterations)
unsigned int index
#define PRECICE_TEST()
Definition Testing.hpp:39
#define PRECICE_TEST_SETUP(...)
Creates and attaches a TestSetup to a Boost test case.
Definition Testing.hpp:29
T begin(T... args)
Performs XML configuration of a participant.
void addConvergenceMeasure(int dataID, bool suffices, bool strict, impl::PtrConvergenceMeasure measure)
Adds a measure to determine the convergence of coupling iterations.
void addDataToReceive(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps)
Adds data to be received on data exchange.
void addDataToSend(const mesh::PtrData &data, mesh::PtrMesh mesh, bool requiresInitialization, bool exchangeSubsteps)
Adds data to be sent on data exchange and possibly be modified during coupling iterations.
void determineInitialDataExchange() override
Determines which data is initialized and therefore has to be exchanged during initialize.
const Eigen::VectorXd & values() const
Returns a const reference to the data values.
const Eigen::VectorXd & previousIteration() const
returns data value from previous iteration
void setSampleAtTime(double time, time::Sample sample)
Add sample at given time to _timeStepsStorage.
int getPreviousIterationSize() const
returns size of previous iteration
mesh::PtrData getData(const std::string &dataName, const std::string &meshName) const
acceleration::PtrAccelerationConfiguration _accelerationConfig
Interface for all coupling schemes.
virtual void initialize()=0
Initializes the coupling scheme and establishes a communication connection to the coupling partner....
virtual ChangedMeshes firstSynchronization(const ChangedMeshes &changes)=0
virtual bool addComputedTime(double timeToAdd)=0
Adds newly computed time. Has to be called before every advance.
virtual double getTime() const =0
Returns the currently computed time of the coupling scheme.
virtual ChangedMeshes secondSynchronization()=0
virtual void markActionFulfilled(Action action)=0
Tells the coupling scheme that the accessor has performed the given action.
virtual bool isActionFulfilled(Action action) const =0
Returns true, if the given action has already been performed by the accessor.
virtual double getNextTimeStepMaxSize() const =0
Returns the maximal size of the next time step to be computed.
virtual bool isTimeWindowComplete() const =0
Returns true, when the accessor can advance to the next time window.
virtual void finalize()=0
Finalizes the coupling and disconnects communication.
virtual bool isCouplingOngoing() const =0
Returns true, when the coupled simulation is still ongoing.
virtual int getTimeWindows() const =0
Returns the currently computed time windows of the coupling scheme.
virtual bool isActionRequired(Action action) const =0
Returns true, if the given action has to be performed by the accessor.
virtual bool hasDataBeenReceived() const =0
Returns true, if data has been exchanged in last call of advance().
Coupling scheme for serial coupling, i.e. staggered execution of two coupled participants.
Measures the convergence from an old data set to a new one.
Configuration for communication channels between solvers.
Performs and provides configuration for Data objects from XML files.
const std::vector< PtrMesh > & meshes() const
Returns all configured meshes.
Container and creator for meshes.
Definition Mesh.hpp:38
Vertex of a mesh.
Definition Vertex.hpp:16
VertexID getID() const
Returns the unique (among vertices of one mesh on one processor) ID of the vertex.
Definition Vertex.hpp:109
Represents an XML tag to be configured automatically.
Definition XMLTag.hpp:28
T empty(T... args)
T end(T... args)
T get(T... args)
contains implementations of coupling schemes for coupled simulations.
std::shared_ptr< DataConfiguration > PtrDataConfiguration
std::shared_ptr< MeshConfiguration > PtrMeshConfiguration
boost::test_tools::predicate_result equals(const std::vector< float > &VectorA, const std::vector< float > &VectorB, float tolerance)
equals to be used in tests. Compares two std::vectors using a given tolerance. Prints both operands o...
Definition Testing.cpp:93
std::string getPathToSources()
Returns the base path to the sources.
Definition Testing.cpp:33
XMLTag getRootTag()
Returns an empty root tag with name "configuration".
Definition XMLTag.cpp:278
std::string configure(XMLTag &tag, const precice::xml::ConfigurationContext &context, std::string_view configurationFilename)
Configures the given configuration from file configurationFilename.
Definition XMLTag.cpp:284
Main namespace of the precice library.
STL namespace.
struct giving access _useOnlyPrimaryCom
Definition M2N.hpp:267
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21