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