|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_matrix_block_h
17 #define dealii_matrix_block_h
39 template <
typename MatrixType>
45 template <
typename MatrixType>
49 template <
typename number>
111 template <
typename MatrixType>
157 operator MatrixType &();
158 operator const MatrixType &()
const;
167 const typename MatrixType::value_type
value);
184 template <
typename number>
186 add(
const std::vector<size_type> &indices,
188 const bool elide_zero_values =
true);
204 template <
typename number>
207 const std::vector<size_type> &col_indices,
209 const bool elide_zero_values =
true);
227 template <
typename number>
230 const std::vector<size_type> &col_indices,
231 const std::vector<number> & values,
232 const bool elide_zero_values =
true);
243 template <
typename number>
248 const number * values,
249 const bool elide_zero_values =
true,
250 const bool col_indices_are_sorted =
false);
257 template <
class VectorType>
266 template <
class VectorType>
275 template <
class VectorType>
284 template <
class VectorType>
301 <<
"Block index " << arg1 <<
" does not match " << arg2);
332 template <
class OTHER_MatrixType>
337 template <
typename number>
353 template <
typename MatrixType>
398 clear(
bool really_clean =
false);
443 template <
typename MatrixType>
519 clear(
bool really_clean =
false);
619 template <
typename MatrixType>
628 template <
typename number>
640 template <
typename MatrixType>
647 template <
typename MatrixType>
654 template <
typename MatrixType>
662 template <
typename MatrixType>
669 template <
typename MatrixType>
676 template <
typename MatrixType>
680 const typename MatrixType::value_type
value)
685 const std::pair<unsigned int, size_type> bi = row_indices.global_to_local(gi);
686 const std::pair<unsigned int, size_type> bj =
687 column_indices.global_to_local(gj);
689 Assert(bi.first == row, ExcBlockIndexMismatch(bi.first, row));
690 Assert(bj.first == column, ExcBlockIndexMismatch(bj.first, column));
696 template <
typename MatrixType>
697 template <
typename number>
700 const std::vector<size_type> &c_indices,
702 const bool elide_zero_values)
710 for (
size_type i = 0; i < row_indices.size(); ++i)
719 template <
typename MatrixType>
720 template <
typename number>
725 const number * values,
732 const std::pair<unsigned int, size_type> bi =
733 row_indices.global_to_local(b_row);
744 Assert(bi.first == row, ExcBlockIndexMismatch(bi.first, row));
748 const std::pair<unsigned int, size_type> bj =
749 column_indices.global_to_local(col_indices[j]);
750 Assert(bj.first == column, ExcBlockIndexMismatch(bj.first, column));
752 matrix.add(bi.second, bj.second, values[j]);
758 template <
typename MatrixType>
759 template <
typename number>
763 const bool elide_zero_values)
771 for (
size_type i = 0; i < indices.size(); ++i)
781 template <
typename MatrixType>
782 template <
typename number>
785 const std::vector<size_type> &col_indices,
786 const std::vector<number> & values,
787 const bool elide_zero_values)
801 template <
typename MatrixType>
802 template <
class VectorType>
810 template <
typename MatrixType>
811 template <
class VectorType>
819 template <
typename MatrixType>
820 template <
class VectorType>
828 template <
typename MatrixType>
829 template <
class VectorType>
837 template <
typename MatrixType>
847 template <
typename MatrixType>
851 const std::string &name)
858 template <
typename MatrixType>
862 for (
size_type i = 0; i < this->size(); ++i)
864 block(i).reinit(sparsity);
869 template <
typename MatrixType>
879 for (
size_type i = 0; i < this->size(); ++i)
886 template <
typename MatrixType>
890 return *this->read<ptr_type>(i);
894 template <
typename MatrixType>
898 return *this->entry<ptr_type>(i);
902 template <
typename MatrixType>
906 return this->entry<ptr_type>(i)->matrix;
913 template <
typename MatrixType>
917 , edge_flux_matrices(f)
921 template <
typename MatrixType>
925 return matrices.size();
929 template <
typename MatrixType>
933 const std::string &name)
937 p[0].column = column;
939 matrices.add(p, name);
942 matrices_in.add(p, name);
943 matrices_out.add(p, name);
945 if (edge_flux_matrices)
947 flux_matrices_up.add(p, name);
948 flux_matrices_down.add(p, name);
953 template <
typename MatrixType>
961 template <
typename MatrixType>
969 template <
typename MatrixType>
977 template <
typename MatrixType>
985 template <
typename MatrixType>
993 template <
typename MatrixType>
1001 template <
typename MatrixType>
1009 template <
typename MatrixType>
1017 template <
typename MatrixType>
1025 template <
typename MatrixType>
1033 template <
typename MatrixType>
1038 for (
size_type i = 0; i < this->size(); ++i)
1048 o[
level].column = col;
1055 template <
typename MatrixType>
1060 for (
size_type i = 0; i < this->size(); ++i)
1070 block_in(i)[
level].row = row;
1071 block_in(i)[
level].column = col;
1073 block_out(i)[
level].row = row;
1074 block_out(i)[
level].column = col;
1081 template <
typename MatrixType>
1086 for (
size_type i = 0; i < this->size(); ++i)
1096 block_up(i)[
level].row = row;
1097 block_up(i)[
level].column = col;
1099 block_down(i)[
level].row = row;
1100 block_down(i)[
level].column = col;
1107 template <
typename MatrixType>
1121 template <
typename MatrixType>
1131 clear_object(matrices);
1132 clear_object(matrices_in);
1133 clear_object(matrices_out);
1134 clear_object(flux_matrices_up);
1135 clear_object(flux_matrices_down);
void Tvmult(VectorType &w, const VectorType &v) const
const value_type & block_up(size_type i) const
MatrixBlock< MatrixType > & operator=(const MatrixBlock< MatrixType > &)=default
AnyData flux_matrices_up
The DG flux from the lower level to a level.
unsigned int max_level() const
void unsubscribe(std::atomic< bool > *const validity, const std::string &identifier="") const
void clear(bool really_clean=false)
unsigned int min_level() const
static ::ExceptionBase & ExcNotImplemented()
const std::string & name(const unsigned int i) const
Name of object at index.
unsigned int size() const
Number of stored data objects.
static ::ExceptionBase & ExcBlockIndexMismatch(size_type arg1, size_type arg2)
void clear_object(AnyData &)
Clear one of the matrix objects.
std::size_t memory_consumption() const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
const value_type & block_out(size_type i) const
std::size_t memory_consumption() const
void subscribe(std::atomic< bool > *const validity, const std::string &identifier="") const
void reinit_edge(const MGLevelObject< BlockSparsityPattern > &sparsity)
const bool edge_matrices
Flag for storing matrices_in and matrices_out.
MGMatrixBlockVector(const bool edge_matrices=false, const bool edge_flux_matrices=false)
static ::ExceptionBase & ExcNotQuadratic()
AnyData flux_matrices_down
The DG flux from a level to the lower level.
std::shared_ptr< value_type > ptr_type
AnyData matrices_out
The matrix from the refinement edge to the interior of a level.
type entry(const std::string &name)
Access to stored data object by name.
typename FullMatrix< number > ::value_type value_type
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void reinit_edge_flux(const MGLevelObject< BlockSparsityPattern > &sparsity)
const bool edge_flux_matrices
Flag for storing flux_matrices_up and flux_matrices_down.
void vmult_add(VectorType &w, const VectorType &v) const
unsigned int global_dof_index
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
const value_type & block(size_type i) const
#define DEAL_II_NAMESPACE_OPEN
void add(const size_type i, const size_type j, const typename MatrixType::value_type value)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void add(size_type row, size_type column, const std::string &name)
void reinit(const BlockSparsityPattern &sparsity)
const value_type & block_down(size_type i) const
@ matrix
Contents is actually a matrix.
AnyData matrices_in
The matrix from the interior of a level to the refinement edge.
void reinit(const BlockSparsityPattern &sparsity)
std::size_t memory_consumption() const
void add(size_type row, size_type column, const std::string &name)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNotInitialized()
AnyData matrices
The level matrices.
SparsityPatternType & block(const size_type row, const size_type column)
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
const BlockIndices & get_column_indices() const
#define Assert(cond, exc)
const value_type & block_in(size_type i) const
void reinit_matrix(const MGLevelObject< BlockSparsityPattern > &sparsity)
void reinit(MatrixBlock<::SparseMatrix< number >> &v, const BlockSparsityPattern &p)
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
const types::global_dof_index invalid_size_type
void Tvmult_add(VectorType &w, const VectorType &v) const
unsigned int size() const
void clear(bool really_clean=false)
BlockIndices column_indices
#define DEAL_II_NAMESPACE_CLOSE
void vmult(VectorType &w, const VectorType &v) const
#define DeclException2(Exception2, type1, type2, outsequence)
const BlockIndices & get_row_indices() const
const value_type & block(size_type i) const
MatrixType & matrix(size_type i)
void add(type entry, const std::string &name)
Add a new data object.