Dune Core Modules (unstable)

A sparse block matrix with compressed row storage. More...

#include <dune/istl/bcrsmatrix.hh>

Classes

class  CreateIterator
 Iterator class for sequential creation of blocks More...
 
class  RealRowIterator
 Iterator access to matrix rows More...
 

Public Types

enum  BuildStage {
  notbuilt =0 , notAllocated =0 , building =1 , rowSizesBuilt =2 ,
  built =3
}
 
enum  BuildMode { row_wise , random , implicit , unknown }
 we support two modes More...
 
using field_type = typename Imp::BlockTraits< B >::field_type
 export the type representing the field
 
typedef B block_type
 export the type representing the components
 
typedef A allocator_type
 export the allocator type
 
typedef A::size_type size_type
 The type for the index access and the size.
 
typedef Imp::CompressedBlockVectorWindow< B, size_typerow_type
 implement row_type with compressed vector
 
typedef ::Dune::CompressionStatistics< size_typeCompressionStatistics
 The type for the statistics object returned by compress()
 
typedef RealRowIterator< row_typeiterator
 The iterator over the (mutable matrix rows.
 
typedef Iterator RowIterator
 rename the iterators for easier access
 
typedef row_type::Iterator ColIterator
 Iterator for the entries of each row.
 
typedef RealRowIterator< const row_typeconst_iterator
 The const iterator over the matrix rows.
 
typedef ConstIterator ConstRowIterator
 rename the const row iterator for easier access
 
typedef row_type::ConstIterator ConstColIterator
 Const iterator to the entries of a row.
 

Public Member Functions

row_typeoperator[] (size_type i)
 random access to the rows
 
const row_typeoperator[] (size_type i) const
 same for read only access
 
Iterator begin ()
 Get iterator to first row.
 
Iterator end ()
 Get iterator to one beyond last row.
 
Iterator beforeEnd ()
 
Iterator beforeBegin ()
 
ConstIterator begin () const
 Get const iterator to first row.
 
ConstIterator end () const
 Get const iterator to one beyond last row.
 
ConstIterator beforeEnd () const
 
ConstIterator beforeBegin () const
 
 BCRSMatrix ()
 an empty matrix
 
 BCRSMatrix (size_type _n, size_type _m, size_type _nnz, BuildMode bm)
 matrix with known number of nonzeroes
 
 BCRSMatrix (size_type _n, size_type _m, BuildMode bm)
 matrix with unknown number of nonzeroes
 
 BCRSMatrix (size_type _n, size_type _m, size_type _avg, double compressionBufferSize, BuildMode bm)
 construct matrix with a known average number of entries per row More...
 
 BCRSMatrix (const BCRSMatrix &Mat)
 copy constructor More...
 
 ~BCRSMatrix ()
 destructor
 
void setBuildMode (BuildMode bm)
 Sets the build mode of the matrix. More...
 
void setSize (size_type rows, size_type columns, size_type nnz=0)
 Set the size of the matrix. More...
 
void setImplicitBuildModeParameters (size_type _avg, double compressionBufferSize)
 Set parameters needed for creation in implicit build mode. More...
 
BCRSMatrixoperator= (const BCRSMatrix &Mat)
 assignment More...
 
BCRSMatrixoperator= (const field_type &k)
 Assignment from a scalar.
 
CreateIterator createbegin ()
 get initial create iterator
 
CreateIterator createend ()
 get create iterator pointing to one after the last block
 
void setrowsize (size_type i, size_type s)
 Set number of indices in row i to s. More...
 
size_type getrowsize (size_type i) const
 get current number of indices in row i
 
void incrementrowsize (size_type i, size_type s=1)
 increment size of row i by s (1 by default)
 
void endrowsizes ()
 indicate that size of all rows is defined
 
void addindex (size_type row, size_type col)
 add index (row,col) to the matrix More...
 
template<typename It >
void setIndicesNoSort (size_type row, It begin, It end)
 Set all column indices for row from the given iterator range. More...
 
template<typename It >
void setIndices (size_type row, It begin, It end)
 Set all column indices for row from the given iterator range. More...
 
void endindices ()
 indicate that all indices are defined, check consistency
 
B & entry (size_type row, size_type col)
 Returns reference to entry (row,col) of the matrix. More...
 
CompressionStatistics compress ()
 Finishes the buildstage in implicit mode. More...
 
BCRSMatrixoperator*= (const field_type &k)
 vector space multiplication with scalar
 
BCRSMatrixoperator/= (const field_type &k)
 vector space division by scalar
 
BCRSMatrixoperator+= (const BCRSMatrix &b)
 Add the entries of another matrix to this one. More...
 
BCRSMatrixoperator-= (const BCRSMatrix &b)
 Subtract the entries of another matrix from this one. More...
 
BCRSMatrixaxpy (field_type alpha, const BCRSMatrix &b)
 Add the scaled entries of another matrix to this one. More...
 
template<class X , class Y >
void mv (const X &x, Y &y) const
 y = A x
 
template<class X , class Y >
void umv (const X &x, Y &y) const
 y += A x
 
template<class X , class Y >
void mmv (const X &x, Y &y) const
 y -= A x
 
template<class X , class Y , class F >
void usmv (F &&alpha, const X &x, Y &y) const
 y += alpha A x
 
template<class X , class Y >
void mtv (const X &x, Y &y) const
 y = A^T x
 
template<class X , class Y >
void umtv (const X &x, Y &y) const
 y += A^T x
 
template<class X , class Y >
void mmtv (const X &x, Y &y) const
 y -= A^T x
 
template<class X , class Y >
void usmtv (const field_type &alpha, const X &x, Y &y) const
 y += alpha A^T x
 
template<class X , class Y >
void umhv (const X &x, Y &y) const
 y += A^H x
 
template<class X , class Y >
void mmhv (const X &x, Y &y) const
 y -= A^H x
 
template<class X , class Y >
void usmhv (const field_type &alpha, const X &x, Y &y) const
 y += alpha A^H x
 
FieldTraits< field_type >::real_type frobenius_norm2 () const
 square of frobenius norm, need for block recursion
 
FieldTraits< field_type >::real_type frobenius_norm () const
 frobenius norm: sqrt(sum over squared values of entries)
 
template<typename ft = field_type, typename std::enable_if<!HasNaN< ft >::value, int >::type = 0>
FieldTraits< ft >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
template<typename ft = field_type, typename std::enable_if<!HasNaN< ft >::value, int >::type = 0>
FieldTraits< ft >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
template<typename ft = field_type, typename std::enable_if< HasNaN< ft >::value, int >::type = 0>
FieldTraits< ft >::real_type infinity_norm () const
 infinity norm (row sum norm, how to generalize for blocks?)
 
template<typename ft = field_type, typename std::enable_if< HasNaN< ft >::value, int >::type = 0>
FieldTraits< ft >::real_type infinity_norm_real () const
 simplified infinity norm (uses Manhattan norm for complex values)
 
size_type N () const
 number of rows (counted in blocks)
 
size_type M () const
 number of columns (counted in blocks)
 
size_type nonzeroes () const
 number of blocks that are stored (the number of blocks that possibly are nonzero)
 
BuildStage buildStage () const
 The current build stage of the matrix.
 
BuildMode buildMode () const
 The currently selected build mode of the matrix.
 
bool exists (size_type i, size_type j) const
 return true if (i,j) is in pattern
 

Protected Member Functions

void setColumnPointers (ConstRowIterator row)
 Copy row sizes from iterator range starting at row and set column index pointers for all rows. More...
 
void setDataPointers ()
 Set data pointers for all rows. More...
 
void copyWindowStructure (const BCRSMatrix &Mat)
 Copy the window structure from another matrix.
 
void deallocate (bool deallocateRows=true)
 deallocate memory of the matrix. More...
 
void allocate (size_type rows, size_type columns, size_type allocationSize, bool allocateRows, bool allocate_data)
 Allocate memory for the matrix structure. More...
 
void implicit_allocate (size_type _n, size_type _m)
 organizes allocation implicit mode calculates correct array size to be allocated and sets the the window pointers to their correct positions for insertion. internally uses allocate() for the real allocation.
 

Detailed Description

template<class B, class A = std::allocator<B>>
class Dune::BCRSMatrix< B, A >

A sparse block matrix with compressed row storage.

Implements a block compressed row storage scheme. The block type B can be any type implementing the matrix interface.

Different ways to build up a compressed row storage matrix are supported:

  1. Row-wise scheme
  2. Random scheme
  3. implicit scheme

Error checking: no error checking is provided normally. Setting the compile time switch DUNE_ISTL_WITH_CHECKING enables error checking.

Details:

  1. Row-wise scheme

Rows are built up in sequential order. Size of the row and the column indices are defined. A row can be used as soon as it is initialized. With respect to memory there are two variants of this scheme: (a) number of non-zeroes known in advance (application finite difference schemes), (b) number of non-zeroes not known in advance (application: Sparse LU, ILU(n)).

...
// third parameter is an optional upper bound for the number
// of nonzeros. If given the matrix will use one array for all values
// as opposed to one for each row.
for(Iter row=B.createbegin(); row!=B.createend(); ++row){
// Add nonzeros for left neighbour, diagonal and right neighbour
if(row.index()>0)
row.insert(row.index()-1);
row.insert(row.index());
if(row.index()<B.N()-1)
row.insert(row.index()+1);
}
// Now the sparsity pattern is fully set up and we can add values
B[0][0]=2;
...
Implementation of the BCRSMatrix class.
Iterator class for sequential creation of blocks
Definition: bcrsmatrix.hh:954
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:2007
A dense n x m matrix.
Definition: fmatrix.hh:117
Implements a matrix constructed from a given type representing a field and compile-time given number ...
  1. Random scheme

For general finite element implementations the number of rows n is known, the number of non-zeroes might also be known (e.g. #edges + #nodes for P2) but the size of a row and the indices of a row can not be defined in sequential order.

...
// initially set row size for each row
B.setrowsize(0,1);
B.setrowsize(3,4);
B.setrowsize(2,1);
B.setrowsize(1,1);
// increase row size for row 2
B.incrementrowsize(2);
// finalize row setup phase
B.endrowsizes();
// add column entries to rows
B.addindex(0,0);
B.addindex(3,1);
B.addindex(2,2);
B.addindex(1,1);
B.addindex(2,0);
B.addindex(3,2);
B.addindex(3,0);
B.addindex(3,3);
// finalize column setup phase
B.endindices();
// set entries using the random access operator
B[0][0] = 1;
B[1][1] = 2;
B[2][0] = 3;
B[2][2] = 4;
B[3][1] = 5;
B[3][2] = 6;
B[3][0] = 7;
B[3][3] = 8;
  1. implicit scheme

With the 'random scheme` described above, the sparsity pattern has to be determined and stored before the matrix is assembled. This requires a dedicated iteration over the grid elements, which can be costly in terms of time. Also, additional memory is needed to store the pattern before it can be given to the 'random' build mode.

On the other hand, often one has good a priori knowledge about the number of entries a row contains on average. The implicit mode tries to make use of that knowledge, and allows the setup of matrix pattern and numerical values together.

Constructing and filling a BCRSMatrix with the 'implicit' mode is performed in two steps: In a setup phase, matrix entries with numerical values can be inserted into the matrix. Then, a compression algorithm is called which defragments and optimizes the memory layout. After this compression step, the matrix is ready to be used, and no further nonzero entries can be added.

To use this mode, either construct a matrix object via

or default-construct the matrix and then call

The parameter _avg specifies the expected number of (block) entries per matrix row.

When the BCRSMatrix object is first constructed with the 'implicit' build mode, two areas for matrix entry storage are allocated:

1) A large continuous chunk of memory that can hold the expected number of entries. In addition, this chunk contains an extra part of memory called the 'compression buffer', located before the memory for the matrix itself. The size of this buffer will be _avg * _n * compressionBufferSize.

2) An associative container indexed by \(i,j\)-pairs, which will hold surplus matrix entries during the setup phase (the 'overflow area'). Its content is merged into the main memory during compression.

