ConicBundle
|
interface for interior point variable and routines specific to primal dual complementarity conditions of a positive semidefinite cone More...
#include <PSCIPBlock.hxx>
Public Member Functions | |
virtual void | clear (CH_Matrix_Classes::Integer dim=0) |
reset all point information to zero for dimension dim, the rest to zero | |
PSCIPBlock (CH_Matrix_Classes::Integer dim=0, CBout *cb=0, int cbinc=-1) | |
default constructor, alos allows to set the order | |
~PSCIPBlock () | |
destructor | |
virtual CH_Matrix_Classes::Integer | get_vecdim () const |
return the length of the svec reprensentation | |
virtual int | center_x (CH_Matrix_Classes::Real val, bool add=false) |
set x to value*"one" to x, or if add==true, add value*"one" to x | |
virtual int | center_z (CH_Matrix_Classes::Real val, bool add=false) |
set z to value*"one" to z, or if add==true, add value*"one" to z | |
virtual int | set_x (const CH_Matrix_Classes::Matrix &vec, CH_Matrix_Classes::Integer startindex, CH_Matrix_Classes::Real &add_center_value) |
set x to the values of vec[startindex+0,+1 ...,+(vecdim-1)] and return in add_center_value a value>=0 that needs to be added to make it feasible | |
virtual int | set_z (const CH_Matrix_Classes::Matrix &vec, CH_Matrix_Classes::Integer startindex, CH_Matrix_Classes::Real &add_center_value) |
set z to the values of vec[startindex+0,+1 ...,+(vecdim-1)] and add sufficient center to make z feasible, return this value>=0 in added_center_value | |
virtual int | vecgetsax (CH_Matrix_Classes::Matrix &vec, CH_Matrix_Classes::Integer startindex, CH_Matrix_Classes::Real a=1., bool add=false) |
on vec[startindex+0,+1 ...,+(vecdim-1)] put or add a * x into vec for a real number a | |
virtual int | vecgetsaz (CH_Matrix_Classes::Matrix &vec, CH_Matrix_Classes::Integer startindex, CH_Matrix_Classes::Real a=1., bool add=false) |
on vec[startindex+0,+1 ...,+(vecdim-1)] put or add a * z into vec for a real number a | |
virtual int | get_mu_info (CH_Matrix_Classes::Integer &mudim, CH_Matrix_Classes::Real &tr_xz, CH_Matrix_Classes::Real &tr_xdzpdxz, CH_Matrix_Classes::Real &tr_dxdz, CH_Matrix_Classes::Real &min_xz, CH_Matrix_Classes::Real &max_xz) const |
add dimensions of the primal-dual pairs to mudim and add the "trace" (the inner product with center) of the respective primal-dual pair products for the current step; update the min and max values of x_i*z_i | |
virtual int | get_nbh_info (CH_Matrix_Classes::Integer mudim, CH_Matrix_Classes::Real tr_xz, CH_Matrix_Classes::Real tr_xdzpdxz, CH_Matrix_Classes::Real tr_dxdz, CH_Matrix_Classes::Real nbh_ubnd, CH_Matrix_Classes::Real &alpha, CH_Matrix_Classes::Real &max_nbh, CH_Matrix_Classes::Real &nrmsqr_xz, CH_Matrix_Classes::Real &nrmsqr_xdzpdxz, CH_Matrix_Classes::Real &nrmsqr_dxdz, CH_Matrix_Classes::Real &ip_xz_xdzpdxz, CH_Matrix_Classes::Real &ip_xz_dxdz, CH_Matrix_Classes::Real &ip_dxdz_xdzpdxz) const |
for limiting the stepsize with respect to the neighborhood this information about norms and inner products of x(.)*z-tr_xz-tr_xz/mudim(.*)1, x.()*dz+dx(.)*z-tr_xdzpdxz/mudim(.*)1, and dx(.)*dz-tr_dxdz/mudim(.)*1 is required, each block adds its contribution to the numbers | |
virtual int | linesearch (CH_Matrix_Classes::Real &alpha) const |
if necessary, reduce alpha to the biggest value so that feasibility is maintained with this step size | |
virtual int | add_muxinv (CH_Matrix_Classes::Matrix &rhs, CH_Matrix_Classes::Integer startindex, CH_Matrix_Classes::Real rhsmu, CH_Matrix_Classes::Real rhscorr, bool minus=false) |
compute the complementarity_rhs=rhsmu*xi-rhscorr*xi*dx*dz (wihtout "-z") for mu=rhsmu and for corrector for factor rhscorr>0., store this and add it to rhs | |
virtual int | set_dx (const CH_Matrix_Classes::Matrix &rhs, CH_Matrix_Classes::Integer startindex) |
extract dx from rhs at startindex and compute at the same time dz (=-sys dx -z +complentarity_rhs); | |
virtual int | set_dx_xizsolverhs (const CH_Matrix_Classes::Matrix &rhs, CH_Matrix_Classes::Integer startindex) |
compute dx=sysinv*rhs and at the same time dz (=-rhs -z +complentarity_rhs); | |
virtual int | apply_xizinv (CH_Matrix_Classes::Matrix &rhs, CH_Matrix_Classes::Integer startindex, bool minus=false) |
compute sysinv*rhs into rhs, possibly with a negative sign | |
virtual int | apply_xiz (CH_Matrix_Classes::Matrix &rhs, CH_Matrix_Classes::Integer startindex, bool minus=false) |
compute sys*rhs into rhs, possibly with a negative sign | |
virtual int | do_step (CH_Matrix_Classes::Real alpha) |
move to (x+alpha*dx, z+alpha*dz) | |
virtual int | add_AxizinvAt (const CH_Matrix_Classes::Matrix &A, CH_Matrix_Classes::Symmatrix &globalsys, bool minus=false, bool Atrans=false) |
add the Schur complement to a big system matrix | |
virtual int | add_xiz (CH_Matrix_Classes::Symmatrix &globalsys, CH_Matrix_Classes::Integer startindex, bool minus=false) |
add (or subract if minus==true) the system matrix to a big system matrix starting at startindex | |
virtual int | get_vecx (CH_Matrix_Classes::Matrix &vecx, CH_Matrix_Classes::Integer startindex) |
return the vector form of x | |
virtual int | get_vecz (CH_Matrix_Classes::Matrix &vecz, CH_Matrix_Classes::Integer startindex) |
return the vector form of z | |
virtual int | get_vecdx (CH_Matrix_Classes::Matrix &vecdx, CH_Matrix_Classes::Integer startindex) |
return the vector form of dx, 1 if not available | |
virtual int | get_vecdz (CH_Matrix_Classes::Matrix &vecdz, CH_Matrix_Classes::Integer startindex) |
return the vector form of dz, 1 if not available | |
Public Member Functions inherited from ConicBundle::InteriorPointBlock | |
InteriorPointBlock (CBout *cb=0, int cbinc=-1) | |
default constructor | |
virtual | ~InteriorPointBlock () |
virtual destructor (implemented in InteriorPointBundleBlock.cxx) | |
Public Member Functions inherited from ConicBundle::CBout | |
virtual void | set_out (std::ostream *out=0, int print_level=1) |
Specifies the output level (out==NULL: no output at all, out!=NULL and level=0: errors and warnings, level>0 increasingly detailed information) More... | |
virtual void | set_cbout (const CBout *cb, int incr=-1) |
Specifies the output level relative to the given CBout class. More... | |
void | clear_cbout () |
reset to default settings (out=0,print_level=1) | |
CBout (const CBout *cb=0, int incr=-1) | |
calls set_cbout | |
CBout (std::ostream *outp, int pl=1) | |
initialize correspondingly | |
CBout (const CBout &cb, int incr=0) | |
copy constructor | |
virtual bool | cb_out (int pl=-1) const |
Returns true if out!=0 and (pl<print_level), pl<0 should be used for WARNINGS and ERRORS only, pl==0 for usual output. | |
std::ostream & | get_out () const |
If cb_out() returned true, this returns the output stream, but it will abort if called with out==0. | |
std::ostream * | get_out_ptr () const |
returns the pointer to the output stream | |
int | get_print_level () const |
returns the print_level | |
virtual int | mfile_data (std::ostream &out) const |
writes problem data to the given outstream | |
Protected Member Functions | |
void | point_changed () |
clear variables that are no longer valid for the current point | |
int | compute_NTscaling () |
computes the NT scaling information for building the system matrix | |
int | compute_Weig_Wvec () |
compute Weig and Wvec with Weig eigenvalues of W and W=Wvec*Wvec' where Wvec = P*Lambda^{.5} for W = P*Lambda*P' | |
Protected Member Functions inherited from ConicBundle::InteriorPointBlock | |
int | pol_le_zero_step (CH_Matrix_Classes::Real &stepsize, CH_Matrix_Classes::Real q0, CH_Matrix_Classes::Real q1, CH_Matrix_Classes::Real q2, CH_Matrix_Classes::Real q3, CH_Matrix_Classes::Real q4, CH_Matrix_Classes::Real abseps=1e-10) const |
for a polynomial q0+alpha*(q1+alpha*(q2+alpha*(q3+alpha*q4))) with q0<=0 find the maximum alpha <= stepsize, that keeps its value <=0. and put stepsize=min(stepsize,alpha) | |
int | minimize_pol_step (CH_Matrix_Classes::Real &stepsize, CH_Matrix_Classes::Real q0, CH_Matrix_Classes::Real q1, CH_Matrix_Classes::Real q2, CH_Matrix_Classes::Real q3, CH_Matrix_Classes::Real q4, CH_Matrix_Classes::Real abseps=1e-10) const |
find the minimizing alpha of a polynomial q0+alpha*(q1+alpha*(q2+alpha*(q3+alpha*q4))) within 0<=alpha<=stepsize and put stepsize=alpha, where q1 is assumed to be <0 or be ==0 becoming negative for small alpa>=0. | |
int | control_nbh_step (CH_Matrix_Classes::Real &stepsize, CH_Matrix_Classes::Real &max_nbh, CH_Matrix_Classes::Real nbh_ubnd, CH_Matrix_Classes::Real mu_xz, CH_Matrix_Classes::Real mu_xdzpdxz, CH_Matrix_Classes::Real mu_dxdz, CH_Matrix_Classes::Real nrmsqr_xz, CH_Matrix_Classes::Real nrmsqr_xdzpdxz, CH_Matrix_Classes::Real nrmsqr_dxdz, CH_Matrix_Classes::Real ip_xz_xdzpdxz, CH_Matrix_Classes::Real ip_xz_dxdz, CH_Matrix_Classes::Real ip_dxdz_xdzpdxz) const |
find a stepsize so that outside the nbh_ubnd some progress is made towards it and inside the upper bound the step size is chosen as large as possible while staying within | |
Protected Attributes | |
CH_Matrix_Classes::Integer | rowdim |
order of the symmetric matrices of the cone, | |
CH_Matrix_Classes::Integer | vecdim |
dimension svec representation | |
CH_Matrix_Classes::Symmatrix | X |
"primal" point X | |
CH_Matrix_Classes::Symmatrix | Z |
"dual" point Z | |
CH_Matrix_Classes::Symmatrix | dX |
current step for X | |
CH_Matrix_Classes::Symmatrix | dZ |
current step for Z | |
CH_Matrix_Classes::Symmatrix | W |
NT scaling matrix. | |
CH_Matrix_Classes::Symmatrix | Winv |
inverse of NT scaling matrix | |
CH_Matrix_Classes::Matrix | G |
Gram Matrix of NT scaling matrix. | |
CH_Matrix_Classes::Matrix | Ginv |
Gram Matrix of inverse NT scaling matrix. | |
CH_Matrix_Classes::Matrix | D |
Diagonal for NT scaling. | |
CH_Matrix_Classes::Matrix | Weig |
eigenvalues of W (if Weig.rowdim()==rowdim(), lazy evaulation) | |
CH_Matrix_Classes::Matrix | Wvec |
W=Wvec*Wvec' where Wvec = P*Lambda^{.5} for W = P*Lambda*P' (if Weig.rowdim()==rowdim(), lazy evaulation) | |
CH_Matrix_Classes::Symmatrix | compl_rhs |
rhs used in solving the complementarity line | |
CH_Matrix_Classes::Real | last_rhs_mu |
the last mu used in rhs computations | |
CH_Matrix_Classes::Real | mu |
in a step mu gets the value of last_rhs_mu | |
CH_Matrix_Classes::Real | old_mu |
in a step old_mu gets the value of mu before this gets last_rhs_mu | |
CH_Matrix_Classes::Real | last_alpha |
last alpha used in do_step() | |
CH_Matrix_Classes::Symmatrix | oldX |
point before x | |
CH_Matrix_Classes::Symmatrix | oldZ |
point before z | |
CH_Matrix_Classes::Matrix | lamX |
eigenvalues of X corresponding to PX (if computed) | |
CH_Matrix_Classes::Matrix | PX |
eigenvectors of X corresponding to lamX (if computed) | |
CH_Matrix_Classes::Symmatrix | tmpsym |
temporary matrices for reuse | |
CH_Matrix_Classes::Symmatrix | tmpsym2 |
temporary matrices for reuse | |
CH_Matrix_Classes::Matrix | tmpvec |
temporary matrices for reuse | |
CH_Matrix_Classes::Matrix | tmpmat |
temporary matrices for reuse | |
interface for interior point variable and routines specific to primal dual complementarity conditions of a positive semidefinite cone
The class implements Nesterov-Todd scaling along Todd, Toh and Tuetuencue.