Class xlifepp::MatrixStorage#

class MatrixStorage#

Inheritence diagram for xlifepp::MatrixStorage:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="xlifepp::ColCsStorage" tooltip="xlifepp::ColCsStorage"] "8" [label="xlifepp::ColDenseStorage" tooltip="xlifepp::ColDenseStorage"] "2" [label="xlifepp::CsStorage" tooltip="xlifepp::CsStorage"] "7" [label="xlifepp::DenseStorage" tooltip="xlifepp::DenseStorage"] "4" [label="xlifepp::DualCsStorage" tooltip="xlifepp::DualCsStorage"] "9" [label="xlifepp::DualDenseStorage" tooltip="xlifepp::DualDenseStorage"] "13" [label="xlifepp::DualSkylineStorage" tooltip="xlifepp::DualSkylineStorage"] "1" [label="xlifepp::MatrixStorage" tooltip="xlifepp::MatrixStorage" fillcolor="#BFBFBF"] "5" [label="xlifepp::RowCsStorage" tooltip="xlifepp::RowCsStorage"] "10" [label="xlifepp::RowDenseStorage" tooltip="xlifepp::RowDenseStorage"] "12" [label="xlifepp::SkylineStorage" tooltip="xlifepp::SkylineStorage"] "6" [label="xlifepp::SymCsStorage" tooltip="xlifepp::SymCsStorage"] "11" [label="xlifepp::SymDenseStorage" tooltip="xlifepp::SymDenseStorage"] "14" [label="xlifepp::SymSkylineStorage" tooltip="xlifepp::SymSkylineStorage"] "3" -> "2" [dir=forward tooltip="public-inheritance"] "8" -> "7" [dir=forward tooltip="public-inheritance"] "2" -> "1" [dir=forward tooltip="public-inheritance"] "7" -> "1" [dir=forward tooltip="public-inheritance"] "4" -> "2" [dir=forward tooltip="public-inheritance"] "9" -> "7" [dir=forward tooltip="public-inheritance"] "13" -> "12" [dir=forward tooltip="public-inheritance"] "5" -> "2" [dir=forward tooltip="public-inheritance"] "10" -> "7" [dir=forward tooltip="public-inheritance"] "12" -> "1" [dir=forward tooltip="public-inheritance"] "6" -> "2" [dir=forward tooltip="public-inheritance"] "11" -> "7" [dir=forward tooltip="public-inheritance"] "14" -> "12" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for xlifepp::MatrixStorage:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="std::vector< T >" tooltip="std::vector< T >"] "2" [label="std::vector< xlifepp::MatrixStorage * >" tooltip="std::vector< xlifepp::MatrixStorage * >"] "1" [label="xlifepp::MatrixStorage" tooltip="xlifepp::MatrixStorage" fillcolor="#BFBFBF"] "2" -> "3" [dir=forward tooltip="template-instance"] "1" -> "2" [dir=forward tooltip="usage"] }

abstract base class of all matrix storage methods

Subclassed by xlifepp::CsStorage, xlifepp::DenseStorage, xlifepp::SkylineStorage

Public Functions

MatrixStorage(const MatrixStorage&)#

copy constructor

copy constructor, NOTE : nbObjectsSharingThis_ is reset to 0

MatrixStorage(StorageType, AccessType at = _noAccess, string_t id = "")#

constructor by storage and access types

MatrixStorage(StorageType, AccessType, number_t, number_t, string_t id = "")#

constructor by storage and access types, number of columns an rows

construct MatrixStorage with size

MatrixStorage(string_t id = "")#

default constructor

virtual ~MatrixStorage()#

virtual destructor

destruct MatrixStorage, error if storage is shared by other matrices

inline AccessType accessType() const#

access to type of access

void add(const MatrixStorage&, const std::vector<number_t> &rowmap = std::vector<number_t>(), const std::vector<number_t> &colmap = std::vector<number_t>())#

add storage in current storage with row/col mapping

add storage in current storage with row/col mapping (generic algorithm, use addRow, addCol, getRow and getCol of children) ms is the matrix storage to add to current storage rowmap(r) gives the row number in current storage of row r of MatrixStorage ms (start index: 1) colmap(c) gives the col number in current storage of col c of MatrixStorage ms (start index: 1) if rowmap (resp.

colmap) is void, it means that row numberings are the same for the both storages !!! warning: ms storage has to be smaller than the current storage

inline virtual void addCol(number_t, const std::set<number_t>&, MatrixPart)#

add row-cols (r,c1), (r,c2), … given by a cols set (r,c>=1)

inline virtual void addIndices(const std::vector<number_t>&, const std::vector<number_t>&)#

add indices (i,j) given by rows (i) and columns (j) (i,j >= 1)

inline virtual void addMatrixMatrix(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&) const#

Matrix + Matrix.

inline virtual void addRow(number_t, const std::set<number_t>&, MatrixPart)#

add row-cols (r,c1), (r,c2), … given by a cols set (r,c>=1)

inline virtual void addSubMatrixIndices(const std::vector<number_t>&, const std::vector<number_t>&)#

add dense submatrix indices in storage

inline StorageBuildType &buildType()#

write access to storage build type

inline StorageBuildType buildType() const#

access to storage build type

inline virtual void clear()#

clear storage vectors

virtual MatrixStorage *clone() const = 0#

create a clone (like a virtual copy constructor)

inline virtual void deleteCols(number_t c1, number_t c2, std::vector<complex_t> &v, SymType s)#

delete cols c1,…, c2 (may be expansive for some storage)

inline virtual void deleteCols(number_t c1, number_t c2, std::vector<Matrix<complex_t>> &v, SymType s)#

delete cols c1,…, c2 (may be expansive for some storage)

inline virtual void deleteCols(number_t c1, number_t c2, std::vector<Matrix<real_t>> &v, SymType s)#

delete cols c1,…, c2 (may be expansive for some storage)

inline virtual void deleteCols(number_t c1, number_t c2, std::vector<real_t> &v, SymType s)#

delete cols c1,…, c2 (may be expansive for some storage)

inline virtual void deleteRows(number_t r1, number_t r2, std::vector<complex_t> &v, SymType s)#

delete rows r1,…, r2 (may be expansive for some storage)

inline virtual void deleteRows(number_t r1, number_t r2, std::vector<Matrix<complex_t>> &v, SymType s)#

delete rows r1,…, r2 (may be expansive for some storage)

inline virtual void deleteRows(number_t r1, number_t r2, std::vector<Matrix<real_t>> &v, SymType s)#

delete rows r1,…, r2 (may be expansive for some storage)

inline virtual void deleteRows(number_t r1, number_t r2, std::vector<real_t> &v, SymType s)#

delete rows r1,…, r2 (may be expansive for some storage)

inline virtual void diagonalMatrixVector(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

diagonal Matrix * Vector (Scalar)

inline number_t diagonalSize() const#

returns number of diagonal entries

inline virtual void diagonalSolver(const std::vector<real_t>&, std::vector<real_t>&, std::vector<real_t>&) const#

Diagonal part factorizations.

MatrixStorage *extract(const std::vector<number_t> &rowIndex, const std::vector<number_t> &colIndex, const string_t &id = "")#

extract storage of submatrix

construct a storage from submatrix of current storage, storage type is not changed

inline virtual void fillSkylineValues(const std::vector<complex_t>&, std::vector<complex_t>&, SymType, MatrixStorage*) const#

fill complex values of current storage as a skyline storage

inline virtual void fillSkylineValues(const std::vector<Matrix<complex_t>>&, std::vector<Matrix<complex_t>>&, SymType, MatrixStorage*) const#

fill values of current storage as a skyline storage

move matrix from from current storage to an other one, return false if not available

inline virtual void fillSkylineValues(const std::vector<Matrix<real_t>>&, std::vector<Matrix<real_t>>&, SymType, MatrixStorage*) const#

fill values of current storage as a skyline storage

inline virtual void fillSkylineValues(const std::vector<real_t>&, std::vector<real_t>&, SymType, MatrixStorage*) const#

fill real values of current storage as a skyline storage

virtual std::vector<std::pair<number_t, number_t>> getCol(SymType, number_t c, number_t r1 = 1, number_t r2 = 0) const#

get (row indices, adress) of col c in set [r1,r2], see child for better algorithm

get (row indices, adress) of col c in set [r1,r2], generic algorithm, see child for better algorithms

virtual std::set<number_t> getCols(number_t r, number_t c1 = 1, number_t c2 = 0) const#

get col indices of row r in set [c1,c2], generic algorithm using position, see children for better algorithms

get col indices of row r in set [c1,c2], generic algorithm using position, see children for better algorithms if c2=0 (default) get all from c1

virtual void getColsV(std::vector<number_t>&, number_t&, number_t r, number_t c1 = 1, number_t c2 = 0) const#

get col indices of row r in set [c1,c2], faster generic algorithm using position, see children for better algorithms

get col indices of row r in set [c1,c2], generic algorithm using position, see children for better algorithms update a vector of numbers with no reallocation, size it before.

nbcol is the number of loaded col indices

virtual std::vector<number_t> getDiag() const#

get diagonal adresses, generic algorithm using position, see children for better algorithms

get diagonal adresses, algorithm using position in non sym/dual storage cases

virtual std::vector<std::pair<number_t, number_t>> getRow(SymType, number_t r, number_t c1 = 1, number_t c2 = 0) const#

get (col indices, adress) of row r in set [c1,c2], see child for better algorithm

virtual std::set<number_t> getRows(number_t c, number_t r1 = 1, number_t r2 = 0) const#

get row indices of col c in set [r1,r2], generic algorithm using position, see children for better algorithms

get row indices of col c in set [r1,r2], generic algorithm using position, see children for better algorithms if r2=0 (default) get all from r1

virtual void getRowsV(std::vector<number_t>&, number_t&, number_t c, number_t r1 = 1, number_t r2 = 0) const#

get row indices of col c in set [r1,r2], faster generic algorithm using position, see children for better algorithms

get row indices of col c in set [r1,r2], generic algorithm using position, see children for better algorithms update a vector of numbers with no reallocation, size it before.

nbrow is the number of loaded row indices

inline virtual void ildlstar(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _selfAdjoint) const#

virtual iLDL* factorizations (Cs child classes only)

inline virtual void ildlt(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _symmetric) const#

virtual incomplete LDLt factorizations (DualCs classes only)

inline virtual void illstar(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _selfAdjoint) const#

virtual iLL* factorizations (Cs child classes only)

inline virtual void illt(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _symmetric) const#

virtual incomplete LLt factorizations (DualCs classes only)

inline virtual void ilu(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _noSymmetry) const#

virtual ILU factorizations ( CS child classes only)

inline virtual bool isDiagonal() const#

true if matrix is diagonal

inline virtual void ldlstar(std::vector<real_t>&, std::vector<real_t>&) const#

virtual LDL* factorizations (Skyline child classes only)

inline virtual void ldlt(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _symmetric) const#

virtual LDLt factorizations (Skyline child classes only)

inline virtual void loadFromFileCoo(std::istream&, std::vector<complex_t>&, SymType, bool)#

load a complex coordinate matrix

inline virtual void loadFromFileCoo(std::istream&, std::vector<Matrix<complex_t>>&, SymType, bool)#

load a coordinate matrix of complex matrices

inline virtual void loadFromFileCoo(std::istream&, std::vector<Matrix<real_t>>&, SymType, bool)#

load a coordinate matrix of real matrices

inline virtual void loadFromFileCoo(std::istream&, std::vector<real_t>&, SymType, bool)#

load a real coordinate matrix

inline virtual void loadFromFileDense(std::istream&, std::vector<complex_t>&, SymType, bool)#

load a complex dense matrix from file

inline virtual void loadFromFileDense(std::istream&, std::vector<Matrix<complex_t>>&, SymType, bool)#

load a dense matrix of complex matrices from file

inline virtual void loadFromFileDense(std::istream&, std::vector<Matrix<real_t>>&, SymType, bool)#

load a dense matrix of real matrices from file

inline virtual void loadFromFileDense(std::istream&, std::vector<real_t>&, SymType, bool)#

load a real dense matrix from file

inline virtual void lowerD1LeftSolver(const std::vector<real_t>&, std::vector<real_t>&, std::vector<real_t>&) const#

Lower + diag(1) part left factorizations.

inline virtual void lowerD1MatrixVector(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

lower Matrix diag 1 (Scalar) * Vector (Scalar)

inline virtual void lowerD1Solver(const std::vector<real_t>&, std::vector<real_t>&, std::vector<real_t>&) const#

Lower + diag(1) part factorizations.

inline virtual void lowerMatrixVector(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

lower Matrix (Scalar) * Vector (Scalar)

virtual number_t lowerPartSize() const = 0#

returns number of entries of lower triangular part

inline virtual void lu(std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _noSymmetry) const#

virtual LU factorizations (Skyline and Dense child classes only)

virtual void multDiagMatrixMatrix(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&) const#

diag Matrix * Matrix (result in current storage), no virtual template

template<typename IteratorA, typename IteratorB, typename IteratorR>
void multDiagMatrixMatrixGeneric(IteratorA, IteratorB, IteratorR) const#

generic diag Matrix * Matrix

template Diag Matrix * Matrix, generic algorithm (see child for better ones) itA: iterator on diag matrix values itB: iterator on matrix values itR: iterator on result matrix values

virtual void multMatrixDiagMatrix(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&) const#

Matrix * diag Matrix (result in current storage), no virtual template.

template<typename IteratorA, typename IteratorB, typename IteratorR>
void multMatrixDiagMatrixGeneric(IteratorA, IteratorB, IteratorR) const#

generic Matrix * diag Matrix

template Diag Matrix * Matrix, generic algorithm (see child for better ones) itA: iterator on matrix values itB: iterator on diag matrix values itR: iterator on result matrix values

virtual void multMatrixMatrix(const std::vector<real_t>&, const MatrixStorage&, const std::vector<real_t>&, std::vector<real_t>&, SymType = _noSymmetry, SymType = _noSymmetry) const#

Matrix * Matrix (result in row dense storage), no virtual template.

template<typename IteratorA, typename IteratorB, typename IteratorR>
void multMatrixMatrixGeneric(IteratorA, const MatrixStorage&, IteratorB, IteratorR, SymType = _noSymmetry, SymType = _noSymmetry) const#

generic Matrix * Matrix

template Matrix * Matrix, generic algorithm (see child for better ones) itmA: iterator on A matrix values itmB: iterator on B matrix values itmR: iterator on result matrix values

inline virtual void multMatrixVector(const std::vector<Matrix<real_t>>&, const std::vector<Vector<real_t>>&, std::vector<Vector<real_t>>&, SymType) const#

Matrix (Matrix) * Vector (Vector)

inline virtual void multMatrixVector(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

Matrix (Scalar) * Vector (Scalar)

inline virtual void multVectorMatrix(const std::vector<Matrix<real_t>>&, const std::vector<Vector<real_t>>&, std::vector<Vector<real_t>>&, SymType) const#

Vector (Vector) * Matrix (Matrix)

inline virtual void multVectorMatrix(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

Vector (Scalar) * Matrix (Scalar)

string_t name() const#

returns type of storage as a string

inline number_t &nbOfColumns()#

returns number of columns

inline number_t nbOfColumns() const#

returns number of columns

inline number_t &nbOfRows()#

returns number of rows

inline number_t nbOfRows() const#

returns number of rows

void noFactorization(const string_t& = "") const#

error handler when no factorization

void noSolver(const string_t& = "") const#

error handler when no solver

inline number_t numberOfObjects() const#

number of objects using this matrix Storage

inline void objectMinus()#

decrease number of objects using this matrix Storage

inline void objectPlus()#

increase number of objects using this matrix Storage

inline virtual number_t pos(number_t i, number_t j, SymType sym = _noSymmetry) const#

returns adress of entry (i,j)

inline virtual void positions(const std::vector<number_t>&, const std::vector<number_t>&, std::vector<number_t>&, bool errorOn = true, SymType = _noSymmetry) const#

access to submatrix positions

inline virtual void print(PrintStream &os = thePrintStream) const#

visualize storage on ostream

virtual void print(std::ostream &os = thePrintStream[0]) const#

visualize storage on ostream

template<typename T>
inline void printCooMatrix(PrintStream &os, const std::vector<T> &m, SymType s = _noSymmetry, real_t tol = theTolerance) const#

output matrix in coordinate format

template<typename T>
void printCooMatrix(std::ostream &os, const std::vector<T> &m, SymType s = _noSymmetry, real_t tol = theTolerance) const#

output matrix in coordinate format

template<typename T>
inline void printDenseMatrix(PrintStream &os, const std::vector<T> &m, SymType s = _noSymmetry) const#

output matrix in dense format

template<typename T>
void printDenseMatrix(std::ostream &os, const std::vector<T> &m, SymType s = _noSymmetry) const#

output matrix in dense format

inline virtual void printEntries(std::ostream&, const std::vector<complex_t>&, number_t vb = 0, const SymType sym = _noSymmetry) const#

output matrix of complex scalars

inline virtual void printEntries(std::ostream&, const std::vector<Matrix<complex_t>>&, number_t vb = 0, const SymType sym = _noSymmetry) const#

output matrix of complex matrices

inline virtual void printEntries(std::ostream&, const std::vector<Matrix<real_t>>&, number_t vb = 0, const SymType sym = _noSymmetry) const#

output matrix of real matrices

inline virtual void printEntries(std::ostream&, const std::vector<real_t>&, number_t vb = 0, const SymType sym = _noSymmetry) const#

output matrix of real scalars

void printHeader(std::ostream&) const#

print header of matrixstorage (name, size, …)

virtual bool sameStorage(const MatrixStorage&) const = 0#

check if two storages have the same structures

virtual std::vector<std::vector<number_t>> scalarColIndices(dimen_t nbr = 1, dimen_t nbc = 1)#

create scalar col indices for each row from current storage and submatrix sizes

create scalar col indices for each row from current storage and submatrix sizes this is a generic tool used by toScalar() member function of storage children it is based on getRow function (virtual function) nbr: row size of submatrix nbc: col size of submatrix if nbc=nbr=1 (scalar case) return the same col indices of current storage See children for optimized algorithms

inline bool &scalarFlag()#

write access to scalarFlag

inline bool scalarFlag() const#

access to scalarFlag

inline virtual void setDiagValue(std::vector<real_t>&, const real_t)#

set value of diagonal

virtual number_t size() const = 0#

returns number of entries

inline virtual std::vector<number_t> skylineColPointer() const#

return skyline col pointer from current storage

inline virtual std::vector<number_t> skylineRowPointer() const#

return skyline row pointer from current storage

inline virtual void sorDiagonalMatrixVector(const std::vector<real_t> &m, const std::vector<real_t> &v, std::vector<real_t> &r, const real_t w) const#

virtual matrix vector products w*Dv for [S]SOR algorithms

template<typename M, typename V, typename R>
void sorDiagonalMatrixVectorG(const std::vector<M> &m, const std::vector<V> &v, std::vector<R> &r, const real_t w) const#

Generic SOR methods, not fast , see child class for best algorithms.

SOR operation w*D*v.

SOR operation w*D*v

inline virtual void sorDiagonalSolver(const std::vector<real_t> &m, const std::vector<real_t> &b, std::vector<real_t> &x, const real_t w) const#

virtual diagonal and triangular linear system solvers [D/w] x = b for [S]SOR algorithms

template<typename M, typename V, typename R>
void sorDiagonalSolverG(const std::vector<M> &m, const std::vector<R> &b, std::vector<V> &x, const real_t w) const#

SOR solver D/w*x = b.

inline virtual void sorLowerMatrixVector(const std::vector<real_t> &m, const std::vector<real_t> &v, std::vector<real_t> &r, const real_t w, const SymType sym = _noSymmetry) const#

virtual matrix vector products [w*D+L]v for [S]SOR algorithms

template<typename M, typename V, typename R>
void sorLowerMatrixVectorG(const std::vector<M> &m, const std::vector<V> &v, std::vector<R> &r, const real_t w) const#

SOR operation [w*D+L]*v.

inline virtual void sorLowerSolver(const std::vector<real_t> &m, const std::vector<real_t> &b, std::vector<real_t> &x, const real_t w) const#

virtual diagonal and triangular linear system solvers [D/w + L] x = b for [S]SOR algorithms

template<typename M, typename V, typename R>
void sorLowerSolverG(const std::vector<M> &m, const std::vector<V> &b, std::vector<R> &x, const real_t w) const#

SOR solver (D/w+L)*x = b.

inline virtual void sorUpperMatrixVector(const std::vector<real_t> &m, const std::vector<real_t> &v, std::vector<real_t> &r, const real_t w, const SymType sym = _noSymmetry) const#

virtual matrix vector products [w*D+U]v for [S]SOR algorithms

template<typename M, typename V, typename R>
void sorUpperMatrixVectorG(const std::vector<M> &m, const std::vector<V> &v, std::vector<R> &r, const real_t w, const SymType sym) const#

SOR operation [w*D+U]*v.

inline virtual void sorUpperSolver(const std::vector<real_t> &m, const std::vector<real_t> &b, std::vector<real_t> &x, const real_t w, const SymType sym = _noSymmetry) const#

virtual diagonal and triangular linear system solvers [D/w+U] x = b for [S]SOR algorithms

template<typename M, typename V, typename R>
void sorUpperSolverG(const std::vector<M> &m, const std::vector<V> &b, std::vector<R> &x, const real_t w, const SymType sym) const#

SOR solver (D/w+U)*x = b.

inline StorageType storageType() const#

returns type of storage as a string

access to type of storage

inline virtual MatrixStorage *toColStorage()#

return new col access storage from current storage

inline virtual MatrixStorage *toDual()#

create a new dual storage from sym storage

inline virtual MatrixStorage *toRowStorage()#

return new row access storage from current storage

inline virtual MatrixStorage *toScalar(dimen_t, dimen_t)#

create a new matrixstorage for scalar storage from non scalar storage

inline virtual void toUmfPack(const std::vector<real_t>&, std::vector<int_t>&, std::vector<int_t>&, std::vector<real_t>&, const SymType sym = _noSymmetry) const#

conversion to umfpack format

inline virtual MatrixStorage *transpose(const std::vector<real_t>&, std::vector<real_t>&) const#

transpose matrix

inline virtual void upperD1MatrixVector(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

upper Matrix diag 1 (Scalar) * Vector (Scalar)

inline virtual void upperD1Solver(const std::vector<real_t>&, std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _noSymmetry) const#

Upper + diag(1) part factorizations.

inline virtual void upperLeftSolver(const std::vector<real_t>&, std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _noSymmetry) const#

upper part left factorizations

inline virtual void upperMatrixVector(const std::vector<real_t>&, const std::vector<real_t>&, std::vector<real_t>&, SymType) const#

upper Matrix (Scalar) * Vector (Scalar)

virtual number_t upperPartSize() const = 0#

returns number of entries of upper triangular part

inline virtual void upperSolver(const std::vector<real_t>&, std::vector<real_t>&, std::vector<real_t>&, const SymType sym = _noSymmetry) const#

upper part factorizations

inline virtual void visual(PrintStream &os = thePrintStream) const#

visualize storage on ostream

virtual void visual(std::ostream &os = thePrintStream[0]) const#

visualize storage on ostream

Public Members

string_t stringId#

string identifier based on characteristic pointers (e.g row/col space pointers)

Public Static Functions

static void clearGlobalVector()#

delete all MatrixStorage objects

static void printAllStorages(std::ostream&)#

print the list of MatrixStorage objects in memory

Public Static Attributes

static std::vector<MatrixStorage*> theMatrixStorages#

list of all matrix storages

list of pointers of all matrix storages

Friends

friend std::ostream &operator<<(std::ostream&, const MatrixStorage&)#

print operator