You can then start filling your matrix by calling entry(size_type row, size_type col), which returns the corresponding matrix entry, creating it on the fly if it does not exist yet. The matrix pattern is hence created implicitly by simply accessing nonzero entries during the initial matrix assembly. Note that new entries are not zero-initialized, though, and hence the first operation on each entry has to be an assignment.

If a row contains more non-zero entries than what was specified in the _avg parameter, the surplus entries are stored in the 'overflow area' during the initial setup phase. After all indices are added, call compress() to trigger the compression step that optimizes the matrix and integrates any entries from the overflow area into the standard BCRS storage. This compression step builds up the final memory layout row by row. It will fail with an exception if the compression buffer is not large enough, which would lead to compressed rows overwriting uncompressed ones. More precisely, if \(\textrm{nnz}_j\) denotes the number of non-zeros in the \(j\)-th row, then the compression algorithm will succeed if the maximal number of non-zeros in the \(i\)-th row is

\[ M_i = \textrm{avg} + A + \sum_{j<i} (\textrm{avg} - \textrm{nnz}_j) \]

for all \(i\), where \( A = \textrm{avg}(n \cdot \textrm{compressionBufferSize}) \) is the total size of the compression buffer determined by the parameters explained above.

The data of the matrix is now located at the beginning of the allocated area, and covers what used to be the compression buffer. In exchange, there is now unused space at the end of the large allocated piece of memory. This will go unused and cannot be freed during the lifetime of the matrix, but it has no negative impact on run-time performance. No matrix entries may be added after the compression step.

