petsc_matrix.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_PETSC_MATRIX_H
21 #define LIBMESH_PETSC_MATRIX_H
22 
23 #include "libmesh/libmesh_common.h"
24 
25 #ifdef LIBMESH_HAVE_PETSC
26 
27 // Local includes
28 #include "libmesh/sparse_matrix.h"
29 #include "libmesh/petsc_macro.h"
30 #include "libmesh/parallel.h"
31 
32 // C++ includes
33 #include <algorithm>
34 
35 // Macro to identify and debug functions which should be called in
36 // parallel on parallel matrices but which may be called in serial on
37 // serial matrices. This macro will only be valid inside non-static
38 // PetscMatrix methods
39 #undef semiparallel_only
40 #ifndef NDEBUG
41 #include <cstring>
42 
43 #define semiparallel_only() do { if (this->initialized()) { const char * mytype; \
44  MatGetType(_mat,&mytype); \
45  if (!strcmp(mytype, MATSEQAIJ)) \
46  parallel_object_only(); } } while (0)
47 #else
48 #define semiparallel_only()
49 #endif
50 
51 
52 // Petsc include files.
53 #include <petscmat.h>
54 
55 
56 
57 namespace libMesh
58 {
59 
60 // Forward Declarations
61 template <typename T> class DenseMatrix;
62 
63 enum PetscMatrixType : int {
64  AIJ=0,
66 
67 
77 template <typename T>
78 class PetscMatrix final : public SparseMatrix<T>
79 {
80 public:
91  explicit
92  PetscMatrix (const Parallel::Communicator & comm_in);
93 
101  explicit
102  PetscMatrix (Mat m,
103  const Parallel::Communicator & comm_in);
104 
110  PetscMatrix (PetscMatrix &&) = delete;
111  PetscMatrix (const PetscMatrix &) = delete;
112  PetscMatrix & operator= (const PetscMatrix &) = delete;
113  PetscMatrix & operator= (PetscMatrix &&) = delete;
114  virtual ~PetscMatrix ();
115 
116  void set_matrix_type(PetscMatrixType mat_type);
117 
118  virtual void init (const numeric_index_type m,
119  const numeric_index_type n,
120  const numeric_index_type m_l,
121  const numeric_index_type n_l,
122  const numeric_index_type nnz=30,
123  const numeric_index_type noz=10,
124  const numeric_index_type blocksize=1) override;
125 
137  void init (const numeric_index_type m,
138  const numeric_index_type n,
139  const numeric_index_type m_l,
140  const numeric_index_type n_l,
141  const std::vector<numeric_index_type> & n_nz,
142  const std::vector<numeric_index_type> & n_oz,
143  const numeric_index_type blocksize=1);
144 
145  virtual void init () override;
146 
153 
154  virtual void clear () override;
155 
156  virtual void zero () override;
157 
158  virtual void zero_rows (std::vector<numeric_index_type> & rows, T diag_value = 0.0) override;
159 
160  virtual void close () override;
161 
162  virtual void flush () override;
163 
164  virtual numeric_index_type m () const override;
165 
166  virtual numeric_index_type n () const override;
167 
168  virtual numeric_index_type row_start () const override;
169 
170  virtual numeric_index_type row_stop () const override;
171 
172  virtual void set (const numeric_index_type i,
173  const numeric_index_type j,
174  const T value) override;
175 
176  virtual void add (const numeric_index_type i,
177  const numeric_index_type j,
178  const T value) override;
179 
180  virtual void add_matrix (const DenseMatrix<T> & dm,
181  const std::vector<numeric_index_type> & rows,
182  const std::vector<numeric_index_type> & cols) override;
183 
184  virtual void add_matrix (const DenseMatrix<T> & dm,
185  const std::vector<numeric_index_type> & dof_indices) override;
186 
187  virtual void add_block_matrix (const DenseMatrix<T> & dm,
188  const std::vector<numeric_index_type> & brows,
189  const std::vector<numeric_index_type> & bcols) override;
190 
191  virtual void add_block_matrix (const DenseMatrix<T> & dm,
192  const std::vector<numeric_index_type> & dof_indices) override
193  { this->add_block_matrix (dm, dof_indices, dof_indices); }
194 
208  virtual void add (const T a, const SparseMatrix<T> & X) override;
209 
210  virtual T operator () (const numeric_index_type i,
211  const numeric_index_type j) const override;
212 
213  virtual Real l1_norm () const override;
214 
215  virtual Real linfty_norm () const override;
216 
217  virtual bool closed() const override;
218 
225  virtual void print_personal(std::ostream & os=libMesh::out) const override;
226 
227  virtual void print_matlab(const std::string & name = "") const override;
228 
229  virtual void get_diagonal (NumericVector<T> & dest) const override;
230 
231  virtual void get_transpose (SparseMatrix<T> & dest) const override;
232 
237  void swap (PetscMatrix<T> &);
238 
247  Mat mat () { libmesh_assert (_mat); return _mat; }
248 
249 protected:
250 
261  virtual void _get_submatrix(SparseMatrix<T> & submatrix,
262  const std::vector<numeric_index_type> & rows,
263  const std::vector<numeric_index_type> & cols,
264  const bool reuse_submatrix) const override;
265 
266 private:
267 
271  Mat _mat;
272 
278 
280 };
281 
282 } // namespace libMesh
283 
284 #endif // #ifdef LIBMESH_HAVE_PETSC
285 #endif // LIBMESH_PETSC_MATRIX_H
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices) override
Definition: petsc_matrix.h:191
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
void set_matrix_type(PetscMatrixType mat_type)
Definition: petsc_matrix.C:121
virtual void print_matlab(const std::string &name="") const override
Definition: petsc_matrix.C:592
virtual void zero() override
Definition: petsc_matrix.C:475
PetscMatrix(const Parallel::Communicator &comm_in)
Definition: petsc_matrix.C:90
virtual numeric_index_type n() const override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: diff_context.h:40
void update_preallocation_and_zero()
Definition: petsc_matrix.C:468
virtual void init() override
Definition: petsc_matrix.C:348
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &brows, const std::vector< numeric_index_type > &bcols) override
Definition: petsc_matrix.C:805
virtual void print_personal(std::ostream &os=libMesh::out) const override
Definition: petsc_matrix.C:671
void swap(PetscMatrix< T > &)
dof_id_type numeric_index_type
Definition: id_types.h:92
PetscMatrix & operator=(const PetscMatrix &)=delete
virtual numeric_index_type row_stop() const override
virtual numeric_index_type m() const override
Definition: petsc_matrix.C:992
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
PetscMatrixType _mat_type
Definition: petsc_matrix.h:279
virtual void get_diagonal(NumericVector< T > &dest) const override
Definition: petsc_matrix.C:908
virtual void clear() override
Definition: petsc_matrix.C:528
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool closed() const override
virtual numeric_index_type row_start() const override
virtual void get_transpose(SparseMatrix< T > &dest) const override
Definition: petsc_matrix.C:921
virtual Real l1_norm() const override
Definition: petsc_matrix.C:546
static const bool value
Definition: xdr_io.C:109
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
Definition: petsc_matrix.C:778
SparseMatrix interface to PETSc Mat.
virtual void _get_submatrix(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols, const bool reuse_submatrix) const override
Definition: petsc_matrix.C:849
virtual void zero_rows(std::vector< numeric_index_type > &rows, T diag_value=0.0) override
Definition: petsc_matrix.C:496
A matrix object used for finite element assembly and numerics.
Definition: dense_matrix.h:54
OStreamProxy out(std::cout)
virtual void flush() override
Definition: petsc_matrix.C:977
virtual void close() override
Definition: petsc_matrix.C:957
virtual Real linfty_norm() const override
Definition: petsc_matrix.C:569