preCICE v3.1.2
Loading...
Searching...
No Matches
WaveformTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
3#include "testing/Testing.hpp"
5#include "time/Time.hpp"
6#include "time/Waveform.hpp"
7
8using namespace precice;
9using namespace precice::time;
10
11BOOST_AUTO_TEST_SUITE(TimeTests)
12BOOST_AUTO_TEST_SUITE(WaveformTests)
13
14BOOST_AUTO_TEST_CASE(testInitialization)
15{
16 PRECICE_TEST(1_rank);
17 const int interpolationDegree = 0;
18 const int valuesSize = 1;
19 Eigen::VectorXd value(valuesSize);
20 Waveform waveform(interpolationDegree);
21 value(0) = 0.0;
22 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
23
25 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
26
27 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
28 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 0.0));
29}
30
31BOOST_AUTO_TEST_CASE(testInitializationVector)
32{
33 PRECICE_TEST(1_rank);
34
35 const int interpolationDegree = 0;
36 const int valuesSize = 3;
37 Eigen::VectorXd value(valuesSize);
38 Waveform waveform(interpolationDegree);
39 value << 0, 0, 0;
40 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
41
43 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
44
45 for (int i = 0; i < valuesSize; i++) {
46 BOOST_TEST(testing::equals(waveform.sample(0.0)(i), 0.0));
47 BOOST_TEST(testing::equals(waveform.sample(1.0)(i), 0.0));
48 }
49}
50
51BOOST_AUTO_TEST_SUITE(InterpolationTests)
52
53BOOST_AUTO_TEST_CASE(testInterpolateDataZerothDegree)
54{
55 PRECICE_TEST(1_rank);
56
58
59 // Test zeroth degree interpolation
60 const int interpolationDegree = 0;
61 const int valuesSize = 1;
62 Eigen::VectorXd value(valuesSize);
63 Waveform waveform(interpolationDegree);
64 value(0) = 0.0;
65 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
66 value(0) = 1.0;
67 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
68
69 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
70 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 2);
71
72 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
73 BOOST_TEST(testing::equals(waveform.sample(0.5)(0), 1.0));
74 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 1.0));
75
76 value(0) = 2.0;
77 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
78
79 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
80 BOOST_TEST(testing::equals(waveform.sample(0.5)(0), 2.0));
81 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
82
83 waveform.timeStepsStorage().move();
84 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
85
86 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
87 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 2.0));
88 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 2.0));
89
90 value(0) = 3.0;
91 waveform.timeStepsStorage().setSampleAtTime(2, time::Sample{1, value});
92
93 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
94 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 3.0));
95 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 3.0));
96}
97
98BOOST_AUTO_TEST_CASE(testInterpolateDataFirstDegree)
99{
100 PRECICE_TEST(1_rank);
101
103
104 // Test first degree interpolation
105 const int interpolationDegree = 1;
106 const int valuesSize = 1;
107 Eigen::VectorXd value(valuesSize);
108 Waveform waveform(interpolationDegree);
109 value(0) = 0.0;
110 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
111 value(0) = 1.0;
112 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
113
114 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
115 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 2);
116
117 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
118 BOOST_TEST(testing::equals(waveform.sample(0.5)(0), 0.5));
119 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 1.0));
120
121 value(0) = 2.0;
122 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
123
124 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
125 BOOST_TEST(testing::equals(waveform.sample(0.5)(0), 1.0));
126 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
127
128 waveform.timeStepsStorage().move();
129 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
130
131 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
132 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 2.0));
133 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 2.0));
134
135 value(0) = 3.0;
136 waveform.timeStepsStorage().setSampleAtTime(2, time::Sample{1, value});
137
138 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
139 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 2.5));
140 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 3.0));
141}
142
143// Remove or modify this feature? Creating a second degree interpolant by using data from previous windows is difficult, because this would require several pieces of data during initialization. What would be useful: Generating a second degree interpolant from multiple samples in a single window (if available). This would go into the least-squares direction
144BOOST_AUTO_TEST_CASE(testInterpolateDataSecondDegree)
145{
146 PRECICE_TEST(1_rank);
147
149
150 // Test second degree interpolation, but there are not enough samples. Therefore, always only first degree.
151 const int interpolationDegree = 2;
152 const int valuesSize = 1;
153 Eigen::VectorXd value(valuesSize);
154 Waveform waveform(interpolationDegree);
155 value(0) = 0.0;
156 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
157 value(0) = 1.0;
158 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
159
160 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
161 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 2);
162
163 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
164 BOOST_TEST(testing::equals(waveform.sample(0.5)(0), 0.5));
165 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 1.0));
166
167 waveform.timeStepsStorage().trim();
168
169 value(0) = 2.0;
170 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
171
172 BOOST_TEST(testing::equals(waveform.sample(0.0)(0), 0.0));
173 BOOST_TEST(testing::equals(waveform.sample(0.5)(0), 1.0));
174 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
175
176 waveform.timeStepsStorage().move();
177 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
178
179 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
180 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 2.0));
181 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 2.0));
182
183 value(0) = 8.0;
184 waveform.timeStepsStorage().setSampleAtTime(2, time::Sample{1, value});
185
186 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
187 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 5.0));
188 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 8.0));
189
190 waveform.timeStepsStorage().trim();
191
192 value(0) = 4.0;
193 waveform.timeStepsStorage().setSampleAtTime(2, time::Sample{1, value});
194
195 BOOST_TEST(testing::equals(waveform.sample(1.0)(0), 2.0));
196 BOOST_TEST(testing::equals(waveform.sample(1.5)(0), 3.0));
197 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 4.0));
198
199 waveform.timeStepsStorage().move();
200 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
201
202 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 4.0));
203 BOOST_TEST(testing::equals(waveform.sample(2.5)(0), 4.0));
204 BOOST_TEST(testing::equals(waveform.sample(3.0)(0), 4.0));
205
206 value(0) = 8.0;
207 waveform.timeStepsStorage().setSampleAtTime(3, time::Sample{1, value});
208
209 BOOST_TEST(testing::equals(waveform.sample(2.0)(0), 4.0));
210 BOOST_TEST(testing::equals(waveform.sample(2.5)(0), 6.0));
211 BOOST_TEST(testing::equals(waveform.sample(3.0)(0), 8.0));
212}
213
214BOOST_AUTO_TEST_CASE(testInterpolateDataFirstDegreeVector)
215{
216 PRECICE_TEST(1_rank);
217
219
220 // Test first degree interpolation
221 const int interpolationDegree = 1;
222 const int valuesSize = 3;
223 Eigen::VectorXd value(valuesSize);
224 Waveform waveform(interpolationDegree);
225 value << 0, 0, 0;
226 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
227 value << 1, 2, 3;
228 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
229
230 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
231 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 2);
232
233 for (int i = 0; i < valuesSize; i++) {
234 BOOST_TEST(testing::equals(waveform.sample(0.0)(i), 0 * value[i]));
235 BOOST_TEST(testing::equals(waveform.sample(0.5)(i), 0.5 * value[i]));
236 BOOST_TEST(testing::equals(waveform.sample(1.0)(i), value[i]));
237 }
238
239 waveform.timeStepsStorage().trim();
240
241 value << 2, 4, 2;
242 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
243
244 for (int i = 0; i < valuesSize; i++) {
245 BOOST_TEST(testing::equals(waveform.sample(0.0)(i), 0 * value[i]));
246 BOOST_TEST(testing::equals(waveform.sample(0.5)(i), 0.5 * value[i]));
247 BOOST_TEST(testing::equals(waveform.sample(1.0)(i), value[i]));
248 }
249
250 waveform.timeStepsStorage().move();
251 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
252
253 for (int i = 0; i < valuesSize; i++) {
254 BOOST_TEST(testing::equals(waveform.sample(1.0)(i), value[i]));
255 BOOST_TEST(testing::equals(waveform.sample(1.5)(i), value[i]));
256 BOOST_TEST(testing::equals(waveform.sample(2.0)(i), value[i]));
257 }
258
259 Eigen::VectorXd value0 = value;
260 value << 1, 2, 3;
261 waveform.timeStepsStorage().setSampleAtTime(2, time::Sample{1, value});
262
263 for (int i = 0; i < valuesSize; i++) {
264 BOOST_TEST(testing::equals(waveform.sample(1.0)(i), value0[i]));
265 BOOST_TEST(testing::equals(waveform.sample(1.5)(i), 0.5 * value0[i] + 0.5 * value[i]));
266 }
267}
268
269BOOST_AUTO_TEST_CASE(testPiecewiseInterpolateDataZerothDegree)
270{
271 PRECICE_TEST(1_rank);
272
274
275 // Test zeroth degree interpolation
276 const int interpolationDegree = 0;
277 const int valuesSize = 1;
278 Eigen::VectorXd value(valuesSize);
279 Waveform waveform(interpolationDegree);
280 value(0) = 0.0;
281 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
282 value(0) = 0.5;
283 waveform.timeStepsStorage().setSampleAtTime(0.5, time::Sample{1, value});
284 value(0) = 1.0;
285 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
286
287 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
288 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 3);
289 BOOST_TEST(testing::equals(waveform.sample(0.00)(0), 0.0));
290 BOOST_TEST(testing::equals(waveform.sample(0.25)(0), 0.5));
291 BOOST_TEST(testing::equals(waveform.sample(0.50)(0), 0.5));
292 BOOST_TEST(testing::equals(waveform.sample(0.75)(0), 1.0));
293 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 1.0));
294
295 waveform.timeStepsStorage().trim();
296
297 value(0) = 1.5;
298 waveform.timeStepsStorage().setSampleAtTime(0.5, time::Sample{1, value});
299
300 value(0) = 2.0;
301 waveform.timeStepsStorage().setSampleAtTime(1.0, time::Sample{1, value});
302 BOOST_TEST(testing::equals(waveform.sample(0.00)(0), 0.0));
303 BOOST_TEST(testing::equals(waveform.sample(0.25)(0), 1.5));
304 BOOST_TEST(testing::equals(waveform.sample(0.50)(0), 1.5));
305 BOOST_TEST(testing::equals(waveform.sample(0.75)(0), 2.0));
306 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.0));
307
308 waveform.timeStepsStorage().move();
309 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
310
311 BOOST_TEST(testing::equals(waveform.sample(0.00)(0), 2.0));
312 BOOST_TEST(testing::equals(waveform.sample(0.25)(0), 2.0));
313 BOOST_TEST(testing::equals(waveform.sample(0.50)(0), 2.0));
314 BOOST_TEST(testing::equals(waveform.sample(0.75)(0), 2.0));
315 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.0));
316
317 value(0) = 3.0;
318 waveform.timeStepsStorage().setSampleAtTime(2.0, time::Sample{1, value});
319 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.0));
320 BOOST_TEST(testing::equals(waveform.sample(1.25)(0), 3.0));
321 BOOST_TEST(testing::equals(waveform.sample(1.50)(0), 3.0));
322 BOOST_TEST(testing::equals(waveform.sample(1.75)(0), 3.0));
323 BOOST_TEST(testing::equals(waveform.sample(2.00)(0), 3.0));
324
325 waveform.timeStepsStorage().trim();
326
327 value(0) = 1.5;
328 waveform.timeStepsStorage().setSampleAtTime(1.5, time::Sample{1, value});
329
330 value(0) = 4.0;
331 waveform.timeStepsStorage().setSampleAtTime(2.0, time::Sample{1, value});
332 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.0));
333 BOOST_TEST(testing::equals(waveform.sample(1.25)(0), 1.5));
334 BOOST_TEST(testing::equals(waveform.sample(1.50)(0), 1.5));
335 BOOST_TEST(testing::equals(waveform.sample(1.75)(0), 4.0));
336 BOOST_TEST(testing::equals(waveform.sample(2.00)(0), 4.0));
337}
338
339BOOST_AUTO_TEST_CASE(testPiecewiseInterpolateDataFirstDegree)
340{
341 PRECICE_TEST(1_rank);
342
344
345 // Test zeroth degree interpolation
346 const int interpolationDegree = 1;
347 const int valuesSize = 1;
348 Eigen::VectorXd value(valuesSize);
349 Waveform waveform(interpolationDegree);
350 value(0) = 0.0;
351 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
352 value(0) = 0.5;
353 waveform.timeStepsStorage().setSampleAtTime(0.5, time::Sample{1, value});
354 value(0) = 1.0;
355 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
356
357 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
358 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 3);
359 BOOST_TEST(testing::equals(waveform.sample(0.00)(0), 0.00));
360 BOOST_TEST(testing::equals(waveform.sample(0.25)(0), 0.25));
361 BOOST_TEST(testing::equals(waveform.sample(0.50)(0), 0.50));
362 BOOST_TEST(testing::equals(waveform.sample(0.75)(0), 0.75));
363 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 1.00));
364
365 value(0) = 1.5;
366 waveform.timeStepsStorage().setSampleAtTime(0.5, time::Sample{1, value});
367
368 value(0) = 2.0;
369 waveform.timeStepsStorage().setSampleAtTime(1.0, time::Sample{1, value});
370 BOOST_TEST(testing::equals(waveform.sample(0.00)(0), 0.00));
371 BOOST_TEST(testing::equals(waveform.sample(0.25)(0), 0.75));
372 BOOST_TEST(testing::equals(waveform.sample(0.50)(0), 1.50));
373 BOOST_TEST(testing::equals(waveform.sample(0.75)(0), 1.75));
374 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.00));
375
376 waveform.timeStepsStorage().move();
377 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 1);
378
379 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.00));
380 BOOST_TEST(testing::equals(waveform.sample(1.25)(0), 2.00));
381 BOOST_TEST(testing::equals(waveform.sample(1.50)(0), 2.00));
382 BOOST_TEST(testing::equals(waveform.sample(1.75)(0), 2.00));
383 BOOST_TEST(testing::equals(waveform.sample(2.00)(0), 2.00));
384
385 value(0) = 3.0;
386 waveform.timeStepsStorage().setSampleAtTime(2.0, time::Sample{1, value});
387 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.00));
388 BOOST_TEST(testing::equals(waveform.sample(1.25)(0), 2.25));
389 BOOST_TEST(testing::equals(waveform.sample(1.50)(0), 2.50));
390 BOOST_TEST(testing::equals(waveform.sample(1.75)(0), 2.75));
391 BOOST_TEST(testing::equals(waveform.sample(2.00)(0), 3.00));
392
393 waveform.timeStepsStorage().trim();
394
395 value(0) = 1.5;
396 waveform.timeStepsStorage().setSampleAtTime(1.5, time::Sample{1, value});
397
398 value(0) = 4.0;
399 waveform.timeStepsStorage().setSampleAtTime(2.0, time::Sample{1, value});
400 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.00));
401 BOOST_TEST(testing::equals(waveform.sample(1.25)(0), 1.75));
402 BOOST_TEST(testing::equals(waveform.sample(1.50)(0), 1.50));
403 BOOST_TEST(testing::equals(waveform.sample(1.75)(0), 2.75));
404 BOOST_TEST(testing::equals(waveform.sample(2.00)(0), 4.00));
405}
406
407BOOST_AUTO_TEST_CASE(testPiecewiseInterpolateDataSecondDegree)
408{
409 PRECICE_TEST(1_rank);
410
412
413 // Test zeroth degree interpolation
414 const int interpolationDegree = 2;
415 const int valuesSize = 1;
416 Waveform waveform(interpolationDegree);
417 Eigen::VectorXd value(valuesSize);
418 value(0) = 0.0;
419 waveform.timeStepsStorage().setSampleAtTime(0, time::Sample{1, value});
420 value(0) = 0.0;
421 waveform.timeStepsStorage().setSampleAtTime(0.5, time::Sample{1, value});
422 value(0) = 2.0;
423 waveform.timeStepsStorage().setSampleAtTime(1, time::Sample{1, value});
424
425 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
426 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 3);
427 BOOST_TEST(testing::equals(waveform.sample(0.00)(0), 0.00));
428 BOOST_TEST(testing::equals(waveform.sample(0.25)(0), -0.25));
429 BOOST_TEST(testing::equals(waveform.sample(0.50)(0), 0.00));
430 BOOST_TEST(testing::equals(waveform.sample(0.75)(0), 0.75));
431 BOOST_TEST(testing::equals(waveform.sample(1.00)(0), 2.00));
432}
433
434BOOST_AUTO_TEST_CASE(testPiecewiseInterpolateDataThirdDegree)
435{
436 PRECICE_TEST(1_rank);
437
439
440 // Test zeroth degree interpolation
441 const int interpolationDegree = 3;
442 const int valuesSize = 1;
443 Waveform waveform(interpolationDegree);
444
445 // linearly increasing values
446 Eigen::VectorXd value(valuesSize);
447 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
448 value(0) = t;
449 waveform.timeStepsStorage().setSampleAtTime(t, time::Sample{1, value});
450 }
451
452 BOOST_TEST(fixture.valuesSize(waveform) == valuesSize);
453 BOOST_TEST(fixture.numberOfStoredSamples(waveform) == 5);
454
455 for (double t : std::vector<double>{0.1, 0.2, 0.3, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0}) {
456 BOOST_TEST(testing::equals(waveform.sample(t)(0), t));
457 }
458
459 waveform.timeStepsStorage().trim();
460
461 // quadratically increasing values
462 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
463 value(0) = t * t;
464 waveform.timeStepsStorage().setSampleAtTime(t, time::Sample{1, value});
465 }
466
467 // interpolates given values
468 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
469 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t));
470 }
471
472 // introduces no approximation error w.r.t function
473 for (double t : std::vector<double>{0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0}) {
474 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t));
475 }
476
477 waveform.timeStepsStorage().trim();
478
479 // cubically increasing values
480 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
481 value(0) = t * t * t;
482 waveform.timeStepsStorage().setSampleAtTime(t, time::Sample{1, value});
483 }
484
485 // interpolates given values
486 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
487 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t * t));
488 }
489
490 // introduces no approximation error w.r.t function
491 for (double t : std::vector<double>{0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0}) {
492 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t * t));
493 }
494
495 waveform.timeStepsStorage().trim();
496
497 // cubically increasing values, but with non-uniform spacing
498 for (double t : std::vector<double>{0, 0.01, 0.1, 0.2, 1}) {
499 value(0) = t * t * t;
500 waveform.timeStepsStorage().setSampleAtTime(t, time::Sample{1, value});
501 }
502
503 // interpolates given values
504 for (double t : std::vector<double>{0, 0.01, 0.1, 0.2, 1}) {
505 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t * t));
506 }
507
508 // introduces no approximation error w.r.t function
509 for (double t : std::vector<double>{0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0}) {
510 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t * t));
511 }
512
513 waveform.timeStepsStorage().trim();
514
515 // quadratically increasing values, but with non-uniform spacing
516 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
517 value(0) = t * t * t * t;
518 waveform.timeStepsStorage().setSampleAtTime(t, time::Sample{1, value});
519 }
520
521 // interpolates given values
522 for (double t : std::vector<double>{0, 0.25, 0.5, 0.75, 1}) {
523 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t * t * t));
524 }
525
526 // introduces approximation error w.r.t function
527 for (double t : std::vector<double>{0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0}) {
528 double tol = 0.015625; // error < h**3 = 0.015625
529 BOOST_TEST(testing::equals(waveform.sample(t)(0), t * t * t * t, tol));
530 }
531}
532
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
#define PRECICE_TEST(...)
Definition Testing.hpp:27
BOOST_AUTO_TEST_CASE(testInitialization)
int valuesSize(time::Waveform &waveform)
int numberOfStoredSamples(time::Waveform &waveform)
void setSampleAtTime(double time, const Sample &sample)
Store Sample at a specific time.
Definition Storage.cpp:27
void trim()
Trims this Storage by deleting all values except values associated with the window start.
Definition Storage.cpp:98
void move()
Move this Storage by deleting all stamples except the one at the end of the window.
Definition Storage.cpp:86
Allows to perform interpolation on samples in storage of given data.
Definition Waveform.hpp:25
Eigen::VectorXd sample(const double time) const
Evaluate waveform at specific point in time. Uses interpolation if necessary.
Definition Waveform.cpp:26
time::Storage & timeStepsStorage()
Returns a reference to the _timeStepsStorage.
Definition Waveform.cpp:16
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
contains the time interpolation logic.
Definition Sample.hpp:6
Main namespace of the precice library.