preCICE
Loading...
Searching...
No Matches
AccelerationSerialTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
18#include "testing/Meshes.hpp"
20#include "testing/Testing.hpp"
22
23using namespace precice;
24using namespace precice::acceleration;
25using namespace precice::acceleration::impl;
26
28
29BOOST_AUTO_TEST_SUITE(AccelerationTests)
30
32 using DataMap = std::map<int, cplscheme::PtrCouplingData>;
33
34 // AccelerationSerialTestsFixture() {}
35};
36
37BOOST_FIXTURE_TEST_SUITE(AccelerationSerialTests, AccelerationSerialTestsFixture)
38
39#ifndef PRECICE_NO_MPI
40
41void testIQNIMVJPP(bool exchangeSubsteps)
42{
44 // use two vectors and see if underrelaxation works
45 double initialRelaxation = 0.01;
46 int maxIterationsUsed = 50;
47 int timeWindowsReused = 6;
48 int reusedTimeWindowsAtRestart = 0;
49 int chunkSize = 0;
50 int filter = Acceleration::QR1FILTER;
51 int restartType = IQNIMVJAcceleration::NO_RESTART;
52 double singularityLimit = 1e-10;
53 double svdTruncationEps = 0.0;
54 bool enforceInitialRelaxation = false;
55 bool alwaysBuildJacobian = false;
56 const double windowStart = 0;
57 const double windowEnd = 1;
58
59 std::vector<int> dataIDs;
60 dataIDs.push_back(0);
61 dataIDs.push_back(1);
62 std::vector<double> factors;
63 factors.resize(2, 1.0);
65 auto dummyMesh = testing::makeDummy2DMesh(4);
66
67 IQNIMVJAcceleration pp(initialRelaxation, enforceInitialRelaxation, maxIterationsUsed,
68 timeWindowsReused, filter, singularityLimit, dataIDs, Acceleration::OnBoundViolation::Ignore, prec, alwaysBuildJacobian,
69 restartType, chunkSize, reusedTimeWindowsAtRestart, svdTruncationEps, !exchangeSubsteps);
70
71 Eigen::VectorXd fcol1;
72
73 mesh::PtrData displacements(new mesh::Data("dvalues", -1, 1));
74 mesh::PtrData forces(new mesh::Data("fvalues", -1, 1));
75
76 // init displacements & forces
77 displacements->emplaceSampleAtTime(windowStart, {1.0, 1.0, 1.0, 1.0});
78 displacements->emplaceSampleAtTime(windowEnd, {1.0, 1.0, 1.0, 1.0});
79 forces->emplaceSampleAtTime(windowStart, {0.2, 0.2, 0.2, 0.2});
80 forces->emplaceSampleAtTime(windowEnd, {0.2, 0.2, 0.2, 0.2});
81
82 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, exchangeSubsteps);
83 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, exchangeSubsteps);
84
85 DataMap data;
86 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
87 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
88 dpcd->storeIteration();
89 fpcd->storeIteration();
90
91 pp.initialize(data);
92
93 displacements->emplaceSampleAtTime(windowEnd, {1.0, 2.0, 3.0, 4.0});
94 forces->emplaceSampleAtTime(windowEnd, {0.1, 0.1, 0.1, 0.1});
95
96 pp.performAcceleration(data, windowStart, windowEnd);
97
98 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(0), 1.00000000000000000000));
99 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(1), 1.01000000000000000888));
100 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(2), 1.02000000000000001776));
101 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(3), 1.03000000000000002665));
102 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(0), 0.199000000000000010214));
103 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(1), 0.199000000000000010214));
104 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(2), 0.199000000000000010214));
105 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(3), 0.199000000000000010214));
106
107 // Update the waveform as well
108 displacements->emplaceSampleAtTime(windowEnd, {10, 10, 10, 10});
109 forces->setSampleAtTime(windowEnd, forces->sample());
110
111 pp.performAcceleration(data, windowStart, windowEnd);
112
113 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(0), -5.63401340929695848558e-01));
114 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(1), 6.10309919173602111186e-01));
115 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(2), 1.78402117927690184729e+00));
116 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(3), 2.95773243938020247157e+00));
117 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(0), 8.28025852497733250157e-02));
118 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(1), 8.28025852497733250157e-02));
119 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(2), 8.28025852497733250157e-02));
120 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(3), 8.28025852497733250157e-02));
121}
122
123PRECICE_TEST_SETUP(1_rank)
124BOOST_AUTO_TEST_CASE(testIQNIMVJPPWithSubsteps)
125{
126 PRECICE_TEST();
127 testIQNIMVJPP(true);
128}
129
130PRECICE_TEST_SETUP(1_rank)
131BOOST_AUTO_TEST_CASE(testIQNIMVJPPWithoutSubsteps)
132{
133 PRECICE_TEST();
134 testIQNIMVJPP(false);
135}
136
137// Test if QR1 works correctly by existence of secondary data
138PRECICE_TEST_SETUP(1_rank)
139BOOST_AUTO_TEST_CASE(testILSQR1WithSecondaryData)
140{
141 PRECICE_TEST();
143 double initialRelaxation = 0.01;
144 int maxIterationsUsed = 50;
145 int timeWindowsReused = 6;
146 int filter = Acceleration::QR1FILTER;
147 double singularityLimit = 1e-3;
148 bool enforceInitialRelaxation = false;
149 const double windowStart = 0;
150 const double windowEnd = 1;
151
152 std::vector<int> dataIDs;
153 std::vector<int> primaryDataIDs;
154 dataIDs.push_back(0);
155 dataIDs.push_back(1);
156 primaryDataIDs.push_back(0);
157 PtrPreconditioner prec(new ResidualSumPreconditioner(-1, false));
158 auto dummyMesh = testing::makeDummy2DMesh(4);
159
160 // only the data with ID 0 is used as primary data for IQNILSAcceleration
161 IQNILSAcceleration pp(initialRelaxation, enforceInitialRelaxation, maxIterationsUsed,
162 timeWindowsReused, filter, singularityLimit, primaryDataIDs, Acceleration::OnBoundViolation::Ignore, prec, false);
163
164 Eigen::VectorXd fcol1;
165
166 mesh::PtrData displacements(new mesh::Data("dvalues", -1, 1));
167 mesh::PtrData forces(new mesh::Data("fvalues", -1, 1));
168
169 // init displacements & forces
170 displacements->emplaceSampleAtTime(windowStart, {1.0, 1.0, 1.0, 1.0});
171 displacements->emplaceSampleAtTime(windowEnd, {1.0, 1.0, 1.0, 1.0});
172 forces->emplaceSampleAtTime(windowStart, {0.2, 0.2, 0.2, 0.2});
173 forces->emplaceSampleAtTime(windowEnd, {0.2, 0.2, 0.2, 0.2});
174
175 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, false);
176 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, false);
177
178 DataMap data;
179 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
180 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
181 dpcd->storeIteration();
182 fpcd->storeIteration();
183
184 pp.initialize(data);
185
186 displacements->emplaceSampleAtTime(windowEnd, {1.0, 2.0, 3.0, 4.0});
187 forces->emplaceSampleAtTime(windowEnd, {0.1, 0.1, 0.1, 0.1});
188
189 pp.performAcceleration(data, windowStart, windowEnd);
190
191 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(0), 1.00000000000000000000));
192 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(1), 1.01000000000000000888));
193 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(2), 1.02000000000000001776));
194 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(3), 1.03000000000000002665));
195 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(0), 0.199000000000000010214));
196 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(1), 0.199000000000000010214));
197 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(2), 0.199000000000000010214));
198 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(3), 0.199000000000000010214));
199
200 // Update the waveform as well
201 displacements->emplaceSampleAtTime(windowEnd, {10, 10, 10, 10});
202 forces->setSampleAtTime(windowEnd, forces->sample());
203
204 // QR1 triggers QRFactorization reset here. Assertions makes sure the reset function get correct data shape also when secondary data is involved.
205 pp.performAcceleration(data, windowStart, windowEnd);
206
207 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(0), -0.565217391304349));
208 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(1), 0.608695652173912));
209 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(2), 1.78260869565217));
210 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(3), 2.95652173913043));
211 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(0), 0.0827826086956522));
212 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(1), 0.0827826086956522));
213 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(2), 0.0827826086956522));
214 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(3), 0.0827826086956522));
215}
216
217void testVIQNPP(bool exchangeSubsteps)
218{
220 // use two vectors and see if underrelaxation works
221
222 double initialRelaxation = 0.01;
223 int maxIterationsUsed = 50;
224 int timeWindowsReused = 6;
226 double singularityLimit = 1e-10;
227 bool enforceInitialRelaxation = false;
228 std::vector<int> dataIDs;
229 const double windowStart = 0;
230 const double windowEnd = 1;
231
232 dataIDs.push_back(0);
233 dataIDs.push_back(1);
234 std::vector<double> factors;
235 factors.resize(2, 1.0);
236 PtrPreconditioner prec(new ConstantPreconditioner(factors));
237
238 std::map<int, double> scalings;
239 scalings.insert(std::make_pair(0, 1.0));
240 scalings.insert(std::make_pair(1, 1.0));
241 auto dummyMesh = testing::makeDummy2DMesh(4);
242
243 IQNILSAcceleration pp(initialRelaxation, enforceInitialRelaxation, maxIterationsUsed,
244 timeWindowsReused, filter, singularityLimit, dataIDs, Acceleration::OnBoundViolation::Ignore, prec, !exchangeSubsteps);
245
246 mesh::PtrData displacements(new mesh::Data("dvalues", -1, 1));
247 mesh::PtrData forces(new mesh::Data("fvalues", -1, 1));
248
249 // init displacements & forces
250 displacements->emplaceSampleAtTime(windowStart, {1.0, 1.0, 1.0, 1.0});
251 displacements->emplaceSampleAtTime(windowEnd, {1.0, 1.0, 1.0, 1.0});
252 forces->emplaceSampleAtTime(windowStart, {0.2, 0.2, 0.2, 0.2});
253 forces->emplaceSampleAtTime(windowEnd, {0.2, 0.2, 0.2, 0.2});
254
255 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, exchangeSubsteps);
256 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, exchangeSubsteps);
257
258 dpcd->storeIteration();
259 fpcd->storeIteration();
260
261 DataMap data;
262 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
263 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
264
265 pp.initialize(data);
266
267 displacements->emplaceSampleAtTime(windowEnd, {1.0, 2.0, 3.0, 4.0});
268 forces->emplaceSampleAtTime(windowEnd, {0.1, 0.1, 0.1, 0.1});
269
270 pp.performAcceleration(data, windowStart, windowEnd);
271
272 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(0), 1.00));
273 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(1), 1.01));
274 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(2), 1.02));
275 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(3), 1.03));
276 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(0), 0.199));
277 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(1), 0.199));
278 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(2), 0.199));
279 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(3), 0.199));
280
281 Eigen::VectorXd newdvalues;
282 utils::append(newdvalues, 10.0);
283 utils::append(newdvalues, 10.0);
284 utils::append(newdvalues, 10.0);
285 utils::append(newdvalues, 10.0);
286
287 displacements->setSampleAtTime(windowEnd, time::Sample(displacements->getDimensions(), newdvalues));
288 forces->setSampleAtTime(windowEnd, forces->sample());
289
290 pp.performAcceleration(data, windowStart, windowEnd);
291
292 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(0), -5.63401340929692295845e-01));
293 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(1), 6.10309919173607440257e-01));
294 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(2), 1.78402117927690717636e+00));
295 BOOST_TEST(testing::equals(data.at(0)->waveform().sample(windowEnd)(3), 2.95773243938020513610e+00));
296 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(0), 8.28025852497733944046e-02));
297 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(1), 8.28025852497733944046e-02));
298 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(2), 8.28025852497733944046e-02));
299 BOOST_TEST(testing::equals(data.at(1)->waveform().sample(windowEnd)(3), 8.28025852497733944046e-02));
300}
301
302PRECICE_TEST_SETUP(1_rank)
303BOOST_AUTO_TEST_CASE(testVIQNPPWithSubsteps)
304{
305 PRECICE_TEST();
306 testVIQNPP(true);
307}
308
309PRECICE_TEST_SETUP(1_rank)
310BOOST_AUTO_TEST_CASE(testVIQNPPWithoutSubsteps)
311{
312 PRECICE_TEST();
313 testVIQNPP(false);
314}
315
316PRECICE_TEST_SETUP(1_rank)
317BOOST_AUTO_TEST_CASE(testConstantUnderrelaxationWithSubsteps)
318{
319 PRECICE_TEST();
320 // use two vectors and see if underrelaxation works
321 double relaxation = 0.4;
322 std::vector<int> dataIDs{0, 1};
323 auto dummyMesh = testing::makeDummy3DMesh(4);
324 const double windowStart = 0;
325 const double windowEnd = 1;
326
327 ConstantRelaxationAcceleration acc(relaxation, dataIDs);
328
329 mesh::PtrData displacements = std::make_shared<mesh::Data>("dvalues", -1, 1);
330 mesh::PtrData forces = std::make_shared<mesh::Data>("fvalues", -1, 1);
331
332 // init displacements & forces
333 displacements->emplaceSampleAtTime(windowStart, {1.0, 2.0, 3.0, 4.0});
334 forces->emplaceSampleAtTime(windowStart, {0.2, 0.2, 0.2, 0.2});
335
336 bool exchangeSubsteps = false;
337
338 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, exchangeSubsteps);
339 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, exchangeSubsteps);
340
341 DataMap data;
342 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
343 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
344 dpcd->storeIteration();
345 fpcd->storeIteration();
346
347 acc.initialize(data);
348
349 displacements->emplaceSampleAtTime(windowEnd, {3.5, 2.0, 2.0, 1.0});
350 forces->emplaceSampleAtTime(windowEnd, {0.1, 0.1, 0.1, 0.1});
351
352 acc.performAcceleration(data, windowStart, windowEnd);
353
354 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 2);
355 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 2);
356 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 2.6);
357 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 2.8);
358 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.16);
359 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.16);
360 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.16);
361 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.16);
362
363 displacements->emplaceSampleAtTime(windowEnd, {10, 10, 10, 10});
364 forces->setSampleAtTime(windowEnd, forces->sample());
365
366 acc.performAcceleration(data, windowStart, windowEnd);
367
368 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 4.6);
369 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 5.2);
370 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 5.8);
371 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 6.4);
372 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.184);
373 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.184);
374 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.184);
375 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.184);
376}
377
378PRECICE_TEST_SETUP(1_rank)
379BOOST_AUTO_TEST_CASE(testAitkenUnderrelaxationWithoutSubsteps)
380{
381 PRECICE_TEST();
382
383 double relaxation = 0.4;
384 std::vector<int> dataIDs{0, 1};
385 std::vector<double> factors{1, 1};
386 auto dummyMesh = testing::makeDummy3DMesh(4);
387 const double windowStart = 0;
388 const double windowEnd = 1;
389
390 PtrPreconditioner prec(new ConstantPreconditioner(factors));
391 AitkenAcceleration acc(relaxation, dataIDs, prec);
392
393 mesh::PtrData displacements = std::make_shared<mesh::Data>("dvalues", -1, 1);
394 mesh::PtrData forces = std::make_shared<mesh::Data>("fvalues", -1, 1);
395
396 // //init displacements & forces
397 displacements->emplaceSampleAtTime(windowStart, {1.0, 2.0, 3.0, 4.0});
398 forces->emplaceSampleAtTime(windowStart, {0.2, 0.2, 0.2, 0.2});
399
400 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, false);
401 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, false);
402
403 DataMap data;
404 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
405 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
406 dpcd->storeIteration();
407 fpcd->storeIteration();
408
409 acc.initialize(data);
410
411 displacements->emplaceSampleAtTime(windowEnd, {3.5, 2.0, 2.0, 1.0});
412 forces->emplaceSampleAtTime(windowEnd, {0.1, 0.1, 0.1, 0.1});
413
414 acc.performAcceleration(data, windowStart, windowEnd);
415
416 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 2);
417 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 2);
418 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 2.6);
419 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 2.8);
420 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.16);
421 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.16);
422 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.16);
423 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.16);
424
425 displacements->emplaceSampleAtTime(windowEnd, {10, 10, 10, 10});
426 forces->setSampleAtTime(windowEnd, forces->sample());
427
428 acc.performAcceleration(data, windowStart, windowEnd);
429
430 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 1.2689851805508461);
431 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 2.2390979382674185);
432 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 3.2092106959839914);
433 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 4.1793234537005644);
434 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.19880451030866292);
435 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.19880451030866292);
436 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.19880451030866292);
437 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.19880451030866292);
438}
439
440PRECICE_TEST_SETUP(1_rank)
441BOOST_AUTO_TEST_CASE(testAitkenUnderrelaxationWithPreconditioner)
442{
443 PRECICE_TEST();
444
445 double relaxation = 0.8;
446 std::vector<int> dataIDs{0, 1, 2, 3};
447 auto dummyMesh = testing::makeDummy2DMesh(2);
448
449 double windowStart = 0;
450 double windowEnd = 1;
451 const double dt = 1;
452
454 AitkenAcceleration acc(relaxation, dataIDs, prec);
455
456 mesh::PtrData data1 = std::make_shared<mesh::Data>("dvalues", -1, 1);
457 mesh::PtrData data2 = std::make_shared<mesh::Data>("fvalues", -1, 1);
458 mesh::PtrData data3 = std::make_shared<mesh::Data>("gvalues", -1, 2);
459 mesh::PtrData data4 = std::make_shared<mesh::Data>("hvalues", -1, 2);
460
461 // init data
462 data1->emplaceSampleAtTime(windowStart, {40, 80});
463 data2->emplaceSampleAtTime(windowStart, {5, 5});
464 data3->emplaceSampleAtTime(windowStart, {1, 2, 3, 4});
465 data4->emplaceSampleAtTime(windowStart, {20, 40, 60, 80});
466
467 cplscheme::PtrCouplingData dpcd = makeCouplingData(data1, dummyMesh, false);
468 cplscheme::PtrCouplingData fpcd = makeCouplingData(data2, dummyMesh, false);
469 cplscheme::PtrCouplingData gpcd = makeCouplingData(data3, dummyMesh, false);
470 cplscheme::PtrCouplingData hpcd = makeCouplingData(data4, dummyMesh, false);
471
472 DataMap data;
473 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
474 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
475 data.insert(std::pair<int, cplscheme::PtrCouplingData>(2, gpcd));
476 data.insert(std::pair<int, cplscheme::PtrCouplingData>(3, hpcd));
477 dpcd->storeIteration();
478 fpcd->storeIteration();
479 gpcd->storeIteration();
480 hpcd->storeIteration();
481
482 acc.initialize(data);
483
484 data1->emplaceSampleAtTime(windowEnd, {1, 7});
485 data2->emplaceSampleAtTime(windowEnd, {10, 10});
486 data3->emplaceSampleAtTime(windowEnd, {10, 11, 12, 13});
487 data4->emplaceSampleAtTime(windowEnd, {40, 60, 80, 100});
488
489 acc.performAcceleration(data, windowStart, windowEnd);
490
491 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 8.8);
492 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 21.6);
493 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 9);
494 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 9);
495 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(0) == 8.2);
496 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(1) == 9.2);
497 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(2) == 10.2);
498 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(0) == 36);
499 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(1) == 56);
500 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(2) == 76);
501 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(3) == 96);
502
503 data1->emplaceSampleAtTime(windowEnd, {2, 14});
504 data2->emplaceSampleAtTime(windowEnd, {8, 8});
505 data3->emplaceSampleAtTime(windowEnd, {13, 14, 15, 16});
506 data4->emplaceSampleAtTime(windowEnd, {41, 61, 81, 90});
507
508 acc.performAcceleration(data, windowStart, windowEnd);
509
510 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == -17.745640722103754);
511 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == -20.295060201548626);
512 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 9.5588663727976648);
513 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 9.5588663727976648);
514 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(0) == 19.235465491190659);
515 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(1) == 20.235465491190659);
516 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(2) == 21.235465491190659);
517 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(0) == 51.912064609583652);
518 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(1) == 71.912064609583652);
519 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(2) == 91.912064609583652);
520 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(3) == 95.196221242658879);
521
522 data1->emplaceSampleAtTime(windowEnd, {2.1, 14.1});
523 data2->emplaceSampleAtTime(windowEnd, {8, 8});
524 data3->emplaceSampleAtTime(windowEnd, {13.05, 14.07, 15.1, 16.1});
525 data4->emplaceSampleAtTime(windowEnd, {42, 60, 81.3, 91});
526
527 acc.iterationsConverged(data, windowStart);
528
529 // move to next window
530 windowStart += dt;
531 windowEnd += dt;
532
533 // move to next window
534 windowStart += dt;
535 windowEnd += dt;
536
537 data1->emplaceSampleAtTime(windowEnd, {3, 16});
538 data2->emplaceSampleAtTime(windowEnd, {7, 7});
539 data3->emplaceSampleAtTime(windowEnd, {18, 19, 20, 21});
540 data4->emplaceSampleAtTime(windowEnd, {50, 70, 90, 110});
541
542 acc.performAcceleration(data, windowStart, windowEnd);
543
544 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 10.4);
545 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 28.8);
546 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 6.6);
547 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 6.6);
548 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(0) == 14.6);
549 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(1) == 15.6);
550 BOOST_TEST(data.at(2)->waveform().sample(windowEnd)(2) == 16.6);
551 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(0) == 44);
552 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(1) == 64);
553 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(2) == 84);
554 BOOST_TEST(data.at(3)->waveform().sample(windowEnd)(3) == 104);
555}
556
557PRECICE_TEST_SETUP(1_rank)
558BOOST_AUTO_TEST_CASE(testConstantUnderrelaxationWithGradientWithSubsteps)
559{
560 PRECICE_TEST();
561 // use two vectors and see if underrelaxation works
562 double relaxation = 0.4;
563 std::vector<int> dataIDs{0, 1};
564 const int dim = 3;
565 auto dummyMesh = testing::makeDummy3DMesh(4);
566 const double windowStart = 0;
567 const double windowEnd = 1;
568
569 ConstantRelaxationAcceleration acc(relaxation, dataIDs);
570
571 mesh::PtrData displacements = std::make_shared<mesh::Data>("dvalues", -1, 1);
572 mesh::PtrData forces = std::make_shared<mesh::Data>("fvalues", -1, 1);
573
574 // init displacements
575 displacements->requireDataGradient();
576 Eigen::MatrixXd displacementGradient(displacements->gradients());
577 displacementGradient.resize(dim, 4);
578 for (unsigned int r = 0; r < dim; ++r) {
579 for (unsigned int c = 0; c < 4; ++c)
580 displacementGradient(r, c) = r + r * c;
581 }
582 displacements->setSampleAtTime(windowStart, time::Sample(displacements->getDimensions(), Eigen::Vector4d{1.0, 2.0, 3.0, 4.0}, displacementGradient));
583 // init forces
584 forces->requireDataGradient();
585 Eigen::MatrixXd forcesGradient(forces->gradients());
586 forcesGradient.resize(dim, 4);
587 forcesGradient.setConstant(-2);
588 forces->setSampleAtTime(windowStart, time::Sample(forces->getDimensions(), Eigen::Vector4d{0.2, 0.2, 0.2, 0.2}, forcesGradient));
589
590 bool exchangeSubsteps = true;
591
592 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, exchangeSubsteps);
593 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, exchangeSubsteps);
594
595 DataMap data;
596 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
597 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
598 dpcd->storeIteration();
599 fpcd->storeIteration();
600
601 acc.initialize(data);
602
603 displacements->setSampleAtTime(windowEnd, time::Sample(displacements->getDimensions(), Eigen::Vector4d{3.5, 2.0, 2.0, 1.0}, Eigen::MatrixXd(displacements->gradients()).setConstant(2.5)));
604 forces->setSampleAtTime(windowEnd, time::Sample(forces->getDimensions(), Eigen::Vector4d{0.1, 0.1, 0.1, 0.1}, Eigen::MatrixXd(forces->gradients()).setConstant(3)));
605
606 acc.performAcceleration(data, windowStart, windowEnd);
607
608 // Test value data
609 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 2);
610 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 2);
611 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 2.6);
612 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 2.8);
613 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.16);
614 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.16);
615 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.16);
616 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.16);
617
618 // Test gradient data
619 BOOST_TEST(data.at(0)->gradients()(0, 0) == 1);
620 BOOST_TEST(data.at(0)->gradients()(0, 1) == 1);
621 BOOST_TEST(data.at(0)->gradients()(0, 2) == 1);
622 BOOST_TEST(data.at(0)->gradients()(1, 0) == 1.6);
623 BOOST_TEST(data.at(0)->gradients()(1, 1) == 2.2);
624 BOOST_TEST(data.at(0)->gradients()(1, 2) == 2.8);
625 BOOST_TEST(data.at(1)->gradients()(0, 0) == 0);
626 BOOST_TEST(data.at(1)->gradients()(0, 1) == 0);
627 BOOST_TEST(data.at(1)->gradients()(0, 2) == 0);
628 BOOST_TEST(data.at(1)->gradients()(1, 0) == 0);
629 BOOST_TEST(data.at(1)->gradients()(1, 1) == 0);
630 BOOST_TEST(data.at(1)->gradients()(1, 2) == 0);
631
632 displacements->setSampleAtTime(windowEnd, time::Sample(displacements->getDimensions(), Eigen::Vector4d{10, 10, 10, 10}, Eigen::MatrixXd(displacements->gradients()).setConstant(4)));
633 forces->setSampleAtTime(windowEnd, forces->sample());
634
635 acc.performAcceleration(data, windowStart, windowEnd);
636
637 // Check that store iteration works properly
638 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 4.6);
639 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 5.2);
640 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 5.8);
641 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 6.4);
642 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.184);
643 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.184);
644 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.184);
645 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.184);
646
647 BOOST_TEST(data.at(0)->gradients()(0, 0) == 1.6);
648 BOOST_TEST(data.at(0)->gradients()(0, 1) == 1.6);
649 BOOST_TEST(data.at(0)->gradients()(0, 2) == 1.6);
650 BOOST_TEST(data.at(0)->gradients()(1, 0) == 2.2);
651 BOOST_TEST(data.at(0)->gradients()(1, 1) == 2.8);
652 BOOST_TEST(data.at(0)->gradients()(1, 2) == 3.4);
653}
654
655PRECICE_TEST_SETUP(1_rank)
656BOOST_AUTO_TEST_CASE(testConstantUnderrelaxationWithoutSubsteps)
657{
658 PRECICE_TEST();
659 // use two vectors and see if underrelaxation works
660 double relaxation = 0.4;
661 std::vector<int> dataIDs{0, 1};
662 auto dummyMesh = testing::makeDummy3DMesh(4);
663 const double windowStart = 0;
664 const double windowEnd = 1;
665
666 ConstantRelaxationAcceleration acc(relaxation, dataIDs);
667
668 mesh::PtrData displacements = std::make_shared<mesh::Data>("dvalues", -1, 1);
669 mesh::PtrData forces = std::make_shared<mesh::Data>("fvalues", -1, 1);
670
671 displacements->emplaceSampleAtTime(windowStart, {1.0, 2.0, 3.0, 4.0});
672 forces->emplaceSampleAtTime(windowStart, {0.2, 0.2, 0.2, 0.2});
673
674 bool exchangeSubsteps = false;
675
676 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, exchangeSubsteps);
677 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, exchangeSubsteps);
678
679 DataMap data;
680 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
681 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
682 dpcd->storeIteration();
683 fpcd->storeIteration();
684
685 acc.initialize(data);
686
687 displacements->emplaceSampleAtTime(windowEnd, {3.5, 2.0, 2.0, 1.0});
688 forces->emplaceSampleAtTime(windowEnd, {0.1, 0.1, 0.1, 0.1});
689
690 acc.performAcceleration(data, windowStart, windowEnd);
691
692 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 2);
693 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 2);
694 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 2.6);
695 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 2.8);
696 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.16);
697 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.16);
698 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.16);
699 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.16);
700
701 displacements->emplaceSampleAtTime(windowEnd, {10, 10, 10, 10});
702
703 forces->setSampleAtTime(windowEnd, forces->sample());
704
705 acc.performAcceleration(data, windowStart, windowEnd);
706
707 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 4.6);
708 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 5.2);
709 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 5.8);
710 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 6.4);
711 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.184);
712 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.184);
713 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.184);
714 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.184);
715}
716
717PRECICE_TEST_SETUP(1_rank)
718BOOST_AUTO_TEST_CASE(testConstantUnderrelaxationWithGradientWithoutSubsteps)
719{
720 PRECICE_TEST();
721 // use two vectors and see if underrelaxation works
722 double relaxation = 0.4;
723 std::vector<int> dataIDs{0, 1};
724 const int dim = 3;
725 auto dummyMesh = testing::makeDummy3DMesh(4);
726 const double windowStart = 0;
727 const double windowEnd = 1;
728
729 ConstantRelaxationAcceleration acc(relaxation, dataIDs);
730
731 mesh::PtrData displacements = std::make_shared<mesh::Data>("dvalues", -1, 1);
732 mesh::PtrData forces = std::make_shared<mesh::Data>("fvalues", -1, 1);
733
734 // init displacements
735 displacements->requireDataGradient();
736 Eigen::MatrixXd displacementGradient(displacements->gradients());
737 displacementGradient.resize(dim, 4);
738 for (unsigned int r = 0; r < dim; ++r) {
739 for (unsigned int c = 0; c < 4; ++c)
740 displacementGradient(r, c) = r + r * c;
741 }
742 displacements->setSampleAtTime(windowStart, time::Sample(displacements->getDimensions(), Eigen::Vector4d{1.0, 2.0, 3.0, 4.0}, displacementGradient));
743 // init forces
744 forces->requireDataGradient();
745 Eigen::MatrixXd forcesGradient(forces->gradients());
746 forcesGradient.resize(dim, 4);
747 forcesGradient.setConstant(-2);
748 forces->setSampleAtTime(windowStart, time::Sample(forces->getDimensions(), Eigen::Vector4d{0.2, 0.2, 0.2, 0.2}, forcesGradient));
749
750 bool exchangeSubsteps = false;
751
752 cplscheme::PtrCouplingData dpcd = makeCouplingData(displacements, dummyMesh, exchangeSubsteps);
753 cplscheme::PtrCouplingData fpcd = makeCouplingData(forces, dummyMesh, exchangeSubsteps);
754
755 DataMap data;
756 data.insert(std::pair<int, cplscheme::PtrCouplingData>(0, dpcd));
757 data.insert(std::pair<int, cplscheme::PtrCouplingData>(1, fpcd));
758 dpcd->storeIteration();
759 fpcd->storeIteration();
760
761 acc.initialize(data);
762
763 displacements->setSampleAtTime(windowEnd, time::Sample(displacements->getDimensions(), Eigen::Vector4d{3.5, 2.0, 2.0, 1.0}, Eigen::MatrixXd(displacements->gradients()).setConstant(2.5)));
764 forces->setSampleAtTime(windowEnd, time::Sample(displacements->getDimensions(), Eigen::Vector4d{0.1, 0.1, 0.1, 0.1}, Eigen::MatrixXd(displacements->gradients()).setConstant(3)));
765
766 acc.performAcceleration(data, windowStart, windowEnd);
767
768 // Test value data
769 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 2);
770 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 2);
771 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 2.6);
772 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 2.8);
773 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.16);
774 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.16);
775 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.16);
776 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.16);
777
778 // Test gradient data
779 BOOST_TEST(data.at(0)->gradients()(0, 0) == 1);
780 BOOST_TEST(data.at(0)->gradients()(0, 1) == 1);
781 BOOST_TEST(data.at(0)->gradients()(0, 2) == 1);
782 BOOST_TEST(data.at(0)->gradients()(1, 0) == 1.6);
783 BOOST_TEST(data.at(0)->gradients()(1, 1) == 2.2);
784 BOOST_TEST(data.at(0)->gradients()(1, 2) == 2.8);
785 BOOST_TEST(data.at(1)->gradients()(0, 0) == 0);
786 BOOST_TEST(data.at(1)->gradients()(0, 1) == 0);
787 BOOST_TEST(data.at(1)->gradients()(0, 2) == 0);
788 BOOST_TEST(data.at(1)->gradients()(1, 0) == 0);
789 BOOST_TEST(data.at(1)->gradients()(1, 1) == 0);
790 BOOST_TEST(data.at(1)->gradients()(1, 2) == 0);
791
792 displacements->setSampleAtTime(windowEnd, time::Sample(displacements->getDimensions(), Eigen::Vector4d{10, 10, 10, 10}, Eigen::MatrixXd(displacements->gradients()).setConstant(4)));
793 forces->setSampleAtTime(windowEnd, forces->sample());
794
795 acc.performAcceleration(data, windowStart, windowEnd);
796
797 // Check that store iteration works properly
798 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(0) == 4.6);
799 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(1) == 5.2);
800 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(2) == 5.8);
801 BOOST_TEST(data.at(0)->waveform().sample(windowEnd)(3) == 6.4);
802 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(0) == 0.184);
803 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(1) == 0.184);
804 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(2) == 0.184);
805 BOOST_TEST(data.at(1)->waveform().sample(windowEnd)(3) == 0.184);
806
807 BOOST_TEST(data.at(0)->gradients()(0, 0) == 1.6);
808 BOOST_TEST(data.at(0)->gradients()(0, 1) == 1.6);
809 BOOST_TEST(data.at(0)->gradients()(0, 2) == 1.6);
810 BOOST_TEST(data.at(0)->gradients()(1, 0) == 2.2);
811 BOOST_TEST(data.at(0)->gradients()(1, 1) == 2.8);
812 BOOST_TEST(data.at(0)->gradients()(1, 2) == 3.4);
813}
814
815#endif // not PRECICE_NO_MPI
816
std::map< int, PtrCouplingData > DataMap
BOOST_AUTO_TEST_CASE(testIQNIMVJPPWithSubsteps)
cplscheme::PtrCouplingData makeCouplingData(mesh::PtrData data, mesh::PtrMesh mesh, bool exchangeSubsteps)
Definition helper.hpp:10
void testIQNIMVJPP(bool exchangeSubsteps)
void testVIQNPP(bool exchangeSubsteps)
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
void performAcceleration(DataMap &cpldata, double windowStart, double windowEnd) final override
void initialize(const DataMap &cpldata) final override
void iterationsConverged(const DataMap &cpldata, double windowStart) final override
void initialize(const DataMap &cplData) final override
Initializes the acceleration.
void performAcceleration(DataMap &cplData, double windowStart, double windowEnd) final override
Performs one acceleration step.
void performAcceleration(DataMap &cplData, double windowStart, double windowEnd) override
Interface quasi-Newton with interface least-squares approximation.
Multi vector quasi-Newton update scheme.
Preconditioner that uses the constant user-defined factors to scale the quasi-Newton system.
Preconditioner that uses the recent residual to scale the quasi-Newton system.
Preconditioner that uses the residuals of all iterations of the current time window summed up to scal...
Describes a set of data values belonging to the vertices of a mesh.
Definition Data.hpp:26
std::shared_ptr< Preconditioner > PtrPreconditioner
contains implementations of acceleration schemes.
std::shared_ptr< CouplingData > PtrCouplingData
std::shared_ptr< Data > PtrData
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
cplscheme::PtrCouplingData makeCouplingData(mesh::PtrData data, mesh::PtrMesh mesh, bool exchangeSubsteps)
Definition helper.hpp:10
auto makeDummy3DMesh(size_t nVertices)
Definition Meshes.hpp:28
auto makeDummy2DMesh(size_t nVertices)
Definition Meshes.hpp:21
void append(Eigen::VectorXd &v, double value)
Main namespace of the precice library.
std::map< int, cplscheme::PtrCouplingData > DataMap