The compress() method returns a value of type Dune::CompressionStatistics, which you can inspect to tune the construction parameters _avg and compressionBufferSize.

Use of copy constructor, assignment operator and matrix vector arithmetic is not supported until the matrix is fully built.

The following sample code constructs a \( 10 \times 10\) matrix, with an expected number of two entries per matrix row. The compression buffer size is set to 0.4. Hence the main chunk of allocated memory will be able to hold 10 * 2 entries in the matrix rows, and 10 * 2 * 0.4 entries in the compression buffer. In total that's 28 entries.

M m(10, 10, 2, 0.4, M::implicit);
// Fill in some arbitrary entries; the order is irrelevant.
// Even operations on these would be possible, you get a reference to the entry!
m.entry(0,0) = 0.;
m.entry(8,0) = 0.;
m.entry(1,8) = 0.; m.entry(1,0) = 0.; m.entry(1,5) = 0.;
m.entry(2,0) = 0.;
m.entry(3,5) = 0.; m.entry(3,0) = 0.; m.entry(3,8) = 0.;
m.entry(4,0) = 0.;
m.entry(9,0) = 0.; m.entry(9,5) = 0.; m.entry(9,8) = 0.;
m.entry(5,0) = 0.; m.entry(5,5) = 0.; m.entry(5,8) = 0.;
m.entry(6,0) = 0.;
m.entry(7,0) = 0.; m.entry(7,5) = 0.; m.entry(7,8) = 0.;

