Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
Classes | Public Member Functions | Private Types | Private Attributes | List of all members
TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d > Class Template Reference

#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 = std::pair< std::bitset<::internal::VectorizedArrayTrait< Number >::width()>, MatrixPairType >
 

Private Attributes

const bool compress_matrices
 
const bool precompute_inverse_diagonal
 
std::vector< MatrixPairTypemass_and_derivative_matrices
 
std::map< MatrixPairTypeWithMask, unsigned int, internal::TensorProductMatrixSymmetricSum::MatrixPairComparator< Number > > cache
 
std::vector< unsigned intindices
 
AlignedVector< Number > mass_matrices
 
AlignedVector< Number > derivative_matrices
 
AlignedVector< Number > eigenvectors
 
AlignedVector< Number > eigenvalues
 
AlignedVector< Number > inverted_eigenvalues
 
std::vector< unsigned intvector_ptr
 
std::vector< unsigned intmatrix_ptr
 
std::vector< unsigned intvector_n_rows_1d
 

Detailed Description

template<int dim, typename Number, int n_rows_1d = -1>
class TensorProductMatrixSymmetricSumCollection< dim, Number, 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.

Member Typedef Documentation

◆ MatrixPairType

template<int dim, typename Number , int n_rows_1d = -1>
using TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::MatrixPairType = std::pair<Table<2, Number>, Table<2, Number> >
private

Definition at line 336 of file tensor_product_matrix.h.

◆ MatrixPairTypeWithMask

template<int dim, typename Number , int n_rows_1d = -1>
using TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::MatrixPairTypeWithMask = std::pair< std::bitset<::internal::VectorizedArrayTrait<Number>::width()>, MatrixPairType>
private

Definition at line 338 of file tensor_product_matrix.h.

Constructor & Destructor Documentation

◆ TensorProductMatrixSymmetricSumCollection()

template<int dim, typename Number , int n_rows_1d = -1>
TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::TensorProductMatrixSymmetricSumCollection ( const AdditionalData additional_data = AdditionalData())

Constructor.

Member Function Documentation

◆ reserve()

template<int dim, typename Number , int n_rows_1d = -1>
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().

◆ insert()

template<int dim, typename Number , int n_rows_1d = -1>
template<typename T >
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.

◆ finalize()

template<int dim, typename Number , int n_rows_1d = -1>
void TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::finalize ( )

Finalize setup. This function computes, e.g., the eigenvalues and the eigenvectors.

◆ apply_inverse()

template<int dim, typename Number , int n_rows_1d = -1>
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.

◆ memory_consumption()

template<int dim, typename Number , int n_rows_1d = -1>
std::size_t TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::memory_consumption ( ) const

Return the memory consumption of this class in bytes.

◆ storage_size()

template<int dim, typename Number , int n_rows_1d = -1>
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.

Member Data Documentation

◆ compress_matrices

template<int dim, typename Number , int n_rows_1d = -1>
const bool TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::compress_matrices
private

Try to compress matrices.

Definition at line 423 of file tensor_product_matrix.h.

◆ precompute_inverse_diagonal

template<int dim, typename Number , int n_rows_1d = -1>
const bool TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::precompute_inverse_diagonal
private

Precompute inverse diagonal.

Definition at line 428 of file tensor_product_matrix.h.

◆ mass_and_derivative_matrices

template<int dim, typename Number , int n_rows_1d = -1>
std::vector<MatrixPairType> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::mass_and_derivative_matrices
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.

◆ cache

template<int dim, typename Number , int n_rows_1d = -1>
std::map< MatrixPairTypeWithMask, unsigned int, internal::TensorProductMatrixSymmetricSum::MatrixPairComparator<Number> > TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::cache
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.

◆ indices

template<int dim, typename Number , int n_rows_1d = -1>
std::vector<unsigned int> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::indices
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.

◆ mass_matrices

template<int dim, typename Number , int n_rows_1d = -1>
AlignedVector<Number> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::mass_matrices
private

Vector of 1d mass matrices.

Definition at line 458 of file tensor_product_matrix.h.

◆ derivative_matrices

template<int dim, typename Number , int n_rows_1d = -1>
AlignedVector<Number> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::derivative_matrices
private

Vector of 1d derivative matrices.

Definition at line 463 of file tensor_product_matrix.h.

◆ eigenvectors

template<int dim, typename Number , int n_rows_1d = -1>
AlignedVector<Number> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::eigenvectors
private

Vector of eigenvectors.

Definition at line 468 of file tensor_product_matrix.h.

◆ eigenvalues

template<int dim, typename Number , int n_rows_1d = -1>
AlignedVector<Number> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::eigenvalues
private

Vector of eigenvalues.

Definition at line 473 of file tensor_product_matrix.h.

◆ inverted_eigenvalues

template<int dim, typename Number , int n_rows_1d = -1>
AlignedVector<Number> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::inverted_eigenvalues
private

Vector of inverted eigenvalues.

Definition at line 478 of file tensor_product_matrix.h.

◆ vector_ptr

template<int dim, typename Number , int n_rows_1d = -1>
std::vector<unsigned int> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::vector_ptr
private

Pointer into mass_matrices, derivative_matrices, and eigenvalues.

Definition at line 483 of file tensor_product_matrix.h.

◆ matrix_ptr

template<int dim, typename Number , int n_rows_1d = -1>
std::vector<unsigned int> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::matrix_ptr
private

Pointer into mass_matrices, derivative_matrices, and eigenvalues.

Definition at line 488 of file tensor_product_matrix.h.

◆ vector_n_rows_1d

template<int dim, typename Number , int n_rows_1d = -1>
std::vector<unsigned int> TensorProductMatrixSymmetricSumCollection< dim, Number, n_rows_1d >::vector_n_rows_1d
private

Number of rows in 1 of each cell.

Definition at line 493 of file tensor_product_matrix.h.


The documentation for this class was generated from the following file: