libMesh
Public Member Functions | Protected Attributes | List of all members
PetscMatrixTest Class Reference
Inheritance diagram for PetscMatrixTest:
[legend]

Public Member Functions

 PetscMatrixTest ()
 
 CPPUNIT_TEST_SUITE (PetscMatrixTest)
 
SPARSEMATRIXTEST CPPUNIT_TEST (testPetscBinaryRead)
 
 CPPUNIT_TEST (testPetscBinaryWrite)
 
 CPPUNIT_TEST (testPetscCopyFromHash)
 
 CPPUNIT_TEST_SUITE_END ()
 
void testPetscBinaryRead ()
 
void testPetscBinaryWrite ()
 
void testPetscHDF5Write ()
 
void testPetscCopyFromHash ()
 
void setUp ()
 
void tearDown ()
 
void setValues ()
 
void testValues ()
 
void testGetAndSet ()
 
void testReadMatlab (const std::string &filename)
 
void testReadMatlab1 ()
 
void testReadMatlab2 ()
 
void testReadMatlab4 ()
 
void testReadHDF5 ()
 
void testTransposeNorms ()
 
void testWriteAndRead ()
 
void testClone ()
 

Protected Attributes

std::string libmesh_suite_name
 
libMesh::Parallel::Communicatormy_comm
 
std::unique_ptr< PetscMatrix< Number > > matrix
 
libMesh::numeric_index_type nonsquare
 
libMesh::numeric_index_type local_m
 
libMesh::numeric_index_type local_n
 
libMesh::numeric_index_type global_m
 
libMesh::numeric_index_type global_n
 
const libMesh::Real _tolerance
 

Detailed Description

Definition at line 11 of file petsc_matrix_test.C.

Constructor & Destructor Documentation

◆ PetscMatrixTest()

PetscMatrixTest::PetscMatrixTest ( )
inline

Definition at line 14 of file petsc_matrix_test.C.

References libMesh::PerfLog::summarized_logs_enabled(), and unitlog.

14  :
17  this->libmesh_suite_name = "SparseMatrixTest";
18  else
19  this->libmesh_suite_name = "PetscMatrixTest";
20  }
bool summarized_logs_enabled()
Definition: perf_log.h:202
libMesh::PerfLog * unitlog
Definition: driver.C:173

Member Function Documentation

◆ CPPUNIT_TEST() [1/3]

SPARSEMATRIXTEST PetscMatrixTest::CPPUNIT_TEST ( testPetscBinaryRead  )

◆ CPPUNIT_TEST() [2/3]

PetscMatrixTest::CPPUNIT_TEST ( testPetscBinaryWrite  )

◆ CPPUNIT_TEST() [3/3]

PetscMatrixTest::CPPUNIT_TEST ( testPetscCopyFromHash  )

◆ CPPUNIT_TEST_SUITE()

PetscMatrixTest::CPPUNIT_TEST_SUITE ( PetscMatrixTest  )

◆ CPPUNIT_TEST_SUITE_END()

PetscMatrixTest::CPPUNIT_TEST_SUITE_END ( )

◆ setUp()

void SparseMatrixTest< PetscMatrix< Number > >::setUp ( )
inlineinherited

Definition at line 39 of file sparse_matrix_test.h.

40  {
41  // By default we'll use the whole communicator in parallel;
42  // Serial-only NumericVector subclasses will need to override
43  // this and set something else first.
44  if (!my_comm)
46 
47  matrix = std::make_unique<DerivedClass>(*my_comm);
48 
49  // Use the same even partitioning that we'll auto-deduce in matrix
50  // read, to make it easier to test parallel reads
51 
52  global_m = global_n = my_comm->size() * 5.75 - 1;
53  if (nonsquare)
54  global_n *= 1.3;
55 
57  row_start = my_comm->rank() * global_m / my_comm->size(),
58  row_stop = (my_comm->rank()+1) * global_m / my_comm->size();
60  col_start = my_comm->rank() * global_n / my_comm->size(),
61  col_stop = (my_comm->rank()+1) * global_n / my_comm->size();
62 
63  local_m = row_stop - row_start;
64  local_n = col_stop - col_start;
65 
66  // Let's just play around with locally dense blocks for now
67  matrix->init(global_m,
68  global_n,
69  local_m,
70  local_n,
71  /*nnz=*/local_n,
72  /*noz=*/5);
73  }
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
libMesh::Parallel::Communicator * my_comm
processor_id_type rank() const
processor_id_type size() const
std::unique_ptr< PetscMatrix< Number > > matrix
uint8_t dof_id_type
Definition: id_types.h:67

