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