https://mooseframework.inl.gov
MatrixEqualityCheck.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "MatrixEqualityCheck.h"
11 
12 #include "MooseUtils.h"
13 #include "libmesh/petsc_matrix.h"
14 
15 registerMooseObject("NavierStokesApp", MatrixEqualityCheck);
16 
19 {
21  params.addClassDescription("Report whether two matrices are the same or not.");
22  params.addRequiredParam<std::string>("mat1", "The petsc binary mat file containing matrix1");
23  params.addRequiredParam<std::string>("mat2", "The petsc binary mat file containing matrix2");
24  params.addParam<Real>(
25  "equivalence_tol", 1e-8, "The relative tolerance for comparing equivalence");
26  return params;
27 }
28 
30  : GeneralPostprocessor(parameters),
31  _equiv_tol(getParam<Real>("equivalence_tol")),
32  _mat1_name(getParam<std::string>("mat1")),
33  _mat2_name(getParam<std::string>("mat2"))
34 {
35 }
36 
37 void
39 {
40  _equiv = true;
41 
44 
45  if ((mat1->row_start() != mat2->row_start()) || (mat1->row_stop() != mat2->row_stop()) ||
46  (mat1->col_start() != mat2->col_start()) || (mat1->col_stop() != mat2->col_stop()))
47  {
48  _equiv = false;
49  return;
50  }
51 
52  for (const auto i : make_range(mat1->row_start(), mat1->row_stop()))
53  for (const auto j : make_range(mat1->col_start(), mat1->col_stop()))
54  {
55  const auto val1 = (*mat1)(i, j);
56  const auto val2 = (*mat2)(i, j);
57  if (!MooseUtils::relativeFuzzyEqual(val1, val2, _equiv_tol) &&
59  {
60  _equiv = false;
61  return;
62  }
63  }
64 }
65 
66 void
68 {
70 
71  if (_petsc_mat1)
72  LibmeshPetscCall(MatDestroy(&_petsc_mat1));
73  if (_petsc_mat2)
74  LibmeshPetscCall(MatDestroy(&_petsc_mat2));
75 }
76 
77 Real
79 {
80  return _equiv;
81 }
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
virtual Real getValue() const override
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
virtual void execute() override
std::unique_ptr< PetscMatrix< Number > > createMatrixFromFile(const libMesh::Parallel::Communicator &comm, Mat &petsc_mat, const std::string &binary_mat_file, unsigned int mat_number_to_load=1)
const Parallel::Communicator & _communicator
void addRequiredParam(const std::string &name, const std::string &doc_string)
static InputParameters validParams()
const std::string & _mat2_name
static InputParameters validParams()
void min(const T &r, T &o, Request &req) const
virtual void finalize() override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool relativeFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
MatrixEqualityCheck(const InputParameters &parameters)
IntRange< T > make_range(T beg, T end)
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
Checks if two matrices are the same by comparing their coefficients.
registerMooseObject("NavierStokesApp", MatrixEqualityCheck)
const std::string & _mat1_name