libMesh
numeric_vector_test.h
Go to the documentation of this file.
1 #ifndef NUMERIC_VECTOR_TEST_H
2 #define NUMERIC_VECTOR_TEST_H
3 
4 // test includes
5 #include "test_comm.h"
6 
7 // libMesh includes
8 #include <libmesh/parallel.h>
9 #include <libmesh/fuzzy_equals.h>
10 
11 #include "libmesh_cppunit.h"
12 
13 #include <memory>
14 
15 #define NUMERICVECTORTEST \
16  CPPUNIT_TEST( testLocalize ); \
17  CPPUNIT_TEST( testLocalizeBase ); \
18  CPPUNIT_TEST( testLocalizeIndices ); \
19  CPPUNIT_TEST( testLocalizeIndicesBase ); \
20  CPPUNIT_TEST( testLocalizeToOne ); \
21  CPPUNIT_TEST( testLocalizeToOneBase ); \
22  CPPUNIT_TEST( testNorms ); \
23  CPPUNIT_TEST( testNormsBase ); \
24  CPPUNIT_TEST( testOperations ); \
25  CPPUNIT_TEST( testOperationsBase ); \
26  CPPUNIT_TEST( testWriteAndRead ); \
27  CPPUNIT_TEST( testWriteAndReadBase );
28 
29 
30 template <class DerivedClass>
31 class NumericVectorTest : public CppUnit::TestCase {
32 
33 protected:
35 
36  std::string libmesh_suite_name;
37 
39 
40 public:
41  void setUp()
42  {
43  // By default we'll use the whole communicator in parallel;
44  // Serial-only NumericVector subclasses will need to override
45  // this and set something else first.
46  if (!my_comm)
48 
49  block_size = 10;
50 
51  // a different size on each processor.
53  static_cast<unsigned int>(my_comm->rank());
54 
55  global_size = 0;
56  for (libMesh::processor_id_type p=0; p<my_comm->size(); p++)
57  global_size += (block_size + static_cast<unsigned int>(p));
58  }
59 
60  void tearDown()
61  {}
62 
63  template <class Base, class Derived>
64  void Operations()
65  {
66  auto v_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
67  Base & v = *v_ptr;
68 
70  first = v.first_local_index(),
71  last = v.last_local_index();
72 
73  for (libMesh::dof_id_type n=first; n != last; n++)
74  v.set (n, static_cast<libMesh::Number>(n+1));
75  v.close();
76 
77  auto v_clone = v.clone();
78  auto & vorig = *v_clone;
79 
80  CPPUNIT_ASSERT(relative_fuzzy_equals(v, vorig));
81 
82  v += v;
83 
84  CPPUNIT_ASSERT(!relative_fuzzy_equals(v, vorig));
85 
86  for (libMesh::dof_id_type n=first; n != last; n++)
87  LIBMESH_ASSERT_NUMBERS_EQUAL
88  (v(n), libMesh::Real(2*n+2),
90 
91  v *= v;
92 
93  for (libMesh::dof_id_type n=first; n != last; n++)
94  LIBMESH_ASSERT_NUMBERS_EQUAL
95  (v(n), libMesh::Real(4*n*n+8*n+4),
97 
98  v -= vorig;
99 
100  for (libMesh::dof_id_type n=first; n != last; n++)
101  LIBMESH_ASSERT_NUMBERS_EQUAL
102  (v(n), libMesh::Real(4*n*n+7*n+3),
104 
105  v /= vorig;
106 
107  for (libMesh::dof_id_type n=first; n != last; n++)
108  LIBMESH_ASSERT_NUMBERS_EQUAL
109  (v(n), libMesh::Real(4*n+3),
111 
112  v.add(1);
113  for (libMesh::dof_id_type n=first; n != last; n++)
114  LIBMESH_ASSERT_NUMBERS_EQUAL
115  (v(n), libMesh::Real(4*n+4),
117 
118  v.add(2, vorig);
119  for (libMesh::dof_id_type n=first; n != last; n++)
120  LIBMESH_ASSERT_NUMBERS_EQUAL
121  (v(n), libMesh::Real(6*n+6),
123 
124  v.scale(1/libMesh::Real(3));
125  for (libMesh::dof_id_type n=first; n != last; n++)
126  LIBMESH_ASSERT_NUMBERS_EQUAL
127  (v(n), libMesh::Real(2*n+2),
129 
130  v = 4;
131  for (libMesh::dof_id_type n=first; n != last; n++)
132  LIBMESH_ASSERT_NUMBERS_EQUAL
133  (v(n), libMesh::Real(4),
135 
136  v = vorig;
137  for (libMesh::dof_id_type n=first; n != last; n++)
138  LIBMESH_ASSERT_NUMBERS_EQUAL
139  (v(n), libMesh::Real(n+1),
141 
142  v.reciprocal();
143  for (libMesh::dof_id_type n=first; n != last; n++)
144  LIBMESH_ASSERT_NUMBERS_EQUAL
145  (v(n), 1/libMesh::Real(n+1),
147 
148  LIBMESH_ASSERT_NUMBERS_EQUAL
149  (v.dot(vorig), libMesh::Real(global_size),
151  }
152 
153  template <class Base, class Derived>
155  {
156  auto v_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
157  Base & v = *v_ptr;
158 
160  first = v.first_local_index(),
161  last = v.last_local_index();
162 
163  for (libMesh::dof_id_type n=first; n != last; n++)
164  v.set (n, static_cast<libMesh::Number>(n+1));
165  v.close();
166 
167  v.print_matlab(libmesh_suite_name+"_vector.m");
168 
169  // Let's make sure we don't have any race conditions; we have
170  // multiple SparseMatrix subclasses that might be trying to read
171  // and write the same file.
172 
174 
175 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
176  // We're not supporting complex reads quite yet
177  return;
178 #endif
179 
180  auto v_new_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
181  Base & v_new = *v_new_ptr;
182 
183  v_new.read_matlab(libmesh_suite_name+"_vector.m");
184 
186 
187  CPPUNIT_ASSERT_EQUAL(v.l1_norm_diff(v_new), libMesh::Real(0));
188  }
189 
190  template <class Base, class Derived>
191  void Norms()
192  {
193  auto v_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
194  Base & v = *v_ptr;
195 
197  first = v.first_local_index(),
198  last = v.last_local_index();
199 
200  for (libMesh::dof_id_type n=first; n != last; n++)
201  v.set (n, -static_cast<libMesh::Number>(n));
202  v.close();
203  for (libMesh::dof_id_type n=first; n != last; n++)
204  v.add (n, -static_cast<libMesh::Number>(n));
205  v.close();
206 
207  const libMesh::Real exact_l1 =
209  const libMesh::Real exact_l2 =
210  std::sqrt(libMesh::Real(global_size-1) *
211  (2*global_size) *
212  (2*global_size-1) / 3);
213  LIBMESH_ASSERT_NUMBERS_EQUAL(v.sum(), -exact_l1,
215  LIBMESH_ASSERT_FP_EQUAL(v.l1_norm(), exact_l1,
217  LIBMESH_ASSERT_FP_EQUAL(v.l2_norm(), exact_l2,
219  LIBMESH_ASSERT_FP_EQUAL(v.linfty_norm(), 2*(global_size-libMesh::Real(1)),
221 
222  auto u_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
223  Base & u = *u_ptr;
224  for (libMesh::dof_id_type n=first; n != last; n++)
225  u.set (n, -static_cast<libMesh::Number>(n * n));
226  u.close();
227 
228  Derived diff_derived(*my_comm, global_size, local_size);
229  Base & diff = diff_derived;
230  diff = u;
231  diff -= v;
232  diff.close();
233 
234  // Use a relative tolerance here; the norms are O(1e5) when our
235  // processor count gets big enough, at which point O(1e-12) is
236  // testing for exact equality...
237  const auto diff_norm = diff.l2_norm();
238  const auto norm_diff = u.l2_norm_diff(v);
239  LIBMESH_ASSERT_FP_EQUAL(diff_norm, norm_diff,
240  (std::abs(diff_norm)+std::abs(norm_diff)) *
242  }
243 
244  template <class Base, class Derived>
245  void Localize(bool to_one=false)
246  {
247  const libMesh::processor_id_type root_pid = 0;
248 
249  {
250  auto v_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
251  Base & v = *v_ptr;
252  std::vector<libMesh::Number> l(global_size);
253 
255  first = v.first_local_index(),
256  last = v.last_local_index();
257 
258  for (libMesh::dof_id_type n=first; n != last; n++)
259  v.set (n, static_cast<libMesh::Number>(n));
260  v.close();
261  for (libMesh::dof_id_type n=first; n != last; n++)
262  v.add (n, static_cast<libMesh::Number>(n));
263  v.close();
264 
265  if (!to_one)
266  v.localize(l);
267  else
268  v.localize_to_one(l,root_pid);
269 
270  if (!to_one || my_comm->rank() == root_pid)
271  // Yes I really mean v.size()
272  for (libMesh::dof_id_type i=0; i<v.size(); i++)
273  LIBMESH_ASSERT_NUMBERS_EQUAL
275 
276  for (libMesh::dof_id_type n=first; n != last; n++)
277  {
278  const auto value = static_cast<libMesh::Number>(n);
279  v.insert (&value, std::vector<libMesh::numeric_index_type>({n}));
280  }
281  v.close();
282 
283  if (!to_one)
284  v.localize(l);
285  else
286  v.localize_to_one(l,root_pid);
287 
288  if (!to_one || my_comm->rank() == root_pid)
289  // Yes I really mean v.size()
290  for (libMesh::dof_id_type i=0; i<v.size(); i++)
291  LIBMESH_ASSERT_NUMBERS_EQUAL
293  }
294  }
295 
296 
297  template <class Base, class Derived>
299  {
300  {
301  auto v_ptr = std::make_unique<Derived>(*my_comm, global_size, local_size);
302  Base & v = *v_ptr;
303 
304  // Let's try pulling the same number of entries from each processor
305  std::vector<libMesh::Number> values(block_size * my_comm->size());
306  std::vector<libMesh::dof_id_type> indices;
307  indices.reserve(block_size * my_comm->size());
308 
310  first = v.first_local_index(),
311  last = v.last_local_index();
312 
313  for (libMesh::dof_id_type n=first; n != last; n++)
314  v.set (n, static_cast<libMesh::Number>(n));
315  v.close();
316 
317  libMesh::dof_id_type end_index = 0;
318  for (libMesh::processor_id_type p=0; p<my_comm->size(); p++)
319  {
320  end_index += block_size + p;
321  for (unsigned int j = 0; j != block_size; ++j)
322  indices.push_back(end_index-j-1);
323  }
324 
325  v.localize(values, indices);
326 
327  end_index = 0;
328  for (libMesh::processor_id_type p=0; p<my_comm->size(); p++)
329  {
330  end_index += block_size + p;
331  for (unsigned int j = 0; j != block_size; ++j)
332  LIBMESH_ASSERT_NUMBERS_EQUAL
333  (values[p*block_size+j], libMesh::Real(end_index-j-1),
335  }
336  }
337  }
338 
340  {
341  LOG_UNIT_TEST;
342 
343  Localize<DerivedClass,DerivedClass>();
344  }
345 
347  {
348  LOG_UNIT_TEST;
349 
350  Localize<libMesh::NumericVector<libMesh::Number>,DerivedClass>();
351  }
352 
354  {
355  LOG_UNIT_TEST;
356 
357  Localize<DerivedClass,DerivedClass >(true);
358  }
359 
361  {
362  LOG_UNIT_TEST;
363 
364  Localize<libMesh::NumericVector<libMesh::Number>,DerivedClass>(true);
365  }
366 
368  {
369  LOG_UNIT_TEST;
370 
371  LocalizeIndices<DerivedClass,DerivedClass >();
372  }
373 
375  {
376  LOG_UNIT_TEST;
377 
378  LocalizeIndices<libMesh::NumericVector<libMesh::Number>,DerivedClass>();
379  }
380 
381  void testNorms()
382  {
383  LOG_UNIT_TEST;
384 
385  Norms<DerivedClass,DerivedClass >();
386  }
387 
389  {
390  LOG_UNIT_TEST;
391 
392  Norms<libMesh::NumericVector<libMesh::Number>,DerivedClass>();
393  }
394 
396  {
397  LOG_UNIT_TEST;
398 
399  Operations<DerivedClass,DerivedClass >();
400  }
401 
403  {
404  LOG_UNIT_TEST;
405 
406  Operations<libMesh::NumericVector<libMesh::Number>,DerivedClass>();
407  }
408 
410  {
411  LOG_UNIT_TEST;
412 
413  WriteAndRead<DerivedClass,DerivedClass >();
414  }
415 
417  {
418  LOG_UNIT_TEST;
419 
420  WriteAndRead<libMesh::NumericVector<libMesh::Number>,DerivedClass>();
421  }
422 };
423 
424 #endif
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
static constexpr Real TOLERANCE
void barrier() const
processor_id_type rank() const
uint8_t processor_id_type
Definition: id_types.h:104
void Localize(bool to_one=false)
processor_id_type size() const
libMesh::Parallel::Communicator * my_comm
unsigned int global_size
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const bool value
Definition: xdr_io.C:54
bool relative_fuzzy_equals(const T &var1, const T2 &var2, const Real tol=TOLERANCE *TOLERANCE)
Function to check whether two variables are equal within a relative tolerance.
Definition: fuzzy_equals.h:78
std::string libmesh_suite_name
uint8_t dof_id_type
Definition: id_types.h:67