Internally the index array now looks like this:

// xxxxxxxx0x800x500x050x050x05
// ........|.|.|.|.|.|.|.|.|.|.

The second row denotes the beginnings of the matrix rows. The eight 'x' on the left are the compression buffer. The overflow area contains the entries (1,5,0.0), (3,8,0.0), (5,8,0.0), (7,8,0.0), and (9,8,0.0). These are entries of rows 1, 3, 5, 7, and 9, which have three entries each, even though only two were anticipated.

//finish building by compressing the array
Statistics about compression achieved in implicit mode.
Definition: bcrsmatrix.hh:88

Internally the index array now looks like this:

// 00580058005800580058xxxxxxxx
// ||..||..||..||..||..........

The compression buffer on the left is gone now, and the matrix has a real CRS layout. The 'x' on the right will be unused for the rest of the matrix' lifetime.

Member Enumeration Documentation

◆ BuildMode

template<class B , class A = std::allocator<B>>
enum Dune::BCRSMatrix::BuildMode

we support two modes

Enumerator
row_wise 

Build in a row-wise manner.

Rows are built up in sequential order. Size of the row and the column indices are defined. A row can be used as soon as it is initialized. With respect to memory there are two variants of this scheme: (a) number of non-zeroes known in advance (application finite difference schemes), (b) number of non-zeroes not known in advance (application: Sparse LU, ILU(n)).

