preCICE v3.1.2
Loading...
Searching...
No Matches
BSplineTest.cpp
Go to the documentation of this file.
1#include "math/Bspline.hpp"
3#include "testing/Testing.hpp"
4
5#include <Eigen/Core>
6
7using namespace precice::testing;
8
11
12BOOST_AUTO_TEST_CASE(TwoPointsLinear)
13{
14 PRECICE_TEST(1_rank);
15 Eigen::Vector2d ts;
16 ts << 1, 2;
17 Eigen::MatrixXd xs(3, 2);
18 xs << 1, 2, 10, 20, 100, 200;
19
20 precice::math::Bspline bspline(ts, xs, 1);
21 // Limits
22 BOOST_TEST(equals(bspline.interpolateAt(1.0), Eigen::Vector3d(1, 10, 100)));
23 BOOST_TEST(equals(bspline.interpolateAt(2.0), Eigen::Vector3d(2, 20, 200)));
24
25 // Midpoint
26 BOOST_TEST(equals(bspline.interpolateAt(1.5), Eigen::Vector3d(1.5, 15, 150)));
27
28 // Quarters
29 BOOST_TEST(equals(bspline.interpolateAt(1.25), Eigen::Vector3d(1.25, 12.5, 125)));
30 BOOST_TEST(equals(bspline.interpolateAt(1.75), Eigen::Vector3d(1.75, 17.5, 175)));
31}
32
33BOOST_AUTO_TEST_CASE(TwoPointsLinearRoundoff)
34{
35 PRECICE_TEST(1_rank);
36 Eigen::Vector2d ts;
37 ts << 1e-6, 2e-6;
38
39 // Evaluate Bspline slightly outside of borders of time window.
40 Eigen::Vector2d teval;
42 // However, cannot distinguish between ts and teval
43 BOOST_ASSERT(equals(ts[0], teval[0]));
44 BOOST_ASSERT(equals(ts[1], teval[1]));
45
46 Eigen::MatrixXd xs(3, 2);
47 xs << 1, 2, 10, 20, 100, 200;
48
49 precice::math::Bspline bspline(ts, xs, 1);
50 // Make sure that evaluating at borders of window (again: with some floating point error within eps) does not introduce observable errors
51 BOOST_TEST(equals(bspline.interpolateAt(teval[0]), Eigen::Vector3d(1, 10, 100)));
52 BOOST_TEST(equals(bspline.interpolateAt(teval[1]), Eigen::Vector3d(2, 20, 200)));
53}
54
55BOOST_AUTO_TEST_CASE(ThreePointsLinear)
56{
57 PRECICE_TEST(1_rank);
58 Eigen::Vector3d ts;
59 ts << 0, 1, 2;
60 Eigen::MatrixXd xs(3, 3);
61 xs << 1, 2, 3, 10, 20, 30, 100, 200, 300;
62 precice::math::Bspline bspline(ts, xs, 1);
63 // Points
64 BOOST_TEST(equals(bspline.interpolateAt(0.0), Eigen::Vector3d(1, 10, 100)));
65 BOOST_TEST(equals(bspline.interpolateAt(1.0), Eigen::Vector3d(2, 20, 200)));
66 BOOST_TEST(equals(bspline.interpolateAt(2.0), Eigen::Vector3d(3, 30, 300)));
67
68 // Midpoints
69 BOOST_TEST(equals(bspline.interpolateAt(0.5), Eigen::Vector3d(1.5, 15, 150)));
70 BOOST_TEST(equals(bspline.interpolateAt(1.5), Eigen::Vector3d(2.5, 25, 250)));
71}
72
73BOOST_AUTO_TEST_CASE(ThreePointsLinearNonEquidistant)
74{
75 Eigen::Vector3d ts;
76 ts << 0, 1, 3;
77 Eigen::MatrixXd xs(3, 3);
78 xs << 1, 2, 3, 10, 20, 30, 100, 200, 300;
79 precice::math::Bspline bspline(ts, xs, 1);
80 // Points
81 BOOST_TEST(equals(bspline.interpolateAt(0.0), Eigen::Vector3d(1, 10, 100)));
82 BOOST_TEST(equals(bspline.interpolateAt(1.0), Eigen::Vector3d(2, 20, 200)));
83 BOOST_TEST(equals(bspline.interpolateAt(3.0), Eigen::Vector3d(3, 30, 300)));
84
85 // Midpoints
86 BOOST_TEST(equals(bspline.interpolateAt(0.5), Eigen::Vector3d(1.5, 15, 150)));
87 BOOST_TEST(equals(bspline.interpolateAt(2.0), Eigen::Vector3d(2.5, 25, 250), 1e-13));
88}
89
90BOOST_AUTO_TEST_CASE(ThreePointsQuadratic)
91{
92 PRECICE_TEST(1_rank);
93 Eigen::Vector3d ts;
94 ts << 0, 1, 2;
95 Eigen::MatrixXd xs(3, 3);
96 xs << 1, 2, 3, 10, 20, 30, 100, 200, 300;
97 precice::math::Bspline bspline(ts, xs, 2);
98 // Points
99 BOOST_TEST(equals(bspline.interpolateAt(0.0), Eigen::Vector3d(1, 10, 100), 1e-13));
100 BOOST_TEST(equals(bspline.interpolateAt(1.0), Eigen::Vector3d(2, 20, 200)));
101 BOOST_TEST(equals(bspline.interpolateAt(2.0), Eigen::Vector3d(3, 30, 300)));
102
103 // Midpoints
104 BOOST_TEST(equals(bspline.interpolateAt(0.5), Eigen::Vector3d(1.5, 15, 150)));
105 BOOST_TEST(equals(bspline.interpolateAt(1.5), Eigen::Vector3d(2.5, 25, 250)));
106}
107
108BOOST_AUTO_TEST_CASE(ThreePointsQuadraticNonEquidistant)
109{
110 PRECICE_TEST(1_rank);
111 Eigen::Vector3d ts;
112 ts << 0, 1, 3;
113 Eigen::MatrixXd xs(3, 3);
114 xs << 1, 2, 3, 10, 20, 30, 100, 200, 300;
115 precice::math::Bspline bspline(ts, xs, 2);
116 // Points
117 BOOST_TEST(equals(bspline.interpolateAt(0.0), Eigen::Vector3d(1, 10, 100), 1e-13));
118 BOOST_TEST(equals(bspline.interpolateAt(1.0), Eigen::Vector3d(2, 20, 200)));
119 BOOST_TEST(equals(bspline.interpolateAt(3.0), Eigen::Vector3d(3, 30, 300)));
120
121 // Midpoints
122 BOOST_TEST(equals(bspline.interpolateAt(0.5), Eigen::Vector3d(37.0 / 24, 370.0 / 24, 3700.0 / 24), 1e-13));
123 BOOST_TEST(equals(bspline.interpolateAt(2.0), Eigen::Vector3d(8.0 / 3, 80.0 / 3, 800.0 / 3), 1e-13));
124}
125
126BOOST_AUTO_TEST_CASE(FloatingPointAccuracy) // see https://github.com/precice/precice/issues/1981
127{
128 PRECICE_TEST(1_rank);
129 Eigen::Vector2d ts;
130 ts << 256.1, 256.2;
131 Eigen::MatrixXd xs(3, 2);
132 xs << 1, 2, 10, 20, 100, 200;
133 precice::math::Bspline bspline(ts, xs, 1);
134 // Points
135 BOOST_TEST(equals(bspline.interpolateAt(256.1), Eigen::Vector3d(1, 10, 100)));
136 BOOST_TEST(equals(bspline.interpolateAt(256.2), Eigen::Vector3d(2, 20, 200)));
137 // 256.1 + 0.1 > 256.2 in floating point numbers!
138 BOOST_TEST(equals(bspline.interpolateAt(256.1 + 0.1), Eigen::Vector3d(2, 20, 200)));
139}
140
141BOOST_AUTO_TEST_SUITE_END() // BSpline
BOOST_AUTO_TEST_CASE(TwoPointsLinear)
Eigen::Vector2d ts
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
#define PRECICE_TEST(...)
Definition Testing.hpp:27
Eigen::VectorXd interpolateAt(double t) const
Samples the B-Spline interpolation.
Definition Bspline.cpp:58
T epsilon(T... args)
contains the testing framework.
Definition helper.hpp:9
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