preCICE v3.1.2
Loading...
Searching...
No Matches
PreconditionerTest.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <stddef.h>
4#include <vector>
12#include "testing/Testing.hpp"
13
14BOOST_AUTO_TEST_SUITE(AccelerationTests)
15
16using namespace precice;
17using namespace precice::acceleration;
18using namespace precice::acceleration::impl;
19
21 Eigen::VectorXd _data;
22 Eigen::VectorXd _res;
23 Eigen::VectorXd _compareDataRes;
24 Eigen::VectorXd _compareDataResSum;
25 Eigen::VectorXd _compareDataResSum2;
26 Eigen::VectorXd _compareDataValue;
27 Eigen::VectorXd _compareDataConstant;
28
30 {
31 _data.resize(8);
32 _data << 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0;
33
34 _res.resize(8);
35 _res << 0.1, 0.1, 0.001, 0.001, 0.001, 0.001, 10.0, 20.0;
36
37 _compareDataRes.resize(8);
38 _compareDataRes << 7.07106781186547372897e+00,
39 1.41421356237309474579e+01,
40 1.50000000000000000000e+03,
41 2.00000000000000000000e+03,
42 2.50000000000000000000e+03,
43 3.00000000000000000000e+03,
44 3.13049516849970566046e-01,
45 3.57770876399966353265e-01;
46
47 _compareDataResSum.resize(8);
48 _compareDataResSum << 7.90585229434499154877e+01,
49 1.58117045886899830975e+02,
50 1.67708453051717078779e+04,
51 2.23611270735622783832e+04,
52 2.79514088419528488885e+04,
53 3.35416906103434157558e+04,
54 3.50007001329973377324e+00,
55 4.00008001519969536020e+00;
56
57 _compareDataResSum2.resize(8);
58 _compareDataResSum2 << 1.58113093108981217938e+02,
59 3.16226186217962435876e+02,
60 4.74339279326943596971e+02,
61 4.00008000319945455914e+00,
62 5.00010000399932064141e+00,
63 6.00012000479918228280e+00,
64 7.00014000559904481236e+00,
65 8.00016000639890734192e+00;
66
67 _compareDataValue.resize(8);
68 _compareDataValue << 4.47213595499957927704e-01,
69 8.94427190999915855407e-01,
70 3.23498319610315276940e-01,
71 4.31331092813753647075e-01,
72 5.39163866017192239255e-01,
73 6.46996639220630553879e-01,
74 6.58504607868518165859e-01,
75 7.52576694706877713514e-01;
76
77 _compareDataConstant.resize(8);
78 _compareDataConstant << 1.00000000000000002082e-03,
79 2.00000000000000004163e-03,
80 1.49999999999999977796e+00,
81 1.99999999999999955591e+00,
82 2.50000000000000044409e+00,
83 2.99999999999999955591e+00,
84 6.99999999999999883585e+05,
85 7.99999999999999650754e+05;
86 }
90};
91
92BOOST_FIXTURE_TEST_SUITE(ResPreconditionerTests, ResPreconditionerFixture)
93
94BOOST_AUTO_TEST_CASE(testResPreconditioner)
95{
96 PRECICE_TEST(1_rank);
98 svs.push_back(2);
99 svs.push_back(4);
100 svs.push_back(2);
101
102 ResidualPreconditioner precond(-1);
103
104 precond.initialize(svs);
105 Eigen::VectorXd backup = _data;
106
107 //should change
108 precond.update(false, _data, _res);
109 BOOST_TEST(precond.requireNewQR());
110 precond.newQRfulfilled();
111 precond.apply(_data);
112 BOOST_TEST(testing::equals(_data, _compareDataRes));
113 precond.revert(_data);
114 BOOST_TEST(testing::equals(_data, backup));
115
116 //should not change weights
117 precond.update(true, _data, _res * 10);
118 BOOST_TEST(not precond.requireNewQR());
119 precond.apply(_data);
120 BOOST_TEST(testing::equals(_data, _compareDataRes));
121 precond.revert(_data);
122 BOOST_TEST(testing::equals(_data, backup));
123}
124
125BOOST_AUTO_TEST_CASE(testResSumPreconditioner)
126{
127 PRECICE_TEST(1_rank);
129 svs.push_back(2);
130 svs.push_back(4);
131 svs.push_back(2);
132
133 ResidualSumPreconditioner precond(-1);
134
135 precond.initialize(svs);
136 Eigen::VectorXd backup = _data;
137
138 //should change, update twice to really test the summation
139 precond.update(false, _data, _res);
140 precond.update(false, _data, _res * 2);
141 BOOST_TEST(precond.requireNewQR());
142 precond.newQRfulfilled();
143 precond.apply(_data);
144 BOOST_TEST(testing::equals(_data, _compareDataResSum));
145
146 precond.revert(_data);
147 BOOST_TEST(testing::equals(_data, backup));
148
149 //should not change weights
150 precond.update(true, _data, _res * 10);
151 BOOST_TEST(not precond.requireNewQR());
152 precond.apply(_data);
153 BOOST_TEST(testing::equals(_data, _compareDataResSum));
154 precond.revert(_data);
155 BOOST_TEST(testing::equals(_data, backup));
156}
157
158BOOST_AUTO_TEST_CASE(testValuePreconditioner)
159{
160 PRECICE_TEST(1_rank);
162 svs.push_back(2);
163 svs.push_back(4);
164 svs.push_back(2);
165
166 ValuePreconditioner precond(-1);
167
168 precond.initialize(svs);
169 Eigen::VectorXd backup = _data;
170
171 //should change, since this is the first time window
172 precond.update(false, _data, _res);
173 BOOST_TEST(precond.requireNewQR());
174 precond.newQRfulfilled();
175 precond.apply(_data);
176 BOOST_TEST(testing::equals(_data, _compareDataValue));
177 precond.revert(_data);
178 BOOST_TEST(testing::equals(_data, backup));
179
180 //now no change
181 precond.update(false, _data, _res);
182 BOOST_TEST(not precond.requireNewQR());
183 precond.apply(_data);
184 BOOST_TEST(testing::equals(_data, _compareDataValue));
185 precond.revert(_data);
186 BOOST_TEST(testing::equals(_data, backup));
187
188 //should change weights
189 precond.update(true, _data * 2, _res);
190 BOOST_TEST(precond.requireNewQR());
191 precond.newQRfulfilled();
192}
193
194BOOST_AUTO_TEST_CASE(testConstPreconditioner)
195{
196 PRECICE_TEST(1_rank);
198 svs.push_back(2);
199 svs.push_back(4);
200 svs.push_back(2);
201
202 std::vector<double> factors;
203 factors.push_back(1e3);
204 factors.push_back(2.0);
205 factors.push_back(1e-5);
206
207 ConstantPreconditioner precond(factors);
208
209 precond.initialize(svs); //new weights already computed here
210 Eigen::VectorXd backup = _data;
211
212 // should have no effect
213 precond.update(false, _data, _res);
214 BOOST_TEST(not precond.requireNewQR());
215 precond.apply(_data);
216 BOOST_TEST(testing::equals(_data, _compareDataConstant));
217 precond.revert(_data);
218 BOOST_TEST(testing::equals(_data, backup));
219
220 //should not change weights
221 precond.update(true, _data, _res);
222 BOOST_TEST(not precond.requireNewQR());
223 precond.apply(_data);
224 BOOST_TEST(testing::equals(_data, _compareDataConstant));
225 precond.revert(_data);
226 BOOST_TEST(testing::equals(_data, backup));
227}
228
229BOOST_AUTO_TEST_CASE(testMultilpleMeshes)
230{
231 PRECICE_TEST(1_rank);
233 svs.push_back(3);
234 svs.push_back(5);
235
236 ResidualSumPreconditioner precond(-1);
237
238 precond.initialize(svs);
239 Eigen::VectorXd backup = _data;
240
241 //should change
242 precond.update(false, _data, _res);
243 BOOST_TEST(precond.requireNewQR());
244 precond.newQRfulfilled();
245 precond.apply(_data);
246 BOOST_TEST(testing::equals(_data, _compareDataResSum2));
247 precond.revert(_data);
248 BOOST_TEST(testing::equals(_data, backup));
249
250 //should not change weights
251 precond.update(true, _data, _res * 10);
252 BOOST_TEST(not precond.requireNewQR());
253 precond.apply(_data);
254 BOOST_TEST(testing::equals(_data, _compareDataResSum2));
255 precond.revert(_data);
256 BOOST_TEST(testing::equals(_data, backup));
257}
258
259#ifndef PRECICE_NO_MPI
260BOOST_AUTO_TEST_CASE(testParallelMatrixScaling)
261{
262 PRECICE_TEST(""_on(4_ranks).setupIntraComm());
263 //setup data
264 int localN = -1;
265 if (context.isPrimary()) {
266 localN = 2;
267 } else if (context.isRank(1)) {
268 localN = 1;
269 } else if (context.isRank(2)) {
270 localN = 0;
271 } else if (context.isRank(3)) {
272 localN = 1;
273 }
274
275 int globalN = 4;
276
277 Eigen::MatrixXd V(localN, 2);
278 Eigen::MatrixXd M(globalN, localN);
279 Eigen::VectorXd x(localN);
280 Eigen::MatrixXd V_back(localN, 2);
281 Eigen::MatrixXd M_back(globalN, localN);
282 Eigen::VectorXd x_back(localN);
283
284 if (context.isPrimary()) {
285 V(0, 0) = 1.0;
286 V(0, 1) = 2.0;
287 V(1, 0) = 3.0;
288 V(1, 1) = 4.0;
289 M(0, 0) = 1.0;
290 M(0, 1) = 2.0;
291 M(1, 0) = 3.0;
292 M(1, 1) = 4.0;
293 M(2, 0) = 1.0;
294 M(2, 1) = 2.0;
295 M(3, 0) = 3.0;
296 M(3, 1) = 4.0;
297 x(0) = 5.0;
298 x(1) = 5.0;
299 } else if (context.isRank(1)) {
300 V(0, 0) = 5.0;
301 V(0, 1) = 6.0;
302 M(0, 0) = 1.0;
303 M(1, 0) = 2.0;
304 M(2, 0) = 3.0;
305 M(3, 0) = 4.0;
306 x(0) = 5.0;
307 } else if (context.isRank(2)) {
308 } else if (context.isRank(3)) {
309 V(0, 0) = 7.0;
310 V(0, 1) = 8.0;
311 M(0, 0) = 1.0;
312 M(1, 0) = 2.0;
313 M(2, 0) = 3.0;
314 M(3, 0) = 4.0;
315 x(0) = 5.0;
316 }
317
318 V_back = V;
319 M_back = M;
320 x_back = x;
321
323 svs.push_back(localN);
324
325 ValuePreconditioner precond(-1);
326 precond.initialize(svs);
327 precond.update(true, x, x);
328 BOOST_TEST(precond.requireNewQR());
329
330 precond.apply(V);
331
332 BOOST_TEST(testing::equals(V, V_back * 0.1));
333
334 precond.revert(V);
335
336 BOOST_TEST(testing::equals(V, V_back));
337}
338#endif
339
BOOST_AUTO_TEST_SUITE(PreProcess)
BOOST_AUTO_TEST_SUITE_END()
BOOST_AUTO_TEST_CASE(testResPreconditioner)
#define PRECICE_TEST(...)
Definition Testing.hpp:27
Preconditioner that uses the constant user-defined factors to scale the quasi-Newton system.
virtual void initialize(std::vector< size_t > &svs)
initialize the preconditioner
void apply(Eigen::MatrixXd &M, bool transpose)
Apply preconditioner to matrix.
void revert(Eigen::MatrixXd &M, bool transpose)
Apply inverse preconditioner to matrix.
void newQRfulfilled()
to tell the preconditioner that QR-decomposition has been recomputed
void update(bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)
Update the scaling after every FSI iteration and require a new QR decomposition (if necessary)
virtual void initialize(std::vector< size_t > &svs)
initialize the preconditioner
bool requireNewQR()
returns true if a QR decomposition from scratch is necessary
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...
virtual void initialize(std::vector< size_t > &svs)
initialize the preconditioner
Preconditioner that uses the values from the previous time window to scale the quasi-Newton system.
contains implementations of acceleration schemes.
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
Main namespace of the precice library.
T push_back(T... args)