DUNE PDELab (git)

Dune::PDELab::ISTL::BCRSMatrixBackend< EntriesPerRow > Struct Template Reference

Backend using (possibly nested) ISTL BCRSMatrices. More...

#include <dune/pdelab/backend/istl/bcrsmatrixbackend.hh>

Public Types

typedef std::size_t size_type
 The size type of the BCRSMatrix.
 
typedef PatternStatistics< size_typeStatistics
 The type of the object holding the statistics generated during pattern construction.
 
template<typename Matrix , typename GFSV , typename GFSU >
using Pattern = typename build_bcrs_pattern_type< typename Matrix::Container, GFSV, GFSU, typename GFSV::Ordering::ContainerAllocationTag >::type
 The type of the pattern object passed to the GridOperator for pattern construction.
 

Public Member Functions

template<typename GridOperator , typename Matrix >
std::vector< StatisticsbuildPattern (const GridOperator &grid_operator, Matrix &matrix) const
 Builds the matrix pattern associated with grid_operator and initializes matrix with it. More...
 
 BCRSMatrixBackend (const EntriesPerRow &entries_per_row)
 Constructs a BCRSMatrixBackend. More...
 

Detailed Description

template<typename EntriesPerRow = std::size_t>
struct Dune::PDELab::ISTL::BCRSMatrixBackend< EntriesPerRow >

Backend using (possibly nested) ISTL BCRSMatrices.

BCRSMatrixBackend is a matrix backend descriptor for ISTL matrices. It expects that both the ansatz and the test function space use ISTL vectors and automatically deduces the correct matrix type from those two vector backends.

The backend uses an accelerated pattern construction scheme, which requires the average number of non-zero entries per matrix row as a priori information. In constrast to the older construction scheme, the improved version never requires more memory than the matrix does after pattern construction and runs a lot faster, as long as it is provided with a reasonable estimate for the number of non-zero entries per row.

Examples
recipe-geometry-grid.cc, recipe-linear-system-assembly.cc, recipe-linear-system-solution-istl.cc, recipe-linear-system-solution-pdelab.cc, and recipe-operator-splitting.cc.

Constructor & Destructor Documentation

◆ BCRSMatrixBackend()

template<typename EntriesPerRow = std::size_t>
Dune::PDELab::ISTL::BCRSMatrixBackend< EntriesPerRow >::BCRSMatrixBackend ( const EntriesPerRow &  entries_per_row)
inline

Constructs a BCRSMatrixBackend.

TODO: Document and flesh out the way this should work for nested matrices (use a nested array as entries_per_row).

Parameters
entries_per_rowThe average number of nonzero entries per row in matrices created with this backend.

Member Function Documentation

◆ buildPattern()

template<typename EntriesPerRow = std::size_t>
template<typename GridOperator , typename Matrix >
std::vector<Statistics> Dune::PDELab::ISTL::BCRSMatrixBackend< EntriesPerRow >::buildPattern ( const GridOperator grid_operator,
Matrix matrix 
) const
inline

Builds the matrix pattern associated with grid_operator and initializes matrix with it.

Returns
a vector with statistics object for all leaf BCRSMatrices in row-major order.

The documentation for this struct was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 3, 22:32, 2024)