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