20 #ifndef LIBMESH_SPARSE_MATRIX_H 21 #define LIBMESH_SPARSE_MATRIX_H 42 template <
typename T>
class SparseMatrix;
43 template <
typename T>
class DenseMatrix;
45 namespace SparsityPattern {
class Graph; }
46 template <
typename T>
class NumericVector;
51 std::ostream & operator << (std::ostream & os, const SparseMatrix<T> & m);
64 class SparseMatrix :
public ReferenceCountedObject<SparseMatrix<T>>,
96 static std::unique_ptr<SparseMatrix<T>>
97 build(
const Parallel::Communicator &
comm,
154 virtual void init () = 0;
159 virtual void clear () = 0;
164 virtual void zero () = 0;
169 virtual void zero_rows (std::vector<numeric_index_type> & rows, T diag_value = 0.0);
175 virtual void close () = 0;
229 const std::vector<numeric_index_type> & rows,
230 const std::vector<numeric_index_type> & cols) = 0;
237 const std::vector<numeric_index_type> & dof_indices) = 0;
246 const std::vector<numeric_index_type> & brows,
247 const std::vector<numeric_index_type> & bcols);
254 const std::vector<numeric_index_type> & dof_indices)
296 virtual bool closed()
const = 0;
331 friend std::ostream & operator << <>(std::ostream & os,
const SparseMatrix<T> &
m);
346 libmesh_not_implemented();
355 const std::vector<numeric_index_type> & rows,
356 const std::vector<numeric_index_type> & cols)
const 371 const std::vector<numeric_index_type> & rows,
372 const std::vector<numeric_index_type> & cols)
const 415 const std::vector<numeric_index_type> & ,
416 const std::vector<numeric_index_type> & ,
419 libmesh_not_implemented();
443 template <
typename T>
444 std::ostream & operator << (std::ostream & os, const SparseMatrix<T> & m)
454 #endif // LIBMESH_SPARSE_MATRIX_H virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &brows, const std::vector< numeric_index_type > &bcols)
static std::unique_ptr< SparseMatrix< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
virtual bool initialized() const
virtual void create_submatrix(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
virtual void print_personal(std::ostream &os=libMesh::out) const =0
SparseMatrix(const Parallel::Communicator &comm)
void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const
virtual ~SparseMatrix()=default
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
const Parallel::Communicator & comm() const
virtual numeric_index_type row_stop() const =0
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
SolverPackage default_solver_package()
Manages the degrees of freedom (DOFs) in a simulation.
virtual void print_matlab(const std::string &="") const
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
dof_id_type numeric_index_type
virtual numeric_index_type m() const =0
virtual void zero_rows(std::vector< numeric_index_type > &rows, T diag_value=0.0)
virtual void get_diagonal(NumericVector< T > &dest) const =0
virtual void reinit_submatrix(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
virtual bool closed() const =0
void attach_dof_map(const DofMap &dof_map)
virtual bool need_full_sparsity_pattern() const
virtual void update_sparsity_pattern(const SparsityPattern::Graph &)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const
virtual numeric_index_type row_start() const =0
virtual void get_transpose(SparseMatrix< T > &dest) const =0
virtual void _get_submatrix(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const
virtual Real linfty_norm() const =0
SparseMatrix & operator=(const SparseMatrix &)=default
virtual Real l1_norm() const =0
A matrix object used for finite element assembly and numerics.
OStreamProxy out(std::cout)
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const =0
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
virtual numeric_index_type n() const =0
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices)