KASKADE 7 development version
|
#include <triplet.hh>
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.
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... | |
MatrixAsTriplet & | operator+= (MatrixAsTriplet const &m) |
Matrix addition: (*this) += m (works also for matrices with non-matching sparsity pattern) More... | |
MatrixAsTriplet & | operator-= (MatrixAsTriplet const &m) |
MatrixAsTriplet & | operator*= (Scalar const &scalar) |
Scaling. More... | |
MatrixAsTriplet & | operator/= (Scalar const &scalar) |
Scaling. More... | |
MatrixAsTriplet & | operator= (MatrixAsTriplet const &m)=default |
Assignment. More... | |
MatrixAsTriplet & | operator= (MatrixAsTriplet &&other)=default |
Move assignment. More... | |
void | scaleRows (std::vector< Scalar >const &scaling) |
MatrixAsTriplet & | transpose () |
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... | |
typedef std::vector<SparseIndexInt>::const_iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::const_index_iterator |
Definition at line 64 of file triplet.hh.
typedef std::vector<Scalar>::const_iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::const_iterator |
Definition at line 62 of file triplet.hh.
typedef std::vector<SparseIndexInt>::iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::index_iterator |
Definition at line 63 of file triplet.hh.
typedef std::vector<Scalar>::iterator Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::iterator |
Definition at line 61 of file triplet.hh.
typedef Scalar Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::value_type |
STL-compliant typedefs.
Definition at line 60 of file triplet.hh.
|
inline |
Constructor.
Definition at line 67 of file triplet.hh.
|
inlineexplicit |
Constructor that allocates memory directly.
Definition at line 70 of file triplet.hh.
|
inlineexplicit |
Constructor from a Dune::BCRSMatrix.
Definition at line 74 of file triplet.hh.
|
inlineexplicit |
Constructor from a NumaBCRSMatrix.
Definition at line 81 of file triplet.hh.
|
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.
|
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.
|
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.
|
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.
|
default |
Copy constructor.
|
default |
Move constructor.
|
inline |
Copy constructor from a triplet matrix with different entry and/or index types.
Definition at line 168 of file triplet.hh.
|
inline |
add a column
Definition at line 536 of file triplet.hh.
|
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().
|
inline |
add to a dense matrix
Definition at line 549 of file triplet.hh.
|
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().
|
inline |
Matrix-vector multiplication: out = (*this) * in.
Definition at line 312 of file triplet.hh.
Referenced by Kaskade::Bridge::DirectInnerSolver< DirectSolver >::ax().
|
inline |
Matrix-vector multiplication: out += alpha * (*this) * in.
?????: Was ist denn nr? number of blocks/rows?? Damit kann man mehrere Vektoren gleichzeitig axpy'en
Definition at line 262 of file triplet.hh.
Referenced by Kaskade::Bridge::KaskadeLinearization< Functional >::d2axpy(), Kaskade::Bridge::NormalStepLinearization< Functional, stateId, adjointId >::d2axpy(), Kaskade::Bridge::TangentialStepLinearization< Functional, stateId, adjointId >::d2axpy(), Kaskade::Bridge::KaskadeLinearization< Functional >::d2taxpy(), Kaskade::Bridge::NormalStepLinearization< Functional, stateId, adjointId >::ddxpy(), Kaskade::Bridge::TangentialStepLinearization< Functional, stateId, adjointId >::ddxpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::mv(), Kaskade::EnergyScalarProduct< X >::operator()(), and Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::usmv().
|
inline |
get iterator on data field
Definition at line 586 of file triplet.hh.
|
inline |
get const iterator on data field
Definition at line 589 of file triplet.hh.
|
inline |
Deletes the data, leaving an empty matrix of size 0 x 0.
Definition at line 182 of file triplet.hh.
|
inline |
Deletes all sub-diagonal entries.
Definition at line 461 of file triplet.hh.
Referenced by Kaskade::HigherOrderRecovery< Grid, Space >::getErrorFunction().
|
inline |
get iterator on end of data field
Definition at line 592 of file triplet.hh.
|
inline |
get const iterator on end of data field
Definition at line 595 of file triplet.hh.
|
inline |
Definition at line 555 of file triplet.hh.
|
inline |
Definition at line 658 of file triplet.hh.
Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::isSymmetric().
|
inline |
Definition at line 642 of file triplet.hh.
|
inline |
Definition at line 562 of file triplet.hh.
|
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.
|
inline |
matrix-vector multiplication
Definition at line 328 of file triplet.hh.
|
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.
|
inline |
Returns number of cols (computes them, if not known)
Definition at line 215 of file triplet.hh.
Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::atxpy(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::ax(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::axpy(), Kaskade::SimpleSparseLinearSystem::cols(), Kaskade::getFactorization(), Kaskade::Bridge::NormalStepLinearization< Functional, stateId, adjointId >::getMatrixBlocks(), Kaskade::Bridge::TangentialStepLinearization< Functional, stateId, adjointId >::getMatrixBlocks(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::M(), Kaskade::DerivativeChecker< Functional, checkD1, SparseInt >::matrixToVTK(), Kaskade::DerivativeCheck::matrixToVTK(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator+=(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::operator-=(), Kaskade::SLAPMatrix< Num, offset >::SLAPMatrix(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toColumns(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toRows(), and Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toVector().
|
inline |
Number of non-zeros. This is the number of entries in the triplet format and can be larger than the number of nonzero entries in the represented sparse matrix.
Definition at line 254 of file triplet.hh.
Referenced by Kaskade::AdditiveSchwarzPreconditioner< Op >::AdditiveSchwarzPreconditioner(), Kaskade::ARMSPreconditioner< Op >::ARMSPreconditioner(), Kaskade::BoomerAMG< Op >::BoomerAMG(), Kaskade::Euclid< Op >::Euclid(), Kaskade::getFactorization(), Kaskade::ILUKPreconditioner< Op >::ILUKPreconditioner(), Kaskade::ILUTPreconditioner< Op >::ILUTPreconditioner(), Kaskade::AdditiveSchwarzPreconditioner< Op >::init(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::lumped(), and Kaskade::PartialDirectPreconditioner< Op >::PartialDirectPreconditioner().
|
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().
|
inline |
Definition at line 426 of file triplet.hh.
|
inline |
Matrix addition: (*this) += m (works also for matrices with non-matching sparsity pattern)
Definition at line 357 of file triplet.hh.
|
inline |
Definition at line 389 of file triplet.hh.
|
inline |
Definition at line 435 of file triplet.hh.
|
default |
Move assignment.
|
default |
Assignment.
|
inline |
Output as Matlab source code.
Definition at line 490 of file triplet.hh.
Referenced by Kaskade::operator<<().
|
inline |
Reserve memory.
Definition at line 242 of file triplet.hh.
|
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().
|
inline |
Definition at line 446 of file triplet.hh.
|
inline |
Shift the row and column indices such that \( A_{00} \) is (structurally) nonzero.
Definition at line 193 of file triplet.hh.
|
inline |
Shifts matrix indices. Can be used together with += to concatenate sparese matrices.
Definition at line 334 of file triplet.hh.
Referenced by Kaskade::Bridge::KaskadeLinearization< Functional >::getDiscreteMatrixBlocks(), Kaskade::Bridge::NormalStepLinearization< Functional, stateId, adjointId >::getMatrixBlocks(), Kaskade::Bridge::TangentialStepLinearization< Functional, stateId, adjointId >::getMatrixBlocks(), Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::getTriplet(), Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::setStartToZero(), and Kaskade::Bridge::DirectInnerSolver< DirectSolver >::solveAdjAndNormal().
|
inline |
get number of entries in sparse matrix
Definition at line 598 of file triplet.hh.
Referenced by Kaskade::MatrixAsTriplet< Scalar, SparseIndexInt >::toBCRS().
|
inline |
copy to Dune::BCRSMatrix
Definition at line 613 of file triplet.hh.
|
inline |
transform into a set of column vectors
Definition at line 500 of file triplet.hh.
|
inline |
write to text file
Definition at line 601 of file triplet.hh.
|
inline |
Definition at line 513 of file triplet.hh.
|
inline |
Definition at line 526 of file triplet.hh.
|
inline |
Transposition.
Definition at line 453 of file triplet.hh.
Referenced by Kaskade::Bridge::KaskadeLinearization< Functional >::d2taxpy(), and Kaskade::HierarchicalBasisErrorEstimator2< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverLA, LinearSolverHA, LinearSolverLM, LinearSolverHM, lumpM, RefinementStrategy >::operator()().
|
inline |
Matrix-vector multiplication.
\( x = x + \alpha A^{T}y \)
Definition at line 296 of file triplet.hh.
|
inline |
scaled matrix-vector multiplication
Definition at line 304 of file triplet.hh.
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().
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().
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().