15#ifndef dealii_matrix_block_h
16#define dealii_matrix_block_h
38template <
typename MatrixType>
44 template <
typename MatrixType>
48 template <
typename number>
109template <
typename MatrixType>
155 operator MatrixType &();
156 operator const MatrixType &()
const;
165 const typename MatrixType::value_type value);
182 template <
typename number>
184 add(
const std::vector<size_type> &indices,
186 const bool elide_zero_values =
true);
202 template <
typename number>
205 const std::vector<size_type> &col_indices,
207 const bool elide_zero_values =
true);
225 template <
typename number>
228 const std::vector<size_type> &col_indices,
229 const std::vector<number> &values,
230 const bool elide_zero_values =
true);
241 template <
typename number>
246 const number *values,
247 const bool elide_zero_values =
true,
248 const bool col_indices_are_sorted =
false);
255 template <
typename VectorType>
257 vmult(VectorType &w,
const VectorType &v)
const;
264 template <
typename VectorType>
273 template <
typename VectorType>
275 Tvmult(VectorType &w,
const VectorType &v)
const;
282 template <
typename VectorType>
299 <<
"Block index " << arg1 <<
" does not match " << arg2);
330 template <
class OTHER_MatrixType>
335 template <
typename number>
350template <
typename MatrixType>
395 clear(
bool really_clean =
false);
439template <
typename MatrixType>
515 clear(
bool really_clean =
false);
615 template <
typename MatrixType>
624 template <
typename number>
636template <
typename MatrixType>
638 : row(
numbers::invalid_size_type)
639 , column(
numbers::invalid_size_type)
643template <
typename MatrixType>
650template <
typename MatrixType>
658template <
typename MatrixType>
665template <
typename MatrixType>
672template <
typename MatrixType>
676 const typename MatrixType::value_type value)
681 const std::pair<unsigned int, size_type> bi = row_indices.global_to_local(gi);
682 const std::pair<unsigned int, size_type> bj =
683 column_indices.global_to_local(gj);
685 Assert(bi.first == row, ExcBlockIndexMismatch(bi.first, row));
686 Assert(bj.first == column, ExcBlockIndexMismatch(bj.first, column));
688 matrix.add(bi.second, bj.second, value);
692template <
typename MatrixType>
693template <
typename number>
696 const std::vector<size_type> &c_indices,
698 const bool elide_zero_values)
706 for (
size_type i = 0; i < row_indices.size(); ++i)
715template <
typename MatrixType>
716template <
typename number>
721 const number *values,
728 const std::pair<unsigned int, size_type> bi =
729 row_indices.global_to_local(b_row);
740 Assert(bi.first == row, ExcBlockIndexMismatch(bi.first, row));
744 const std::pair<unsigned int, size_type> bj =
745 column_indices.global_to_local(col_indices[j]);
746 Assert(bj.first == column, ExcBlockIndexMismatch(bj.first, column));
748 matrix.add(bi.second, bj.second, values[j]);
754template <
typename MatrixType>
755template <
typename number>
759 const bool elide_zero_values)
767 for (
size_type i = 0; i < indices.size(); ++i)
777template <
typename MatrixType>
778template <
typename number>
781 const std::vector<size_type> &col_indices,
782 const std::vector<number> &values,
783 const bool elide_zero_values)
797template <
typename MatrixType>
798template <
typename VectorType>
806template <
typename MatrixType>
807template <
typename VectorType>
811 matrix.vmult_add(w, v);
815template <
typename MatrixType>
816template <
typename VectorType>
824template <
typename MatrixType>
825template <
typename VectorType>
829 matrix.Tvmult_add(w, v);
833template <
typename MatrixType>
843template <
typename MatrixType>
847 const std::string &name)
854template <
typename MatrixType>
858 for (
size_type i = 0; i < this->size(); ++i)
860 block(i).reinit(sparsity);
865template <
typename MatrixType>
875 for (
size_type i = 0; i < this->size(); ++i)
882template <
typename MatrixType>
886 return *this->read<ptr_type>(i);
890template <
typename MatrixType>
894 return *this->entry<ptr_type>(i);
898template <
typename MatrixType>
902 return this->entry<ptr_type>(i)->matrix;
909template <
typename MatrixType>
913 , edge_flux_matrices(f)
917template <
typename MatrixType>
921 return matrices.size();
925template <
typename MatrixType>
929 const std::string &name)
933 p[0].column = column;
935 matrices.add(p, name);
938 matrices_in.add(p, name);
939 matrices_out.add(p, name);
941 if (edge_flux_matrices)
943 flux_matrices_up.add(p, name);
944 flux_matrices_down.add(p, name);
949template <
typename MatrixType>
957template <
typename MatrixType>
965template <
typename MatrixType>
973template <
typename MatrixType>
981template <
typename MatrixType>
989template <
typename MatrixType>
997template <
typename MatrixType>
1005template <
typename MatrixType>
1013template <
typename MatrixType>
1021template <
typename MatrixType>
1029template <
typename MatrixType>
1034 for (
size_type i = 0; i < this->size(); ++i)
1044 o[
level].column = col;
1051template <
typename MatrixType>
1056 for (
size_type i = 0; i < this->size(); ++i)
1066 block_in(i)[
level].row = row;
1067 block_in(i)[
level].column = col;
1069 block_out(i)[
level].row = row;
1070 block_out(i)[
level].column = col;
1077template <
typename MatrixType>
1082 for (
size_type i = 0; i < this->size(); ++i)
1092 block_up(i)[
level].row = row;
1093 block_up(i)[
level].column = col;
1095 block_down(i)[
level].row = row;
1096 block_down(i)[
level].column = col;
1103template <
typename MatrixType>
1112 o[
level].matrix.clear();
1117template <
typename MatrixType>
1127 clear_object(matrices);
1128 clear_object(matrices_in);
1129 clear_object(matrices_out);
1130 clear_object(flux_matrices_up);
1131 clear_object(flux_matrices_down);
const std::string & name(const unsigned int i) const
Name of object at index.
type entry(const std::string &name)
Access to stored data object by name.
void add(type entry, const std::string &name)
Add a new data object.
unsigned int size() const
Number of stored data objects.
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
SparsityPatternType & block(const size_type row, const size_type column)
const BlockIndices & get_column_indices() const
const BlockIndices & get_row_indices() const
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&...args)
unsigned int max_level() const
unsigned int min_level() const
AnyData flux_matrices_up
The DG flux from the lower level to a level.
const value_type & block_down(size_type i) const
void clear_object(AnyData &)
Clear one of the matrix objects.
const value_type & block(size_type i) const
void clear(bool really_clean=false)
void reinit_edge_flux(const MGLevelObject< BlockSparsityPattern > &sparsity)
const bool edge_flux_matrices
Flag for storing flux_matrices_up and flux_matrices_down.
AnyData matrices_out
The matrix from the refinement edge to the interior of a level.
void reinit_matrix(const MGLevelObject< BlockSparsityPattern > &sparsity)
const value_type & block_in(size_type i) const
MGMatrixBlockVector(const bool edge_matrices=false, const bool edge_flux_matrices=false)
const value_type & block_out(size_type i) const
std::size_t memory_consumption() const
AnyData flux_matrices_down
The DG flux from a level to the lower level.
AnyData matrices
The level matrices.
void reinit_edge(const MGLevelObject< BlockSparsityPattern > &sparsity)
AnyData matrices_in
The matrix from the interior of a level to the refinement edge.
const bool edge_matrices
Flag for storing matrices_in and matrices_out.
void add(size_type row, size_type column, const std::string &name)
unsigned int size() const
const value_type & block_up(size_type i) const
const std::string & name(const unsigned int i) const
std::shared_ptr< value_type > ptr_type
void add(size_type row, size_type column, const std::string &name)
const value_type & block(size_type i) const
void clear(bool really_clean=false)
std::size_t memory_consumption() const
void reinit(const BlockSparsityPattern &sparsity)
MatrixType & matrix(size_type i)
MatrixBlock< MatrixType > & operator=(const MatrixBlock< MatrixType > &)=default
void vmult_add(VectorType &w, const VectorType &v) const
void reinit(const BlockSparsityPattern &sparsity)
BlockIndices column_indices
MatrixBlock(size_type i, size_type j)
void add(const size_type i, const size_type j, const typename MatrixType::value_type value)
std::size_t memory_consumption() const
void add(const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true)
void Tvmult_add(VectorType &w, const VectorType &v) const
void add(const std::vector< size_type > &indices, const FullMatrix< number > &full_matrix, const bool elide_zero_values=true)
void add(const size_type row_index, const std::vector< size_type > &col_indices, const std::vector< number > &values, const bool elide_zero_values=true)
void vmult(VectorType &w, const VectorType &v) const
void Tvmult(VectorType &w, const VectorType &v) const
typename MatrixType::value_type value_type
MatrixBlock(const MatrixBlock< MatrixType > &M)=default
void add(const size_type row, const size_type n_cols, const size_type *col_indices, const number *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
friend void internal::reinit(MatrixBlock<::SparseMatrix< number > > &v, const BlockSparsityPattern &p)
void unsubscribe(std::atomic< bool > *const validity, const std::string &identifier="") const
void subscribe(std::atomic< bool > *const validity, const std::string &identifier="") const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcBlockIndexMismatch(size_type arg1, size_type arg2)
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcNotQuadratic()
#define DEAL_II_NOT_IMPLEMENTED()
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
unsigned int global_dof_index