Class xlifepp::RowDenseStorage#

class RowDenseStorage : public xlifepp::DenseStorage#

Inheritence diagram for xlifepp::RowDenseStorage:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="xlifepp::DenseStorage" tooltip="xlifepp::DenseStorage"] "3" [label="xlifepp::MatrixStorage" tooltip="xlifepp::MatrixStorage"] "1" [label="xlifepp::RowDenseStorage" tooltip="xlifepp::RowDenseStorage" fillcolor="#BFBFBF"] "2" -> "3" [dir=forward tooltip="public-inheritance"] "1" -> "2" [dir=forward tooltip="public-inheritance"] }

Collaboration diagram for xlifepp::RowDenseStorage:

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

handles dense storage of matrix stored row by row

Public Functions

RowDenseStorage(number_t, number_t, string_t id = "RowDenseStorage")#

constructor by access type, number of columns and rows

RowDenseStorage(number_t, string_t id = "RowDenseStorage")#

constructor by access type, number of columns and rows (square matrix)

RowDenseStorage(string_t id = "RowDenseStorage")#

default constructor

inline virtual ~RowDenseStorage()#

virtual destructor

template<typename M1, typename M2, typename R>
void addMatrixMatrix(const std::vector<M1>&, const std::vector<M2>&, std::vector<R>&) const#

template Matrix + Matrix addition and specializations (see below)

Add two matrices.

templated row dense Matrix + Matrix

Parameters:
  • m1 – vector values_ of first matrix

  • m2 – vector values_ of second matrix

  • r – vector values_ of result matrix

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

Matrix + Matrix.

inline virtual RowDenseStorage *clone() const#

create a clone (virtual copy constructor, covariant)

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

diag Matrix * Vector (Scalar)

template<typename M, typename V, typename X>
void diagonalSolver(const std::vector<M> &m, std::vector<V> &v, std::vector<X> &x) const#

Diagonal linear system solver: D x = b.

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

Specializations of diagonal linear solvers D x = v.

inline virtual void gaussSolver(std::vector<real_t> &A, std::vector<std::vector<real_t>> &bs) const#

specializations of gaussSolver and LU factorization thats call generic functions

template<typename T>
void gaussSolver(std::vector<T> &A, std::vector<std::vector<T>> &b) const#

Gauss elimination with pivoting (omp); n rhs.

template<typename T>
void gaussSolver(std::vector<T> &A, std::vector<T> &b) const#

Gauss elimination with pivoting (omp); 1 rhs.

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

get (row indices, adress) of col c in set [r1,r2]

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

get (col indices, adress) of row r in set [c1,c2]

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

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

template<typename M, typename V, typename X>
void lowerD1Solver(const std::vector<M> &m, std::vector<V> &v, std::vector<X> &x) const#

Template diagonal and triangular part solvers (after factorization) Lower triangular with unit diagonal linear system solver: (I+L) x = b.

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

Specializations of lower triangular part with unit diagonal linear solvers (I + L) x = v */.

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

lower Matrix (Scalar) * Vector (Scalar)

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

virtual LU factorizations (Skyline and Dense child classes only)

template<typename T>
void lu(std::vector<T> &A, std::vector<T> &LU) const#

LU factorization with no pivoting (omp)

LU factorization with no pivoting.

template<typename T>
void lu(std::vector<T> &A, std::vector<T> &LU, std::vector<number_t> &p) const#

LU factorization with row pivoting (omp)

LU factorization with real row pivoting LU is a row dense matrix where its rows have been permuted according to the row permutation vector (LU=PA) be care: solving A*X=B is equivalent to solving P*A*X=P*B <=> L*U*X=PB, the product by P has to be done before calling triangular part solvers! product B=A*X is equivalent to product B=inv(P)*L*U*X, the product by inv(P) has to be done after calling triangular part products! product B=X*A is equivalent to product B=X*inv(P)*L*U=tr(tr(U)*tr(L)*P*tr(X)), the product by P has to be done before calling triangular part products! this products are not taken into account by storage class !

