20 #ifndef LIBMESH_TRILINOS_EPETRA_MATRIX_H 21 #define LIBMESH_TRILINOS_EPETRA_MATRIX_H 25 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA 29 #include <Epetra_FECrsMatrix.h> 30 #include <Epetra_Map.h> 31 #include <Epetra_MpiComm.h> 34 #ifdef LIBMESH_TRILINOS_HAVE_EPETRAEXT 35 # include <EpetraExt_MatrixMatrix.h> 50 template <
typename T>
class DenseMatrix;
121 virtual void init ()
override;
123 virtual void clear ()
override;
125 virtual void zero ()
override;
127 virtual void close ()
override;
139 const T
value)
override;
143 const T
value)
override;
146 const std::vector<numeric_index_type> & rows,
147 const std::vector<numeric_index_type> & cols)
override;
150 const std::vector<numeric_index_type> & dof_indices)
override;
171 virtual bool closed()
const override;
192 Epetra_FECrsMatrix *
mat () { libmesh_assert(
_mat);
return _mat; }
194 const Epetra_FECrsMatrix *
mat ()
const { libmesh_assert(
_mat);
return _mat; }
229 #endif // LIBMESH_TRILINOS_HAVE_EPETRA 230 #endif // LIBMESH_TRILINOS_EPETRA_MATRIX_H
virtual numeric_index_type row_stop() const override
virtual void get_transpose(SparseMatrix< T > &dest) const override
bool _destroy_mat_on_exit
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
const Parallel::Communicator & comm() const
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
virtual numeric_index_type n() const override
virtual void init() override
dof_id_type numeric_index_type
virtual numeric_index_type m() const override
virtual void clear() override
EpetraMatrix & operator=(const EpetraMatrix &)=delete
virtual numeric_index_type row_start() const override
virtual void close() override
virtual void print_personal(std::ostream &os=libMesh::out) const override
EpetraMatrix(const Parallel::Communicator &comm)
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
virtual bool closed() const override
virtual void zero() override
Epetra_FECrsMatrix * _mat
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
Epetra_FECrsMatrix * mat()
A matrix object used for finite element assembly and numerics.
void swap(EpetraMatrix< T > &)
OStreamProxy out(std::cout)
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
virtual void update_sparsity_pattern(const SparsityPattern::Graph &) override
const Epetra_FECrsMatrix * mat() const
virtual bool need_full_sparsity_pattern() const override