CH_Matrix_Classes::Matrix Class Reference
[Matrix (dense, real, m by n)]

Matrix class for real values of type Real More...

#include <matrix.hxx>

Inheritance diagram for CH_Matrix_Classes::Matrix:

CH_Matrix_Classes::Memarrayuser ConicBundle::PrimalMatrix

List of all members.

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
 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)
 generate a matrix of size nr x nc initializing the elements from the (one dimensional) array dp with increment incr
 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)
Matrixinit (const Matrix &A, Real d=1., int atrans=0)
 initialize to *this=A*d where A may be transposed
Matrixinit (const Indexmatrix &A, Real d=1.)
 initialize to *this=A*d
Matrixinit (const Sparsemat &A, Real d=1.)
 initialize to *this=A*d
Matrixinit (const Symmatrix &S, Real d=1.)
 initialize to *this=A*d
Matrixinit (const Sparsesym &, Real d=1.)
 initialize to *this=A*d
Matrixinit (const Realrange &)
 initialize *this to a column vector holding the elements of Realrange
Matrixinit (Integer nr, Integer nc, Real d)
 intialize *this to a matrix of size nr x nc initializing all elements to the value d
Matrixinit (Integer nr, Integer nc, const Real *dp, Integer incr=1)
 generate a matrix of size nr x nc initializing the elements from the (one dimensional) array dp with increment incr
Matrixinit (const std::vector< Real > &vec)
 use std::vector<Real> to initialize this to a column vector of the same size and content
void newsize (Integer nr, Integer nc)
 resize the matrix to nr x nc elements but WITHOUT initializing the memory
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)
Realoperator() (Integer i, Integer j)
 returns reference to element (i,j) of the matrix (rowindex i, columnindex j)
Realoperator() (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()
Realoperator[] (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->rowdim(), with row i a copy of row vec(i) of *this
Matrixtriu (Integer d=0)
 keeps everything above and including diagonal d, everything below is set to zero, returns *this
Matrixtril (Integer d=0)
 keeps everything below and including diagonal d, everything above is set to zero, returns *this
Matrixsubassign (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()
Matrixsubassign (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
Matrixdelete_rows (const Indexmatrix &ind)
 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
Matrixdelete_cols (const Indexmatrix &ind)
 all colmuns indexed by vector ind are deleted, no column should appear twice in ind, remaining columns are moved up keeping their order, returns *this
Matrixinsert_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
Matrixinsert_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
Matrixreduce_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
Matrixconcat_right (const Matrix &A)
 concats matrix A to the right of *this, A or *this may be the 0x0 matrix initally, returns *this
Matrixconcat_below (const Matrix &A)
 concats matrix A to the bottom of *this, A or *this may be the 0x0 matrix initally, returns *this
Matrixconcat_right (Real d)
 concat value d at the bottom of *this, *this must be a column vector or the 0x0 matrix, returns *this
Matrixconcat_below (Real d)
 concat value d at the right of *this, *this must be a row vector or the 0x0 matrix, returns *this
Realget_store ()
 returns the current address of the internal value array; use cautiously, do not use delete!
const Realget_store () const
 returns the current address of the internal value array; use cautiously!
BLAS-like Routines (Members)
Matrixxeya (const Matrix &A, Real d=1., int atrans=0)
 sets *this=d*A where A may be transposed and returns *this
Matrixxpeya (const Matrix &A, Real d=1.)
 sets *this+=d*A and returns *this
Matrixxeya (const Indexmatrix &A, Real d=1.)
 sets *this=d*A and returns *this
Matrixxpeya (const Indexmatrix &A, Real d=1.)
 sets *this+=d*A and returns *this
Usual Arithmetic Operators (Members)
Matrixoperator= (const Matrix &A)
Matrixoperator*= (const Matrix &s)
Matrixoperator+= (const Matrix &v)
Matrixoperator-= (const Matrix &v)
Matrixoperator%= (const Matrix &A)
 ATTENTION: this is redefined as the Hadamard product, (*this)(i,j)=(*this)(i,j)*A(i,j) for all i,j.
Matrixoperator/= (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
Matrixoperator*= (Real d)
Matrixoperator/= (Real d)
 ATTENTION: d is NOT checked for 0.
Matrixoperator+= (Real d)
 sets (*this)(i,j)+=d for all i,j
Matrixoperator-= (Real d)
 sets (*this)(i,j)-=d for all i,j
Matrixtranspose ()
 transposes itself (cheap for vectors, expensive for matrices)
Connections to other Classes (Members)
Matrixxeya (const Symmatrix &A, Real d=1.)
 sets *this=d*A and returns *this
Matrixxpeya (const Symmatrix &A, Real d=1.)
 sets *this+=d*A and returns *this
Matrixoperator= (const Symmatrix &S)
Matrixoperator*= (const Symmatrix &S)
Matrixoperator+= (const Symmatrix &S)
Matrixoperator-= (const Symmatrix &S)
Matrixxeya (const Sparsesym &A, Real d=1.)
 sets *this=d*A and returns *this
Matrixxpeya (const Sparsesym &A, Real d=1.)
 sets *this+=d*A and returns *this
Matrixoperator= (const Sparsesym &)
Matrixoperator*= (const Sparsesym &S)
Matrixoperator+= (const Sparsesym &S)
Matrixoperator-= (const Sparsesym &S)
Matrixxeya (const Sparsemat &A, Real d=1.)
 sets *this=d*A and returns *this
Matrixxpeya (const Sparsemat &A, Real d=1.)
 sets *this+=d*A and returns *this
Matrixoperator= (const Sparsemat &A)
Matrixoperator*= (const Sparsemat &A)
Matrixoperator+= (const Sparsemat &A)
Matrixoperator-= (const Sparsemat &A)
Elementwise Operations (Members)
Matrixrand (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
Matrixshuffle (CH_Tools::GB_rand *random_generator=0)
 shuffle the elements randomly (does not change dimensions)
Matrixinv (void)
 sets (*this)(i,j)=1./(*this)(i,j) for all i,j and returns *this
Matrixsqrt (void)
 sets (*this)(i,j)=sqrt((*this)(i,j)) for all i,j and returns *this
Matrixsign (Real tol=1e-12)
 sets (*this)(i,j)=sign((*this)(i,j),tol) for all i,j using sign(double,double) and returns *this
Matrixfloor (void)
 sets (*this)(i,j)=floor((*this)(i,j)) for all i,j and returns *this
Matrixceil (void)
 sets (*this)(i,j)=ceil((*this)(i,j)) for all i,j and returns *this
Matrixrint (void)
 sets (*this)(i,j)=rint((*this)(i,j)) for all i,j and returns *this
Matrixround (void)
 sets (*this)(i,j)=round((*this)(i,j)) for all i,j and returns *this
Matrixabs (void)
 sets (*this)(i,j)=abs((*this)(i,j)) for all i,j and returns *this
Numerical Methods (Members)
Matrixscale_rows (const Matrix &vec)
 scales each row i of (*this) by vec(i), i.e., (*this)=diag(vec)*(*this), and returns (*this)
Matrixscale_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); currently it always returns 0
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=1e-10) 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 (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=1e-10) 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.
int QR_concat_right (const Matrix &A, Indexmatrix &piv, Integer r, Real tol=1e-10)
int ls (Matrix &rhs, Real tol=1e-10)
 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
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
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.

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
Realm
 pointer to store, order is columnwise (a11,a21,...,anr1,a12,a22,.....)
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=0)
 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=0)
 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)
Matrixxbpeya (Matrix &x, const Matrix &y, Real alpha=1., Real beta=0., int ytrans=0)
 returns x= alpha*y+beta*x, where y may be transposed (ytrans=1); if beta==0. then x is initialized to the correct size
Matrixxeyapzb (Matrix &x, const Matrix &y, const Matrix &z, Real alpha=1., Real beta=1.)
 returns x= alpha*y+beta*z; x is initialized to the correct size
Matrixgenmult (const Matrix &A, const Matrix &B, Matrix &C, Real alpha=1., Real beta=0., int atrans=0, int btrans=0)
 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)