template<typename M, typename V, typename R>
void multMatrixVector(const std::vector<M>&, const std::vector<V>&, std::vector<R>&) const#

templated row dense Matrix x Vector

template<typename M, typename V, typename R>
void multMatrixVector(const std::vector<M> &m, V *vp, R *rp) const#

templated row dense Matrix x Vector (pointer form)

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

Matrix (Matrix) * Vector (Vector)

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

Matrix (Scalar) * Vector (Scalar)

inline virtual void multMatrixVector(const std::vector<real_t> &m, real_t *vp, real_t *rp, SymType sym) const#

Matrix * Vector (pointer form)

template<typename M, typename V, typename R>
void multVectorMatrix(const std::vector<M>&, const std::vector<V>&, std::vector<R>&) const#

templated Vector x row dense Matrix

template<typename M, typename V, typename R>
void multVectorMatrix(const std::vector<M> &m, V *vp, R *rp) const#

templated row dense Vector x Matrix (pointer form)

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

Matrix (Matrix) * Vector (Vector)

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

Vector (Scalar) * Matrix (Scalar)

inline virtual void multVectorMatrix(const std::vector<real_t> &m, real_t *vp, real_t *rp, SymType sym) const#

Vector * Matrix (pointer form)

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

overloaded pos returns adress of entry (i,j)

inline virtual number_t pos_(number_t i, number_t j, SymType s = _noSymmetry) const#

overloaded pos fast returns adress of entry (i,j)

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

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

output row dense matrix of complex scalars

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

output row dense matrix of complex matrices

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

output row dense matrix of real matrices

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

output row dense matrix of real scalars

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

output row dense matrix of complex vectors

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

output row dense matrix of real vectors

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

set value of diagonal

template<typename T>
void setDiagValueRowDense(std::vector<T> &m, const T k)#

Set value of Diagonal.

virtual RowDenseStorage *toScalar(dimen_t, dimen_t)#

create a new scalar RowDense storage from current RowDense storage and submatrix sizes

template<typename M, typename OrdinalType>
void toUmfPack(const std::vector<M> &values, std::vector<OrdinalType> &colPointer, std::vector<OrdinalType> &rowIndex, std::vector<M> &mat) const#

conversion to umfpack

template<typename M1, typename Idx>
void toUmfPack(const std::vector<M1> &m1, std::vector<Idx> &colPtUmf, std::vector<Idx> &rowIdxUmf, std::vector<M1> &resultUmf) const#

Extract and convert matrix storage to UMFPack format (Matlab sparse matrix)

Parameters:
  • m1 – vector values_ current matrix

  • colPtUmf – vector column Pointer of UMFPack format

  • rowIdxUmf – vector row Index of UMFPack format

  • resultUmf – vector values of UMFPack format

inline virtual void toUmfPack(const std::vector<real_t> &values, std::vector<int_t> &colPointer, std::vector<int_t> &rowIndex, std::vector<real_t> &mat, const SymType sym) const#

specializations of umfpack conversion

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

transpose matrix

template<typename T>
MatrixStorage *transpose(const std::vector<T> &m, std::vector<T> &mt) const#

transpose matrix

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

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

template<typename M, typename V, typename X>
void upperD1Solver(const std::vector<M> &m, std::vector<V> &v, std::vector<X> &x) const#

Upper triangular with unit diagonal linear system solver: (I+U) x = b.

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

Specializations of upper triangular part with unit diagonal linear solvers (I + U) x = v.

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

upper Matrix (Scalar) * Vector (Scalar)

template<typename M, typename V, typename X>
void upperSolver(const std::vector<M> &m, std::vector<V> &v, std::vector<X> &x) const#

upper triangular linear system solver: (D+U) x = b

inline virtual void upperSolver(const std::vector<real_t> &m, std::vector<real_t> &v, std::vector<real_t> &x, const SymType sym) const#

specializations of upper triangular part linear solvers (D + U) x = v