Go to the documentation of this file.
18 #include "libmesh/diagonal_matrix.h"
19 #include "libmesh/numeric_vector.h"
20 #include "libmesh/dense_matrix.h"
21 #include "libmesh/dof_map.h"
22 #include "libmesh/libmesh_common.h"
62 _diagonal->init(m, m_l);
71 _diagonal->init(this->_dof_map->n_dofs());
78 _diagonal->init(other, fast);
102 template <
typename T>
109 template <
typename T>
113 return _diagonal->size();
116 template <
typename T>
120 return _diagonal->size();
123 template <
typename T>
127 return _diagonal->first_local_index();
130 template <
typename T>
134 return _diagonal->last_local_index();
137 template <
typename T>
142 _diagonal->set(i,
value);
145 template <
typename T>
150 _diagonal->add(i,
value);
153 template <
typename T>
156 const std::vector<numeric_index_type> & rows,
157 const std::vector<numeric_index_type> & cols)
162 for (decltype(m) i = 0; i < m; ++i)
163 for (decltype(n) j = 0; j < n; ++j)
165 auto global_i = rows[i];
166 auto global_j = cols[j];
167 if (global_i == global_j)
168 _diagonal->add(global_i, dm(i, j));
172 template <
typename T>
175 const std::vector<numeric_index_type> & dof_indices)
177 _diagonal->add_vector(dm.
diagonal(), dof_indices);
180 template <
typename T>
184 auto x_diagonal = _diagonal->zero_clone();
186 _diagonal->add(a, *x_diagonal);
189 template <
typename T>
194 return (*_diagonal)(i);
199 template <
typename T>
203 return _diagonal->l1_norm();
206 template <
typename T>
210 return _diagonal->linfty_norm();
213 template <
typename T>
217 return _diagonal->closed();
220 template <
typename T>
224 _diagonal->print(os);
227 template <
typename T>
234 template <
typename T>
240 *diagonal_dest = *_diagonal;
242 libmesh_error_msg(
"DenseMatrix<T>::get_transpose currently only accepts another DenseMatrix<T> "
246 template <
typename T>
250 for (
auto row : rows)
251 _diagonal->set(row, val);
254 template <
typename T>
DiagonalMatrix(const Parallel::Communicator &comm)
Constructor; initializes the matrix to be empty, without any structure, i.e.
DenseVector< T > diagonal() const
Return the matrix diagonal.
void zero_rows(std::vector< numeric_index_type > &rows, T val=0) override
Sets all row entries to 0 then puts diag_value in the diagonal entry.
numeric_index_type row_stop() const override
DiagonalMatrix & operator=(DiagonalMatrix &&)=default
The libMesh namespace provides an interface to certain functionality in the library.
const NumericVector< T > & diagonal() const
virtual void get_diagonal(NumericVector< T > &dest) const =0
Copies the diagonal part of the matrix into dest.
numeric_index_type m() const override
Defines a dense matrix for use in Finite Element-type computations.
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.
Diagonal matrix class whose underlying storage is a vector.
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
void set(const numeric_index_type i, const numeric_index_type j, const T value) override
Set the element (i,j) to value.
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
bool closed() const override
void get_transpose(SparseMatrix< T > &dest) const override
Copies the transpose of the matrix into dest, which may be *this.
void get_diagonal(NumericVector< T > &dest) const override
Copies the diagonal part of the matrix into dest.
Real l1_norm() const override
Real linfty_norm() const override
numeric_index_type n() const override
dof_id_type numeric_index_type
void init() override
Initialize this matrix using the sparsity structure computed by dof_map.
std::unique_ptr< NumericVector< T > > _diagonal
Underlying diagonal matrix storage.
void zero() override
Set all entries to 0.
void close() override
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
void add(const numeric_index_type i, const numeric_index_type j, const T value) override
Add value to the element (i,j).
T operator()(const numeric_index_type i, const numeric_index_type j) const override
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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
numeric_index_type row_start() const override
void clear() override
Restores the SparseMatrix<T> to a pristine state.