ConicBundle
|
Matrix class for real values of type Real More...
#include <matrix.hxx>
Public Member Functions | |
Constructors, Destructor, and Initialization (Members) | |
Matrix () | |
empty matrix | |
Matrix (const Matrix &, Real d=1., int atrans=0) | |
copy constructor, *this=d*A | |
Matrix (const Realrange &) | |
generate a column vector holding the elements of this Realrange | |
Matrix (Integer nr, Integer nc) | |
generate a matrix of size nr x nc but WITHOUT initializing the memory More... | |
Matrix (Integer nr, Integer nc, Real d) | |
generate a matrix of size nr x nc initializing all elements to the value d | |
Matrix (Integer nr, Integer nc, const Real *dp, Integer incr=1, Real d=1.) | |
generate a matrix of size nr x nc initializing the elements from the (one dimensional) array dp with increment incr and scaled by d | |
Matrix (const std::vector< Real > &vec) | |
copy a std::vector<Real> to a column vector of the same size and content | |
~Matrix () | |
void | set_init (bool) |
after external initialization, call matrix.set_init(true) (not needed if CONICBUNDLE_DEBUG is undefined) | |
bool | get_init () const |
returns true if the matrix has been declared initialized (not needed if CONICBUNDLE_DEBUG is undefined) | |
Matrix & | init (const Matrix &A, Real d=1., int atrans=0) |
initialize to *this=A*d where A may be transposed | |
Matrix & | init (const Indexmatrix &A, Real d=1.) |
initialize to *this=A*d | |
Matrix & | init (const Sparsemat &A, Real d=1.) |
initialize to *this=A*d | |
Matrix & | init (const Symmatrix &S, Real d=1.) |
initialize to *this=A*d | |
Matrix & | init (const Sparsesym &, Real d=1.) |
initialize to *this=A*d | |
Matrix & | init (const Realrange &) |
initialize *this to a column vector holding the elements of Realrange | |
Matrix & | init (Integer nr, Integer nc, Real d) |
intialize *this to a matrix of size nr x nc initializing all elements to the value d | |
Matrix & | init (Integer nr, Integer nc, const Real *dp, Integer incr=1, Real d=1.) |
generate a matrix of size nr x nc initializing the elements from the (one dimensional) array dp with increment incr and scaled by d | |
Matrix & | init (const std::vector< Real > &vec) |
use std::vector<Real> to initialize this to a column vector of the same size and content | |
Matrix & | init_diag (int nr, Real d=1.) |
initialize to a diagonal nr x nr matrix with constant diagonal value d | |
Matrix & | init_diag (const Matrix &vec, Real d=1.) |
initialize to a diagonal matrix with diagonal given by vec | |
Matrix & | init_diag (const Indexmatrix &vec, Real d=1.) |
initialize to a diagonal matrix with diagonal given by vec | |
void | newsize (Integer nr, Integer nc) |
resize the matrix to nr x nc elements but WITHOUT initializing the memory More... | |
Conversions from other Matrix Classes (Members) | |
Matrix (const Indexmatrix &A, Real d=1.) | |
(*this)=d*A | |
Matrix (const Sparsemat &A, Real d=1.) | |
(*this)=d*A | |
Matrix (const Symmatrix &S, Real d=1.) | |
(*this)=d*A | |
Matrix (const Sparsesym &, Real d=1.) | |
(*this)=d*A | |
Size and Type Information (Members) | |
void | dim (Integer &_nr, Integer &_nc) const |
returns the number of rows in _nr and the number of columns in _nc | |
Integer | dim () const |
returns the dimension rows * columns when the matrix is regarded as a vector | |
Integer | rowdim () const |
returns the row dimension | |
Integer | coldim () const |
returns the column dimension | |
Mtype | get_mtype () const |
returns the type of the matrix, MTmatrix | |
Indexing and Submatrices (Members) | |
Real & | operator() (Integer i, Integer j) |
returns reference to element (i,j) of the matrix (rowindex i, columnindex j) | |
Real & | operator() (Integer i) |
returns reference to element (i) of the matrix if regarded as vector of stacked columns [element (irowdim, i/rowdim)] | |
Real | operator() (Integer i, Integer j) const |
returns value of element (i,j) of the matrix (rowindex i, columnindex j) | |
Real | operator() (Integer i) const |
returns value of element (i) of the matrix if regarded as vector of stacked columns [element (irowdim, i/rowdim)] | |
Matrix | operator() (const Indexmatrix &vecrow, const Indexmatrix &veccol) const |
returns a new submatrix as indexed by vecrow and veccol, A(i,j)=(*this)(vecrow(i),veccol(j)) for 0<=i<vecrow.dim(), 0<=j<veccol.dim() | |
Matrix | operator() (const Indexmatrix &A) const |
returns a new matrix B of the same shape as A with B(i,j)=(*this)(A(i),A(j)) for 0<=i<A.rowdim(), 0<=j<A.coldim() | |
Real & | operator[] (Integer i) |
returns reference to element (i) of the matrix if regarded as vector of stacked columns [element (irowdim, i/rowdim)] | |
Real | operator[] (Integer i) const |
returns value of element (i) of the matrix if regarded as vector of stacked columns [element (irowdim, i/rowdim)] | |
Matrix | col (Integer i) const |
returns column i copied to a new matrix | |
Matrix | row (Integer i) const |
returns row i copied to a new matrix | |
Matrix | cols (const Indexmatrix &vec) const |
returns a matrix of size this->rowdim() x vec.dim(), with column i a copy of column vec(i) of *this | |
Matrix | rows (const Indexmatrix &vec) const |
returns a matrix of size vec.dim() x this->coldim(), with row i a copy of row vec(i) of *this | |
Matrix & | swap_rowsij (Integer i, Integer j) |
row i of this matrix is swapped with row j | |
Matrix & | swap_colsij (Integer i, Integer j) |
column i of this matrix is swapped with column j | |
Matrix & | pivot_permute_rows (const Indexmatrix &piv, bool inverse=false) |
for i=0 to rowdim row i of this matrix is swapped with row piv(i); for inverse =true the inverse permutation is generated | |
Matrix & | pivot_permute_cols (const Indexmatrix &piv, bool inverse=false) |
for j=0 to coldim column j of this matrix is swapped with column piv(j); for inverse =true the inverse permutation is generated | |
Matrix & | triu (Integer d=0) |
keeps everything above and including diagonal d, everything below is set to zero, returns *this | |
Matrix & | tril (Integer d=0) |
keeps everything below and including diagonal d, everything above is set to zero, returns *this | |
Matrix & | subassign (const Indexmatrix &vecrow, const Indexmatrix &veccol, const Matrix &A) |
assigns A to a submatrix of *this, (*this)(vecrow(i),veccol(j))=A(i,j) for 0<=i<vecrow.dim(), 0<=j<veccol.dim() | |
Matrix & | subassign (const Indexmatrix &vec, const Matrix &A) |
assigns vector A to a subvector of *this, (*this)(vec(i))=A(i) for 0<=i<vec.dim(), *this, vec, and A may be rectangular matrices, their dimesions are not changed, returns *this | |
Matrix & | delete_rows (const Indexmatrix &ind, bool sorted_increasingly=false) |
all rows indexed by vector ind are deleted, no row should appear twice in ind, remaining rows are moved up keeping their order, returns *this | |
Matrix & | delete_cols (const Indexmatrix &ind, bool sorted_increasingly=false) |
all colmuns indexed by vector ind are deleted, no column should appear twice in ind, remaining columns are moved left keeping their order, returns *this | |
Matrix & | insert_row (Integer i, const Matrix &v) |
insert the row vector v before row i, 0<=i<= row dimension, for i==row dimension the row is appended below; appending to a 0x0 matrix is allowed, returns *this | |
Matrix & | insert_col (Integer i, const Matrix &v) |
insert a column before column i, 0<=i<= column dimension, for i==column dimension the column is appended at the right; appending to a 0x0 matrix is allowed, returns *this | |
Matrix & | reduce_length (Integer n) |
(*this) is set to a column vector of length min{max{0,n},dim()}; usually used to truncate a vector, returns *this | |
Matrix & | concat_right (const Matrix &A, int Atrans=0) |
concats matrix A (or its tranpose) to the right of *this, A or *this may be the 0x0 matrix initally, returns *this | |
Matrix & | concat_below (const Matrix &A) |
concats matrix A to the bottom of *this, A or *this may be the 0x0 matrix initally, returns *this | |
Matrix & | concat_right (Real d) |
concat value d at the bottom of *this, *this must be a column vector or the 0x0 matrix, returns *this | |
Matrix & | concat_below (Real d) |
concat value d at the right of *this, *this must be a row vector or the 0x0 matrix, returns *this | |
Matrix & | enlarge_right (Integer addnc) |
enlarge the matrix by addnc>=0 columns without intializaton of the new columns, returns *this (marked as not initialized if nr>0) | |
Matrix & | enlarge_below (Integer addnr) |
enlarge the matrix by addnr>=0 rows without intializaton of the new rows, returns *this (marked as not initialized if addnr>0) | |
Matrix & | enlarge_right (Integer addnc, Real d) |
enlarge the matrix by addnc>=0 columns intializing the new columns by value d, returns *this | |
Matrix & | enlarge_below (Integer addnr, Real d) |
enlarge the matrix by addnr>=0 rows intializing the new rows by value d, returns *this | |
Matrix & | enlarge_right (Integer addnc, const Real *dp, Real d=1.) |
enlarge the matrix by addnc>=0 columns intializing the new columns by the values pointed to by dp times d, returns *this | |
Matrix & | enlarge_below (Integer addnr, const Real *dp, Real d=1.) |
enlarge the matrix by addnr>=0 rows intializing the new rows by the values pointed to by dp times d, returns *this | |
Real * | get_store () |
returns the current address of the internal value array; use cautiously, do not use delete! | |
const Real * | get_store () const |
returns the current address of the internal value array; use cautiously! | |
BLAS-like Routines (Members) | |
Matrix & | xeya (const Matrix &A, Real d=1., int atrans=0) |
sets *this=d*A where A may be transposed and returns *this | |
Matrix & | xpeya (const Matrix &A, Real d=1.) |
sets *this+=d*A and returns *this | |
Matrix & | xeya (const Indexmatrix &A, Real d=1.) |
sets *this=d*A and returns *this | |
Matrix & | xpeya (const Indexmatrix &A, Real d=1.) |
sets *this+=d*A and returns *this | |
Usual Arithmetic Operators (Members) | |
Matrix & | operator= (const Matrix &A) |
Matrix & | operator*= (const Matrix &s) |
Matrix & | operator+= (const Matrix &v) |
Matrix & | operator-= (const Matrix &v) |
Matrix & | operator%= (const Matrix &A) |
ATTENTION: this is redefined as the Hadamard product, (*this)(i,j)=(*this)(i,j)*A(i,j) for all i,j. | |
Matrix & | operator/= (const Matrix &A) |
ATTENTION: this is redefined to act componentwise without checking for zeros, (*this)(i,j)=(*this)(i,j)/A(i,j) for all i,j. | |
Matrix | operator- () const |
Matrix & | operator*= (Real d) |
Matrix & | operator/= (Real d) |
ATTENTION: d is NOT checked for 0. | |
Matrix & | operator+= (Real d) |
sets (*this)(i,j)+=d for all i,j | |
Matrix & | operator-= (Real d) |
sets (*this)(i,j)-=d for all i,j | |
Matrix & | transpose () |
transposes itself (cheap for vectors, expensive for matrices) | |
Connections to other Classes (Members) | |
Matrix & | xeya (const Symmatrix &A, Real d=1.) |
sets *this=d*A and returns *this | |
Matrix & | xpeya (const Symmatrix &A, Real d=1.) |
sets *this+=d*A and returns *this | |
Matrix & | operator= (const Symmatrix &S) |
Matrix & | operator*= (const Symmatrix &S) |
Matrix & | operator+= (const Symmatrix &S) |
Matrix & | operator-= (const Symmatrix &S) |
Matrix & | xeya (const Sparsesym &A, Real d=1.) |
sets *this=d*A and returns *this | |
Matrix & | xpeya (const Sparsesym &A, Real d=1.) |
sets *this+=d*A and returns *this | |
Matrix & | operator= (const Sparsesym &) |
Matrix & | operator*= (const Sparsesym &S) |
Matrix & | operator+= (const Sparsesym &S) |
Matrix & | operator-= (const Sparsesym &S) |
Matrix & | xeya (const Sparsemat &A, Real d=1.) |
sets *this=d*A and returns *this | |
Matrix & | xpeya (const Sparsemat &A, Real d=1.) |
sets *this+=d*A and returns *this | |
Matrix & | operator= (const Sparsemat &A) |
Matrix & | operator*= (const Sparsemat &A) |
Matrix & | operator+= (const Sparsemat &A) |
Matrix & | operator-= (const Sparsemat &A) |
Elementwise Operations (Members) | |
Matrix & | rand (Integer nr, Integer nc, CH_Tools::GB_rand *random_generator=0) |
resize *this to an nr x nc matrix and assign to (i,j) a random number uniformly from [0,1] for all i,j | |
Matrix & | rand_normal (Integer nr, Integer nc, Real mean=0., Real variance=1., int generator_type=0) |
resize *this to an nr x nc matrix and assign to (i,j) a random number from the normal distribution with given mean and variance (generators: 0 std, 1 mt, 2 mt64) | |
Matrix & | shuffle (CH_Tools::GB_rand *random_generator=0) |
shuffle the elements randomly (does not change dimensions) | |
Matrix & | inv (void) |
sets (*this)(i,j)=1./(*this)(i,j) for all i,j and returns *this | |
Matrix & | sqrt (void) |
sets (*this)(i,j)=sqrt((*this)(i,j)) for all i,j and returns *this | |
Matrix & | sqr (void) |
sets (*this)(i,j)=sqr((*this)(i,j)) for all i,j and returns *this | |
Matrix & | sign (Real tol=1e-12) |
sets (*this)(i,j)=sign((*this)(i,j),tol) for all i,j using sign(double,double) and returns *this | |
Matrix & | floor (void) |
sets (*this)(i,j)=floor((*this)(i,j)) for all i,j and returns *this | |
Matrix & | ceil (void) |
sets (*this)(i,j)=ceil((*this)(i,j)) for all i,j and returns *this | |
Matrix & | rint (void) |
sets (*this)(i,j)=rint((*this)(i,j)) for all i,j and returns *this | |
Matrix & | round (void) |
sets (*this)(i,j)=round((*this)(i,j)) for all i,j and returns *this | |
Matrix & | abs (void) |
sets (*this)(i,j)=abs((*this)(i,j)) for all i,j and returns *this | |
Numerical Methods (Members) | |
bool | contains_nan () |
returns true if any matrix element returns true on std::isnan | |
Matrix & | scale_rows (const Matrix &vec) |
scales each row i of (this) by vec(i), i.e., (*this)=diag(vec)(*this), and returns (*this) | |
Matrix & | scale_cols (const Matrix &vec) |
scales each column i of (*this) by vec(i), i.e., (*this)=(*this)*diag(vec), and returns (*this) | |
int | triu_solve (Matrix &rhs, Real tol=1e-10) |
solves (*this)*x=rhs for x by back substitution regarding (*this) as an upper triangle matrix and stores x in rhs. Returns 0 on success, otherwise i+1 if abs(*this)(i,i)<tol and the remaining row of rhs is nonzero. | |
int | tril_solve (Matrix &rhs, Real tol=1e-10) |
solves (*this)*x=rhs for x by forward substitution regarding (*this) as an upper triangle matrix and stores x in rhs. Returns 0 on success, otherwise i+1 if abs(*this)(i,i)<tol and the reduced row of rhs is nonzero. | |
int | QR_factor (Real tol=1e-10) |
computes a Householder QR_factorization overwriting (*this); returns 0 on success, otherwise column index +1 if the norm of the column is below tol when reaching it. | |
int | QR_factor (Matrix &Q, Real tol=1e-10) |
computes a Householder QR_factorization computing the Q matrix explicitly and setting (*this)=R; it always returns 0 | |
int | QR_factor (Matrix &Q, Matrix &R, Real tol) const |
computes a Householder QR_factorization computing matrices Q and R explicitly and leaving (*this) unchanged; it always returns 0 | |
int | QR_factor (Indexmatrix &piv, Real tol=1e-10) |
computes a Householder QR_factorization with pivoting and overwriting (*this); the pivoting permutation is stored in piv; returns the rank | |
int | QR_factor_relpiv (Indexmatrix &piv, Real tol=1e-10) |
computes a Householder QR_factorization with pivoting on those with tolerable norm and tolerable reduction in norm and overwriting (*this); the pivoting permutation is stored in piv; returns the rank | |
int | QR_factor (Matrix &Q, Indexmatrix &piv, Real tol=1e-10) |
Computes a Householder QR_factorization with pivoting computing the Q matrix explicitly and setting (*this)=R; the pivoting permutation is stored in piv; returns the rank. | |
int | QR_factor (Matrix &Q, Matrix &R, Indexmatrix &piv, Real tol) const |
computes a Householder QR_factorization with pivoting computing matrices Q and R explicitly and leaving (*this) unchanged; the pivoting permutation is stored in piv; returns the rank | |
int | Qt_times (Matrix &A, Integer r) const |
computes A=transpose(Q)*A, assuming a housholder Q is coded in the first r columns of the lower triangle of (*this); it always returns 0 | |
int | Q_times (Matrix &A, Integer r) const |
computes A=Q*A, assuming a housholder Q is coded in the first r columns of the lower triangle of (*this); it always returns 0 | |
int | times_Q (Matrix &A, Integer r) const |
computes A=A*Q, assuming a housholder Q is coded in the first r columns of the lower triangle of (*this); it always returns 0 | |
int | QR_solve (Matrix &rhs, Real tol=1e-10) |
solves (*this)*x=rhs by factorizing and overwriting (*this); rhs is overwritten with the solution. Returns 0 on success, otherwise i+1 if in the backsolve abs(*this)(i,i)<tol and the reduced row of rhs is nonzero. More... | |
int | QR_concat_right (const Matrix &A, Indexmatrix &piv, Integer r, Real tol=1e-10) |
int | ls (Matrix &rhs, Real tol) |
computes a least squares solution by QR_solve, overwriting (*this). rhs is overwritten with the solution. In fact, the full code is return this->QR_solve(rhs,tol); | |
int | nnls (Matrix &rhs, Matrix *dual=0, Real tol=1e-10) const |
computes a nonnegative least squares solution; rhs is overwritten by the solution; if dual!=0, the dual variables are stored there; returns 0 on success, 1 on failure More... | |
Find (Members) | |
Indexmatrix | find (Real tol=1e-10) const |
returns an Indexmatrix ind so that (*this)(ind(i)) 0<=i<ind.dim() runs through all nonzero elements | |
Indexmatrix | find_number (Real num=0., Real tol=1e-10) const |
returns an Indexmatrix ind so that (*this)(ind(i)) 0<=i<ind.dim() runs through all elements of value num | |
Input, Output (Members) | |
void | display (std::ostream &out, int precision=0, int width=0, int screenwidth=0) const |
displays a matrix in a pretty way for bounded screen widths; for variables of value zero default values are used. More... | |
void | mfile_output (std::ostream &out, int precision=16, int width=0) const |
outputs a matrix A in the format "[ A(0,1) ... A(0,nc-1)\n ... A(nr-1,nc-1)];\n" so that it can be read e.g. by octave as an m-file More... | |
Private Types | |
enum | AuxTask { none, qr } |
the auxilliary storage aux_m may be used for different pruposes which are listed here; this is used to indicate the latest purpose in aux_task More... | |
Private Member Functions | |
void | init_to_zero () |
initialize the matrix to a 0x0 matrix without storage | |
Private Attributes | |
Integer | mem_dim |
amount of memory currently allocated | |
Integer | nr |
number of rows | |
Integer | nc |
number of columns | |
Real * | m |
pointer to store, order is columnwise (a11,a21,...,anr1,a12,a22,.....) | |
AuxTask | aux_task |
the current purpose of the auxilliary information | |
Integer | aux_dim |
amount of memory currently allocated in auxm | |
Real * | aux_m |
pointer to supplementary store for repeated use in some routines like qr | |
bool | is_init |
flag whether memory is initialized, it is only used if CONICBUNDLE_DEBUG is defined | |
Static Private Attributes | |
static const Mtype | mtype |
used for MatrixError templates (runtime type information was not yet existing) | |
Friends | |
class | Indexmatrix |
class | Symmatrix |
class | Sparsemat |
class | Sparsesym |
Indexing and Submatrices (Friends) | |
Matrix | diag (const Matrix &A) |
returns a column vector v consisting of the elements v(i)=(*this)(i,i), 0<=i<min(row dimension,column dimension) | |
Matrix | triu (const Matrix &A, Integer i) |
retuns a matrix that keeps the upper triangle of A starting with diagonal d, i.e., (i,j)=A(i,j) for 0<=i<row dimension, max(0,i+d)<=j<column dimension, and sets (i,j)=0 otherwise | |
Matrix | tril (const Matrix &A, Integer i) |
retuns a matrix that keeps the lower triangle of A starting with diagonal d, i.e., (i,j)=A(i,j) for 0<=i<row dimension, 0<=j<min(i+d+1,column dimension), and sets (i,j)=0 otherwise | |
Matrix | concat_right (const Matrix &A, const Matrix &B) |
returns a new matrix [A, B], i.e., it concats matrices A and B rowwise; A or B may be a 0x0 matrix | |
Matrix | concat_below (const Matrix &A, const Matrix &B) |
returns a bew matrix [A; B], i.e., it concats matrices A and B columnwise; A or B may be a 0x0 matrix | |
void | swap (Matrix &A, Matrix &B) |
swap the content of the two matrices A and B (involves no copying) | |
BLAS-like Routines (Friends) | |
Matrix & | xbpeya (Matrix &x, const Matrix &y, Real alpha, Real beta, int ytrans) |
returns x= alpha*y+beta*x, where y may be transposed (ytrans=1); if beta==0. then x is initialized to the correct size | |
Matrix & | xeyapzb (Matrix &x, const Matrix &y, const Matrix &z, Real alpha, Real beta) |
returns x= alpha*y+beta*z; x is initialized to the correct size | |
Matrix & | genmult (const Matrix &A, const Matrix &B, Matrix &C, Real alpha, Real beta, int atrans, int btrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size | |
Usual Arithmetic Operators (Friends) | |
Matrix | operator* (const Matrix &A, const Matrix &B) |
returns Matrix equal to A*B | |
Matrix | operator+ (const Matrix &A, const Matrix &B) |
returns Matrix equal to A+B | |
Matrix | operator- (const Matrix &A, const Matrix &B) |
returns Matrix equal to A-B | |
Matrix | operator% (const Matrix &A, const Matrix &B) |
ATTENTION: this is redefined as the Hadamard product, C(i,j)=A(i,j)*B(i,j) for all i,j. | |
Matrix | operator/ (const Matrix &A, const Matrix &B) |
ATTENTION: this is redefined to act componentwise without checking for zeros, C(i,j)=A(i,j)/B(i,j) for all i,j. | |
Matrix | operator* (const Matrix &A, Real d) |
returns Matrix equal to A*d | |
Matrix | operator* (Real d, const Matrix &A) |
returns Matrix equal to d*A | |
Matrix | operator/ (const Matrix &A, Real d) |
ATTENTION: d is NOT checked for 0. | |
Matrix | operator+ (const Matrix &A, Real d) |
returns (i,j)=A(i,j)+d for all i,j | |
Matrix | operator+ (Real d, const Matrix &A) |
returns (i,j)=A(i,j)+d for all i,j | |
Matrix | operator- (const Matrix &A, Real d) |
returns (i,j)=A(i,j)-d for all i,j | |
Matrix | operator- (Real d, const Matrix &A) |
returns (i,j)=d-A(i,j) for all i,j | |
Matrix | transpose (const Matrix &A) |
returns a transposed matrix of A | |
Connections to other Classes (Friends) | |
std::vector< double > & | assign (std::vector< double > &vec, const Matrix &A) |
interpret A as a vector and copy it to a std::vector<double> which is also returned | |
Matrix & | genmult (const Symmatrix &A, const Matrix &B, Matrix &C, Real alpha, Real beta, int btrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size More... | |
Matrix & | genmult (const Matrix &A, const Symmatrix &B, Matrix &C, Real alpha, Real beta, int atrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size More... | |
Matrix & | genmult (const Sparsesym &A, const Matrix &B, Matrix &C, Real alpha, Real beta, int btrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size | |
Matrix & | genmult (const Matrix &A, const Sparsesym &B, Matrix &C, Real alpha, Real beta, int atrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size | |
Matrix & | genmult (const Sparsemat &A, const Matrix &B, Matrix &C, Real alpha, Real beta, int atrans, int btrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size | |
Matrix & | genmult (const Sparsemat &A, const Matrix &B, Integer colB, Matrix &C, Integer colC, Real alpha, Real beta, int atrans, int btrans) |
returns C.col(colC)=beta*C.col(colC)+alpha*A*B.col(colB), where A and B may be transposed first; C must not be equal to A and B; if beta==0. then C is initialized, but the size of C must be correct already | |
Matrix & | genmult (const Matrix &A, const Sparsemat &B, Matrix &C, Real alpha, Real beta, int atrans, int btrans) |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size | |
Elementwise Operations (Friends) | |
Matrix | rand (Integer nr, Integer nc, CH_Tools::GB_rand *random_generator) |
return a nr x nc matrix with (i,j) assigned a random number uniformly from [0,1] for all i,j | |
Matrix | inv (const Matrix &A) |
returns a matrix with elements (i,j)=1./((*this)(i,j)) for all i,j; ATTENTION: no check for division by zero | |
Matrix | sqrt (const Matrix &A) |
returns a matrix with elements (i,j)=sqrt((*this)(i,j)) for all i,j | |
Matrix | sqr (const Matrix &A) |
returns a matrix with elements (i,j)=sqr((*this)(i,j)) for all i,j | |
Matrix | sign (const Matrix &A, Real tol) |
returns a matrix with elements (i,j)=sign((*this)(i,j)) for all i,j using sign(double,double) | |
Matrix | floor (const Matrix &A) |
returns a matrix with elements (i,j)=floor((*this)(i,j)) for all i,j | |
Matrix | ceil (const Matrix &A) |
returns a matrix with elements (i,j)=ceil((*this)(i,j)) for all i,j | |
Matrix | rint (const Matrix &A) |
returns a matrix with elements (i,j)=rint((*this)(i,j)) for all i,j | |
Matrix | round (const Matrix &A) |
returns a matrix with elements (i,j)=round((*this)(i,j)) for all i,j | |
Matrix | abs (const Matrix &A) |
returns a matrix with elements (i,j)=abs((*this)(i,j)) for all i,j | |
Numerical Methods (Friends) | |
Real | trace (const Matrix &A) |
returns the sum of the diagonal elements A(i,i) over all i | |
Real | ip (const Matrix &A, const Matrix &B) |
returns the usual inner product of A and B, i.e., the sum of A(i,j)*B(i,j) over all i,j | |
Real | ip_min_max (const Matrix &A, const Matrix &B, Real &minval, Real &maxval) |
returns in addition to the usual inner product of A and B the minimum and maximum value of A(i,j)*B(i,j) over all (i,j) | |
Real | colip (const Matrix &A, Integer j, const Matrix *scaling) |
returns the squared Frobenius norm of column j of A, i.e., the sum of A(i,j)*A(i,j) over all i with possibly (if scaling!=0) each term i multiplied by (*scaling)(i) | |
Real | rowip (const Matrix &A, Integer i, const Matrix *scaling) |
returns the squared Frobenius norm of row i of A, i.e., the sum of A(i,j)*A(i,j) over all j with possibly (if scaling!=0) each term j multiplied by (*scaling)(j) | |
Matrix | colsip (const Matrix &A) |
returns the column vector of the squared Frobenius norm of all columns j of A, i.e., the sum of A(i,j)*A(i,j) over all i for each j | |
Matrix | rowsip (const Matrix &A) |
returns the row vector of the squared Frobenius norm of all rowd i of A, i.e., the sum of A(i,j)*A(i,j) over all j for each i | |
Real | norm2 (const Matrix &A) |
returns the Frobenius norm of A, i.e., the square root of the sum of A(i,j)*A(i,j) over all i,j | |
Real | normDsquared (const Matrix &A, const Matrix &d, int atrans, int dinv) |
returns trace(A^TDA)=|A|^2_D with D=Diag(d). A may be transposed, D may be inverted but there is no check for division by zero | |
Matrix | sumrows (const Matrix &A) |
returns a row vector holding the sum over all rows, i.e., (1 1 ... 1)*A | |
Matrix | sumcols (const Matrix &A) |
returns a column vector holding the sum over all columns, i.e., A*(1 1 ... 1)^T | |
Real | sum (const Matrix &A) |
returns the sum over all elements of A, i.e., (1 1 ... 1)*A*(1 1 ... 1)^T | |
Matrix | house (const Matrix &A, Integer i, Integer j, Real tol) |
returns the Householder vector of size A.rowdim() for the subcolumn A(i:A.rowdim(),j) | |
int | rowhouse (Matrix &A, const Matrix &v, Integer i, Integer j) |
Housholder pre-multiplication of A with Householder vector v; the first nonzero of v is index i, the multplication is applied to all columns of A with index >=j; always returns 0. | |
int | colhouse (Matrix &A, const Matrix &v, Integer i, Integer j) |
Housholder post-multiplication of A with Householder vector v; the first nonzero of v is index i, the multplication is applied to all rows of A with index >=j; always returns 0. | |
int | QR_factor (const Matrix &A, Matrix &Q, Matrix &R, Real tol) |
computes a Householder QR factorization of A and outputs Q and R leaving A unchanged; always returns 0 | |
int | QR_factor (const Matrix &A, Matrix &Q, Matrix &R, Indexmatrix &piv, Real tol) |
computes a Householder QR factorization of A with pivating. It outputs Q, R, and the pivoting permuation in piv; returns the rank of A | |
Comparisons, Max, Min, Sort, Find (Friends) | |
Matrix | operator< (const Matrix &A, const Matrix &B) |
returns a matrix having elements (i,j)=Real(A(i,j)<B(i,j)) for all i,j | |
Matrix | operator> (const Matrix &A, const Matrix &B) |
returns a matrix having elements (i,j)=Real(A(i,j)>B(i,j)) for all i,j | |
Matrix | operator<= (const Matrix &A, const Matrix &B) |
returns a matrix having elements (i,j)=Real(A(i,j)<=B(i,j)) for all i,j | |
Matrix | operator>= (const Matrix &A, const Matrix &B) |
returns a matrix having elements (i,j)=Real(A(i,j)>=B(i,j)) for all i,j | |
Matrix | operator== (const Matrix &A, const Matrix &B) |
returns a matrix having elements (i,j)=Real(A(i,j)==B(i,j)) for all i,j | |
Matrix | operator!= (const Matrix &A, const Matrix &B) |
returns a matrix having elements (i,j)=Real(A(i,j)!=B(i,j)) for all i,j | |
Matrix | operator< (const Matrix &A, Real d) |
returns a matrix having elements (i,j)=Real(A(i,j)<d) for all i,j | |
Matrix | operator> (const Matrix &A, Real d) |
returns a matrix having elements (i,j)=Real(A(i,j)>d) for all i,j | |
Matrix | operator<= (const Matrix &A, Real d) |
returns a matrix having elements (i,j)=Real(A(i,j)<=d) for all i,j | |
Matrix | operator>= (const Matrix &A, Real d) |
returns a matrix having elements (i,j)=Real(A(i,j)>=d) for all i,j | |
Matrix | operator== (const Matrix &A, Real d) |
returns a matrix having elements (i,j)=Real(A(i,j)==d) for all i,j | |
Matrix | operator!= (const Matrix &A, Real d) |
returns a matrix having elements (i,j)=Real(A(i,j)!=d) for all i,j | |
Matrix | operator< (Real d, const Matrix &A) |
returns a matrix having elements (i,j)=Real(d<A(i,j)) for all i,j | |
Matrix | operator> (Real d, const Matrix &A) |
returns a matrix having elements (i,j)=Real(d>A(i,j)) for all i,j | |
Matrix | operator<= (Real d, const Matrix &A) |
returns a matrix having elements (i,j)=Real(d<=A(i,j)) for all i,j | |
Matrix | operator>= (Real d, const Matrix &A) |
returns a matrix having elements (i,j)=Real(d>=A(i,j)) for all i,j | |
Matrix | operator== (Real d, const Matrix &A) |
returns a matrix having elements (i,j)=Real(d==A(i,j)) for all i,j | |
Matrix | operator!= (Real d, const Matrix &A) |
returns a matrix having elements (i,j)=Real(d!=A(i,j)) for all i,j | |
bool | equal (const Matrix &A, const Matrix &B) |
returns true if both matrices have the same size and the same elements | |
Matrix | minrows (const Matrix &A) |
returns a row vector holding in each column the minimum over all rows in this column | |
Matrix | mincols (const Matrix &A) |
returns a column vector holding in each row the minimum over all columns in this row | |
Real | min (const Matrix &A, Integer *iindex, Integer *jindex) |
returns the minimum value over all elements of the matrix | |
Matrix | maxrows (const Matrix &A) |
returns a row vector holding in each column the maximum over all rows in this column | |
Matrix | maxcols (const Matrix &A) |
returns a column vector holding in each row the maximum over all columns in this row | |
Real | max (const Matrix &A, Integer *iindex, Integer *jindex) |
returns the maximum value over all elements of the matrix | |
Indexmatrix | sortindex (const Matrix &vec, bool nondecreasing) |
returns an Indexmatrix ind so that vec(ind(0))<=vec(ind(1))<=...<=vec(ind(vec.dim()-1)) (vec may be rectangular, set nondecreasing=false for opposite order) | |
void | sortindex (const Matrix &vec, Indexmatrix &ind, bool nondecreasing) |
sets ind so that vec(ind(0))<=vec(ind(1))<=...<=vec(ind(vec.dim()-1)) (vec may be rectangular, set nondecreasing=false for opposite order) | |
Indexmatrix | find (const Matrix &A, Real tol) |
returns an Indexmatrix ind so that A(ind(i)) 0<=i<ind.dim() runs through all nonzero elements with abs(A(j))>tol | |
Indexmatrix | find_number (const Matrix &A, Real num, Real tol) |
returns an Indexmatrix ind so that A(ind(i)) 0<=i<ind.dim() runs through all elements of A having value num, i.e., abs(A(j)-num)<tol | |
Input, Output (Friends) | |
std::ostream & | operator<< (std::ostream &o, const Matrix &v) |
output format (nr and nc are Integer values, all others Real values): nr nc \n A(1,1) A(1,2) ... A(1,nc) \n A(2,1) ... A(nr,nc) \n | |
std::istream & | operator>> (std::istream &i, Matrix &v) |
input format (nr and nc are Integer values, all others Real values): nr nc \n A(1,1) A(1,2) ... A(1,nc) \n A(2,1) ... A(nr,nc) \n | |
Additional Inherited Members | |
Protected Member Functions inherited from CH_Matrix_Classes::Memarrayuser | |
Memarrayuser () | |
if memarray is NULL, then a new Memarray is generated. In any case the number of users of the Memarray is incremented | |
virtual | ~Memarrayuser () |
the number of users is decremented and the Memarray memory manager is destructed, if the number is zero. | |
Static Protected Attributes inherited from CH_Matrix_Classes::Memarrayuser | |
static Memarray * | memarray |
pointer to common memory manager for all Memarrayusers, instantiated in memarray.cxx | |
Matrix class for real values of type Real
Internally a matrix of size nr x nc is stored in a one dimensional array of Real variables, the elements are arranged in columnwise order (a11,a21,...,anr1,a12,a22,...).
Any matrix element can be indexed by (i,j), which internally refers to m[i+j*nr], or directly by the one dimensional index (i+j*nr). The latter view directly corresponds to the vec() operator often used in the linear algebra literature, i.e., the matrix is transformed to a vector by stacking the columns on top of each other.
|
private |
the auxilliary storage aux_m may be used for different pruposes which are listed here; this is used to indicate the latest purpose in aux_task
Enumerator | |
---|---|
none | aux_m currently not in use |
qr | aux_m is now / was last used by a QR factorization routine which stores the beta multipliers to be used in QR_mult type routines |
generate a matrix of size nr x nc but WITHOUT initializing the memory
If initializing the memory externally and CONICBUNDLE_DEBUG is defined, please use set_init() via matrix.set_init(true) in order to avoid warnings concerning improper initialization
References init_to_zero(), and newsize().
void CH_Matrix_Classes::Matrix::display | ( | std::ostream & | out, |
int | precision = 0 , |
||
int | width = 0 , |
||
int | screenwidth = 0 |
||
) | const |
displays a matrix in a pretty way for bounded screen widths; for variables of value zero default values are used.
out | output stream |
precision | number of most significant digits, default=4 |
width | field width, default = precision+6 |
screenwidth | maximum number of characters in one output line, default = 80 |
Referenced by ConicBundle::CMgramdense::display(), ConicBundle::CMlowranksd::display(), and ConicBundle::CMlowrankdd::display().
void CH_Matrix_Classes::Matrix::mfile_output | ( | std::ostream & | out, |
int | precision = 16 , |
||
int | width = 0 |
||
) | const |
outputs a matrix A in the format "[ A(0,1) ... A(0,nc-1)\n ... A(nr-1,nc-1)];\n" so that it can be read e.g. by octave as an m-file
out | output stream |
precision | number of most significant digits, default=16 |
width | field width, default = precision+6 |
resize the matrix to nr x nc elements but WITHOUT initializing the memory
If initializing the memory externally and CONICBUNDLE_DEBUG is defined, please use set_init() via matrix.set_init(true) in order to avoid warnings concerning improper initialization
Referenced by ConicBundle::UQPSolver::init_size(), and Matrix().
computes a nonnegative least squares solution; rhs is overwritten by the solution; if dual!=0, the dual variables are stored there; returns 0 on success, 1 on failure
Computes the least squares solution of min ||Ax-b|| s.t. x >=0;
The KKT system A'*A*x - A'*b - l = 0; x >=0, l>=0, x'*l=0 is solved solved by interior point method with QR-solution of the extended system.
The current implementation is based on [P. Matsoms, "Sparse Linear Least Squares Problems in Optimization", Comput. Opt. and Appl., 7, 89-110 (1997)] but is only a quick and rather sloppy implementation of it ...
int CH_Matrix_Classes::Matrix::QR_concat_right | ( | const Matrix & | A, |
Indexmatrix & | piv, | ||
Integer | r, | ||
Real | tol = 1e-10 |
||
) |
extend the current Householder QR-factorization stored in this by appending and factorizing the columns of A yielding a new QR_factorization
[in] | A | contains the addtionial columns to be factorized |
[in,out] | piv | contains on input, the permution vector returned by QR_factor with pivoting, and on output the new entire permutation |
[in] | r | gives the initial rank of *this as returned by QR_factor with pviaton |
[in] | tol | gives the tolerance for regarding a vector as having norm zero |
solves (*this)*x=rhs by factorizing and overwriting (*this); rhs is overwritten with the solution. Returns 0 on success, otherwise i+1 if in the backsolve abs(*this)(i,i)<tol and the reduced row of rhs is nonzero.
To avoid overwriting (*this), use the appropriate version of QR_factor and triu_solve.
Referenced by ls().
|
friend |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size
returns C=beta*C+alpha*A*B, where B may be transposed; C must not be equal to B; if beta==0. then C is initialized to the correct size
|
friend |
returns C=beta*C+alpha*A*B, where A and B may be transposed; C must not be equal to A and B; if beta==0. then C is initialized to the correct size
returns C=beta*C+alpha*A*B, where A may be transposed; C must not be equal to A; if beta==0. then C is initialized to the correct size