preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
WatchPointTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <istream>
4#include <iterator>
5#include <memory>
6#include <string>
7#include <vector>
9#include "mesh/Data.hpp"
10#include "mesh/Mesh.hpp"
14#include "testing/Testing.hpp"
15#include "utils/assertion.hpp"
16
17namespace precice::mesh {
18class Vertex;
19} // namespace precice::mesh
20
21using namespace precice;
23
24BOOST_AUTO_TEST_SUITE(PreciceTests)
25
26namespace {
27std::vector<double> readDoublesFromTXTFile(std::string const &filename, int skip = 0)
28{
29 std::ifstream is{filename};
30 if (skip > 0) {
32 while (skip--) {
33 is >> ignore;
34 }
35 }
37}
38} // namespace
39
40void testWatchPoint(const TestContext &context,
41 bool withEdge,
42 std::vector<double> watchPosition,
43 std::vector<double> expected)
44{
45 using namespace mesh;
46 using Eigen::VectorXd;
47 // Setup geometry
48 std::string name("rectangle");
49 PtrMesh mesh(new Mesh(name, 2, testing::nextMeshID()));
50
51 if (context.size > 1) {
52 if (context.isPrimary()) {
53 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector2d(0.0, 0.0));
54 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector2d(0.0, 1.0));
55 if (withEdge) {
56 mesh->createEdge(v1, v2);
57 }
58 } else {
59 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector2d(0.0, 1.0));
60 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector2d(0.0, 2.0));
61 if (withEdge) {
62 mesh->createEdge(v2, v3);
63 }
64 }
65 } else {
66 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector2d(0.0, 0.0));
67 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector2d(0.0, 1.0));
68 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector2d(0.0, 2.0));
69 if (withEdge) {
70 mesh->createEdge(v1, v2);
71 mesh->createEdge(v2, v3);
72 }
73 }
74
75 using precice::testing::operator""_dataID;
76 PtrData doubleData = mesh->createData("DoubleData", 1, 0_dataID);
77 PtrData vectorData = mesh->createData("VectorData", 2, 1_dataID);
78
79 if (context.size > 1) {
80 if (context.isPrimary()) {
81 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{1.0, 2.0}}));
82 vectorData->setSampleAtTime(0, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{1.0, 2.0, 3.0, 4.0}}));
83 } else {
84 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{2.0, 3.0}}));
85 vectorData->setSampleAtTime(0, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{3.0, 4.0, 5.0, 6.0}}));
86 }
87 } else {
88 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{1., 2., 3.}}));
89 vectorData->setSampleAtTime(0, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{1., 2., 3., 4., 5., 6.}}));
90 }
91
92 std::string filename0("precice-WatchPointTest-timeseries-0.log");
93
94 bool isWatchpointClosest = false;
95 // this scope forces the filestreams to be closed
96 {
97 // Create watchpoints
98 Eigen::Vector2d pointToWatch0(watchPosition[0], watchPosition[1]);
99 impl::WatchPoint watchpoint0(pointToWatch0, mesh, filename0);
100
101 // Initialize
102 watchpoint0.initialize();
103 isWatchpointClosest = watchpoint0.isClosest();
104 // Write output
105 watchpoint0.exportPointData(0.0);
106
107 // Change data (next timestep)
108 if (context.size > 1) {
109 if (context.isPrimary()) {
110 doubleData->setSampleAtTime(1, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{2., 3.}}));
111 vectorData->setSampleAtTime(1, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{2., 3., 4., 5.}}));
112 } else {
113 doubleData->setSampleAtTime(1, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{3., 4.}}));
114 vectorData->setSampleAtTime(1, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{4., 5., 6., 7.}}));
115 }
116 } else {
117 doubleData->setSampleAtTime(1, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{2., 3., 4.}}));
118 vectorData->setSampleAtTime(1, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{2., 3., 4., 5., 6., 7.}}));
119 }
120
121 // Write output again
122 watchpoint0.exportPointData(1.0);
123
124 doubleData->moveToNextWindow();
125 vectorData->moveToNextWindow();
126
127 // Write output again to check if the data stays the same
128 watchpoint0.exportPointData(2.0);
129 }
130
131 // File Format: Time Coordinate0 Coordinate1 DoubleData VectorData0 VectorData1
132 if (isWatchpointClosest) {
133 BOOST_TEST_CONTEXT("Validating watchpoint0")
134 {
135 auto result = readDoublesFromTXTFile(filename0, 6);
136 BOOST_TEST(result.size() == expected.size());
137 for (size_t i = 0; i < result.size(); ++i) {
138 BOOST_TEST_CONTEXT("entry index: " << i)
139 {
140 using testing::equals;
141 BOOST_TEST(equals(result.at(i), expected.at(i)));
142 }
143 }
144 }
145 }
146}
147
148BOOST_AUTO_TEST_SUITE(WatchPoint)
149
150PRECICE_TEST_SETUP(1_rank)
151BOOST_AUTO_TEST_CASE(TimeSeriesNoEdgeSerialPoint1)
152{
153 PRECICE_TEST();
154 bool withEdge = false;
155 auto watchPointPosition = std::vector<double>{-0.5, 0.6};
156 auto expected = std::vector<double>{
157 0.0, 0.0, 1.0, 2.0, 3.0, 4.0,
158 1.0, 0.0, 1.0, 3.0, 4.0, 5.0,
159 2.0, 0.0, 1.0, 3.0, 4.0, 5.0};
160 testWatchPoint(context, withEdge, watchPointPosition, expected);
161}
162
163PRECICE_TEST_SETUP(1_rank)
164BOOST_AUTO_TEST_CASE(TimeSeriesWithEdgeSerialPoint1)
165{
166 PRECICE_TEST();
167 bool withEdge = true;
168 auto watchPointPosition = std::vector<double>{-0.5, 0.6};
169 auto expected = std::vector<double>{
170 0.0, 0.0, 0.6, 1.6, 2.2, 3.2,
171 1.0, 0.0, 0.6, 2.6, 3.2, 4.2,
172 2.0, 0.0, 0.6, 2.6, 3.2, 4.2};
173 testWatchPoint(context, withEdge, watchPointPosition, expected);
174}
175
176PRECICE_TEST_SETUP(1_rank)
177BOOST_AUTO_TEST_CASE(TimeSeriesNoEdgeSerialPoint2)
178{
179 PRECICE_TEST();
180 bool withEdge = false;
181 auto watchPointPosition = std::vector<double>{0.0, 1.6};
182 auto expected = std::vector<double>{
183 0.0, 0.0, 2.0, 3.0, 5.0, 6.0,
184 1.0, 0.0, 2.0, 4.0, 6.0, 7.0,
185 2.0, 0.0, 2.0, 4.0, 6.0, 7.0};
186 testWatchPoint(context, withEdge, watchPointPosition, expected);
187}
188
189PRECICE_TEST_SETUP(1_rank)
190BOOST_AUTO_TEST_CASE(TimeSeriesWithEdgeSerialPoint2)
191{
192 PRECICE_TEST();
193 bool withEdge = true;
194 auto watchPointPosition = std::vector<double>{0.0, 1.6};
195 auto expected = std::vector<double>{
196 0.0, 0.0, 1.6, 2.6, 4.2, 5.2,
197 1.0, 0.0, 1.6, 3.6, 5.2, 6.2,
198 2.0, 0.0, 1.6, 3.6, 5.2, 6.2};
199 testWatchPoint(context, withEdge, watchPointPosition, expected);
200}
201
202PRECICE_TEST_SETUP(""_on(2_ranks).setupIntraComm(), Require::Events)
203BOOST_AUTO_TEST_CASE(TimeSeriesNoEdgeParallelPoint1)
204{
205 PRECICE_TEST();
206 bool withEdge = false;
207 auto watchPointPosition = std::vector<double>{-0.5, 0.6};
208 auto expected = std::vector<double>{
209 0.0, 0.0, 1.0, 2.0, 3.0, 4.0,
210 1.0, 0.0, 1.0, 3.0, 4.0, 5.0,
211 2.0, 0.0, 1.0, 3.0, 4.0, 5.0};
212 testWatchPoint(context, withEdge, watchPointPosition, expected);
213}
214
215PRECICE_TEST_SETUP(""_on(2_ranks).setupIntraComm(), Require::Events)
216BOOST_AUTO_TEST_CASE(TimeSeriesWithEdgeParallelPoint1)
217{
218 PRECICE_TEST();
219 bool withEdge = true;
220 auto watchPointPosition = std::vector<double>{-0.5, 0.6};
221 auto expected = std::vector<double>{
222 0.0, 0.0, 0.6, 1.6, 2.2, 3.2,
223 1.0, 0.0, 0.6, 2.6, 3.2, 4.2,
224 2.0, 0.0, 0.6, 2.6, 3.2, 4.2};
225 testWatchPoint(context, withEdge, watchPointPosition, expected);
226}
227
228PRECICE_TEST_SETUP(""_on(2_ranks).setupIntraComm(), Require::Events)
229BOOST_AUTO_TEST_CASE(TimeSeriesNoEdgeParallelPoint2)
230{
231 PRECICE_TEST();
232 bool withEdge = false;
233 auto watchPointPosition = std::vector<double>{0.0, 1.6};
234 auto expected = std::vector<double>{
235 0.0, 0.0, 2.0, 3.0, 5.0, 6.0,
236 1.0, 0.0, 2.0, 4.0, 6.0, 7.0,
237 2.0, 0.0, 2.0, 4.0, 6.0, 7.0};
238 testWatchPoint(context, withEdge, watchPointPosition, expected);
239}
240
241PRECICE_TEST_SETUP(""_on(2_ranks).setupIntraComm(), Require::Events)
242BOOST_AUTO_TEST_CASE(TimeSeriesWithEdgeParallelPoint2)
243{
244 PRECICE_TEST();
245 bool withEdge = true;
246 auto watchPointPosition = std::vector<double>{0.0, 1.6};
247 auto expected = std::vector<double>{
248 0.0, 0.0, 1.6, 2.6, 4.2, 5.2,
249 1.0, 0.0, 1.6, 3.6, 5.2, 6.2,
250 2.0, 0.0, 1.6, 3.6, 5.2, 6.2};
251 testWatchPoint(context, withEdge, watchPointPosition, expected);
252}
253
254PRECICE_TEST_SETUP(1_rank)
256{
257 PRECICE_TEST();
258 using namespace mesh;
259 using Eigen::VectorXd;
260 // Setup geometry
261 std::string name("rectangle");
262 PtrMesh mesh(new Mesh(name, 2, testing::nextMeshID()));
263
264 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector2d(0.0, 0.0));
265 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector2d(0.0, 1.0));
266 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector2d(1.0, 2.0));
267 mesh->createEdge(v1, v2);
268
269 PtrData doubleData = mesh->createData("DoubleData", 1, 0_dataID);
270 PtrData vectorData = mesh->createData("VectorData", 2, 1_dataID);
271
272 // v1, v2 carry data 1
273 // v2 and later v3 carry data 2
274 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{1., 1., 2.}}));
275 vectorData->setSampleAtTime(0, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{1., 1., 1., 1., 2., 2.}}));
276
277 std::string filename0("precice-WatchPointTest-reinit-0.log");
278 std::string filename1("precice-WatchPointTest-reinit-1.log");
279
280 // this scope forces the filestreams to be closed
281 {
282 // Create watchpoints
283 Eigen::Vector2d pointToWatch0(0.1, 0.5);
284 impl::WatchPoint watchpoint0(pointToWatch0, mesh, filename0);
285 Eigen::Vector2d pointToWatch1(1.1, 0.5);
286 impl::WatchPoint watchpoint1(pointToWatch1, mesh, filename1);
287
288 // Initialize
289 watchpoint0.initialize();
290 watchpoint1.initialize();
291
292 // Write output
293 watchpoint0.exportPointData(0.0);
294 watchpoint1.exportPointData(0.0);
295
296 // Change Mesh and data
297 mesh::Vertex &v4 = mesh->createVertex(Eigen::Vector2d(1.0, 0.0));
298 mesh->createEdge(v3, v4);
299 mesh->index().clear();
300 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{1., 1., 2., 2.}}));
301 vectorData->setSampleAtTime(0, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{1., 1., 1., 1., 2., 2., 2., 2.}}));
302
303 // Re-Initialize
304 watchpoint0.initialize();
305 watchpoint1.initialize();
306
307 // Write output
308 watchpoint0.exportPointData(1.0);
309 watchpoint1.exportPointData(1.0);
310 }
311
312 // File Format: Time Coordinate0 Coordinate1 DoubleData VectorData0 VectorData1
313 BOOST_TEST_CONTEXT("Validating watchpoint0")
314 {
315 auto result = readDoublesFromTXTFile(filename0, 6);
316 auto expected = std::vector<double>{
317 0.0, 0.0, 0.5, 1.0, 1.0, 1.0,
318 1.0, 0.0, 0.5, 1.0, 1.0, 1.0};
319 BOOST_TEST(result.size() == expected.size());
320 for (size_t i = 0; i < result.size(); ++i) {
321 BOOST_TEST_CONTEXT("entry index: " << i)
322 {
323 using testing::equals;
324 BOOST_TEST(equals(result.at(i), expected.at(i)));
325 }
326 }
327 }
328
329 BOOST_TEST_CONTEXT("Validating watchpoint1")
330 {
331 auto result = readDoublesFromTXTFile(filename1, 6);
332 auto expected = std::vector<double>{
333 0.0, 0.0, 0.5, 1.0, 1.0, 1.0,
334 1.0, 1.0, 0.5, 2.0, 2.0, 2.0};
335 BOOST_TEST(result.size() == expected.size());
336 for (size_t i = 0; i < result.size(); ++i) {
337 BOOST_TEST_CONTEXT("entry index: " << i)
338 {
339 using testing::equals;
340 BOOST_TEST(equals(result.at(i), expected.at(i)));
341 }
342 }
343 }
344}
345
346PRECICE_TEST_SETUP(1_rank)
347BOOST_AUTO_TEST_CASE(VolumetricInterpolation2D)
348{
349 PRECICE_TEST();
350 using namespace mesh;
351 using Eigen::VectorXd;
352 // Setup geometry
353 std::string name("triangle");
354 PtrMesh mesh(new Mesh(name, 2, testing::nextMeshID()));
355
356 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector2d(0.0, 0.0));
357 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector2d(1.0, 0.0));
358 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector2d(0.0, 1.0));
359 mesh->createEdge(v1, v2);
360 mesh->createEdge(v2, v3);
361 mesh->createEdge(v3, v1);
362 mesh->createTriangle(v1, v2, v3);
363
364 PtrData doubleData = mesh->createData("DoubleData", 1, 0_dataID);
365 PtrData vectorData = mesh->createData("VectorData", 2, 1_dataID);
366
367 // Data is (1,1,2) for the scalar, and same for each vector component.
368 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{1., 1., 2., 2.}}));
369 vectorData->setSampleAtTime(0, time::Sample(vectorData->getDimensions(), Eigen::VectorXd{{1., 1., 1., 1., 2., 2., 2., 2.}}));
370
371 std::string filename0("precice-WatchPointTest-volumetric2d-0.log");
372 std::string filename1("precice-WatchPointTest-volumetric2d-1.log");
373
374 // this scope forces the filestreams to be closed
375 {
376 // Create watchpoints
377 Eigen::Vector2d pointToWatch0(0.1, 0.5);
378 impl::WatchPoint watchpoint0(pointToWatch0, mesh, filename0);
379 Eigen::Vector2d pointToWatch1(0.25, 0.25);
380 impl::WatchPoint watchpoint1(pointToWatch1, mesh, filename1);
381
382 // Initialize
383 watchpoint0.initialize();
384 watchpoint1.initialize();
385
386 // Write output
387 watchpoint0.exportPointData(0.0);
388 watchpoint1.exportPointData(0.0);
389 }
390
391 // File Format: Time Coordinate0 Coordinate1 DoubleData VectorData0 VectorData1
392 BOOST_TEST_CONTEXT("Validating watchpoint0")
393 {
394 auto result = readDoublesFromTXTFile(filename0, 6);
395 auto expected = std::vector<double>{
396 0.0, 0.1, 0.5, 1.5, 1.5, 1.5};
397 BOOST_TEST(result.size() == expected.size());
398 for (size_t i = 0; i < result.size(); ++i) {
399 BOOST_TEST_CONTEXT("entry index: " << i)
400 {
401 using testing::equals;
402 BOOST_TEST(equals(result.at(i), expected.at(i)));
403 }
404 }
405 }
406
407 BOOST_TEST_CONTEXT("Validating watchpoint1")
408 {
409 auto result = readDoublesFromTXTFile(filename1, 6);
410 auto expected = std::vector<double>{
411 0.0, 0.25, 0.25, 1.25, 1.25, 1.25};
412 BOOST_TEST(result.size() == expected.size());
413 for (size_t i = 0; i < result.size(); ++i) {
414 BOOST_TEST_CONTEXT("entry index: " << i)
415 {
416 using testing::equals;
417 BOOST_TEST(equals(result.at(i), expected.at(i)));
418 }
419 }
420 }
421}
422
423PRECICE_TEST_SETUP(1_rank)
424BOOST_AUTO_TEST_CASE(VolumetricInterpolation3D)
425{
426 PRECICE_TEST();
427 using namespace mesh;
428 using Eigen::VectorXd;
429 // Setup geometry
430 std::string name("triangle");
431 PtrMesh mesh(new Mesh(name, 3, testing::nextMeshID()));
432
433 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
434 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
435 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
436 mesh::Vertex &v4 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 1.0));
437
438 mesh->createTetrahedron(v1, v2, v3, v4);
439
440 PtrData doubleData = mesh->createData("DoubleData", 1, 0_dataID);
441
442 // Data is (1,1,2) for the scalar, and same for each vector component.
443 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd{{1., 2., 3., 4.}}));
444
445 std::string filename0("precice-WatchPointTest-volumetric3d-0.log");
446 std::string filename1("precice-WatchPointTest-volumetric3d-1.log");
447
448 // this scope forces the filestreams to be closed
449 {
450 // Create watchpoints
451 Eigen::Vector3d pointToWatch0(0.1, 0.5, 0.2);
452 impl::WatchPoint watchpoint0(pointToWatch0, mesh, filename0);
453 Eigen::Vector3d pointToWatch1(0.25, 0.25, 0.25);
454 impl::WatchPoint watchpoint1(pointToWatch1, mesh, filename1);
455
456 // Initialize
457 watchpoint0.initialize();
458 watchpoint1.initialize();
459
460 // Write output
461 watchpoint0.exportPointData(0.0);
462 watchpoint1.exportPointData(0.0);
463 }
464
465 // File Format: Time Coordinate0 Coordinate1 Coordinate2 DoubleData
466 BOOST_TEST_CONTEXT("Validating watchpoint0")
467 {
468 auto result = readDoublesFromTXTFile(filename0, 5);
469 auto expected = std::vector<double>{
470 0.0, 0.1, 0.5, 0.2, 2.7};
471 BOOST_TEST(result.size() == expected.size());
472 for (size_t i = 0; i < result.size(); ++i) {
473 BOOST_TEST_CONTEXT("entry index: " << i)
474 {
475 using testing::equals;
476 BOOST_TEST(equals(result.at(i), expected.at(i)));
477 }
478 }
479 }
480
481 BOOST_TEST_CONTEXT("Validating watchpoint1")
482 {
483 auto result = readDoublesFromTXTFile(filename1, 5);
484 auto expected = std::vector<double>{
485 0.0, 0.25, 0.25, 0.25, 2.5};
486 BOOST_TEST(result.size() == expected.size());
487 for (size_t i = 0; i < result.size(); ++i) {
488 BOOST_TEST_CONTEXT("entry index: " << i)
489 {
490 using testing::equals;
491 BOOST_TEST(equals(result.at(i), expected.at(i)));
492 }
493 }
494 }
495}
496
497PRECICE_TEST_SETUP(""_on(4_ranks).setupIntraComm(), Require::Events)
498BOOST_AUTO_TEST_CASE(VolumetricParallel)
499{
500 PRECICE_TEST();
501 using namespace mesh;
502 using Eigen::VectorXd;
503 // Setup geometry
504 std::string name("splitTetra");
505 PtrMesh mesh(new Mesh(name, 3, testing::nextMeshID()));
506
507 switch (context.rank) {
508 case 0: {
509 mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
510 } break;
511 case 1: {
512 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
513 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
514 mesh->createEdge(v1, v2);
515 } break;
516 case 2: {
517 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
518 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
519 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 1.0));
520 mesh->createTriangle(v1, v2, v3);
521 } break;
522
523 case 3: {
524 mesh::Vertex &v1 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 0.0));
525 mesh::Vertex &v2 = mesh->createVertex(Eigen::Vector3d(1.0, 0.0, 0.0));
526 mesh::Vertex &v3 = mesh->createVertex(Eigen::Vector3d(0.0, 1.0, 0.0));
527 mesh::Vertex &v4 = mesh->createVertex(Eigen::Vector3d(0.0, 0.0, 1.0));
528
529 mesh->createTetrahedron(v1, v2, v3, v4);
530 } break;
531 }
532
533 PtrData doubleData = mesh->createData("DoubleData", 1, 0_dataID);
534
535 // Data is (1, 2, 3, 4) on the tetra, other ranks agree on their subset
536 doubleData->setSampleAtTime(0, time::Sample(doubleData->getDimensions(), Eigen::VectorXd(context.rank + 1).setLinSpaced(1., context.rank + 1)));
537
538 std::string filename0("precice-WatchPointTest-volumetricParallel-0.log");
539 bool isClosest;
540
541 // this scope forces the filestreams to be closed
542 {
543 // Create watchpoints
544 Eigen::Vector3d pointToWatch0(0.1, 0.2, 0.3);
545 // Barycentric coordinates are (0.4, 0.1, 0.2, 0.3)
546 // Thus expected value inside tetra is 0.4 + 0.1*2 + 0.2*3 + 0.3*4 = 2.4
547 impl::WatchPoint watchpoint0(pointToWatch0, mesh, filename0);
548
549 watchpoint0.initialize();
550 watchpoint0.exportPointData(0.0);
551 isClosest = watchpoint0.isClosest();
552 }
553
554 // File Format: Time Coordinate0 Coordinate1 Coordinate2 DoubleData
555 // Closest rank does the check so we know the file is closed.
556 if (isClosest) {
557 BOOST_TEST_CONTEXT("Validating watchpoint0")
558 {
559 auto result = readDoublesFromTXTFile(filename0, 5);
560 auto expected = std::vector<double>{
561 0.0, 0.1, 0.2, 0.3, 2.4};
562 BOOST_TEST(result.size() = expected.size());
563 for (size_t i = 0; i < result.size(); ++i) {
564 BOOST_TEST_CONTEXT("entry index: " << i)
565 {
566 using testing::equals;
567 BOOST_TEST(equals(result.at(i), expected.at(i)));
568 }
569 }
570 }
571 }
572}
573
575
576BOOST_AUTO_TEST_SUITE_END() // Precice
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
std::string name
#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
void testWatchPoint(const TestContext &context, bool withEdge, std::vector< double > watchPosition, std::vector< double > expected)
BOOST_AUTO_TEST_CASE(TimeSeriesNoEdgeSerialPoint1)
Observes and exports coordinates of a point on the geometry.
void exportPointData(double time)
Writes one line with data of the watchpoint into the output file.
Vertex of a mesh.
Definition Vertex.hpp:16
int size
the size of the Communicator of the current participant
provides Mesh, Data and primitives.
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
void ignore(Arguments &&...)
Definition ignore.hpp:10
Main namespace of the precice library.