libMesh
sparse_matrix_test.h
Go to the documentation of this file.
1 #ifndef SPARSE_MATRIX_TEST_H
2 #define SPARSE_MATRIX_TEST_H
3 
4 // test includes
5 #include "test_comm.h"
6 
7 // libMesh includes
8 #include <libmesh/dense_matrix.h>
9 #include <libmesh/sparse_matrix.h>
10 #include <libmesh/fuzzy_equals.h>
11 
12 #include "libmesh_cppunit.h"
13 
14 // C++ includes
15 #include <vector>
16 #ifdef LIBMESH_HAVE_CXX11_THREAD
17 #include <thread>
18 #include <algorithm>
19 #endif
20 
21 
22 #define SPARSEMATRIXTEST \
23  CPPUNIT_TEST(testGetAndSet); \
24  CPPUNIT_TEST(testReadMatlab1); \
25  CPPUNIT_TEST(testReadMatlab2); \
26  CPPUNIT_TEST(testReadMatlab4); \
27  CPPUNIT_TEST(testWriteAndRead); \
28  CPPUNIT_TEST(testClone);
29 
30 
31 
32 template <class DerivedClass>
33 class SparseMatrixTest : public CppUnit::TestCase
34 {
35 public:
36 
37  void setUp()
38  {
39  // By default we'll use the whole communicator in parallel;
40  // Serial-only NumericVector subclasses will need to override
41  // this and set something else first.
42  if (!my_comm)
44 
45  matrix = std::make_unique<DerivedClass>(*my_comm);
46 
47  // Use the same even partitioning that we'll auto-deduce in matrix
48  // read, to make it easier to test parallel reads
49 
50  global_m = global_n = my_comm->size() * 5.75 - 1;
51  if (nonsquare)
52  global_n *= 1.3;
53 
55  row_start = my_comm->rank() * global_m / my_comm->size(),
56  row_stop = (my_comm->rank()+1) * global_m / my_comm->size();
58  col_start = my_comm->rank() * global_n / my_comm->size(),
59  col_stop = (my_comm->rank()+1) * global_n / my_comm->size();
60 
61  local_m = row_stop - row_start;
62  local_n = col_stop - col_start;
63 
64  // Let's just play around with locally dense blocks for now
65  matrix->init(global_m,
66  global_n,
67  local_m,
68  local_n,
69  /*nnz=*/local_n,
70  /*noz=*/5);
71  }
72 
73 
74  void tearDown() {}
75 
76 
77  void setValues()
78  {
79  std::vector<libMesh::numeric_index_type> rows(local_m);
80  std::iota(rows.begin(), rows.end(), matrix->row_start());
81  std::vector<libMesh::numeric_index_type> cols(local_n);
82  std::iota(cols.begin(), cols.end(), matrix->col_start());
83 
85 
86  for (auto i : libMesh::make_range(local_m))
87  for (auto j : libMesh::make_range(local_n))
88  local(i, j) = (i + 1) * (j + 1) * (my_comm->rank() + 1);
89 
90  matrix->zero();
91 
92  matrix->add_matrix(local, rows, cols);
93  matrix->close();
94  }
95 
96 
97  void testValues()
98  {
99  auto functor = [this]()
100  {
101  std::vector<libMesh::numeric_index_type> cols_to_get;
102  std::vector<libMesh::Number> values;
103  const libMesh::numeric_index_type col_start =
104  matrix->col_start();
106  libMesh::make_range(matrix->row_start(),
107  matrix->row_stop()))
108  {
109  matrix->get_row(i, cols_to_get, values);
110  for (libMesh::numeric_index_type col_j :
111  libMesh::index_range(cols_to_get))
112  {
113  CPPUNIT_ASSERT_EQUAL(cols_to_get[col_j], col_start + col_j);
114  LIBMESH_ASSERT_NUMBERS_EQUAL
115  ((i - matrix->row_start() + 1) * (col_j + 1) * (my_comm->rank() + 1),
116  values[col_j], _tolerance);
117  }
118  }
119  };
120 
121 #ifdef LIBMESH_HAVE_CXX11_THREAD
122  auto num_threads = std::min(unsigned(2),
123  std::max(
124  std::thread::hardware_concurrency(),
125  unsigned(1)));
126  std::vector<std::thread> threads(num_threads);
127  for (unsigned int thread = 0; thread < num_threads; ++thread)
128  threads[thread] = std::thread(functor);
129  std::for_each(threads.begin(), threads.end(),
130  [](std::thread & x){x.join();});
131 #else
132  functor();
133 #endif
134  }
135 
136 
138  {
139  LOG_UNIT_TEST;
140 
141  setValues();
142 
143  testValues();
144  }
145 
146  void testReadMatlab(const std::string & filename)
147  {
148  // Laspack doesn't handle non-square matrices)
149  if (matrix->solver_package() == libMesh::LASPACK_SOLVERS)
150  return;
151 
152  matrix->clear();
153 
154  auto matrix2 = std::make_unique<DerivedClass>(*my_comm);
155 
156  matrix->read(filename);
157 
158 #ifndef LIBMESH_HAVE_GZSTREAM
159  return;
160 #endif
161 
162  matrix2->read(std::string(filename)+".gz");
163 
164  // We need some more SparseMatrix operators, but not today
165  CPPUNIT_ASSERT(matrix->l1_norm() == matrix2->l1_norm());
166  CPPUNIT_ASSERT(matrix->linfty_norm() == matrix2->linfty_norm());
167  }
168 
169 
171  {
172  LOG_UNIT_TEST;
173  testReadMatlab("matrices/geom_1_extraction_op.m");
174  }
175 
176 
178  {
179  LOG_UNIT_TEST;
180  testReadMatlab("matrices/geom_2_extraction_op.m");
181  }
182 
183 
185  {
186  LOG_UNIT_TEST;
187  testReadMatlab("matrices/geom_4_extraction_op.m");
188  }
189 
190 
192  {
193  LOG_UNIT_TEST;
194 
195  setValues();
196 
197  // If we're working with serial matrices then just print one of
198  // them so they don't step on the others' toes.
199  //
200  // Use a very short filename, because we had a bug with that and
201  // we want to test it.
202  if (matrix->n_processors() > 1 ||
203  TestCommWorld->rank() == 0)
204  matrix->print_matlab("M.m");
205 
206  matrix->clear();
207 
208  // Let's make sure we don't have any race conditions; we have
209  // multiple SparseMatrix subclasses that might be trying to read
210  // and write the same file.
211 
213 
214 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
215  // We're not supporting complex reads quite yet
216  return;
217 #endif
218 
219  matrix->read("M.m");
220 
222 
223  testValues();
224  }
225 
226  void testClone()
227  {
228  LOG_UNIT_TEST;
229 
230  setValues();
231 
232  // Matrix must be closed before it can be cloned.
233  matrix->close();
234 
235  {
236  // Create copy, test that it can go out of scope
237  auto copy = matrix->clone();
238 
239  // Check that matrices have the same local/global sizes
240  CPPUNIT_ASSERT_EQUAL(copy->m(), matrix->m());
241  CPPUNIT_ASSERT_EQUAL(copy->n(), matrix->n());
242  CPPUNIT_ASSERT_EQUAL(copy->local_m(), matrix->local_m());
243  CPPUNIT_ASSERT_EQUAL(copy->row_start(), matrix->row_start());
244  CPPUNIT_ASSERT_EQUAL(copy->row_stop(), matrix->row_stop());
245 
246  // Check that copy has same values as original
247  LIBMESH_ASSERT_FP_EQUAL(copy->l1_norm(), matrix->l1_norm(), _tolerance);
248  CPPUNIT_ASSERT(relative_fuzzy_equals(*matrix, *copy));
249  copy->scale(2);
250  CPPUNIT_ASSERT(!relative_fuzzy_equals(*matrix, *copy));
251  }
252 
253  {
254  // Create zero copy
255  auto zero_copy = matrix->zero_clone();
256 
257  // Check that matrices have the same local/global sizes
258  CPPUNIT_ASSERT_EQUAL(zero_copy->m(), matrix->m());
259  CPPUNIT_ASSERT_EQUAL(zero_copy->n(), matrix->n());
260  CPPUNIT_ASSERT_EQUAL(zero_copy->local_m(), matrix->local_m());
261  CPPUNIT_ASSERT_EQUAL(zero_copy->row_start(), matrix->row_start());
262  CPPUNIT_ASSERT_EQUAL(zero_copy->row_stop(), matrix->row_stop());
263 
264  // Check that zero_copy has same values as original
265  LIBMESH_ASSERT_FP_EQUAL(0.0, zero_copy->l1_norm(), _tolerance);
266  }
267  }
268 
269 protected:
270 
271  std::string libmesh_suite_name;
272 
274 
275  std::unique_ptr<DerivedClass> matrix;
276 
278  local_m, local_n,
280 
282 };
283 
284 #endif // SPARSE_MATRIX_TEST_H
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
static constexpr Real TOLERANCE
void barrier() const
libMesh::Parallel::Communicator * my_comm
processor_id_type rank() const
libMesh::numeric_index_type local_n
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota was created back when std::iota was just an SGI STL extension.
Definition: utility.h:229
processor_id_type size() const
dof_id_type numeric_index_type
Definition: id_types.h:99
void testReadMatlab(const std::string &filename)
std::unique_ptr< DerivedClass > matrix
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
libMesh::numeric_index_type global_m
std::string libmesh_suite_name
const libMesh::Real _tolerance
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
libMesh::numeric_index_type nonsquare
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
libMesh::numeric_index_type global_n
Defines a dense matrix for use in Finite Element-type computations.
Definition: dof_map.h:75
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
libMesh::numeric_index_type local_m
uint8_t dof_id_type
Definition: id_types.h:67