◆ setValues()

void SparseMatrixTest< PetscMatrix< Number > >::setValues ( )
inlineinherited

Definition at line 79 of file sparse_matrix_test.h.

80  {
81  std::vector<libMesh::numeric_index_type> rows(local_m);
82  std::iota(rows.begin(), rows.end(), matrix->row_start());
83  std::vector<libMesh::numeric_index_type> cols(local_n);
84  std::iota(cols.begin(), cols.end(), matrix->col_start());
85 
87 
88  for (auto i : libMesh::make_range(local_m))
89  for (auto j : libMesh::make_range(local_n))
90  local(i, j) = (i + 1) * (j + 1) * (my_comm->rank() + 1);
91 
92  matrix->zero();
93 
94  matrix->add_matrix(local, rows, cols);
95  matrix->close();
96  }
libMesh::Parallel::Communicator * my_comm
processor_id_type rank() const
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
std::unique_ptr< PetscMatrix< Number > > matrix
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
Defines a dense matrix for use in Finite Element-type computations.
Definition: dof_map.h:75

◆ tearDown()

void SparseMatrixTest< PetscMatrix< Number > >::tearDown ( )
inlineinherited

Definition at line 76 of file sparse_matrix_test.h.

76 {}

◆ testClone()

void SparseMatrixTest< PetscMatrix< Number > >::testClone ( )
inlineinherited

Definition at line 259 of file sparse_matrix_test.h.

260  {
261  LOG_UNIT_TEST;
262 
263  setValues();
264 
265  // Matrix must be closed before it can be cloned.
266  matrix->close();
267 
268  {
269  // Create copy, test that it can go out of scope
270  auto copy = matrix->clone();
271 
272  // Check that matrices have the same local/global sizes
273  CPPUNIT_ASSERT_EQUAL(copy->m(), matrix->m());
274  CPPUNIT_ASSERT_EQUAL(copy->n(), matrix->n());
275  CPPUNIT_ASSERT_EQUAL(copy->local_m(), matrix->local_m());
276  CPPUNIT_ASSERT_EQUAL(copy->row_start(), matrix->row_start());
277  CPPUNIT_ASSERT_EQUAL(copy->row_stop(), matrix->row_stop());
278 
279  // Check that copy has same values as original
280  LIBMESH_ASSERT_FP_EQUAL(copy->l1_norm(), matrix->l1_norm(), _tolerance);
281  CPPUNIT_ASSERT(relative_fuzzy_equals(*matrix, *copy));
282  copy->scale(2);
283  CPPUNIT_ASSERT(!relative_fuzzy_equals(*matrix, *copy));
284  }
285 
286  {
287  // Create zero copy
288  auto zero_copy = matrix->zero_clone();
289 
290  // Check that matrices have the same local/global sizes
291  CPPUNIT_ASSERT_EQUAL(zero_copy->m(), matrix->m());
292  CPPUNIT_ASSERT_EQUAL(zero_copy->n(), matrix->n());
293  CPPUNIT_ASSERT_EQUAL(zero_copy->local_m(), matrix->local_m());
294  CPPUNIT_ASSERT_EQUAL(zero_copy->row_start(), matrix->row_start());
295  CPPUNIT_ASSERT_EQUAL(zero_copy->row_stop(), matrix->row_stop());
296 
297  // Check that zero_copy has same values as original
298  LIBMESH_ASSERT_FP_EQUAL(0.0, zero_copy->l1_norm(), _tolerance);
299  }
300  }
std::unique_ptr< PetscMatrix< Number > > matrix
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

◆ testGetAndSet()

void SparseMatrixTest< PetscMatrix< Number > >::testGetAndSet ( )
inlineinherited

Definition at line 139 of file sparse_matrix_test.h.

140  {
141  LOG_UNIT_TEST;
142 
143  setValues();
144 
145  testValues();
146  }

◆ testPetscBinaryRead()

void PetscMatrixTest::testPetscBinaryRead ( )
inline

Definition at line 42 of file petsc_matrix_test.C.

References libMesh::TOLERANCE.

