libMesh
trilinos_epetra_matrix.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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_TRILINOS_EPETRA_MATRIX_H
21 #define LIBMESH_TRILINOS_EPETRA_MATRIX_H
22 
23 #include "libmesh/libmesh_common.h"
24 
25 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
26 
27 // Trilinos includes
28 #include "libmesh/ignore_warnings.h"
29 #include <Epetra_FECrsMatrix.h>
30 #include <Epetra_Map.h>
31 #include <Epetra_MpiComm.h>
32 
33 // The EpetraExt interface is only needed for EpetraMatrix::add()
34 #ifdef LIBMESH_TRILINOS_HAVE_EPETRAEXT
35 # include <EpetraExt_MatrixMatrix.h>
36 #endif
37 #include "libmesh/restore_warnings.h"
38 
39 // Local includes
40 #include "libmesh/sparse_matrix.h"
41 
42 // C++ includes
43 #include <algorithm>
44 #include <cstddef>
45 
46 namespace libMesh
47 {
48 
49 // Forward Declarations
50 template <typename T> class DenseMatrix;
51 
52 
53 
62 template <typename T>
63 class EpetraMatrix final : public SparseMatrix<T>
64 {
65 public:
77 
86  EpetraMatrix (Epetra_FECrsMatrix * m,
88 
94  EpetraMatrix (EpetraMatrix &&) = delete;
95  EpetraMatrix (const EpetraMatrix &) = delete;
96  EpetraMatrix & operator= (const EpetraMatrix &) = delete;
97  EpetraMatrix & operator= (EpetraMatrix &&) = delete;
98  virtual ~EpetraMatrix ();
99  virtual SparseMatrix<T> & operator=(const SparseMatrix<T> &) override
100  {
101  libmesh_not_implemented();
102  return *this;
103  }
104 
105  virtual SolverPackage solver_package() override
106  {
107  return TRILINOS_SOLVERS;
108  }
109 
113  virtual bool need_full_sparsity_pattern () const override
114  { return true; }
115 
121  virtual void update_sparsity_pattern (const SparsityPattern::Graph &) override;
122 
123  virtual void init (const numeric_index_type m,
124  const numeric_index_type n,
125  const numeric_index_type m_l,
126  const numeric_index_type n_l,
127  const numeric_index_type nnz=30,
128  const numeric_index_type noz=10,
129  const numeric_index_type blocksize=1) override;
130 
131  virtual void init (ParallelType = PARALLEL) override;
132 
136  virtual void clear () noexcept override;
137 
138  virtual void zero () override;
139 
140  virtual std::unique_ptr<SparseMatrix<T>> zero_clone () const override;
141 
142  virtual std::unique_ptr<SparseMatrix<T>> clone () const override;
143 
144  virtual void close () override;
145 
146  virtual numeric_index_type m () const override;
147 
148  virtual numeric_index_type n () const override;
149 
150  virtual numeric_index_type row_start () const override;
151 
152  virtual numeric_index_type row_stop () const override;
153 
154  virtual numeric_index_type col_start () const override;
155 
156  virtual numeric_index_type col_stop () const override;
157 
158  virtual void set (const numeric_index_type i,
159  const numeric_index_type j,
160  const T value) override;
161 
162  virtual void add (const numeric_index_type i,
163  const numeric_index_type j,
164  const T value) override;
165 
166  virtual void add_matrix (const DenseMatrix<T> & dm,
167  const std::vector<numeric_index_type> & rows,
168  const std::vector<numeric_index_type> & cols) override;
169 
170  virtual void add_matrix (const DenseMatrix<T> & dm,
171  const std::vector<numeric_index_type> & dof_indices) override;
172 
183  virtual void add (const T a, const SparseMatrix<T> & X) override;
184 
185  virtual T operator () (const numeric_index_type i,
186  const numeric_index_type j) const override;
187 
188  virtual Real l1_norm () const override;
189 
190  virtual Real linfty_norm () const override;
191 
192  virtual bool closed() const override;
193 
194  virtual void print_personal(std::ostream & os=libMesh::out) const override;
195 
196  virtual void get_diagonal (NumericVector<T> & dest) const override;
197 
198  virtual void get_transpose (SparseMatrix<T> & dest) const override;
199 
200  virtual void get_row(numeric_index_type i,
201  std::vector<numeric_index_type> & indices,
202  std::vector<T> & values) const override;
203 
207  void swap (EpetraMatrix<T> &);
208 
217  Epetra_FECrsMatrix * mat () { libmesh_assert(_mat); return _mat; }
218 
219  const Epetra_FECrsMatrix * mat () const { libmesh_assert(_mat); return _mat; }
220 
221 
222 private:
223 
227  Epetra_FECrsMatrix * _mat;
228 
232  Epetra_Map * _map;
233 
237  Epetra_CrsGraph * _graph;
238 
244 
250 };
251 
252 } // namespace libMesh
253 
254 #endif // LIBMESH_TRILINOS_HAVE_EPETRA
255 #endif // LIBMESH_TRILINOS_EPETRA_MATRIX_H
virtual SparseMatrix< T > & operator=(const SparseMatrix< T > &) override
This looks like a copy assignment operator, but note that, unlike normal copy assignment operators...
virtual std::unique_ptr< SparseMatrix< T > > clone() const override
virtual numeric_index_type row_stop() const override
virtual SolverPackage solver_package() override
virtual void get_transpose(SparseMatrix< T > &dest) const override
Copies the transpose of the matrix into dest, which may be *this.
virtual void clear() noexcept override
clear() is called from the destructor, so it should not throw.
bool _destroy_mat_on_exit
This boolean value should only be set to false for the constructor which takes an Epetra_FECrsMatrix ...
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
const Parallel::Communicator & comm() const
bool _use_transpose
Epetra has no GetUseTranspose so we need to keep track of whether we&#39;re transposed manually...
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const override
The libMesh namespace provides an interface to certain functionality in the library.
virtual void init(const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type nnz=30, const numeric_index_type noz=10, const numeric_index_type blocksize=1) override
Initialize SparseMatrix with the specified sizes.
virtual numeric_index_type col_stop() const override
virtual void get_row(numeric_index_type i, std::vector< numeric_index_type > &indices, std::vector< T > &values) const override
Get a row from the matrix.
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
Add the full matrix dm to the SparseMatrix.
virtual numeric_index_type n() const override
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
Set the element (i,j) to value.
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
virtual numeric_index_type col_start() const override
dof_id_type numeric_index_type
Definition: id_types.h:99
Epetra_CrsGraph * _graph
Holds the sparsity pattern.
virtual numeric_index_type m() const override
libmesh_assert(ctx)
EpetraMatrix & operator=(const EpetraMatrix &)=delete
virtual numeric_index_type row_start() const override
virtual void close() override
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
virtual void print_personal(std::ostream &os=libMesh::out) const override
Print the contents of the matrix to the screen in a package-personalized style, if available...
EpetraMatrix(const Parallel::Communicator &comm)
Constructor; initializes the matrix to be empty, without any structure, i.e.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real l1_norm() const override
virtual Real linfty_norm() const override
virtual void get_diagonal(NumericVector< T > &dest) const override
Copies the diagonal part of the matrix into dest.
OStreamProxy out
virtual bool closed() const override
virtual void zero() override
Set all entries to 0.
static const bool value
Definition: xdr_io.C:54
Epetra_FECrsMatrix * _mat
Actual Epetra datatype to hold matrix entries.
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
Add value to the element (i,j).
Epetra_FECrsMatrix * mat()
SolverPackage
Defines an enum for various linear solver packages.
Defines a dense matrix for use in Finite Element-type computations.
Definition: dof_map.h:75
void swap(EpetraMatrix< T > &)
Swaps the internal data pointers, no actual values are swapped.
Epetra_Map * _map
Holds the distributed Map.
This class provides a nice interface to the Epetra data structures for parallel, sparse matrices...
virtual void update_sparsity_pattern(const SparsityPattern::Graph &) override
Updates the matrix sparsity pattern.
const Epetra_FECrsMatrix * mat() const
virtual bool need_full_sparsity_pattern() const override
The EpetraMatrix needs the full sparsity pattern.
ParallelType
Defines an enum for parallel data structure types.