KASKADE 7 development version
Classes | Public Types | Public Member Functions | Public Attributes | List of all members
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt > Class Template Reference

More...

#include <triplet.hh>

Detailed Description

template<class Scalar, class SparseIndexInt = int>
class Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >

Sparse Matrix in triplet format.

A simple sparse matrix format consisting of scalar-valued (row,col,value) entries. Use this as the least common denominator of sparse linear algebra interfaces.

Multiple entries with same row and col are semantically summed up.

Template Parameters
Scalar(double or float, field_type)
SparseIndexInt(int,long)

long has to be used as SparseIndexInt if the number of nonzeros of the sparse matrix exceeds INT_MAX 2,147,483,647 or if the umfpack direct solver is used and the internal workspace exceeds this number, what will happen much earlier.

Definition at line 56 of file triplet.hh.

Public Types

typedef Scalar value_type
 STL-compliant typedefs. More...
 
typedef std::vector< Scalar >::iterator iterator
 
typedef std::vector< Scalar >::const_iterator const_iterator
 
typedef std::vector< SparseIndexInt >::iterator index_iterator
 
typedef std::vector< SparseIndexInt >::const_iterator const_index_iterator
 

Public Member Functions

 MatrixAsTriplet ()
 Constructor. More...
 
 MatrixAsTriplet (size_t s)
 Constructor that allocates memory directly. More...
 
template<int n, int m, class Allocator >
 MatrixAsTriplet (Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Allocator > const &other)
 Constructor from a Dune::BCRSMatrix. More...
 
template<int n, int m, class Index >
 MatrixAsTriplet (NumaBCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Index > const &other)
 Constructor from a NumaBCRSMatrix. More...
 
 MatrixAsTriplet (std::vector< SparseIndexInt > &&row, std::vector< SparseIndexInt > &&col, std::vector< Scalar > &&value)
 Move constructor taking raw data. More...
 
 MatrixAsTriplet (SparseIndexInt nrows, SparseIndexInt ncols, Scalar const &value)
 Constructs a matrix of size (nrows,ncols), filled with value. If value==0.0, then it is an empty matrix. More...
 
 MatrixAsTriplet (SparseIndexInt nrows, SparseIndexInt ncols, SparseIndexInt shift, Scalar const &value)
 Constructs a diagonal matrix of size (nrows,ncols). The diagonal is shifted down shift rows (if shift is negative, it is shifted up) More...
 
 MatrixAsTriplet (SparseIndexInt nrows, SparseIndexInt ncols, SparseIndexInt shift, std::vector< Scalar > const &value)
 Constructs a diagonal matrix of size (nrows,ncols). The diagonal is shifted down shift rows (if shift is negative, it is shifted up) More...
 
 MatrixAsTriplet (MatrixAsTriplet const &other)=default
 Copy constructor. More...
 
 MatrixAsTriplet (MatrixAsTriplet &&other)=default
 Move constructor. More...
 
template<class OtherScalar , class OtherIndex >
 MatrixAsTriplet (MatrixAsTriplet< OtherScalar, OtherIndex > const &A)
 Copy constructor from a triplet matrix with different entry and/or index types. More...
 
void clear ()
 Deletes the data, leaving an empty matrix of size 0 x 0. More...
 
void setStartToZero ()
 Shift the row and column indices such that \( A_{00} \) is (structurally) nonzero. More...
 
SparseIndexInt nrows () const
 Returns number of rows (computes them, if not known) More...
 
SparseIndexInt N () const
 Returns number of rows (computes them, if not known) This is a Dune ISTL compatibility method. More...
 
SparseIndexInt ncols () const
 Returns number of cols (computes them, if not known) More...
 
SparseIndexInt M () const
 Returns number of columns (computes them, if not known) This is a Dune ISTL compatibility method. More...
 
void resize (size_t s)
 Resizes the memory. More...
 
void reserve (size_t s)
 Reserve memory. More...
 
size_t nnz () const
 
template<class X , class Y >
void axpy (Y &out, X const &in, Scalar alpha=1.0, int nvectors=1) const
 Matrix-vector multiplication: out += alpha * (*this) * in. More...
 
template<class X , class Y >
void atxpy (Y &y, X const &x, Scalar alpha=1.0) const
 Matrix-vector multiplication. More...
 