43  {
44  auto mat_to_read = std::make_unique<PetscMatrix<Number>>(*my_comm);
45 
46  // Petsc binary formats depend on sizeof(PetscInt) and
47  // sizeof(PetscScalar)
48 #if LIBMESH_DOF_ID_BYTES == 4 && LIBMESH_DEFAULT_DOUBLE_PRECISION
49  mat_to_read->read("matrices/geom_1_extraction_op.petsc32");
50 #elif LIBMESH_DOF_ID_BYTES == 8 && LIBMESH_DEFAULT_DOUBLE_PRECISION
51  mat_to_read->read("matrices/geom_1_extraction_op.petsc64");
52 #else
53  return;
54 #endif
55 
56  CPPUNIT_ASSERT_EQUAL(mat_to_read->m(), dof_id_type(27));
57  CPPUNIT_ASSERT_EQUAL(mat_to_read->n(), dof_id_type(27));
58 
59  // Our read_matlab partitioning doesn't necessarily match PETSc's
60  // MatLoad partitioning, so we'll partition a matrix to compare
61  // against manually. These particular files have bandwidth 8,
62  // so we'll ask for sufficient n_nz and n_oz to handle that
63  // regardless of partitioning.
64  auto mat_ascii_format = std::make_unique<PetscMatrix<Number>>(*my_comm);
65  mat_ascii_format->init(mat_to_read->m(), mat_to_read->n(),
66  mat_to_read->local_m(), mat_to_read->local_n(),
67  8, 7);
68 
69  mat_ascii_format->read_matlab("matrices/geom_1_extraction_op.m");
70 
71  mat_ascii_format->add(-1, *mat_to_read);
72  CPPUNIT_ASSERT_LESS(TOLERANCE, mat_ascii_format->l1_norm());
73  }
static constexpr Real TOLERANCE
uint8_t dof_id_type
Definition: id_types.h:67

◆ testPetscBinaryWrite()

void PetscMatrixTest::testPetscBinaryWrite ( )
inline

Definition at line 76 of file petsc_matrix_test.C.

References libMesh::TOLERANCE.

77  {
78  auto mat_to_read = std::make_unique<PetscMatrix<Number>>(*my_comm);
79  mat_to_read->read_matlab("matrices/geom_1_extraction_op.m");
80  mat_to_read->print_petsc_binary("geom_1_extraction_op.petsc");
81 
82  // Our read_matlab partitioning doesn't necessarily match PETSc's
83  // MatLoad partitioning, so we'll partition a matrix to compare
84  // against manually. These particular files have bandwidth 8,
85  // so we'll ask for sufficient n_nz and n_oz to handle that
86  // regardless of partitioning.
87  auto mat_reread = std::make_unique<PetscMatrix<Number>>(*my_comm);
88  mat_reread->init(mat_to_read->m(), mat_to_read->n(),
89  mat_to_read->local_m(), mat_to_read->local_n(),
90  8, 7);
91 
92  mat_reread->read_petsc_binary("geom_1_extraction_op.petsc");
93 
94  mat_to_read->add(-1, *mat_reread);
95  CPPUNIT_ASSERT_LESS(TOLERANCE, mat_to_read->l1_norm());
96  }
static constexpr Real TOLERANCE

◆ testPetscCopyFromHash()

void PetscMatrixTest::testPetscCopyFromHash ( )
inline

Definition at line 113 of file petsc_matrix_test.C.

References libMesh::PetscMatrix< T >::add(), libMesh::PetscMatrix< T >::copy_from_hash(), libMesh::PetscMatrix< T >::finish_initialization(), libMesh::PetscMatrix< T >::init_without_preallocation(), libMesh::make_range(), and libMesh::SparseMatrix< T >::use_hash_table().

114  {
116  const numeric_index_type M = my_comm->size();
117  const numeric_index_type m = 1;
118  const numeric_index_type blocksize = 1;
119  mat.use_hash_table(true);
120  mat.init_without_preallocation(M, M, m, m, blocksize);
121  mat.finish_initialization();
122 
123  // We'll just write a dense row
124  for (const auto j : make_range(my_comm->size()))
125  mat.add(my_comm->rank(), j, my_comm->rank() * j);
126 
127  auto copy = mat.copy_from_hash();
128  for (const auto j : make_range(my_comm->size()))
129  CPPUNIT_ASSERT_EQUAL((*copy)(my_comm->rank(), j), static_cast<Number>(my_comm->rank() * j));
130  }
libMesh::Parallel::Communicator * my_comm
processor_id_type rank() const
processor_id_type size() const
dof_id_type numeric_index_type
Definition: id_types.h:99
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
This class provides a nice interface to the PETSc C-based AIJ data structures for parallel...
Definition: petsc_matrix.h:61

◆ testPetscHDF5Write()

void PetscMatrixTest::testPetscHDF5Write ( )
inline

Definition at line 99 of file petsc_matrix_test.C.

References libMesh::TOLERANCE.

