ConicBundle
|
If in Lagrangean relaxation primal solutions are in the form of a real vector or, more generally a matrix, then an approximate primal solution can be generated by supplying primal information of this form for each epsilon subgradient within ConicBundle::MatrixFunctionOracle::evaluate(). More...
#include <MatrixCBSolver.hxx>
Public Member Functions | |
PrimalMatrix () | |
empty matrix | |
PrimalMatrix (CH_Matrix_Classes::Integer nr, CH_Matrix_Classes::Integer nc) | |
generate a matrix of size nr x nc but WITHOUT initializing the memory More... | |
PrimalMatrix (CH_Matrix_Classes::Integer r, CH_Matrix_Classes::Integer c, CH_Matrix_Classes::Real d) | |
generate a matrix of size nr x nc initializing all elements to the value d | |
PrimalMatrix (const PrimalMatrix &pm) | |
copy constructor, *this=pm | |
PrimalMatrix (const CH_Matrix_Classes::Matrix &pm) | |
copy constructor, *this=pm | |
PrimalMatrix & | operator= (const CH_Matrix_Classes::Matrix &pd) |
copy operator, *this=pm | |
PrimalData * | clone_primal_data () const |
produces a new PrimalMatrix that is a copy of itself; the caller has to delete the returned object at some point | |
int | aggregate_primal_data (const PrimalData &it, double itsfactor) |
multiply *this Matrix with myfactor and add itsfactor*it (it must dynamic_cast to a PrimalMatrix) | |
int | scale_primal_data (double myfactor) |
multiply/scale *this with a nonnegative myfactor | |
Public Member Functions inherited from CH_Matrix_Classes::Matrix | |
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... | |
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 | |
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 | |
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! | |
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 | |
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) | |
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) |
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 | |
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... | |
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 | |
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... | |
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 | |
If in Lagrangean relaxation primal solutions are in the form of a real vector or, more generally a matrix, then an approximate primal solution can be generated by supplying primal information of this form for each epsilon subgradient within ConicBundle::MatrixFunctionOracle::evaluate().
In many applications, e.g. in Lagrangean relaxation, the convex minimization problem arises as the dual of a convex primal maximization problem. In this case one is typically interested in obtaining a primal approximate solution in addition to the dual solution. Under reasonable conditions this is possible if the primal solutions that give rise to the subgradients are aggregated along with the subgradients within the bundle algorithm. If the primal data can be represented as a CH_Matrix_Classes::Matrix then the user has to supply in the oracle for each sugradient the corresponding primal data in a PrimalMatrix and the algorithm will do the rest. Observe that a PrimalMatrix can be used exactly in the same way as a CH_Matrix_Classes::Matrix and that they are assignable among each other.
The primal data has to be supplied within ConicBundle::MatrixFunctionOracle::Evaluate() and can be retrieved via the methods ConicBundle::MatrixCBSolver::get_approximate_primal() and ConicBundle::MatrixCBSolver::get_center_primal()
|
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