Reference documentation for deal.II version 9.6.0
|
#include <deal.II/lac/tensor_product_matrix.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
TensorProductMatrixSymmetricSumCollection (const AdditionalData &additional_data=AdditionalData()) | |
void | reserve (const unsigned int size) |
template<typename T > | |
void | insert (const unsigned int index, const T &Ms, const T &Ks) |
void | finalize () |
void | apply_inverse (const unsigned int index, const ArrayView< Number > &dst_in, const ArrayView< const Number > &src_in) const |
std::size_t | memory_consumption () const |
std::size_t | storage_size () const |
Private Types | |
using | MatrixPairType = std::pair<Table<2, Number>, Table<2, Number>> |
using | MatrixPairTypeWithMask |
Private Attributes | |
const bool | compress_matrices |
const bool | precompute_inverse_diagonal |
std::vector< MatrixPairType > | mass_and_derivative_matrices |
std::map< MatrixPairTypeWithMask, unsigned int, internal::TensorProductMatrixSymmetricSum::MatrixPairComparator< Number > > | cache |
std::vector< unsigned int > | indices |
AlignedVector< Number > | mass_matrices |
AlignedVector< Number > | derivative_matrices |
AlignedVector< Number > | eigenvectors |
AlignedVector< Number > | eigenvalues |
AlignedVector< Number > | inverted_eigenvalues |
std::vector< unsigned int > | vector_ptr |
std::vector< unsigned int > | matrix_ptr |
std::vector< unsigned int > | vector_n_rows_1d |
A class similar to TensorProductMatrixSymmetricSum.
The class TensorProductMatrixSymmetricSum stores a 1d mass matrix, 1d stiffness matrix, eigenvalues and eigenvectors for each direction. If one uses one TensorProductMatrixSymmetricSum instance for, e.g., each cell, these quantities are stored for each cell. There is no possibility to reuse quantities between TensorProductMatrixSymmetricSum instances even if the values of the internal data structures might be the same. This class targets the case of many TensorProductMatrixSymmetricSum instances, where some of them might possibly share the underlying 1d matrices and hence re-use the same data.
This class is flexible and allows to interpret the parameter index
arbitrarily. In the case of an element-centric patch smoother, the index might correspond to the cell index and, in the case of a vertex patch smoother, the index might correspond to the vertex index defining the vertex patch. If n_rows_1d
is set to -1, the sizes of the mass matrices and the stiffness matrices can differ between cells, which might be useful if the class is used in the context of hp-refinement to construct a smoother.
Definition at line 334 of file tensor_product_matrix.h.
|
private |
Definition at line 336 of file tensor_product_matrix.h.
|
private |
Definition at line 338 of file tensor_product_matrix.h.
TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::TensorProductMatrixSymmetricSumCollection | ( | const AdditionalData & | additional_data = AdditionalData() | ) |
Constructor.
void TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::reserve | ( | const unsigned int | size | ) |
Allocate memory. The parameter specifies
the maximum value of the index used in invert() and apply_inverse().
void TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::insert | ( | const unsigned int | index, |
const T & | Ms, | ||
const T & | Ks ) |
For a given index
, attach the mass matrices Ms
and stiffness matrices Ks
to the stored data, looking out for compression possibilities.
void TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::finalize | ( | ) |
Finalize setup. This function computes, e.g., the eigenvalues and the eigenvectors.
void TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::apply_inverse | ( | const unsigned int | index, |
const ArrayView< Number > & | dst_in, | ||
const ArrayView< const Number > & | src_in ) const |
Apply the inverse matrix for a given index
.
std::size_t TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
std::size_t TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::storage_size | ( | ) | const |
Return the number of 1d matrices of each type stored internally. In the case that no compression could be performed, its value is the parameter passed to the function reserve() times the number of dimension. If compression could be performed, the value returned is less. In the optimal case (uniform Cartesian mesh), the value is one.
|
private |
Try to compress matrices.
Definition at line 423 of file tensor_product_matrix.h.
|
private |
Precompute inverse diagonal.
Definition at line 428 of file tensor_product_matrix.h.
|
private |
Container used to collect 1d matrices if no compression is requested. The memory is freed during finalize().
Definition at line 434 of file tensor_product_matrix.h.
|
private |
Container used during setup to determine the unique 1d matrices. The memory is freed during finalize().
Definition at line 444 of file tensor_product_matrix.h.
|
private |
Map from index to the storage position within mass_matrices, derivative_matrices, eigenvectors, and eigenvalues. If compression was not successful, this vector is empty, since the vectors can be access directly with the given index.
Definition at line 453 of file tensor_product_matrix.h.
|
private |
Vector of 1d mass matrices.
Definition at line 458 of file tensor_product_matrix.h.
|
private |
Vector of 1d derivative matrices.
Definition at line 463 of file tensor_product_matrix.h.
|
private |
Vector of eigenvectors.
Definition at line 468 of file tensor_product_matrix.h.
|
private |
Vector of eigenvalues.
Definition at line 473 of file tensor_product_matrix.h.
|
private |
Vector of inverted eigenvalues.
Definition at line 478 of file tensor_product_matrix.h.
|
private |
Pointer into mass_matrices, derivative_matrices, and eigenvalues.
Definition at line 483 of file tensor_product_matrix.h.
|
private |
Pointer into mass_matrices, derivative_matrices, and eigenvalues.
Definition at line 488 of file tensor_product_matrix.h.
|
private |
Number of rows in 1 of each cell.
Definition at line 493 of file tensor_product_matrix.h.