template<class X , class Y >
void usmtv (Scalar const alpha, X const &x, Y &y) const
 Matrix-vector multiplication. More...
 
template<class X , class Y >
void usmv (Scalar const alpha, X const &x, Y &y) const
 scaled matrix-vector multiplication More...
 
template<class X , class Y >
void ax (Y &out, X const &in) const
 Matrix-vector multiplication: out = (*this) * in. More...
 
template<class X , class Y >
void mv (X const &x, Y &y) const
 matrix-vector multiplication More...
 
void shiftIndices (SparseIndexInt row, SparseIndexInt col)
 Shifts matrix indices. Can be used together with += to concatenate sparese matrices. More...
 
void addEntry (SparseIndexInt row, SparseIndexInt col, Scalar const &value)
 Adds the given value to the specified matrix entry. The entry is created in case it has been structurally zero before. More...
 
MatrixAsTripletoperator+= (MatrixAsTriplet const &m)
 Matrix addition: (*this) += m (works also for matrices with non-matching sparsity pattern) More...
 
MatrixAsTripletoperator-= (MatrixAsTriplet const &m)
 
MatrixAsTripletoperator*= (Scalar const &scalar)
 Scaling. More...
 
MatrixAsTripletoperator/= (Scalar const &scalar)
 Scaling. More...
 
MatrixAsTripletoperator= (MatrixAsTriplet const &m)=default
 Assignment. More...
 
MatrixAsTripletoperator= (MatrixAsTriplet &&other)=default
 Move assignment. More...
 
void scaleRows (std::vector< Scalar >const &scaling)
 
MatrixAsTriplettranspose ()
 Transposition. More...
 
void deleteLowerTriangle ()
 Deletes all sub-diagonal entries. More...
 
std::ostream & print (std::ostream &out=std::cout, double threshold=0.0) const
 Output as Matlab source code. More...
 
void toColumns (std::vector< std::vector< Scalar > > &colsvector) const
 transform into a set of column vectors More...
 
void toRows (std::vector< std::vector< Scalar > > &rows) const
 
void toVector (std::vector< Scalar > &colsvector)
 
void addColumn (std::vector< Scalar > &colsvector, size_t position)
 add a column More...
 
template<class Mat >
void addToMatrix (Mat &mat)
 add to a dense matrix More...
 
void erase (SparseIndexInt i)
 
MatrixAsTriplet lumped ()
 
iterator begin ()
 get iterator on data field More...
 
const_iterator begin () const
 get const iterator on data field More...
 
iterator end ()
 get iterator on end of data field More...
 
const_iterator end () const
 get const iterator on end of data field More...
 
size_t size () const
 get number of entries in sparse matrix More...
 
void toFile (const char *filename, unsigned int precision=10) const
 write to text file More...
 
template<int bcrsN = 1>
std::unique_ptr< Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, bcrsN, bcrsN > > > toBCRS () const
 copy to Dune::BCRSMatrix More...
 
bool isSymmetric () const
 
std::pair< size_t, bool > findEntry (SparseIndexInt row, SparseIndexInt col) const
 

Public Attributes

std::vector< SparseIndexInt > ridx
 row indices More...
 
std::vector< SparseIndexInt > cidx
 column indices More...
 
std::vector< Scalar > data
 data More...
 

Member Typedef Documentation

◆ const_index_iterator

template<class Scalar , class SparseIndexInt = int>
typedef std::vector<SparseIndexInt>::const_iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::const_index_iterator

Definition at line 64 of file triplet.hh.

◆ const_iterator

template<class Scalar , class SparseIndexInt = int>
typedef std::vector<Scalar>::const_iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::const_iterator

Definition at line 62 of file triplet.hh.

◆ index_iterator

template<class Scalar , class SparseIndexInt = int>
typedef std::vector<SparseIndexInt>::iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::index_iterator

Definition at line 63 of file triplet.hh.

◆ iterator

template<class Scalar , class SparseIndexInt = int>
typedef std::vector<Scalar>::iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::iterator

Definition at line 61 of file triplet.hh.

◆ value_type

template<class Scalar , class SparseIndexInt = int>
typedef Scalar Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::value_type

STL-compliant typedefs.