Matrix operator+ (const Matrix &A, const Matrix &B)
Matrix operator- (const Matrix &A, const Matrix &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)
Matrix operator* (Real d, const Matrix &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)
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
Matrixgenmult (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
Matrixgenmult (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
Matrixgenmult (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
Matrixgenmult (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
Matrixgenmult (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
Matrixgenmult (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=0)
 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)=abs((*this)(i,j)) for all i,j
Matrix sqrt (const Matrix &A)
 returns a matrix with elements (i,j)=abs((*this)(i,j)) for all i,j
Matrix sign (const Matrix &A, Real tol=1e-12)
 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 colip (const Matrix &A, Integer j)
 returns the squared Frobenius norm of column j of A, i.e., the sum of A(i,j)*A(i,j) over all i
Real rowip (const Matrix &A, Integer i)
 returns the squared Frobenius norm of row i of A, i.e., the sum of A(i,j)*A(i,j) over all j
Matrix colsip (const Matrix &A)
 returns the row 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 column 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=0, int dinv=0)
 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=0, Integer j=0, Real tol=1e-10)
 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=0, Integer j=0)
 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=0, Integer j=0)
 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=1e-10)
 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=1e-10)
 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
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=0, Integer *jindex=0)
 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=0, Integer *jindex=0)
 returns the maximum value over all elements of the matrix
Indexmatrix sortindex (const Matrix &vec)
 returns an Indexmatrix ind so that vec(ind(0))<=vec(ind(1))<=...<=vec(ind(vec.dim()-1)) (vec may be rectangular)
void sortindex (const Matrix &vec, Indexmatrix &ind)
 sets ind so that vec(ind(0))<=vec(ind(1))<=...<=vec(ind(vec.dim()-1)) (vec may be rectangular)
Indexmatrix find (const Matrix &A, Real tol=1e-10)
 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=0., Real tol=1e-10)
 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


Detailed Description

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.


Constructor & Destructor Documentation

CH_Matrix_Classes::Matrix::Matrix ( Integer  nr,
Integer  nc 
) [inline]

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().


Member Function Documentation

void CH_Matrix_Classes::Matrix::newsize ( Integer  nr,
Integer  nc 
)

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 init(), and Matrix().

int CH_Matrix_Classes::Matrix::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.

To avoid overwriting (*this), use the appropriate version of QR_factor and triu_solve.

Referenced by ls().

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

Parameters:
[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
Returns:
the rank of the new QR-facotrization

int CH_Matrix_Classes::Matrix::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

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 ...

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.

Parameters:
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


The documentation for this class was generated from the following files:

Generated on Tue May 3 16:52:54 2011 for ConicBundle by  doxygen 1.5.6