Nektar++
testNekVector.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: testNekVector.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: Test code for NekVector
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38#include <boost/test/tools/floating_point_comparison.hpp>
39#include <boost/test/unit_test.hpp>
40
41#include <iostream>
42
43using std::cout;
44using std::endl;
45
46namespace Nektar::UnitTests
47{
48using namespace Nektar;
49
51{
52public:
54 {
55 }
56 explicit VectorTestClass(int data) : m_dataValue(data)
57 {
58 }
59 VectorTestClass(const VectorTestClass &in) = default;
61 {
62 }
63
65 {
67 return *this;
68 }
69
70 int value() const
71 {
72 return m_dataValue;
73 }
74
75private:
77};
78
79bool operator==(const VectorTestClass &lhs, const VectorTestClass &rhs)
80{
81 return lhs.value() == rhs.value();
82}
83
84bool operator!=(const VectorTestClass &lhs, const VectorTestClass &rhs)
85{
86 return !(lhs == rhs);
87}
88
89struct TenD
90{
91 static const unsigned int Value = 10;
92};
93
94BOOST_AUTO_TEST_CASE(TestNekVectorConstruction)
95{
96 // Test the constructors on a numeric type.
97 {
98 NekVector<double> p1(10, 2.7);
99 BOOST_CHECK(p1.GetDimension() == 10);
100 for (unsigned int i = 0; i < p1.GetDimension(); ++i)
101 {
102 BOOST_CHECK(p1(i) == 2.7);
103 BOOST_CHECK(p1[i] == 2.7);
104 }
105
106 NekVector<double> p2(1, 0.0);
107 BOOST_CHECK(p2.GetDimension() == 1);
108 BOOST_CHECK(p2(0) == 0.0);
109
110 NekVector<double> p3(p1);
111 BOOST_CHECK(p3.GetDimension() == 10);
112 for (unsigned int i = 0; i < p3.GetDimension(); ++i)
113 {
114 BOOST_CHECK(p3(i) == 2.7);
115 BOOST_CHECK(p3[i] == 2.7);
116 }
117
118 p2 = p3;
119 BOOST_CHECK(p2.GetDimension() == 10);
120 for (unsigned int i = 0; i < p2.GetDimension(); ++i)
121 {
122 BOOST_CHECK_EQUAL(p2(i), 2.7);
123 BOOST_CHECK_EQUAL(p2[i], 2.7);
124 }
125 }
126}
127
128BOOST_AUTO_TEST_CASE(TestNekVectorOperators)
129{
130 NekVector<double> v1(3, 1.0);
131 v1(0) = 1.1;
132 v1(1) = 1.2;
133 v1(2) = 1.3;
134
135 BOOST_CHECK(v1(0) == 1.1);
136 BOOST_CHECK(v1(1) == 1.2);
137 BOOST_CHECK(v1(2) == 1.3);
138
139 v1[0] = 1.4;
140 v1[1] = 1.5;
141 v1[2] = 1.6;
142
143 BOOST_CHECK(v1[0] == 1.4);
144 BOOST_CHECK(v1[1] == 1.5);
145 BOOST_CHECK(v1[2] == 1.6);
146
147 v1.x() = 1.7;
148 v1.y() = 1.8;
149 v1.z() = 1.9;
150
151 BOOST_CHECK(v1[0] == 1.7);
152 BOOST_CHECK(v1[1] == 1.8);
153 BOOST_CHECK(v1[2] == 1.9);
154
155 BOOST_CHECK(v1.x() == 1.7);
156 BOOST_CHECK(v1.y() == 1.8);
157 BOOST_CHECK(v1.z() == 1.9);
158
159 NekVector<double> v2(3, 1.0);
160 v2.x() = 1.7;
161 v2.y() = 1.8;
162 v2.z() = 1.9;
163
164 BOOST_CHECK(v1 == v2);
165 BOOST_CHECK(!(v1 != v2));
166
167 NekVector<double> v3(4, 2.0);
168 BOOST_CHECK(v3 != v1);
169
170 NekVector<double> v4(3, 0.0);
171 BOOST_CHECK(v4 != v1);
172}
173
174BOOST_AUTO_TEST_CASE(TestNekVectorArithmetic)
175{
176}
177
179{
180 double vals[] = {1, -2, 3};
181 NekVector<double> v(3, vals);
182
183 double epsilon = 1e-11;
184 BOOST_CHECK_EQUAL(v.L1Norm(), 1.0 + 2.0 + 3.0);
185 BOOST_CHECK_CLOSE(v.L2Norm(), sqrt(1.0 + 4.0 + 9.0), epsilon);
186 BOOST_CHECK_EQUAL(v.InfinityNorm(), 3);
187}
188
189BOOST_AUTO_TEST_CASE(TestMatrixVectorMultiply)
190{
191 {
192 // double matrix_buf[] = {1.0, 2.0, 3.0,
193 // 4.0, 5.0, 6.0,
194 // 7.0, 8.0, 9.0};
195 double matrix_buf[] = {1.0, 4.0, 7.0, 2.0, 5.0, 8.0, 3.0, 6.0, 9.0};
196 double vector_buf[] = {20.0, 30.0, 40.0};
197
198 NekMatrix<double> m(3, 3, matrix_buf);
199 NekVector<double> v(3, vector_buf);
200 NekVector<double> result = m * v;
201
202 double result_buf[] = {200.0, 470.0, 740.0};
203 NekVector<double> expected_result(3, result_buf);
204 BOOST_CHECK_EQUAL(result, expected_result);
205 BOOST_CHECK_EQUAL(3u, result.GetDimension());
206 }
207
208 {
209 // double matrix_buf[] = {1.0, 2.0, 3.0,
210 // 4.0, 5.0, 6.0,
211 // 7.0, 8.0, 9.0};
212 double matrix_buf[] = {1.0, 4.0, 7.0, 2.0, 5.0, 8.0, 3.0, 6.0, 9.0};
213 double vector_buf[] = {20.0, 30.0, 40.0};
214
215 std::shared_ptr<NekMatrix<double>> m(
216 new NekMatrix<double>(3, 3, matrix_buf));
217 NekMatrix<NekMatrix<double>, ScaledMatrixTag> s(2.0, m);
218 NekVector<double> v(3, vector_buf);
219 NekVector<double> result = s * v;
220
221 double result_buf[] = {400.0, 940.0, 1480.0};
222 NekVector<double> expected_result(3, result_buf);
223 BOOST_CHECK_EQUAL(result, expected_result);
224 BOOST_CHECK_EQUAL(3u, result.GetDimension());
225 }
226
227 {
228 // double m1_buf[] = {1.0, 2.0, 3.0, 4.0};
229 // double m2_buf[] = {5.0, 6.0, 7.0, 8.0};
230 // double m3_buf[] = {9.0, 10.0, 11.0, 12.0};
231 double m1_buf[] = {1.0, 3.0, 2.0, 4.0};
232 double m2_buf[] = {5.0, 7.0, 6.0, 8.0};
233 double m3_buf[] = {9.0, 11.0, 10.0, 12.0};
234 double vector_buf[] = {20.0, 30.0};
235
236 std::shared_ptr<NekMatrix<double>> m1(
237 new NekMatrix<double>(2, 2, m1_buf));
238 std::shared_ptr<NekMatrix<double>> m2(
239 new NekMatrix<double>(2, 2, m2_buf));
240 std::shared_ptr<NekMatrix<double>> m3(
241 new NekMatrix<double>(2, 2, m3_buf));
242
243 NekMatrix<NekMatrix<double>, BlockMatrixTag> b(3, 1, 2, 2);
244 b.SetBlock(0, 0, m1);
245 b.SetBlock(1, 0, m2);
246 b.SetBlock(2, 0, m3);
247
248 NekVector<double> v(2, vector_buf);
249 NekVector<double> result = b * v;
250
251 double result_buf[] = {80.0, 180.0, 280.0, 380.0, 480.0, 580.0};
252 NekVector<double> expected_result(6, result_buf);
253 BOOST_CHECK_EQUAL(result, expected_result);
254 BOOST_CHECK_EQUAL(6u, result.GetDimension());
255 }
256
257 {
258 double m_buf[] = {1, 4, 7, 2, 5, 8, 3, 6, 9};
259 double v_buf[] = {10, 11, 12};
260 double expected_result_buf[] = {68, 167, 266};
261 double expected_transpose_result_buf[] = {138, 171, 204};
262
263 NekMatrix<double> m(3, 3, m_buf);
264 NekVector<double> v_variable(3, v_buf);
265 NekVector<double> v_constant(3, v_buf);
266 NekVector<double> expected_result(3, expected_result_buf);
267 NekVector<double> expected_transpose_result(
268 3, expected_transpose_result_buf);
269
270 NekVector<double> constantResult = m * v_constant;
271 NekVector<double> variableResult = m * v_variable;
272 BOOST_CHECK_EQUAL(variableResult, expected_result);
273 BOOST_CHECK_EQUAL(constantResult, expected_result);
274
275 m.Transpose();
276 constantResult = m * v_constant;
277 variableResult = m * v_variable;
278 BOOST_CHECK_EQUAL(variableResult, expected_transpose_result);
279 BOOST_CHECK_EQUAL(constantResult, expected_transpose_result);
280
281 m.Transpose();
282 constantResult = m * v_constant;
283 variableResult = m * v_variable;
284 BOOST_CHECK_EQUAL(variableResult, expected_result);
285 BOOST_CHECK_EQUAL(constantResult, expected_result);
286
287 NekMatrix<double> transposed = Transpose(m);
288 constantResult = transposed * v_constant;
289 variableResult = transposed * v_variable;
290 BOOST_CHECK_EQUAL(variableResult, expected_transpose_result);
291 BOOST_CHECK_EQUAL(constantResult, expected_transpose_result);
292 }
293
294 {
295 double m_buf[] = {1, 4, 2, 5, 3, 6};
296 double v_non_transposed_buf[] = {10, 11, 12};
297 double v_transposed_buf[] = {20, 21};
298 double expected_result_buf[] = {68, 167};
299 double expected_transpose_result_buf[] = {104, 145, 186};
300
301 NekMatrix<double> m(2, 3, m_buf);
302 NekVector<double> v_non_transposed_variable(3, v_non_transposed_buf);
303 NekVector<double> v_non_transposed_constant(3, v_non_transposed_buf);
304 NekVector<double> v_transposed_variable(2, v_transposed_buf);
305 NekVector<double> v_transposed_constant(2, v_transposed_buf);
306 NekVector<double> expected_result(2, expected_result_buf);
307 NekVector<double> expected_transpose_result(
308 3, expected_transpose_result_buf);
309
310 NekVector<double> variableResult = m * v_non_transposed_variable;
311 NekVector<double> constantResult = m * v_non_transposed_constant;
312 BOOST_CHECK_EQUAL(variableResult, expected_result);
313 BOOST_CHECK_EQUAL(constantResult, expected_result);
314
315 m.Transpose();
316 variableResult = m * v_transposed_variable;
317 constantResult = m * v_transposed_constant;
318 BOOST_CHECK_EQUAL(variableResult, expected_transpose_result);
319 BOOST_CHECK_EQUAL(constantResult, expected_transpose_result);
320
321 m.Transpose();
322 variableResult = m * v_non_transposed_variable;
323 constantResult = m * v_non_transposed_constant;
324 BOOST_CHECK_EQUAL(variableResult, expected_result);
325 BOOST_CHECK_EQUAL(constantResult, expected_result);
326
327 NekMatrix<double> transposed = Transpose(m);
328 variableResult = transposed * v_transposed_variable;
329 constantResult = transposed * v_transposed_constant;
330 BOOST_CHECK_EQUAL(variableResult, expected_transpose_result);
331 BOOST_CHECK_EQUAL(constantResult, expected_transpose_result);
332 }
333}
334
335BOOST_AUTO_TEST_CASE(TestVectorConstructorsWithSizeArguments)
336{
337 {
338 double buf[] = {1.0, 2.0, 3.0, 4.0};
339 Array<OneD, double> a(4, buf);
340 NekVector<double> b(3, a);
341
342 BOOST_CHECK_EQUAL(b.GetRows(), 3u);
343 BOOST_CHECK_EQUAL(b.GetDimension(), 3u);
344 BOOST_CHECK_EQUAL(b[0], 1.0);
345 BOOST_CHECK_EQUAL(b[1], 2.0);
346 BOOST_CHECK_EQUAL(b[2], 3.0);
347
349 BOOST_CHECK_EQUAL(*iter, 1.0);
350 ++iter;
351 BOOST_CHECK_EQUAL(*iter, 2.0);
352 ++iter;
353 BOOST_CHECK_EQUAL(*iter, 3.0);
354 ++iter;
355 BOOST_CHECK(iter == b.end());
356 }
357
358 {
359 double buf[] = {1.0, 2.0, 3.0, 4.0};
361 NekVector<double> b(3, a);
362
363 BOOST_CHECK_EQUAL(b.GetRows(), 3u);
364 BOOST_CHECK_EQUAL(b.GetDimension(), 3u);
365 BOOST_CHECK_EQUAL(b[0], 1.0);
366 BOOST_CHECK_EQUAL(b[1], 2.0);
367 BOOST_CHECK_EQUAL(b[2], 3.0);
368
370 BOOST_CHECK_EQUAL(*iter, 1.0);
371 ++iter;
372 BOOST_CHECK_EQUAL(*iter, 2.0);
373 ++iter;
374 BOOST_CHECK_EQUAL(*iter, 3.0);
375 ++iter;
376 BOOST_CHECK(iter == b.end());
377 }
378
379 {
380 double buf[] = {1.0, 2.0, 3.0, 4.0};
381 Array<OneD, double> a(4, buf);
383
384 BOOST_CHECK_EQUAL(b.GetRows(), 3u);
385 BOOST_CHECK_EQUAL(b.GetDimension(), 3u);
386 BOOST_CHECK_EQUAL(b[0], 1.0);
387 BOOST_CHECK_EQUAL(b[1], 2.0);
388 BOOST_CHECK_EQUAL(b[2], 3.0);
389
391 BOOST_CHECK_EQUAL(*iter, 1.0);
392 ++iter;
393 BOOST_CHECK_EQUAL(*iter, 2.0);
394 ++iter;
395 BOOST_CHECK_EQUAL(*iter, 3.0);
396 ++iter;
397 BOOST_CHECK(iter == b.end());
398
399 b[0] = 5.0;
400 BOOST_CHECK_EQUAL(b[0], a[0]);
401 }
402}
403
404} // namespace Nektar::UnitTests
boost::call_traits< DataType >::reference x()
Definition: NekVector.cpp:275
DataType L1Norm() const
Definition: NekVector.cpp:426
unsigned int GetDimension() const
Returns the number of dimensions for the point.
Definition: NekVector.cpp:201
unsigned int GetRows() const
Definition: NekVector.cpp:206
boost::call_traits< DataType >::reference y()
Definition: NekVector.cpp:282
iterator begin()
Definition: NekVector.cpp:235
DataType * iterator
Definition: NekVector.hpp:109
iterator end()
Definition: NekVector.cpp:241
DataType InfinityNorm() const
Definition: NekVector.cpp:436
DataType L2Norm() const
Definition: NekVector.cpp:431
boost::call_traits< DataType >::reference z()
Definition: NekVector.cpp:289
VectorTestClass operator=(const VectorTestClass &rhs)
VectorTestClass(const VectorTestClass &in)=default
bool operator!=(const PointTestClass &lhs, const PointTestClass &rhs)
bool operator==(const PointTestClass &lhs, const PointTestClass &rhs)
BOOST_AUTO_TEST_CASE(TestGaussInt)
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
static const unsigned int Value