Definition at line 60 of file triplet.hh.

Constructor & Destructor Documentation

◆ MatrixAsTriplet() [1/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( )
inline

Constructor.

Definition at line 67 of file triplet.hh.

◆ MatrixAsTriplet() [2/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( size_t  s)
inlineexplicit

Constructor that allocates memory directly.

Definition at line 70 of file triplet.hh.

◆ MatrixAsTriplet() [3/11]

template<class Scalar , class SparseIndexInt = int>
template<int n, int m, class Allocator >
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Allocator > const &  other)
inlineexplicit

Constructor from a Dune::BCRSMatrix.

Definition at line 74 of file triplet.hh.

◆ MatrixAsTriplet() [4/11]

template<class Scalar , class SparseIndexInt = int>
template<int n, int m, class Index >
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( NumaBCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Index > const &  other)
inlineexplicit

Constructor from a NumaBCRSMatrix.

Definition at line 81 of file triplet.hh.

◆ MatrixAsTriplet() [5/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( std::vector< SparseIndexInt > &&  row,
std::vector< SparseIndexInt > &&  col,
std::vector< Scalar > &&  value 
)
inline

Move constructor taking raw data.

The size of row, col, and value has to be identical. All entries in row and col need to be nonnegative. The supplied containers are empty after construction.

Definition at line 94 of file triplet.hh.

◆ MatrixAsTriplet() [6/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( SparseIndexInt  nrows,
SparseIndexInt  ncols,
Scalar const &  value 
)
inline

Constructs a matrix of size (nrows,ncols), filled with value. If value==0.0, then it is an empty matrix.

Definition at line 107 of file triplet.hh.

◆ MatrixAsTriplet() [7/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( SparseIndexInt  nrows,
SparseIndexInt  ncols,
SparseIndexInt  shift,
Scalar const &  value 
)
inline

Constructs a diagonal matrix of size (nrows,ncols). The diagonal is shifted down shift rows (if shift is negative, it is shifted up)

Definition at line 124 of file triplet.hh.

◆ MatrixAsTriplet() [8/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( SparseIndexInt  nrows,
SparseIndexInt  ncols,
SparseIndexInt  shift,
std::vector< Scalar > const &  value 
)
inline

Constructs a diagonal matrix of size (nrows,ncols). The diagonal is shifted down shift rows (if shift is negative, it is shifted up)

Definition at line 143 of file triplet.hh.

◆ MatrixAsTriplet() [9/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( MatrixAsTriplet< Scalar, SparseIndexInt > const &  other)
default

Copy constructor.

◆ MatrixAsTriplet() [10/11]

template<class Scalar , class SparseIndexInt = int>
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( MatrixAsTriplet< Scalar, SparseIndexInt > &&  other)
default

Move constructor.

◆ MatrixAsTriplet() [11/11]

template<class Scalar , class SparseIndexInt = int>
template<class OtherScalar , class OtherIndex >
Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet ( MatrixAsTriplet< OtherScalar, OtherIndex > const &  A)
inline

Copy constructor from a triplet matrix with different entry and/or index types.

Definition at line 168 of file triplet.hh.

Member Function Documentation

◆ addColumn()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addColumn ( std::vector< Scalar > &  colsvector,
size_t  position 
)
inline

add a column

Definition at line 536 of file triplet.hh.

◆ addEntry()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addEntry ( SparseIndexInt  row,
SparseIndexInt  col,
Scalar const &  value 
)
inline

Adds the given value to the specified matrix entry. The entry is created in case it has been structurally zero before.

Definition at line 349 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped().

◆ addToMatrix()

template<class Scalar , class SparseIndexInt = int>
template<class Mat >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addToMatrix ( Mat &  mat)
inline

add to a dense matrix

Definition at line 549 of file triplet.hh.

◆ atxpy()

template<class Scalar , class SparseIndexInt = int>
template<class X , class Y >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::atxpy ( Y &  y,
X const &  x,
Scalar  alpha = 1.0 
) const
inline

Matrix-vector multiplication.

\( x = x + \alpha A^{T}y \)

Definition at line 283 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::usmtv().

◆ ax()

template<class Scalar , class SparseIndexInt = int>
template<class X , class Y >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ax ( Y &  out,
X const &  in 
) const
inline

Matrix-vector multiplication: out = (*this) * in.

Definition at line 312 of file triplet.hh.

Referenced by Kaskade::Bridge::DirectInnerSolver< DirectSolver >::ax().

◆ axpy()

template<class Scalar , class SparseIndexInt = int>
template<class X , class Y >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::axpy ( Y &  out,
X const &  in,
Scalar  alpha = 1.0,
int  nvectors = 1 
) const
inline

◆ begin() [1/2]

template<class Scalar , class SparseIndexInt = int>
iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::begin ( )
inline

get iterator on data field

Definition at line 586 of file triplet.hh.

◆ begin() [2/2]

template<class Scalar , class SparseIndexInt = int>
const_iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::begin ( ) const
inline

get const iterator on data field

Definition at line 589 of file triplet.hh.

◆ clear()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::clear ( )
inline

Deletes the data, leaving an empty matrix of size 0 x 0.

Definition at line 182 of file triplet.hh.

◆ deleteLowerTriangle()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::deleteLowerTriangle ( )
inline

Deletes all sub-diagonal entries.

Definition at line 461 of file triplet.hh.

Referenced by Kaskade::HigherOrderRecovery< Grid, Space >::getErrorFunction().

◆ end() [1/2]

template<class Scalar , class SparseIndexInt = int>
iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::end ( )
inline

get iterator on end of data field

Definition at line 592 of file triplet.hh.

◆ end() [2/2]

template<class Scalar , class SparseIndexInt = int>
const_iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::end ( ) const
inline

get const iterator on end of data field

Definition at line 595 of file triplet.hh.

◆ erase()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::erase ( SparseIndexInt  i)
inline

Definition at line 555 of file triplet.hh.

◆ findEntry()

template<class Scalar , class SparseIndexInt = int>
std::pair< size_t, bool > Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::findEntry ( SparseIndexInt  row,
SparseIndexInt  col 
) const
inline
Returns
std::pair<size_t,bool> retVal. retVal.first: index. retVal.second: true, if entry was found

Definition at line 658 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::isSymmetric().

◆ isSymmetric()

template<class Scalar , class SparseIndexInt = int>
bool Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::isSymmetric ( ) const
inline

Definition at line 642 of file triplet.hh.

◆ lumped()

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped ( )
inline

Definition at line 562 of file triplet.hh.

◆ M()

template<class Scalar , class SparseIndexInt = int>
SparseIndexInt Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::M ( ) const
inline

Returns number of columns (computes them, if not known) This is a Dune ISTL compatibility method.

Definition at line 225 of file triplet.hh.

◆ mv()

template<class Scalar , class SparseIndexInt = int>
template<class X , class Y >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::mv ( X const &  x,
Y &  y 
) const
inline

matrix-vector multiplication

Definition at line 328 of file triplet.hh.

◆ N()

template<class Scalar , class SparseIndexInt = int>
SparseIndexInt Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::N ( ) const
inline

Returns number of rows (computes them, if not known) This is a Dune ISTL compatibility method.

Definition at line 212 of file triplet.hh.

◆ ncols()

template<class Scalar , class SparseIndexInt = int>
SparseIndexInt Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ncols ( ) const
inline

◆ nnz()

template<class Scalar , class SparseIndexInt = int>
size_t Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::nnz ( ) const
inline

◆ nrows()

template<class Scalar , class SparseIndexInt = int>
SparseIndexInt Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::nrows ( ) const
inline

Returns number of rows (computes them, if not known)

Definition at line 202 of file triplet.hh.

Referenced by Kaskade::AdditiveSchwarzPreconditioner< Op >::AdditiveSchwarzPreconditioner(), Kaskade::ARMSPreconditioner< Op >::ARMSPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::atxpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ax(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::axpy(), Kaskade::BlockLUFactorization< Factorization >::BlockLUFactorization(), Kaskade::BoomerAMG< Op >::BoomerAMG(), Kaskade::Euclid< Op >::Euclid(), Kaskade::getFactorization(), Kaskade::Bridge::NormalStepLinearization< Functional, stateId, adjointId >::getMatrixBlocks(), Kaskade::Bridge::TangentialStepLinearization< Functional, stateId, adjointId >::getMatrixBlocks(), Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::getTriplet(), Kaskade::ILUKPreconditioner< Op >::ILUKPreconditioner(), Kaskade::ILUTPreconditioner< Op >::ILUTPreconditioner(), Kaskade::AdditiveSchwarzPreconditioner< Op >::init(), Kaskade::SchurPreconditionerDetail::InvertLumpedMatrix< Scalar_ >::InvertLumpedMatrix(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped(), Kaskade::DerivativeChecker< Functional, checkD1, SparseInt >::matrixToVTK(), Kaskade::DerivativeCheck::matrixToVTK(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::N(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator+=(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator-=(), Kaskade::SimpleSparseLinearSystem::rows(), Kaskade::SLAPMatrix< Num, offset >::SLAPMatrix(), Kaskade::DirectLinearSolver< Scalar, SparseInt >::solve(), Kaskade::Bridge::DirectInnerSolver< DirectSolver >::solveAdjAndNormal(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toColumns(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toRows(), and Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toVector().

◆ operator*=()

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator*= ( Scalar const &  scalar)
inline

Scaling.

Definition at line 426 of file triplet.hh.

◆ operator+=()

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator+= ( MatrixAsTriplet< Scalar, SparseIndexInt > const &  m)
inline

Matrix addition: (*this) += m (works also for matrices with non-matching sparsity pattern)

Definition at line 357 of file triplet.hh.

◆ operator-=()

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator-= ( MatrixAsTriplet< Scalar, SparseIndexInt > const &  m)
inline

Definition at line 389 of file triplet.hh.

◆ operator/=()

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator/= ( Scalar const &  scalar)
inline

Scaling.

Definition at line 435 of file triplet.hh.

◆ operator=() [1/2]

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator= ( MatrixAsTriplet< Scalar, SparseIndexInt > &&  other)
default

Move assignment.

◆ operator=() [2/2]

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator= ( MatrixAsTriplet< Scalar, SparseIndexInt > const &  m)
default

Assignment.

◆ print()

template<class Scalar , class SparseIndexInt = int>
std::ostream & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::print ( std::ostream &  out = std::cout,
double  threshold = 0.0 
) const
inline

Output as Matlab source code.

Definition at line 490 of file triplet.hh.

Referenced by Kaskade::operator<<().

◆ reserve()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::reserve ( size_t  s)
inline

Reserve memory.

Definition at line 242 of file triplet.hh.

◆ resize()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::resize ( size_t  s)
inline

Resizes the memory.

Existing values are preserved if s >= nnz(), and the new entries are undefined. For s<nnz(), all the values are undefined.

Definition at line 234 of file triplet.hh.

Referenced by Kaskade::BlockLUFactorization< Factorization >::BlockLUFactorization().

◆ scaleRows()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::scaleRows ( std::vector< Scalar >const &  scaling)
inline

Definition at line 446 of file triplet.hh.

◆ setStartToZero()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::setStartToZero ( )
inline

Shift the row and column indices such that \( A_{00} \) is (structurally) nonzero.

Definition at line 193 of file triplet.hh.

◆ shiftIndices()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::shiftIndices ( SparseIndexInt  row,
SparseIndexInt  col 
)
inline

◆ size()

template<class Scalar , class SparseIndexInt = int>
size_t Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::size ( ) const
inline

get number of entries in sparse matrix

Definition at line 598 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS().

◆ toBCRS()

template<class Scalar , class SparseIndexInt = int>
template<int bcrsN = 1>
std::unique_ptr< Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, bcrsN, bcrsN > > > Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS ( ) const
inline

copy to Dune::BCRSMatrix

Definition at line 613 of file triplet.hh.

◆ toColumns()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toColumns ( std::vector< std::vector< Scalar > > &  colsvector) const
inline

transform into a set of column vectors

Definition at line 500 of file triplet.hh.

◆ toFile()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toFile ( const char *  filename,
unsigned int  precision = 10 
) const
inline

write to text file

Definition at line 601 of file triplet.hh.

◆ toRows()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toRows ( std::vector< std::vector< Scalar > > &  rows) const
inline

Definition at line 513 of file triplet.hh.

◆ toVector()

template<class Scalar , class SparseIndexInt = int>
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toVector ( std::vector< Scalar > &  colsvector)
inline

Definition at line 526 of file triplet.hh.

◆ transpose()

template<class Scalar , class SparseIndexInt = int>
MatrixAsTriplet & Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::transpose ( )
inline

◆ usmtv()

template<class Scalar , class SparseIndexInt = int>
template<class X , class Y >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::usmtv ( Scalar const  alpha,
X const &  x,
Y &  y 
) const
inline

◆ usmv()

template<class Scalar , class SparseIndexInt = int>
template<class X , class Y >
void Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::usmv ( Scalar const  alpha,
X const &  x,
Y &  y 
) const
inline

Member Data Documentation

◆ cidx

template<class Scalar , class SparseIndexInt = int>
std::vector<SparseIndexInt> Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::cidx

column indices

Definition at line 671 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addColumn(), Kaskade::AdditiveSchwarzPreconditioner< Op >::AdditiveSchwarzPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addToMatrix(), Kaskade::SchurPreconditionerDetail::JacobiIteration< Scalar_ >::apply(), Kaskade::AssemblyDetail::Fill< MatrixAsTriplet< Scalar, SparseIndexInt > >::apply(), Kaskade::ARMSPreconditioner< Op >::ARMSPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::atxpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ax(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::axpy(), Kaskade::BlockLUFactorization< Factorization >::BlockLUFactorization(), Kaskade::BoomerAMG< Op >::BoomerAMG(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::clear(), Kaskade::ThreadedMatrixDetail::CRSChunk< Entry, Index >::CRSChunk(), Kaskade::ThreadedMatrixDetail::CRSChunkPattern< Index >::CRSChunkPattern(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::deleteLowerTriangle(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::erase(), Kaskade::Euclid< Op >::Euclid(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::findEntry(), Kaskade::HigherOrderRecovery< Grid, Space >::getErrorFunction(), Kaskade::ThreadedMatrixDetail::getRowCount(), Kaskade::ILUKPreconditioner< Op >::ILUKPreconditioner(), Kaskade::ILUTPreconditioner< Op >::ILUTPreconditioner(), Kaskade::AdditiveSchwarzPreconditioner< Op >::init(), Kaskade::SchurPreconditionerDetail::InvertLumpedMatrix< Scalar_ >::InvertLumpedMatrix(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::isSymmetric(), Kaskade::SchurPreconditionerDetail::JacobiIteration< Scalar_ >::JacobiIteration(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet(), Kaskade::DerivativeChecker< Functional, checkD1, SparseInt >::matrixToVTK(), Kaskade::DerivativeCheck::matrixToVTK(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ncols(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator+=(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator-=(), Kaskade::PartialDirectPreconditioner< Op >::PartialDirectPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::print(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::reserve(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::resize(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::setStartToZero(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::shiftIndices(), Kaskade::SimpleSparseLinearSystem::SimpleSparseLinearSystem(), Kaskade::SLAPMatrix< Num, offset >::SLAPMatrix(), Kaskade::DirectLinearSolver< Scalar, SparseInt >::solve(), Kaskade::Bridge::DirectInnerSolver< DirectSolver >::solveAdjAndNormal(), Kaskade::Limex< Eq >::step(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toColumns(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toFile(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toRows(), Kaskade::SLAPMatrix< Num, offset >::toTriplet(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toVector(), and Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::transpose().

◆ data

template<class Scalar , class SparseIndexInt = int>
std::vector<Scalar> Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::data

data

Definition at line 673 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addColumn(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addToMatrix(), Kaskade::SchurPreconditionerDetail::JacobiIteration< Scalar_ >::apply(), Kaskade::AssemblyDetail::Fill< MatrixAsTriplet< Scalar, SparseIndexInt > >::apply(), Kaskade::ARMSPreconditioner< Op >::ARMSPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::atxpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ax(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::axpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::begin(), Kaskade::BlockLUFactorization< Factorization >::BlockLUFactorization(), Kaskade::BoomerAMG< Op >::BoomerAMG(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::clear(), Kaskade::ThreadedMatrixDetail::CRSChunk< Entry, Index >::CRSChunk(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::deleteLowerTriangle(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::end(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::erase(), Kaskade::Euclid< Op >::Euclid(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::findEntry(), Kaskade::HigherOrderRecovery< Grid, Space >::getErrorFunction(), Kaskade::ILUKPreconditioner< Op >::ILUKPreconditioner(), Kaskade::ILUTPreconditioner< Op >::ILUTPreconditioner(), Kaskade::AdditiveSchwarzPreconditioner< Op >::init(), Kaskade::SchurPreconditionerDetail::InvertLumpedMatrix< Scalar_ >::InvertLumpedMatrix(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::isSymmetric(), Kaskade::SchurPreconditionerDetail::JacobiIteration< Scalar_ >::JacobiIteration(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet(), Kaskade::DerivativeChecker< Functional, checkD1, SparseInt >::matrixToVTK(), Kaskade::DerivativeCheck::matrixToVTK(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator*=(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator+=(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator-=(), Kaskade::PartialDirectPreconditioner< Op >::PartialDirectPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::print(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::reserve(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::resize(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::scaleRows(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::shiftIndices(), Kaskade::SimpleSparseLinearSystem::SimpleSparseLinearSystem(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::size(), Kaskade::SLAPMatrix< Num, offset >::SLAPMatrix(), Kaskade::DirectLinearSolver< Scalar, SparseInt >::solve(), Kaskade::Bridge::DirectInnerSolver< DirectSolver >::solveAdjAndNormal(), Kaskade::Limex< Eq >::step(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toColumns(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toFile(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toRows(), Kaskade::SLAPMatrix< Num, offset >::toTriplet(), and Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toVector().

◆ ridx

template<class Scalar , class SparseIndexInt = int>
std::vector<SparseIndexInt> Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ridx

row indices

Definition at line 669 of file triplet.hh.

Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addColumn(), Kaskade::AdditiveSchwarzPreconditioner< Op >::AdditiveSchwarzPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::addToMatrix(), Kaskade::SchurPreconditionerDetail::JacobiIteration< Scalar_ >::apply(), Kaskade::AssemblyDetail::Fill< MatrixAsTriplet< Scalar, SparseIndexInt > >::apply(), Kaskade::ARMSPreconditioner< Op >::ARMSPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::atxpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ax(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::axpy(), Kaskade::BlockLUFactorization< Factorization >::BlockLUFactorization(), Kaskade::BoomerAMG< Op >::BoomerAMG(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::clear(), Kaskade::ThreadedMatrixDetail::CRSChunk< Entry, Index >::CRSChunk(), Kaskade::ThreadedMatrixDetail::CRSChunkPattern< Index >::CRSChunkPattern(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::deleteLowerTriangle(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::erase(), Kaskade::Euclid< Op >::Euclid(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::findEntry(), Kaskade::HigherOrderRecovery< Grid, Space >::getErrorFunction(), Kaskade::ThreadedMatrixDetail::getRowCount(), Kaskade::ILUKPreconditioner< Op >::ILUKPreconditioner(), Kaskade::ILUTPreconditioner< Op >::ILUTPreconditioner(), Kaskade::AdditiveSchwarzPreconditioner< Op >::init(), Kaskade::SchurPreconditionerDetail::InvertLumpedMatrix< Scalar_ >::InvertLumpedMatrix(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::isSymmetric(), Kaskade::SchurPreconditionerDetail::JacobiIteration< Scalar_ >::JacobiIteration(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::MatrixAsTriplet(), Kaskade::DerivativeChecker< Functional, checkD1, SparseInt >::matrixToVTK(), Kaskade::DerivativeCheck::matrixToVTK(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::nnz(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::nrows(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator+=(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator-=(), Kaskade::PartialDirectPreconditioner< Op >::PartialDirectPreconditioner(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::print(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::reserve(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::resize(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::scaleRows(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::setStartToZero(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::shiftIndices(), Kaskade::SimpleSparseLinearSystem::SimpleSparseLinearSystem(), Kaskade::SLAPMatrix< Num, offset >::SLAPMatrix(), Kaskade::DirectLinearSolver< Scalar, SparseInt >::solve(), Kaskade::Bridge::DirectInnerSolver< DirectSolver >::solveAdjAndNormal(), Kaskade::Limex< Eq >::step(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toColumns(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toFile(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toRows(), Kaskade::SLAPMatrix< Num, offset >::toTriplet(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toVector(), and Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::transpose().


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