random 

Build entries randomly.

For general finite element implementations the number of rows n is known, the number of non-zeroes might also be known (e.g. #edges + #nodes for P2) but the size of a row and the indices of a row cannot be defined in sequential order.

implicit 

Build entries randomly with an educated guess for the number of entries per row.

Allows random order generation as in random mode, but row sizes do not need to be given first. Instead an average number of non-zeroes per row is passed to the constructor. Matrix setup is finished with compress(), full data access during build stage is possible.

unknown 

Build mode not set!

◆ BuildStage

template<class B , class A = std::allocator<B>>
enum Dune::BCRSMatrix::BuildStage
Enumerator
notbuilt 

Matrix is not built at all, no memory has been allocated, build mode and size can still be set.

notAllocated 

Matrix is not built at all, no memory has been allocated, build mode and size can still be set.

building 

Matrix is currently being built, some memory has been allocated, build mode and size are fixed.

rowSizesBuilt 

The row sizes of the matrix are known.

 Only used in random mode.
built 

The matrix structure is fully built.

Constructor & Destructor Documentation

◆ BCRSMatrix() [1/2]

template<class B , class A = std::allocator<B>>
Dune::BCRSMatrix< B, A >::BCRSMatrix ( size_type  _n,
size_type  _m,
size_type  _avg,
double  compressionBufferSize,
BuildMode  bm 
)
inline

construct matrix with a known average number of entries per row

Constructs a matrix in implicit buildmode.

Parameters
_nnumber of rows of the matrix
_mnumber of columns of the matrix
_avgexpected average number of entries per row
compressionBufferSizefraction of _n*_avg which is expected to be needed for elements that exceed _avg entries per row.
bmMatrix build mode

References DUNE_THROW, Dune::BCRSMatrix< B, A >::implicit, and Dune::BCRSMatrix< B, A >::implicit_allocate().

◆ BCRSMatrix() [2/2]

template<class B , class A = std::allocator<B>>
Dune::BCRSMatrix< B, A >::BCRSMatrix ( const BCRSMatrix< B, A > &  Mat)
inline

Member Function Documentation

◆ addindex()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::addindex ( size_type  row,
size_type  col 
)
inline

add index (row,col) to the matrix

This method can only be used when building the BCRSMatrix in random mode.

addindex adds a new column entry to the row. If this column entry already exists, nothing is done.

Don't call addindex after the setup phase is finished (after endindices is called).

Referenced by Dune::BDMatrix< B, A >::setSize(), and Dune::BTDMatrix< B, A >::setSize().

◆ allocate()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::allocate ( size_type  rows,
size_type  columns,
size_type  allocationSize,
bool  allocateRows,
bool  allocate_data 
)
inlineprotected

Allocate memory for the matrix structure.

Sets the number of rows and columns of the matrix and allocates the memory needed for the storage of the matrix entries.

Warning
After calling this methods on an already allocated (and probably setup matrix) results in all the structure and data being lost. Please call deallocate() before calling allocate in this case.
Parameters
rowsThe number of rows the matrix should contain.
columnsthe number of columns the matrix should contain.
allocationSizeThe number of nonzero entries the matrix should hold (if omitted defaults to 0).
allocateRowsWhether we have to allocate the row pointers, too. (Defaults to true)
allocate_dataWhether data array needs to be allocated

References Dune::BCRSMatrix< B, A >::building, DUNE_THROW, and Dune::size().

Referenced by Dune::BCRSMatrix< B, A >::BCRSMatrix(), Dune::BCRSMatrix< B, A >::endrowsizes(), Dune::BCRSMatrix< B, A >::implicit_allocate(), Dune::BCRSMatrix< B, A >::operator=(), and Dune::BCRSMatrix< B, A >::setSize().

◆ axpy()

template<class B , class A = std::allocator<B>>
BCRSMatrix & Dune::BCRSMatrix< B, A >::axpy ( field_type  alpha,
const BCRSMatrix< B, A > &  b 
)
inline

Add the scaled entries of another matrix to this one.

Matrix axpy operation: *this += alpha * b

Parameters
alphaScaling factor.
bThe matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.

References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, DUNE_THROW, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

◆ beforeBegin() [1/2]

template<class B , class A = std::allocator<B>>
Iterator Dune::BCRSMatrix< B, A >::beforeBegin ( )
inline
Returns
an iterator that is positioned before the first row of the matrix.

◆ beforeBegin() [2/2]

template<class B , class A = std::allocator<B>>
ConstIterator Dune::BCRSMatrix< B, A >::beforeBegin ( ) const
inline
Returns
an iterator that is positioned before the first row of the matrix.

◆ beforeEnd() [1/2]

template<class B , class A = std::allocator<B>>
Iterator Dune::BCRSMatrix< B, A >::beforeEnd ( )
inline
Returns
an iterator that is positioned before the end iterator of the rows, i.e. at the last row.

◆ beforeEnd() [2/2]

template<class B , class A = std::allocator<B>>
ConstIterator Dune::BCRSMatrix< B, A >::beforeEnd ( ) const
inline
Returns
an iterator that is positioned before the end iterator of the rows. i.e. at the last row.

◆ compress()

template<class B , class A = std::allocator<B>>
CompressionStatistics Dune::BCRSMatrix< B, A >::compress ( )
inline

Finishes the buildstage in implicit mode.

Performs compression of index and data arrays with linear complexity in the number of nonzeroes.

After calling this method, the matrix is in the built state and no more entries can be added.

Returns
An object with some statistics about the compression for future optimization.

References Dune::CompressionStatistics< size_type >::avg, Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::building, Dune::BCRSMatrix< B, A >::built, DUNE_THROW, Dune::BCRSMatrix< B, A >::implicit, Dune::CompressionStatistics< size_type >::maximum, Dune::CompressionStatistics< size_type >::mem_ratio, Dune::BCRSMatrix< B, A >::notAllocated, Dune::CompressionStatistics< size_type >::overflow_total, and Dune::size().

◆ deallocate()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::deallocate ( bool  deallocateRows = true)
inlineprotected

deallocate memory of the matrix.

Parameters
deallocateRowsWhether we have to deallocate the row pointers, too. If false they will not be touched. (Defaults to true).

References Dune::BCRSMatrix< B, A >::notAllocated.

Referenced by Dune::BCRSMatrix< B, A >::operator=(), Dune::BCRSMatrix< B, A >::setSize(), and Dune::BCRSMatrix< B, A >::~BCRSMatrix().

◆ entry()

template<class B , class A = std::allocator<B>>
B & Dune::BCRSMatrix< B, A >::entry ( size_type  row,
size_type  col 
)
inline

Returns reference to entry (row,col) of the matrix.

This method can only be used when the matrix is in implicit building mode.

A reference to entry (row, col) of the matrix is returned. If entry (row, col) is accessed for the first time, it is created on the fly.

This method can only be used while building the matrix, after compression operator[] gives a much better performance.

◆ operator+=()

template<class B , class A = std::allocator<B>>
BCRSMatrix & Dune::BCRSMatrix< B, A >::operator+= ( const BCRSMatrix< B, A > &  b)
inline

Add the entries of another matrix to this one.

Parameters
bThe matrix to add to this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.

References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, DUNE_THROW, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

◆ operator-=()

template<class B , class A = std::allocator<B>>
BCRSMatrix & Dune::BCRSMatrix< B, A >::operator-= ( const BCRSMatrix< B, A > &  b)
inline

Subtract the entries of another matrix from this one.

Parameters
bThe matrix to subtract from this one. Its sparsity pattern has to be subset of the sparsity pattern of this matrix.

References Dune::BCRSMatrix< B, A >::begin(), Dune::BCRSMatrix< B, A >::built, DUNE_THROW, Dune::BCRSMatrix< B, A >::end(), Dune::BCRSMatrix< B, A >::M(), and Dune::BCRSMatrix< B, A >::N().

◆ operator=()

template<class B , class A = std::allocator<B>>
BCRSMatrix & Dune::BCRSMatrix< B, A >::operator= ( const BCRSMatrix< B, A > &  Mat)
inline

◆ setBuildMode()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setBuildMode ( BuildMode  bm)
inline

◆ setColumnPointers()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setColumnPointers ( ConstRowIterator  row)
inlineprotected

Copy row sizes from iterator range starting at row and set column index pointers for all rows.

This method does not modify the data pointers, as those are set only after building the pattern (to allow for a delayed allocation).

Referenced by Dune::BCRSMatrix< B, A >::endrowsizes().

◆ setDataPointers()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setDataPointers ( )
inlineprotected

Set data pointers for all rows.

This method assumes that column pointers and row sizes have been correctly set up by a prior call to setColumnPointers().

Referenced by Dune::BCRSMatrix< B, A >::endindices(), and Dune::BCRSMatrix< B, A >::CreateIterator::operator++().

◆ setImplicitBuildModeParameters()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setImplicitBuildModeParameters ( size_type  _avg,
double  compressionBufferSize 
)
inline

Set parameters needed for creation in implicit build mode.

Use this method before setSize() to define storage behaviour of a matrix in implicit build mode

Parameters
_avgexpected average number of entries per row
compressionBufferSizefraction of _n*_avg which is expected to be needed for elements that exceed _avg entries per row.

References DUNE_THROW, and Dune::BCRSMatrix< B, A >::notAllocated.

◆ setIndices()

template<class B , class A = std::allocator<B>>
template<typename It >
void Dune::BCRSMatrix< B, A >::setIndices ( size_type  row,
It  begin,
It  end 
)
inline

Set all column indices for row from the given iterator range.

The iterator range has to be of the same length as the previously set row size. The entries in the iterator range do not have to be in any particular order, but must not contain duplicate values. This method will insert the indices and sort them afterwards.

Calling this method overwrites any previously set column indices!

◆ setIndicesNoSort()

template<class B , class A = std::allocator<B>>
template<typename It >
void Dune::BCRSMatrix< B, A >::setIndicesNoSort ( size_type  row,
It  begin,
It  end 
)
inline

Set all column indices for row from the given iterator range.

The iterator range has to be of the same length as the previously set row size. The entries in the iterator range must be sorted and must not contain duplicate values. This method will insert the indices in the given order.

Calling this method overwrites any previously set column indices!

◆ setrowsize()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setrowsize ( size_type  i,
size_type  s 
)
inline

Set number of indices in row i to s.

The number s may actually be larger than the true number of nonzero entries in row i. In that case, the extra memory goes wasted. You will receive run-time warnings about this, sent to the Dune::dwarn stream.

References Dune::BCRSMatrix< B, A >::building, DUNE_THROW, and Dune::BCRSMatrix< B, A >::random.

Referenced by Dune::BDMatrix< B, A >::setSize().

◆ setSize()

template<class B , class A = std::allocator<B>>
void Dune::BCRSMatrix< B, A >::setSize ( size_type  rows,
size_type  columns,
size_type  nnz = 0 
)
inline

Set the size of the matrix.

Sets the number of rows and columns of the matrix and allocates the memory needed for the storage of the matrix entries.

Warning
After calling this methods on an already allocated (and probably setup matrix) results in all the structure and data being deleted. I.~e. one has to setup the matrix again.
Parameters
rowsThe number of rows the matrix should contain.
columnsthe number of columns the matrix should contain.
nnzThe number of nonzero entries the matrix should hold (if omitted defaults to 0). Must be omitted in implicit mode.

References Dune::BCRSMatrix< B, A >::allocate(), Dune::BCRSMatrix< B, A >::deallocate(), DUNE_THROW, Dune::BCRSMatrix< B, A >::implicit, and Dune::BCRSMatrix< B, A >::implicit_allocate().

Referenced by Dune::BDMatrix< B, A >::setSize(), and Dune::BTDMatrix< B, A >::setSize().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 20, 23:30, 2024)