100  {
101  auto mat_to_read = std::make_unique<PetscMatrix<Number>>(*my_comm);
102  mat_to_read->read_matlab("matrices/geom_1_extraction_op.m");
103  mat_to_read->print_petsc_hdf5("geom_1_extraction_op.hdf5");
104 
105  auto mat_reread = std::make_unique<PetscMatrix<Number>>(*my_comm);
106  mat_reread->read_petsc_hdf5("geom_1_extraction_op.hdf5");
107 
108  mat_to_read->add(-1, *mat_reread);
109  CPPUNIT_ASSERT_LESS(TOLERANCE, mat_to_read->l1_norm());
110  }
static constexpr Real TOLERANCE

◆ testReadHDF5()

void SparseMatrixTest< PetscMatrix< Number > >::testReadHDF5 ( )
inlineinherited

Definition at line 194 of file sparse_matrix_test.h.

195  {
196 #ifdef LIBMESH_HAVE_HDF5
197  LOG_UNIT_TEST;
198 
199  matrix->clear();
200  auto matrix2 = std::make_unique<DerivedClass>(*my_comm);
201  matrix->read("matrices/geom_1_extraction_op.m");
202  matrix2->read("matrices/geom_1_extraction_op.h5");
203 
204  // We need some more SparseMatrix operators, but not today
205  CPPUNIT_ASSERT(matrix->l1_norm() == matrix2->l1_norm());
206  CPPUNIT_ASSERT(matrix->linfty_norm() == matrix2->linfty_norm());
207 #endif // LIBMESH_HAVE_HDF5
208  }
std::unique_ptr< PetscMatrix< Number > > matrix

◆ testReadMatlab()

void SparseMatrixTest< PetscMatrix< Number > >::testReadMatlab ( const std::string &  filename)
inlineinherited

Definition at line 149 of file sparse_matrix_test.h.

150  {
151  // Laspack doesn't handle non-square matrices)
152  if (matrix->solver_package() == libMesh::LASPACK_SOLVERS)
153  return;
154 
155  matrix->clear();
156 
157  auto matrix2 = std::make_unique<DerivedClass>(*my_comm);
158 
159  matrix->read(filename);
160 
161 #ifndef LIBMESH_HAVE_GZSTREAM
162  return;
163 #endif
164 
165  matrix2->read(std::string(filename)+".gz");
166 
167  // We need some more SparseMatrix operators, but not today
168  CPPUNIT_ASSERT(matrix->l1_norm() == matrix2->l1_norm());
169  CPPUNIT_ASSERT(matrix->linfty_norm() == matrix2->linfty_norm());
170  }
std::unique_ptr< PetscMatrix< Number > > matrix

◆ testReadMatlab1()

void SparseMatrixTest< PetscMatrix< Number > >::testReadMatlab1 ( )
inlineinherited

Definition at line 173 of file sparse_matrix_test.h.

174  {
175  LOG_UNIT_TEST;
176  testReadMatlab("matrices/geom_1_extraction_op.m");
177  }
void testReadMatlab(const std::string &filename)

◆ testReadMatlab2()

void SparseMatrixTest< PetscMatrix< Number > >::testReadMatlab2 ( )
inlineinherited

Definition at line 180 of file sparse_matrix_test.h.

181  {
182  LOG_UNIT_TEST;
183  testReadMatlab("matrices/geom_2_extraction_op.m");
184  }
void testReadMatlab(const std::string &filename)

◆ testReadMatlab4()

void SparseMatrixTest< PetscMatrix< Number > >::testReadMatlab4 ( )
inlineinherited

Definition at line 187 of file sparse_matrix_test.h.

188  {
189  LOG_UNIT_TEST;
190  testReadMatlab("matrices/geom_4_extraction_op.m");
191  }
void testReadMatlab(const std::string &filename)

◆ testTransposeNorms()

void SparseMatrixTest< PetscMatrix< Number > >::testTransposeNorms ( )
inlineinherited

Definition at line 211 of file sparse_matrix_test.h.

212  {
213  LOG_UNIT_TEST;
214 
215  setValues();
216 
217  auto matrix2 = std::make_unique<DerivedClass>(*my_comm);
218  matrix->get_transpose(*matrix2);
219  LIBMESH_ASSERT_FP_EQUAL(matrix->l1_norm(), matrix2->linfty_norm(), _tolerance);
220  LIBMESH_ASSERT_FP_EQUAL(matrix2->l1_norm(), matrix->linfty_norm(), _tolerance);
221  }
std::unique_ptr< PetscMatrix< Number > > matrix

◆ testValues()

void SparseMatrixTest< PetscMatrix< Number > >::testValues ( )
inlineinherited

