Nektar++
testNekMatrix.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: testNekMatrix.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Tests NekMatrix functionality.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36#include <UnitTests/util.h>
37#include <boost/test/tools/floating_point_comparison.hpp>
38#include <boost/test/unit_test.hpp>
39#include <iostream>
40
41namespace Nektar
42{
43namespace MatrixUnitTests
44{
45BOOST_AUTO_TEST_CASE(TestFullMatrixInversion)
46{
47 {
48 double buf[] = {1.0, 3.0, 2.0, 4.0};
49 NekMatrix<double> m(2, 2, buf);
50 m.Invert();
51
52 BOOST_CHECK_CLOSE(*m(0, 0), -2.0, .00001);
53 BOOST_CHECK_CLOSE(*m(0, 1), 1.0, .00001);
54 BOOST_CHECK_CLOSE(*m(1, 0), 3.0 / 2.0, .00001);
55 BOOST_CHECK_CLOSE(*m(1, 1), -1.0 / 2.0, .00001);
56 }
57
58 {
59 double buf[] = {1.7, 5.6, 9.1, -3.4, -2.0, 4.5, 7.8, 8.2, -2.5,
60 5.0, 9.0, 3.6, 7.3, .98, 1.0, -12.6, 2.9, 6.4,
61 .02, -4.0, -1, 1.7, 5.6, 2.98, 5.0};
62 NekMatrix<double> m(5, 5, buf);
63 m.Invert();
64
65 BOOST_CHECK_CLOSE(*m(0, 0), 0.0005010995978, .000001);
66 BOOST_CHECK_CLOSE(*m(0, 1), -0.4704403712, .000001);
67 BOOST_CHECK_CLOSE(*m(0, 2), 0.2719063614, .000001);
68 BOOST_CHECK_CLOSE(*m(0, 3), -0.3941557805, .000001);
69 BOOST_CHECK_CLOSE(*m(0, 4), 0.09043166650, .000001);
70
71 BOOST_CHECK_CLOSE(*m(1, 0), -0.01087166322, .000001);
72 BOOST_CHECK_CLOSE(*m(1, 1), 0.3242048735, .000001);
73 BOOST_CHECK_CLOSE(*m(1, 2), -0.1605116333, .000001);
74 BOOST_CHECK_CLOSE(*m(1, 3), 0.09133673974, .000001);
75 BOOST_CHECK_CLOSE(*m(1, 4), 0.01293234281, .000001);
76
77 BOOST_CHECK_CLOSE(*m(2, 0), 0.06465598241, .000001);
78 BOOST_CHECK_CLOSE(*m(2, 1), 0.2593895742, .000001);
79 BOOST_CHECK_CLOSE(*m(2, 2), -0.09057233795, .000001);
80 BOOST_CHECK_CLOSE(*m(2, 3), 0.3106562459, .000001);
81 BOOST_CHECK_CLOSE(*m(2, 4), -0.1589713628, .000001);
82
83 BOOST_CHECK_CLOSE(*m(3, 0), -0.03464982687, .000001);
84 BOOST_CHECK_CLOSE(*m(3, 1), 0.2655177914, .000001);
85 BOOST_CHECK_CLOSE(*m(3, 2), -0.1016866584, .000001);
86 BOOST_CHECK_CLOSE(*m(3, 3), 0.2125363445, .000001);
87 BOOST_CHECK_CLOSE(*m(3, 4), -0.1099886183, .000001);
88
89 BOOST_CHECK_CLOSE(*m(4, 0), -0.02957895492, .000001);
90 BOOST_CHECK_CLOSE(*m(4, 1), -0.3518447037, .000001);
91 BOOST_CHECK_CLOSE(*m(4, 2), 0.2060393188, .000001);
92 BOOST_CHECK_CLOSE(*m(4, 3), -0.1411012255, .000001);
93 BOOST_CHECK_CLOSE(*m(4, 4), 0.1670437017, .000001);
94 }
95
96 {
97 }
98}
99
100BOOST_AUTO_TEST_CASE(TestDiagonalMatrixInversion)
101{
102 double buf[] = {1.0, 2.0, 3.0, 4.0};
103 NekMatrix<double> m(4, 4, buf, eDIAGONAL);
104 m.Invert();
105
106 BOOST_CHECK_EQUAL(m(0, 0), 1.0 / 1.0);
107 BOOST_CHECK_EQUAL(m(1, 1), 1.0 / 2.0);
108 BOOST_CHECK_EQUAL(m(2, 2), 1.0 / 3.0);
109 BOOST_CHECK_EQUAL(m(3, 3), 1.0 / 4.0);
110}
111
112BOOST_AUTO_TEST_CASE(TestNekMatrixConstruction)
113{
114 {
115 std::shared_ptr<Nektar::Matrix<double>> a(
116 new Nektar::NekMatrix<double>(3, 4));
117 std::shared_ptr<Nektar::NekMatrix<double>> b(
118 new Nektar::NekMatrix<double>(5, 6));
119
120 BOOST_CHECK_EQUAL(a->GetRows(), 3u);
121 BOOST_CHECK_EQUAL(a->GetColumns(), 4u);
122 BOOST_CHECK_EQUAL(a->GetSize()[0], 3u);
123 BOOST_CHECK_EQUAL(a->GetSize()[1], 4u);
124
125 BOOST_CHECK_EQUAL(b->GetRows(), 5u);
126 BOOST_CHECK_EQUAL(b->GetColumns(), 6u);
127 BOOST_CHECK_EQUAL(b->GetSize()[0], 5u);
128 BOOST_CHECK_EQUAL(b->GetSize()[1], 6u);
129
130 BOOST_CHECK_EQUAL(a->GetStorageType(), eFULL);
131 BOOST_CHECK_EQUAL(b->GetStorageType(), eFULL);
132
133 BOOST_CHECK_EQUAL(a->GetStorageSize(), 3 * 4);
134 BOOST_CHECK_EQUAL(b->GetStorageSize(), 5 * 6);
135 }
136
137 {
138 std::shared_ptr<Nektar::Matrix<double>> a(
140 std::shared_ptr<Nektar::NekMatrix<double>> b(
142
143 BOOST_CHECK_EQUAL(a->GetRows(), 3);
144 BOOST_CHECK_EQUAL(a->GetColumns(), 3);
145 BOOST_CHECK_EQUAL(a->GetSize()[0], 3);
146 BOOST_CHECK_EQUAL(a->GetSize()[1], 3);
147
148 BOOST_CHECK_EQUAL(b->GetRows(), 5);
149 BOOST_CHECK_EQUAL(b->GetColumns(), 5);
150 BOOST_CHECK_EQUAL(b->GetSize()[0], 5);
151 BOOST_CHECK_EQUAL(b->GetSize()[1], 5);
152
153 BOOST_CHECK_EQUAL(a->GetStorageType(), eDIAGONAL);
154 BOOST_CHECK_EQUAL(b->GetStorageType(), eDIAGONAL);
155
156 BOOST_CHECK_EQUAL(a->GetStorageSize(), 3);
157 BOOST_CHECK_EQUAL(b->GetStorageSize(), 5);
158 }
159}
160
161BOOST_AUTO_TEST_CASE(TestFullNekMatrixGetValue)
162{
164 double data[] = {1.0, 89.0, 0.0, 45.12, 2.0, -12.3,
165 892.2532, 76.12, 3.0, -56.7, 211.22, 45.23};
166
167 Nektar::NekMatrix<double> m(4, 3, data);
168
169 BOOST_CHECK_EQUAL(m(0, 0), 1.0);
170 BOOST_CHECK_EQUAL(m(0, 1), 2.0);
171 BOOST_CHECK_EQUAL(m(0, 2), 3.0);
172
173 BOOST_CHECK_EQUAL(m(1, 0), 89.0);
174 BOOST_CHECK_EQUAL(m(1, 1), -12.3);
175 BOOST_CHECK_EQUAL(m(1, 2), -56.7);
176
177 BOOST_CHECK_EQUAL(m(2, 0), 0.0);
178 BOOST_CHECK_EQUAL(m(2, 1), 892.2532);
179 BOOST_CHECK_EQUAL(m(2, 2), 211.22);
180
181 BOOST_CHECK_EQUAL(m(3, 0), 45.12);
182 BOOST_CHECK_EQUAL(m(3, 1), 76.12);
183 BOOST_CHECK_EQUAL(m(3, 2), 45.23);
184
185#ifdef NEKTAR_FULLDEBUG
186 BOOST_CHECK_THROW(m(0, 3), ErrorUtil::NekError);
187 BOOST_CHECK_THROW(m(4, 1), ErrorUtil::NekError);
188 BOOST_CHECK_THROW(m(4, 3), ErrorUtil::NekError);
189#endif // NEKTAR_FULLDEBUG
190
191 std::shared_ptr<Nektar::Matrix<double>> m1(
192 new Nektar::NekMatrix<double>(4, 3, data));
193 std::shared_ptr<Nektar::NekMatrix<double>> m2(
194 new Nektar::NekMatrix<double>(4, 3, data));
195
196 BOOST_CHECK_EQUAL((*m1)(0, 0), 1.0);
197 BOOST_CHECK_EQUAL((*m1)(0, 1), 2.0);
198 BOOST_CHECK_EQUAL((*m1)(0, 2), 3.0);
199
200 BOOST_CHECK_EQUAL((*m1)(1, 0), 89.0);
201 BOOST_CHECK_EQUAL((*m1)(1, 1), -12.3);
202 BOOST_CHECK_EQUAL((*m1)(1, 2), -56.7);
203
204 BOOST_CHECK_EQUAL((*m1)(2, 0), 0.0);
205 BOOST_CHECK_EQUAL((*m1)(2, 1), 892.2532);
206 BOOST_CHECK_EQUAL((*m1)(2, 2), 211.22);
207
208 BOOST_CHECK_EQUAL((*m1)(3, 0), 45.12);
209 BOOST_CHECK_EQUAL((*m1)(3, 1), 76.12);
210 BOOST_CHECK_EQUAL((*m1)(3, 2), 45.23);
211
212#ifdef NEKTAR_FULLDEBUG
213 BOOST_CHECK_THROW((*m1)(0, 3), ErrorUtil::NekError);
214 BOOST_CHECK_THROW((*m1)(4, 1), ErrorUtil::NekError);
215 BOOST_CHECK_THROW((*m1)(4, 3), ErrorUtil::NekError);
216#endif // NEKTAR_FULLDEBUG
217
218 BOOST_CHECK_EQUAL((*m2)(0, 0), 1.0);
219 BOOST_CHECK_EQUAL((*m2)(0, 1), 2.0);
220 BOOST_CHECK_EQUAL((*m2)(0, 2), 3.0);
221
222 BOOST_CHECK_EQUAL((*m2)(1, 0), 89.0);
223 BOOST_CHECK_EQUAL((*m2)(1, 1), -12.3);
224 BOOST_CHECK_EQUAL((*m2)(1, 2), -56.7);
225
226 BOOST_CHECK_EQUAL((*m2)(2, 0), 0.0);
227 BOOST_CHECK_EQUAL((*m2)(2, 1), 892.2532);
228 BOOST_CHECK_EQUAL((*m2)(2, 2), 211.22);
229
230 BOOST_CHECK_EQUAL((*m2)(3, 0), 45.12);
231 BOOST_CHECK_EQUAL((*m2)(3, 1), 76.12);
232 BOOST_CHECK_EQUAL((*m2)(3, 2), 45.23);
233
234#ifdef NEKTAR_FULLDEBUG
235 BOOST_CHECK_THROW((*m2)(0, 3), ErrorUtil::NekError);
236 BOOST_CHECK_THROW((*m2)(4, 1), ErrorUtil::NekError);
237 BOOST_CHECK_THROW((*m2)(4, 3), ErrorUtil::NekError);
238#endif // NEKTAR_FULLDEBUG
239}
240
241BOOST_AUTO_TEST_CASE(TestDiagonalMatrixGetValue)
242{
244 double data[] = {8.9, 3.4, 5.7};
245 Nektar::NekMatrix<double> m1(3, 3, data, eDIAGONAL);
246 std::shared_ptr<Nektar::Matrix<double>> m2(
247 new Nektar::NekMatrix<double>(3, 3, data, eDIAGONAL));
248 std::shared_ptr<Nektar::NekMatrix<double>> m3(
249 new Nektar::NekMatrix<double>(3, 3, data, eDIAGONAL));
250
251 BOOST_CHECK_EQUAL(m1(0, 0), 8.9);
252 BOOST_CHECK_EQUAL((*m2)(0, 0), 8.9);
253 BOOST_CHECK_EQUAL((*m3)(0, 0), 8.9);
254
255 BOOST_CHECK_EQUAL(m1(0, 1), 0.0);
256 BOOST_CHECK_EQUAL((*m2)(0, 1), 0.0);
257 BOOST_CHECK_EQUAL((*m3)(0, 1), 0.0);
258
259 BOOST_CHECK_EQUAL(m1(0, 2), 0.0);
260 BOOST_CHECK_EQUAL((*m2)(0, 2), 0.0);
261 BOOST_CHECK_EQUAL((*m3)(0, 2), 0.0);
262
263 BOOST_CHECK_EQUAL(m1(1, 0), 0.0);
264 BOOST_CHECK_EQUAL((*m2)(1, 0), 0.0);
265 BOOST_CHECK_EQUAL((*m3)(1, 0), 0.0);
266
267 BOOST_CHECK_EQUAL(m1(1, 1), 3.4);
268 BOOST_CHECK_EQUAL((*m2)(1, 1), 3.4);
269 BOOST_CHECK_EQUAL((*m3)(1, 1), 3.4);
270
271 BOOST_CHECK_EQUAL(m1(1, 2), 0.0);
272 BOOST_CHECK_EQUAL((*m2)(1, 2), 0.0);
273 BOOST_CHECK_EQUAL((*m3)(1, 2), 0.0);
274
275 BOOST_CHECK_EQUAL(m1(2, 0), 0.0);
276 BOOST_CHECK_EQUAL((*m2)(2, 0), 0.0);
277 BOOST_CHECK_EQUAL((*m3)(2, 0), 0.0);
278
279 BOOST_CHECK_EQUAL(m1(2, 1), 0.0);
280 BOOST_CHECK_EQUAL((*m2)(2, 1), 0.0);
281 BOOST_CHECK_EQUAL((*m3)(2, 1), 0.0);
282
283 BOOST_CHECK_EQUAL(m1(2, 2), 5.7);
284 BOOST_CHECK_EQUAL((*m2)(2, 2), 5.7);
285 BOOST_CHECK_EQUAL((*m3)(2, 2), 5.7);
286
287#ifdef NEKTAR_FULLDEBUG
288 BOOST_CHECK_THROW(m1(0, 4), ErrorUtil::NekError);
289 BOOST_CHECK_THROW((*m2)(0, 4), ErrorUtil::NekError);
290 BOOST_CHECK_THROW((*m3)(0, 4), ErrorUtil::NekError);
291
292 BOOST_CHECK_THROW(m1(4, 0), ErrorUtil::NekError);
293 BOOST_CHECK_THROW((*m2)(4, 0), ErrorUtil::NekError);
294 BOOST_CHECK_THROW((*m3)(4, 0), ErrorUtil::NekError);
295
296 BOOST_CHECK_THROW(m1(4, 4), ErrorUtil::NekError);
297 BOOST_CHECK_THROW((*m2)(4, 4), ErrorUtil::NekError);
298 BOOST_CHECK_THROW((*m3)(4, 4), ErrorUtil::NekError);
299#endif // NEKTAR_FULLDEBUG
300}
301
302BOOST_AUTO_TEST_CASE(TestFullNekMatrixSetValue)
303{
304 double data[] = {1.0, 3.0, 5.0, 2.0, 4.0, 6.0};
305 Nektar::NekMatrix<double> m1(3, 2, data);
306 std::shared_ptr<Nektar::Matrix<double>> m2(
307 new Nektar::NekMatrix<double>(3, 2, data));
308 std::shared_ptr<Nektar::NekMatrix<double>> m3(
309 new Nektar::NekMatrix<double>(3, 2, data));
310
311 m1.SetValue(0, 0, -1.0);
312 m2->SetValue(1, 1, -2.0);
313 m3->SetValue(2, 1, -3.0);
314
315 BOOST_CHECK_EQUAL(m1(0, 0), -1.0);
316 BOOST_CHECK_EQUAL(m1(0, 1), 2.0);
317 BOOST_CHECK_EQUAL(m1(1, 0), 3.0);
318 BOOST_CHECK_EQUAL(m1(1, 1), 4.0);
319 BOOST_CHECK_EQUAL(m1(2, 0), 5.0);
320 BOOST_CHECK_EQUAL(m1(2, 1), 6.0);
321
322 BOOST_CHECK_EQUAL((*m2)(0, 0), 1.0);
323 BOOST_CHECK_EQUAL((*m2)(0, 1), 2.0);
324 BOOST_CHECK_EQUAL((*m2)(1, 0), 3.0);
325 BOOST_CHECK_EQUAL((*m2)(1, 1), -2.0);
326 BOOST_CHECK_EQUAL((*m2)(2, 0), 5.0);
327 BOOST_CHECK_EQUAL((*m2)(2, 1), 6.0);
328
329 BOOST_CHECK_EQUAL((*m3)(0, 0), 1.0);
330 BOOST_CHECK_EQUAL((*m3)(0, 1), 2.0);
331 BOOST_CHECK_EQUAL((*m3)(1, 0), 3.0);
332 BOOST_CHECK_EQUAL((*m3)(1, 1), 4.0);
333 BOOST_CHECK_EQUAL((*m3)(2, 0), 5.0);
334 BOOST_CHECK_EQUAL((*m3)(2, 1), -3.0);
335}
336
337BOOST_AUTO_TEST_CASE(TestDiagonalNekMatrixSetValue)
338{
340 double data[] = {8.9, 3.4, 5.7};
341 Nektar::NekMatrix<double> m1(3, 3, data, eDIAGONAL);
342 std::shared_ptr<Nektar::Matrix<double>> m2(
343 new Nektar::NekMatrix<double>(3, 3, data, eDIAGONAL));
344 std::shared_ptr<Nektar::NekMatrix<double>> m3(
345 new Nektar::NekMatrix<double>(3, 3, data, eDIAGONAL));
346
347 m1.SetValue(0, 0, 1.0);
348 m2->SetValue(1, 1, 2.0);
349 m3->SetValue(2, 2, 3.0);
350
351 BOOST_CHECK_EQUAL(m1.GetStorageType(), eDIAGONAL);
352 BOOST_CHECK_EQUAL(m2->GetStorageType(), eDIAGONAL);
353 BOOST_CHECK_EQUAL(m3->GetStorageType(), eDIAGONAL);
354
355 BOOST_CHECK_EQUAL(m1(0, 0), 1.0);
356 BOOST_CHECK_EQUAL(m1(1, 1), 3.4);
357 BOOST_CHECK_EQUAL(m1(2, 2), 5.7);
358
359 BOOST_CHECK_EQUAL((*m2)(0, 0), 8.9);
360 BOOST_CHECK_EQUAL((*m2)(1, 1), 2.0);
361 BOOST_CHECK_EQUAL((*m2)(2, 2), 5.7);
362
363 BOOST_CHECK_EQUAL((*m3)(0, 0), 8.9);
364 BOOST_CHECK_EQUAL((*m3)(1, 1), 3.4);
365 BOOST_CHECK_EQUAL((*m3)(2, 2), 3.0);
366
367#ifdef NEKTAR_FULLDEBUG
368 BOOST_CHECK_THROW(m1.SetValue(0, 3, 2.0), ErrorUtil::NekError);
369 BOOST_CHECK_THROW(m1.SetValue(3, 0, 2.0), ErrorUtil::NekError);
370 BOOST_CHECK_THROW(m1.SetValue(3, 3, 2.0), ErrorUtil::NekError);
371
372 BOOST_CHECK_THROW(m2->SetValue(0, 3, 2.0), ErrorUtil::NekError);
373 BOOST_CHECK_THROW(m2->SetValue(3, 0, 2.0), ErrorUtil::NekError);
374 BOOST_CHECK_THROW(m2->SetValue(3, 3, 2.0), ErrorUtil::NekError);
375
376 BOOST_CHECK_THROW(m3->SetValue(0, 3, 2.0), ErrorUtil::NekError);
377 BOOST_CHECK_THROW(m3->SetValue(3, 0, 2.0), ErrorUtil::NekError);
378 BOOST_CHECK_THROW(m3->SetValue(3, 3, 2.0), ErrorUtil::NekError);
379#endif // NEKTAR_FULLDEBUG
380}
381
382BOOST_AUTO_TEST_CASE(TestFullFullMatrixAddition)
383{
384 double lhs_data[] = {1.0, 4.0, 2.0, 5.0, 3.0, 6.0};
385 double rhs_data[] = {10.0, 13.0, 11.0, 14.0, 12.0, 15.0};
386
387 Nektar::NekMatrix<double> lhs(2, 3, lhs_data);
388 Nektar::NekMatrix<double> rhs(2, 3, rhs_data);
389 Nektar::NekMatrix<double> result = lhs + rhs;
390
391 BOOST_CHECK_EQUAL(result.GetRows(), 2);
392 BOOST_CHECK_EQUAL(result.GetColumns(), 3);
393 BOOST_CHECK_EQUAL(result.GetStorageType(), eFULL);
394 BOOST_CHECK_EQUAL(result(0, 0), 11.0);
395 BOOST_CHECK_EQUAL(result(0, 1), 13.0);
396 BOOST_CHECK_EQUAL(result(0, 2), 15.0);
397 BOOST_CHECK_EQUAL(result(1, 0), 17.0);
398 BOOST_CHECK_EQUAL(result(1, 1), 19.0);
399 BOOST_CHECK_EQUAL(result(1, 2), 21.0);
400}
401} // namespace MatrixUnitTests
402
403namespace UnitTests
404{
405
406using namespace Nektar;
407
408BOOST_AUTO_TEST_CASE(testNekMatrixConstruction)
409{
410 // Basic, dense matrix construction.
411 {
412 double buf[] = {1.0, 4.0, 7.0, 10.0, 2.0, 5.0,
413 8.0, 11.0, 3.0, 6.0, 9.0, 12.0};
414 NekMatrix<double> dynamic_matrix(4, 3, buf);
415
416 BOOST_CHECK(dynamic_matrix.GetRows() == 4);
417 BOOST_CHECK(dynamic_matrix.GetColumns() == 3);
418
419 for (unsigned int i = 0; i < 4; ++i)
420 {
421 for (unsigned int j = 0; j < 3; ++j)
422 {
423 BOOST_CHECK(dynamic_matrix(i, j) == buf[4 * j + i]);
424 }
425 }
426 }
427}
428
429BOOST_AUTO_TEST_CASE(testNekMatrixBasicMath)
430{
431 // Addition tests.
432 {
433 double buf[] = {1.0, 4.0, 7.0, 2.0, 5.0, 8.0, 3.0, 6.0, 9.0};
434 NekMatrix<double> m1(3, 3, buf);
435 NekMatrix<double> m2(3, 3, buf);
436 NekMatrix<double> m3 = m1 + m2;
437
438 for (unsigned int i = 0; i < 3; ++i)
439 {
440 for (unsigned int j = 0; j < 3; ++j)
441 {
442 BOOST_CHECK(m3(i, j) == buf[3 * j + i] + buf[3 * j + i]);
443 }
444 }
445
446 NekMatrix<double> m4(3, 3, buf);
447 NekMatrix<double> m5(3, 3, buf);
448 NekMatrix<double> m6 = m4 + m5;
449
450 for (unsigned int i = 0; i < 3; ++i)
451 {
452 for (unsigned int j = 0; j < 3; ++j)
453 {
454 BOOST_CHECK(m6(i, j) == buf[3 * j + i] + buf[3 * j + i]);
455 }
456 }
457 }
458 // Multiply
459
460 {
461 double buf1[] = {1, 4, 7, 2, 5, 8, 3, 6, 9};
462 double buf2[] = {10, 15, 19, 11, 16, 20, 12, 17, 21};
463 NekMatrix<double> lhs(3, 3, buf1);
464 NekMatrix<double> rhs(3, 3, buf2);
465
466 NekMatrix<double> result = lhs * rhs;
467
468 BOOST_CHECK(result.GetRows() == 3);
469 BOOST_CHECK(result.GetColumns() == 3);
470
471 double epsilon = 1e-12;
472 BOOST_CHECK_CLOSE(*result(0, 0), 97.0, epsilon);
473 BOOST_CHECK_CLOSE(*result(0, 1), 103.0, epsilon);
474 BOOST_CHECK_CLOSE(*result(0, 2), 109.0, epsilon);
475
476 BOOST_CHECK_CLOSE(*result(1, 0), 229.0, epsilon);
477 BOOST_CHECK_CLOSE(*result(1, 1), 244.0, epsilon);
478 BOOST_CHECK_CLOSE(*result(1, 2), 259.0, epsilon);
479
480 BOOST_CHECK_CLOSE(*result(2, 0), 361.0, epsilon);
481 BOOST_CHECK_CLOSE(*result(2, 1), 385.0, epsilon);
482 BOOST_CHECK_CLOSE(*result(2, 2), 409.0, epsilon);
483 }
484
485 {
486 double buf1[] = {1, 4, 7, 2, 5, 8, 3, 6, 9};
487 double buf2[] = {10, 15, 19, 11, 16, 20, 12, 17, 21, 14, 18, 22};
488 NekMatrix<double> lhs(3, 3, buf1);
489 NekMatrix<double> rhs(3, 4, buf2);
490
491 NekMatrix<double> result = lhs * rhs;
492
493 BOOST_CHECK(result.GetRows() == 3);
494 BOOST_CHECK(result.GetColumns() == 4);
495
496 double epsilon = 1e-12;
497 BOOST_CHECK_CLOSE(*result(0, 0), 97.0, epsilon);
498 BOOST_CHECK_CLOSE(*result(0, 1), 103.0, epsilon);
499 BOOST_CHECK_CLOSE(*result(0, 2), 109.0, epsilon);
500 BOOST_CHECK_CLOSE(*result(0, 3), 116.0, epsilon);
501
502 BOOST_CHECK_CLOSE(*result(1, 0), 229.0, epsilon);
503 BOOST_CHECK_CLOSE(*result(1, 1), 244.0, epsilon);
504 BOOST_CHECK_CLOSE(*result(1, 2), 259.0, epsilon);
505 BOOST_CHECK_CLOSE(*result(1, 3), 278.0, epsilon);
506
507 BOOST_CHECK_CLOSE(*result(2, 0), 361.0, epsilon);
508 BOOST_CHECK_CLOSE(*result(2, 1), 385.0, epsilon);
509 BOOST_CHECK_CLOSE(*result(2, 2), 409.0, epsilon);
510 BOOST_CHECK_CLOSE(*result(2, 3), 440.0, epsilon);
511 }
512
513 //
514 // // Transpose
515 //
516 // // Determinant.
517 //
518 // // Invert/Check for singularity.
519 //
520 // // Eigenvalues/vectors
521 //
522 // // Condition number wrt various norms.
523 //
524 // // Various norm computations.
525 //
526 // // LU Decomposition? More appropriate in LinAlg?
527}
528
529} // namespace UnitTests
530} // namespace Nektar
BOOST_AUTO_TEST_CASE(TestFullMatrixInversion)
void RedirectCerrIfNeeded()
Definition: util.cpp:41
BOOST_AUTO_TEST_CASE(TestGaussInt)