qpOASES
3.2.1
An Implementation of the Online Active Set Strategy
|
Interfaces matrix-vector operations tailored to symmetric sparse matrices. More...
#include <Matrices.hpp>
Public Member Functions | |
SymSparseMat () | |
SymSparseMat (int_t nr, int_t nc, sparse_int_t *r, sparse_int_t *c, real_t *v) | |
SymSparseMat (int_t nr, int_t nc, int_t ld, const real_t *const v) | |
virtual | ~SymSparseMat () |
virtual Matrix * | duplicate () const |
virtual SymmetricMatrix * | duplicateSym () const |
virtual returnValue | bilinear (const Indexlist *const icols, int_t xN, const real_t *x, int_t xLD, real_t *y, int_t yLD) const |
virtual void | free ()=0 |
virtual real_t | diag (int_t i) const =0 |
virtual BooleanType | isDiag () const =0 |
virtual real_t | getNorm (int_t type=2) const =0 |
virtual real_t | getRowNorm (int_t rNum, int_t type=2) const =0 |
virtual returnValue | getRowNorm (real_t *norm, int_t type=2) const =0 |
virtual returnValue | getRow (int_t rNum, const Indexlist *const icols, real_t alpha, real_t *row) const =0 |
virtual returnValue | getCol (int_t cNum, const Indexlist *const irows, real_t alpha, real_t *col) const =0 |
virtual returnValue | getSparseSubmatrix (const Indexlist *const irows, const Indexlist *const icols, int_t rowoffset, int_t coloffset, int_t &numNonzeros, int_t *irn, int_t *jcn, real_t *avals, BooleanType only_lower_triangular=BT_FALSE) const |
virtual returnValue | getSparseSubmatrix (const Indexlist *const irows, int_t idx_icol, int_t rowoffset, int_t coloffset, int_t &numNonzeros, int_t *irn, int_t *jcn, real_t *avals, BooleanType only_lower_triangular=BT_FALSE) const |
virtual returnValue | getSparseSubmatrix (int_t idx_row, const Indexlist *const icols, int_t rowoffset, int_t coloffset, int_t &numNonzeros, int_t *irn, int_t *jcn, real_t *avals, BooleanType only_lower_triangular=BT_FALSE) const |
virtual returnValue | getSparseSubmatrix (int_t irowsLength, const int_t *const irowsNumber, int_t icolsLength, const int_t *const icolsNumber, int_t rowoffset, int_t coloffset, int_t &numNonzeros, int_t *irn, int_t *jcn, real_t *avals, BooleanType only_lower_triangular=BT_FALSE) const =0 |
virtual returnValue | times (int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD) const =0 |
virtual returnValue | times (const Indexlist *const irows, const Indexlist *const icols, int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD, BooleanType yCompr=BT_TRUE) const =0 |
virtual returnValue | transTimes (int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD) const =0 |
virtual returnValue | transTimes (const Indexlist *const irows, const Indexlist *const icols, int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD) const =0 |
virtual returnValue | addToDiag (real_t alpha)=0 |
virtual real_t * | full () const =0 |
virtual returnValue | print (const char *name=0) const =0 |
virtual returnValue | writeToFile (FILE *output_file, const char *prefix) const =0 |
BooleanType | needToFreeMemory () const |
void | doFreeMemory () |
void | doNotFreeMemory () |
virtual void | free () |
virtual void | setVal (const real_t *newVal) |
virtual real_t | diag (int_t i) const |
virtual BooleanType | isDiag () const |
virtual real_t | getNorm (int_t type=2) const |
virtual real_t | getRowNorm (int_t rNum, int_t type=2) const |
virtual returnValue | getRowNorm (real_t *norm, int_t type=2) const |
virtual returnValue | getRow (int_t rNum, const Indexlist *const icols, real_t alpha, real_t *row) const |
virtual returnValue | getCol (int_t cNum, const Indexlist *const irows, real_t alpha, real_t *col) const |
virtual returnValue | getSparseSubmatrix (int_t irowsLength, const int_t *const irowsNumber, int_t icolsLength, const int_t *const icolsNumber, int_t rowoffset, int_t coloffset, int_t &numNonzeros, int_t *irn, int_t *jcn, real_t *avals, BooleanType only_lower_triangular=BT_FALSE) const |
virtual returnValue | times (int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD) const |
virtual returnValue | times (const Indexlist *const irows, const Indexlist *const icols, int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD, BooleanType yCompr=BT_TRUE) const |
virtual returnValue | transTimes (int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD) const |
virtual returnValue | transTimes (const Indexlist *const irows, const Indexlist *const icols, int_t xN, real_t alpha, const real_t *x, int_t xLD, real_t beta, real_t *y, int_t yLD) const |
virtual returnValue | addToDiag (real_t alpha) |
sparse_int_t * | createDiagInfo () |
virtual real_t * | full () const |
virtual returnValue | print (const char *name=0) const |
virtual returnValue | writeToFile (FILE *output_file, const char *prefix) const |
Protected Attributes | |
BooleanType | freeMemory |
int_t | nRows |
int_t | nCols |
sparse_int_t * | ir |
sparse_int_t * | jc |
sparse_int_t * | jd |
real_t * | val |
Symmetric sparse matrix class (column compressed format).
|
inline |
Default constructor.
|
inline |
Constructor with arguments.
nr | Number of rows. |
nc | Number of columns. |
r | Row indices (length). |
c | Indices to first entry of columns (nCols+1). |
v | Vector of entries (length). |
Constructor from dense matrix.
nr | Number of rows. |
nc | Number of columns. |
ld | Leading dimension. |
v | Row major stored matrix elements. |
|
inlinevirtual |
Destructor.
References Matrix::duplicate(), END_NAMESPACE_QPOASES, and real_t.
|
pure virtualinherited |
Adds given offset to diagonal of matrix.
alpha | Diagonal offset. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), QProblemB::regulariseHessian(), and Matrix::~Matrix().
|
virtualinherited |
Adds given offset to diagonal of matrix.
alpha | Diagonal offset. |
Implements Matrix.
References BT_FALSE, SparseMatrix::ir, isZero(), SparseMatrix::jd, SparseMatrix::nCols, SparseMatrix::nRows, RET_DIAGONAL_NOT_INITIALISED, RET_NO_DIAGONAL_AVAILABLE, SUCCESSFUL_RETURN, THROWERROR, and SparseMatrix::val.
|
virtual |
Compute bilinear form y = x'*H*x using submatrix given by index list.
icols | Index list specifying columns of x. |
xN | Number of vectors to multiply. |
x | Input vector to be multiplied (uncompressed). |
xLD | Leading dimension of input x. |
y | Output vector of results (compressed). |
yLD | Leading dimension of output y. |
Implements SymmetricMatrix.
References END_NAMESPACE_QPOASES, Indexlist::iSort, SparseMatrixRow::jd, Indexlist::length, Indexlist::number, RET_DIAGONAL_NOT_INITIALISED, SUCCESSFUL_RETURN, THROWERROR, and SparseMatrixRow::val.
Referenced by QProblem::computeProjectedCholesky().
|
inherited |
Create jd field from ir and jc.
References SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::jd, and SparseMatrix::nCols.
Referenced by QProblemB::createDiagSparseMat(), and solveOqpBenchmark().
Returns i-th diagonal entry.
i | Index. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), QProblemB::determineHessianType(), QProblem::removeBound(), QProblemB::removeBound(), and Matrix::~Matrix().
Returns i-th diagonal entry.
i | Index. |
Implements Matrix.
References INFTY, SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::jd, RET_DIAGONAL_NOT_INITIALISED, THROWERROR, and SparseMatrix::val.
|
inlineinherited |
Enables de-allocation of internal memory.
References BT_TRUE, and Matrix::freeMemory.
Referenced by QProblemB::createDiagSparseMat(), DenseMatrix::duplicate(), SparseMatrix::duplicate(), SparseMatrixRow::duplicate(), SymDenseMat::duplicateSym(), duplicateSym(), QProblemB::setupQPdataFromFile(), QProblem::setupQPdataFromFile(), solveOqpBenchmark(), SparseMatrix::SparseMatrix(), and SparseMatrixRow::SparseMatrixRow().
|
inlineinherited |
Disables de-allocation of internal memory.
Referenced by SparseMatrix::free(), SparseMatrixRow::free(), Matrix::Matrix(), SparseMatrix::SparseMatrix(), and SparseMatrixRow::SparseMatrixRow().
|
virtual |
Returns a deep-copy of the Matrix object.
Reimplemented from SparseMatrix.
|
virtual |
Returns a deep-copy of the SymmetricMatrix object.
Implements SymmetricMatrix.
References Matrix::doFreeMemory(), SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::jd, SparseMatrixRow::jd, SparseMatrix::nCols, SparseMatrixRow::nCols, SparseMatrix::nRows, SparseMatrixRow::nRows, real_t, SparseMatrix::val, and SparseMatrixRow::val.
|
pure virtualinherited |
Frees all internal memory.
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), DenseMatrix::~DenseMatrix(), and Matrix::~Matrix().
|
virtualinherited |
Frees all internal memory.
Implements Matrix.
References Matrix::doNotFreeMemory(), SparseMatrix::ir, SparseMatrix::jc, and SparseMatrix::val.
Referenced by SparseMatrix::~SparseMatrix().
|
pure virtualinherited |
Allocates and creates dense matrix array in row major format.
Note: Calling function has to free allocated memory!
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), SolutionAnalysis::getKktViolation(), QProblem::writeQpDataIntoMatFile(), and Matrix::~Matrix().
|
virtualinherited |
Allocates and creates dense matrix array in row major format.
Note: Calling function has to free allocated memory!
Implements Matrix.
References SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::nCols, SparseMatrix::nRows, real_t, and SparseMatrix::val.
Referenced by SparseMatrix::print().
|
pure virtualinherited |
Retrieve indexed entries of matrix column multiplied by alpha.
cNum | Column number. |
irows | Index list specifying rows. |
alpha | Scalar factor. |
col | Output column vector. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by QProblemB::computeCholesky(), QProblem::computeProjectedCholesky(), DenseMatrix::DenseMatrix(), QProblem::removeBound(), and Matrix::~Matrix().
|
virtualinherited |
Retrieve indexed entries of matrix column multiplied by alpha.
cNum | Column number. |
irows | Index list specifying rows. |
alpha | Scalar factor. |
col | Output column vector. |
Implements Matrix.
References BT_TRUE, SparseMatrix::ir, isEqual(), Indexlist::iSort, SparseMatrix::jc, Indexlist::number, SUCCESSFUL_RETURN, and SparseMatrix::val.
Get the N-norm of the matrix
type | Norm type, 1: one-norm, 2: Euclidean norm. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), DenseMatrix::getNorm(), DenseMatrix::getRowNorm(), QProblemB::regulariseHessian(), and Matrix::~Matrix().
Get the N-norm of the matrix
type | Norm type, 1: one-norm, 2: Euclidean norm. |
Implements Matrix.
References SparseMatrix::jc, SparseMatrix::nCols, REFER_NAMESPACE_QPOASES, and SparseMatrix::val.
|
pure virtualinherited |
Retrieve indexed entries of matrix row multiplied by alpha.
rNum | Row number. |
icols | Index list specifying columns. |
alpha | Scalar factor. |
row | Output row vector. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by QProblem::addConstraint(), QProblem::addConstraint_checkLI(), SQProblemSchur::addConstraint_checkLISchur(), QProblem::addConstraint_ensureLI(), DenseMatrix::DenseMatrix(), QProblemB::removeBound(), and Matrix::~Matrix().
|
virtualinherited |
Retrieve indexed entries of matrix row multiplied by alpha.
rNum | Row number. |
icols | Index list specifying columns. |
alpha | Scalar factor. |
row | Output row vector. |
Implements Matrix.
References BT_TRUE, SparseMatrix::ir, isEqual(), Indexlist::iSort, SparseMatrix::jc, Indexlist::length, SparseMatrix::nCols, Indexlist::number, SUCCESSFUL_RETURN, and SparseMatrix::val.
Get the N-norm of a row
rNum | Row number. |
type | Norm type, 1: one-norm, 2: Euclidean norm. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), QProblem::setA(), and Matrix::~Matrix().
|
pure virtualinherited |
Get the N-norm of all rows
norm | Norm of each row. |
type | Norm type, 1: one-norm, 2: Euclidean norm. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Get the N-norm of a row
rNum | Row number. |
type | Norm type, 1: one-norm, 2: Euclidean norm. |
Implements Matrix.
References getAbs(), getSqrt(), INFTY, SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::nCols, real_t, RET_INVALID_ARGUMENTS, THROWERROR, and SparseMatrix::val.
|
virtualinherited |
Get the N-norm of all rows
norm | Norm of each row. |
type | Norm type, 1: one-norm, 2: Euclidean norm. |
Implements Matrix.
References getAbs(), getSqrt(), SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::nCols, SparseMatrix::nRows, RET_INVALID_ARGUMENTS, SUCCESSFUL_RETURN, and SparseMatrix::val.
|
virtualinherited |
Retrieve entries of submatrix in Harwell-Boeing sparse format. If irn, jcn, and avals are null, this only counts the number of nonzeros. Otherwise, numNonzeros containts the size of irn, jcn, and avals on entry, and the written number of entries on return.
irows | Index list specifying rows. |
icols | Index list specifying columns. |
rowoffset | Offset for row entries. |
coloffset | Offset for row entries. |
numNonzeros | Number of nonzeros in submatrix. |
irn | Row position of entries (as position in irows) plus rowoffset. |
jcn | Column position of entries (as position in irows) plus coloffset. |
avals | Numerical values of the entries. |
only_lower_triangular | if true, only the lower triangular portion is returned. This can only be true for symmetric matrices and if irows==jcols. |
References Indexlist::getLength(), and Indexlist::getNumberArray().
Referenced by DenseMatrix::DenseMatrix(), Matrix::getSparseSubmatrix(), and Matrix::~Matrix().
|
virtualinherited |
Retrieve entries of submatrix in Harwell-Boeing sparse format. If irn, jcn, and avals are null, this only counts the number of nonzeros. Otherwise, numNonzeros containts the size of irn, jcn, and avals on entry, and the written number of entries on return. This version retrieves one column.
irows | Index list specifying rows. |
idx_icol | Index list specifying columns. |
rowoffset | Offset for row entries. |
coloffset | Offset for row entries. |
numNonzeros | Number of nonzeros in submatrix. |
irn | Row position of entries (as position in irows) plus rowoffset. |
jcn | Column position of entries (as position in irows) plus coloffset. |
avals | Numerical values of the entries. |
only_lower_triangular | if true, only the lower triangular portion is returned. This can only be true for symmetric matrices and if irows==jcols. |
References Indexlist::getLength(), Indexlist::getNumberArray(), and Matrix::getSparseSubmatrix().
|
virtualinherited |
Retrieve entries of submatrix in Harwell-Boeing sparse format. If irn, jcn, and avals are null, this only counts the number of nonzeros. Otherwise, numNonzeros containts the size of irn, jcn, and avals on entry, and the written number of entries on return. This version retrieves one row.
idx_row | Row number. |
icols | Index list specifying columns. |
rowoffset | Offset for row entries. |
coloffset | Offset for row entries. |
numNonzeros | Number of nonzeros in submatrix. |
irn | Row position of entries (as position in irows) plus rowoffset. |
jcn | Column position of entries (as position in irows) plus coloffset. |
avals | Numerical values of the entries. |
only_lower_triangular | if true, only the lower triangular portion is returned. This can only be true for symmetric matrices and if irows==jcols. |
References Indexlist::getLength(), Indexlist::getNumberArray(), and Matrix::getSparseSubmatrix().
|
pure virtualinherited |
Retrieve entries of submatrix in Harwell-Boeing sparse format. If irn, jcn, and avals are null, this only counts the number of nonzeros. Otherwise, numNonzeros containts the size of irn, jcn, and avals on entry, and the written number of entries on return.
irowsLength | Number of rows. |
irowsNumber | Array with row numbers. |
icolsLength | Number of columns. |
icolsNumber | Array with column numbers. |
rowoffset | Offset for row entries. |
coloffset | Offset for row entries. |
numNonzeros | Number of nonzeros in submatrix. |
irn | Row position of entries (as position in irows) plus rowoffset. |
jcn | Column position of entries (as position in irows) plus coloffset. |
avals | Numerical values of the entries. |
only_lower_triangular | if true, only the lower triangular portion is returned. This can only be true for symmetric matrices and if irows==jcols. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
|
virtualinherited |
Retrieve entries of submatrix in Harwell-Boeing sparse format. If irn, jcn, and avals are null, this only counts the number of nonzeros. Otherwise, numNonzeros containts the size of irn, jcn, and avals on entry, and the written number of entries on return.
irowsLength | Number of rows. |
irowsNumber | Array with row numbers. |
icolsLength | Number of columns. |
icolsNumber | Array with column numbers. |
rowoffset | Offset for row entries. |
coloffset | Offset for row entries. |
numNonzeros | Number of nonzeros in submatrix. |
irn | Row position of entries (as position in irows) plus rowoffset. |
jcn | Column position of entries (as position in irows) plus coloffset. |
avals | Numerical values of the entries. |
only_lower_triangular | if true, only the lower triangular portion is returned. This can only be true for symmetric matrices and if irows==jcols. |
Implements Matrix.
References BT_FALSE, SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::nRows, RET_INVALID_ARGUMENTS, SUCCESSFUL_RETURN, THROWERROR, and SparseMatrix::val.
|
pure virtualinherited |
Checks whether matrix is square and diagonal.
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), QProblemB::determineHessianType(), and Matrix::~Matrix().
|
virtualinherited |
Checks whether matrix is square and diagonal.
Implements Matrix.
References BT_FALSE, BT_TRUE, SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::nCols, and SparseMatrix::nRows.
|
inlineinherited |
Returns whether internal memory needs to be de-allocated.
References Matrix::freeMemory.
Referenced by DenseMatrix::duplicate(), SymDenseMat::duplicateSym(), DenseMatrix::~DenseMatrix(), SparseMatrix::~SparseMatrix(), and SparseMatrixRow::~SparseMatrixRow().
|
pure virtualinherited |
Prints matrix to screen.
name | Name of matrix. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), DenseMatrix::print(), and Matrix::~Matrix().
|
virtualinherited |
Prints matrix to screen.
name | Name of matrix. |
Implements Matrix.
References SparseMatrix::full(), SparseMatrix::nCols, SparseMatrix::nRows, real_t, and REFER_NAMESPACE_QPOASES.
|
virtualinherited |
Sets value array.
Thanks to Frank Chuang.
newVal | ... |
References SparseMatrix::jc, SparseMatrix::nCols, and SparseMatrix::val.
|
pure virtualinherited |
Evaluate Y=alpha*A*X + beta*Y.
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), QProblem::determineStepDirection(), QProblemB::determineStepDirection(), QProblem::ensureNonzeroCurvature(), QProblemB::getObjVal(), QProblem::performStep(), QProblem::printIteration(), QProblemB::printIteration(), QProblem::removeBound(), QProblem::removeConstraint(), QProblem::setA(), QProblem::setupAuxiliaryQP(), QProblem::setupAuxiliaryQPgradient(), QProblemB::setupAuxiliaryQPgradient(), QProblem::setupAuxiliaryQPsolution(), and Matrix::~Matrix().
|
pure virtualinherited |
Evaluate matrix vector product with submatrix given by Indexlist.
irows | Index list specifying rows. |
icols | Index list specifying columns. |
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
yCompr | Compressed storage for y. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
|
virtualinherited |
Evaluate Y=alpha*A*X + beta*Y.
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
Implements Matrix.
References BT_FALSE, BT_TRUE, SparseMatrix::ir, isEqual(), isZero(), SparseMatrix::jc, SparseMatrix::nCols, SparseMatrix::nRows, SUCCESSFUL_RETURN, and SparseMatrix::val.
|
virtualinherited |
Evaluate matrix vector product with submatrix given by Indexlist.
irows | Index list specifying rows. |
icols | Index list specifying columns. |
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
yCompr | Compressed storage for y. |
Implements Matrix.
References BT_FALSE, BT_TRUE, getAbs(), getMax(), SparseMatrix::ir, isEqual(), Indexlist::iSort, isZero(), SparseMatrix::jc, Indexlist::length, SparseMatrix::nCols, SparseMatrix::nRows, Indexlist::number, real_t, SUCCESSFUL_RETURN, and SparseMatrix::val.
|
pure virtualinherited |
Evaluate Y=alpha*A'*X + beta*Y.
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by QProblem::addBound_ensureLI(), QProblem::addConstraint_ensureLI(), DenseMatrix::DenseMatrix(), QProblem::determineStepDirection(), QProblem::printIteration(), QProblem::setupAuxiliaryQPgradient(), and Matrix::~Matrix().
|
pure virtualinherited |
Evaluate matrix transpose vector product.
irows | Index list specifying rows. |
icols | Index list specifying columns. |
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
|
virtualinherited |
Evaluate Y=alpha*A'*X + beta*Y.
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
Implements Matrix.
References BT_FALSE, BT_TRUE, SparseMatrix::ir, isEqual(), isZero(), SparseMatrix::jc, SparseMatrix::nCols, SUCCESSFUL_RETURN, and SparseMatrix::val.
|
virtualinherited |
Evaluate matrix transpose vector product.
irows | Index list specifying rows. |
icols | Index list specifying columns. |
xN | Number of vectors to multiply. |
alpha | Scalar factor for matrix vector product. |
x | Input vector to be multiplied. |
xLD | Leading dimension of input x. |
beta | Scalar factor for y. |
y | Output vector of results. |
yLD | Leading dimension of output y. |
Implements Matrix.
References BT_FALSE, BT_TRUE, SparseMatrix::ir, isEqual(), Indexlist::iSort, isZero(), SparseMatrix::jc, Indexlist::length, SparseMatrix::nRows, Indexlist::number, real_t, SUCCESSFUL_RETURN, and SparseMatrix::val.
|
pure virtualinherited |
Write matrix to file.
Implemented in SparseMatrixRow, SparseMatrix, and DenseMatrix.
Referenced by DenseMatrix::DenseMatrix(), and Matrix::~Matrix().
|
virtualinherited |
Write matrix to file.
Implements Matrix.
References SparseMatrix::ir, SparseMatrix::jc, SparseMatrix::nCols, SUCCESSFUL_RETURN, and SparseMatrix::val.
|
protectedinherited |
Indicating whether internal memory needs to be de-allocated.
Referenced by Matrix::doFreeMemory(), and Matrix::needToFreeMemory().
|
protectedinherited |
Row indices (length).
Referenced by SparseMatrix::addToDiag(), SparseMatrix::createDiagInfo(), SparseMatrix::diag(), SparseMatrix::duplicate(), duplicateSym(), SparseMatrix::free(), SparseMatrix::full(), SparseMatrix::getCol(), SparseMatrix::getRow(), SparseMatrix::getRowNorm(), SparseMatrix::getSparseSubmatrix(), SparseMatrix::isDiag(), SparseMatrix::SparseMatrix(), SparseMatrix::times(), SparseMatrix::transTimes(), and SparseMatrix::writeToFile().
|
protectedinherited |
Indices to first entry of columns (nCols+1).
Referenced by SparseMatrix::createDiagInfo(), SparseMatrix::diag(), SparseMatrix::duplicate(), duplicateSym(), SparseMatrix::free(), SparseMatrix::full(), SparseMatrix::getCol(), SparseMatrix::getNorm(), SparseMatrix::getRow(), SparseMatrix::getRowNorm(), SparseMatrix::getSparseSubmatrix(), SparseMatrix::isDiag(), SparseMatrix::setVal(), SparseMatrix::SparseMatrix(), SparseMatrix::times(), SparseMatrix::transTimes(), and SparseMatrix::writeToFile().
|
protectedinherited |
Indices to first entry of lower triangle (including diagonal) (nCols).
Referenced by SparseMatrix::addToDiag(), SparseMatrix::createDiagInfo(), SparseMatrix::diag(), SparseMatrix::duplicate(), duplicateSym(), and SparseMatrix::~SparseMatrix().
|
protectedinherited |
Number of columns.
Referenced by SparseMatrix::addToDiag(), SparseMatrix::createDiagInfo(), SparseMatrix::duplicate(), duplicateSym(), SparseMatrix::full(), SparseMatrix::getNorm(), SparseMatrix::getRow(), SparseMatrix::getRowNorm(), SparseMatrix::isDiag(), SparseMatrix::print(), SparseMatrix::setVal(), SparseMatrix::SparseMatrix(), SparseMatrix::times(), SparseMatrix::transTimes(), and SparseMatrix::writeToFile().
|
protectedinherited |
Number of rows.
Referenced by SparseMatrix::addToDiag(), SparseMatrix::duplicate(), duplicateSym(), SparseMatrix::full(), SparseMatrix::getRowNorm(), SparseMatrix::getSparseSubmatrix(), SparseMatrix::isDiag(), SparseMatrix::print(), SparseMatrix::SparseMatrix(), SparseMatrix::times(), and SparseMatrix::transTimes().
|
protectedinherited |
Vector of entries (length).
Referenced by SparseMatrix::addToDiag(), SparseMatrix::diag(), SparseMatrix::duplicate(), duplicateSym(), SparseMatrix::free(), SparseMatrix::full(), SparseMatrix::getCol(), SparseMatrix::getNorm(), SparseMatrix::getRow(), SparseMatrix::getRowNorm(), SparseMatrix::getSparseSubmatrix(), SparseMatrix::setVal(), SparseMatrix::SparseMatrix(), SparseMatrix::times(), SparseMatrix::transTimes(), and SparseMatrix::writeToFile().