Definition at line 99 of file sparse_matrix_test.h.

100  {
101  auto functor = [this]()
102  {
103  std::vector<libMesh::numeric_index_type> cols_to_get;
104  std::vector<libMesh::Number> values;
105  const libMesh::numeric_index_type col_start =
106  matrix->col_start();
108  libMesh::make_range(matrix->row_start(),
109  matrix->row_stop()))
110  {
111  matrix->get_row(i, cols_to_get, values);
112  for (libMesh::numeric_index_type col_j :
113  libMesh::index_range(cols_to_get))
114  {
115  CPPUNIT_ASSERT_EQUAL(cols_to_get[col_j], col_start + col_j);
116  LIBMESH_ASSERT_NUMBERS_EQUAL
117  ((i - matrix->row_start() + 1) * (col_j + 1) * (my_comm->rank() + 1),
118  values[col_j], _tolerance);
119  }
120  }
121  };
122 
123 #ifdef LIBMESH_HAVE_CXX11_THREAD
124  auto num_threads = std::min(unsigned(2),
125  std::max(
126  std::thread::hardware_concurrency(),
127  unsigned(1)));
128  std::vector<std::thread> threads(num_threads);
129  for (unsigned int thread = 0; thread < num_threads; ++thread)
130  threads[thread] = std::thread(functor);
131  std::for_each(threads.begin(), threads.end(),
132  [](std::thread & x){x.join();});
133 #else
134  functor();
135 #endif
136  }
libMesh::Parallel::Communicator * my_comm
processor_id_type rank() const
dof_id_type numeric_index_type
Definition: id_types.h:99
std::unique_ptr< PetscMatrix< Number > > matrix
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
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

◆ testWriteAndRead()

void SparseMatrixTest< PetscMatrix< Number > >::testWriteAndRead ( )
inlineinherited

Definition at line 224 of file sparse_matrix_test.h.

225  {
226  LOG_UNIT_TEST;
227 
228  setValues();
229 
230  // If we're working with serial matrices then just print one of
231  // them so they don't step on the others' toes.
232  //
233  // Use a very short filename, because we had a bug with that and
234  // we want to test it.
235  if (matrix->n_processors() > 1 ||
236  TestCommWorld->rank() == 0)
237  matrix->print_matlab("M.m");
238 
239  matrix->clear();
240 
241  // Let's make sure we don't have any race conditions; we have
242  // multiple SparseMatrix subclasses that might be trying to read
243  // and write the same file.
244 
246 
247 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
248  // We're not supporting complex reads quite yet
249  return;
250 #endif
251 
252  matrix->read("M.m");
253 
255 
256  testValues();
257  }
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:171
void barrier() const
processor_id_type rank() const
std::unique_ptr< PetscMatrix< Number > > matrix

Member Data Documentation

◆ _tolerance

const libMesh::Real SparseMatrixTest< PetscMatrix< Number > >::_tolerance
protectedinherited

Definition at line 314 of file sparse_matrix_test.h.

◆ global_m

libMesh::numeric_index_type SparseMatrixTest< PetscMatrix< Number > >::global_m
protectedinherited

Definition at line 310 of file sparse_matrix_test.h.

◆ global_n

libMesh::numeric_index_type SparseMatrixTest< PetscMatrix< Number > >::global_n
protectedinherited

Definition at line 310 of file sparse_matrix_test.h.

◆ libmesh_suite_name

std::string SparseMatrixTest< PetscMatrix< Number > >::libmesh_suite_name
protectedinherited

Definition at line 304 of file sparse_matrix_test.h.

◆ local_m

libMesh::numeric_index_type SparseMatrixTest< PetscMatrix< Number > >::local_m
protectedinherited

Definition at line 310 of file sparse_matrix_test.h.

◆ local_n

libMesh::numeric_index_type SparseMatrixTest< PetscMatrix< Number > >::local_n
protectedinherited

Definition at line 310 of file sparse_matrix_test.h.

◆ matrix

std::unique_ptr<PetscMatrix< Number > > SparseMatrixTest< PetscMatrix< Number > >::matrix
protectedinherited

Definition at line 308 of file sparse_matrix_test.h.

◆ my_comm

libMesh::Parallel::Communicator* SparseMatrixTest< PetscMatrix< Number > >::my_comm
protectedinherited

Definition at line 306 of file sparse_matrix_test.h.

◆ nonsquare

libMesh::numeric_index_type SparseMatrixTest< PetscMatrix< Number > >::nonsquare
protectedinherited

Definition at line 310 of file sparse_matrix_test.h.


The documentation for this class was generated from the following file: