NDOSolver / FiOracle
Interfaces and Solvers for NonDifferentiable Optimization
Loading...
Searching...
No Matches
OSIMPSolver Class Reference

#include <OSIMPSolver.h>

Inheritance diagram for OSIMPSolver:
MPSolver

Public Types

Public Types
enum  StabFun { unset = 0 , none , boxstep , quadratic }
 Public enum for handling the Bundle stabilization term [see below].
 
enum  OsiAlg {
  kAuto = 0 , kPrim , kDual , kNet ,
  kBar , kSif , kCon
}
 
enum  OsiRed { rNo = 0 , rPrim , rDual , rBoth }
 
- Public Types inherited from MPSolver
enum  MPParam { kMaxTme = 0 , kOptEps , kFsbEps , kLastMPParam }
 
enum  MPStatus {
  kOK = 0 , kUnbndd , kUnfsbl , kStppd ,
  kError
}
 

Public Member Functions

Constructor
 OSIMPSolver (std::istream *iStrm=nullptr)
 
Other initializations

@ -------------------------------------------------------------------—

void SetDim (cIndex MxBSz=0, FiOracle *Oracle=nullptr, const bool UsAvSt=false) override
 
void Sett (cHpNum tt=1) override
 
void SetPar (const int wp, cHpNum value) override
 
void SetThreads (int nthreads) override
 
void SetLowerBound (cHpNum LwBnd=- Inf< HpNum >(), cIndex wFi=Inf< Index >()) override
 
void SetMPLog (ostream *outs=0, const char lvl=0) override
 
void SetAlgo (const OsiAlg algo=kAuto, const OsiRed reduc=rNo)
 
OSIMPSolver-specific methods
void SetOsi (OsiSolverInterface *osi=nullptr)
 Provides OSIMPSolver with the actual OsiSolverInterface.
 
OsiSolverInterface * GetOsi (void) const
 gives back the current OsiSolverInterface
 
void SetStabType (const StabFun sf=none)
 Sets the stabilizing term in the Master Problem.
 
Solving the problem
MPStatus SolveMP (void) override
 
Reading results
HpNum ReadFiBLambda (cIndex wFi=Inf< Index >()) override
 
HpNum ReadDt (cHpNum tt=1) override
 
HpNum ReadSigma (cIndex wFi=Inf< Index >()) override
 
HpNum ReadDStart (cHpNum tt=1) override
 
cLMRow Readd (bool Fulld=false) override
 
void ReadZ (LMRow tz, cIndex_Set &I, Index &D, cIndex wFi=Inf< Index >()) override
 
cHpRow ReadMult (cIndex_Set &I, Index &D, cIndex wFi=Inf< Index >(), const bool IncldCnst=true) override
 
HpNum ReadLBMult (cIndex wFi=Inf< Index >()) override
 
cHpRow ReadDualEasy (cIndex wFi) override
 returns the dual variables of the constraints of an easy component
 
cHpRow ReadReducedCostsEasy (cIndex wFi) override
 returns the reduced costs of the variables of an easy component
 
HpNum ReadGid (cIndex Nm=Inf< Index >()) override
 
void MakeLambda1 (cHpRow Lmbd, HpRow Lmbd1, cHpNum Tau) override
 
void SensitAnals (HpNum &lp, HpNum &cp) override
 
Reading the data of the problem
Index BSize (cIndex wFi=Inf< Index >()) override
 returns the current number of items (of either type) in the bundle.
 
Index BCSize (cIndex wFi=Inf< Index >()) override
 returns the current number of constraints in the bundle.
 
Index MaxName (cIndex wFi=Inf< Index >()) override
 
Index WComponent (cIndex i) override
 
bool IsSubG (cIndex i) override
 
Index NumNNVars (void) override
 
Index NumBxdVars (void) override
 
bool IsNN (cIndex i) override
 
void CheckIdentical (const bool Chk=true) override
 
cHpRow ReadLinErr (void) override
 
HpNum ReadLowerBound (cIndex wFi=Inf< Index >()) override
 
HpNum EpsilonD (void) override
 
bool FiBLambdaIsExact (cIndex wFi)
 
Adding / removing / changing data
SgRow GetItem (cIndex wFi=Inf< Index >()) override
 
void SetItemBse (cIndex_Set SGBse=nullptr, cIndex SGBDm=0) override
 
Index CheckSubG (cHpNum DFi, cHpNum Tau, HpNum &Ai, HpNum &ScPri) override
 
Index CheckCnst (HpNum &Ai, HpNum &ScPri, cHpRow CrrPnt) override
 
bool ChangesMPSol (void) override
 
void SetItem (cIndex Nm=Inf< Index >()) override
 
void SubstItem (cIndex Nm) override
 
void RmvItem (cIndex i) override
 
void RmvItems (void) override
 
void SetActvSt (cIndex_Set AVrs=nullptr, cIndex AVDm=0) override
 
void AddActvSt (cIndex_Set Addd, cIndex AdDm, cIndex_Set AVrs) override
 
void RmvActvSt (cIndex_Set Rmvd, cIndex RmDm, cIndex_Set AVrs) override
 
void AddVars (cIndex NNwVrs) override
 
void RmvVars (cIndex_Set whch=nullptr, Index hwmny=0) override
 
void ChgAlfa (cHpRow DeltaAlfa) override
 
void ChgAlfa (cHpRow NewAlfa, cIndex wFi) override
 
void ChgAlfa (cIndex i, cHpNum Ai) override
 
void ChangeCurrPoint (cLMRow DLambda, cHpRow DFi) override
 
void ChangeCurrPoint (cHpNum Tau, cHpRow DFi) override
 
void ChgSubG (cIndex strt, Index stp, cIndex wFi) override
 
void ChgCosts (Index wFi, cLMRow Lambda) override
 
void ChgRLHS (Index wFi) override
 
void ChgLUBD (Index wFi) override
 
Destructor
virtual ~OSIMPSolver ()
 
- Public Member Functions inherited from MPSolver
virtual void SetMPLog (std::ostream *outs=0, const char lvl=0)
 
virtual void SetMPTime (const bool TimeIt=true)
 
 MPSolver (void)
 
void MPTime (double &t_us, double &t_ss)
 
double MPTime (void)
 
Index GetMaxBSize (void)
 returns the maximum bundle size
 
Index GetMaxSGLen (void)
 returns the maximum subgradient length
 
Index GetCrrSGLen (void)
 returns the current subgradient length
 
Index GetNrFi (void)
 returns the number of components of Fi()
 
virtual ~MPSolver ()
 destructor: deletes the timer. it is virtual, as it must be
 

Private Member Functions

void cleanup (void)
 
void tUpdatePrices (cIndex strt=0, Index stp=Inf< Index >())
 
void ptUpdatePrices (cIndex strt=0, Index stp=Inf< Index >())
 
void ptUpdatePricesInPlace (void)
 
void UpdateRhoCol (void)
 
void switchToQP (void)
 
void activate (Index i)
 
void deactivate (Index i)
 
void resizeHP (Index i)
 
void resizeI (Index i)
 
void CheckDS (void)
 
Index CheckBCopy (void)
 

Private Attributes

OsiSolverInterface * osiSlvr
 pointer to the OsiSolver
 
FiOracleFIO
 pointer to the FiOracle
 
OsiAlg algorithm
 algorithm used in the solver
 
bool useactiveset
 true if the active set are taken in account
 
bool first
 
bool checkID
 
Bool_Vec weasy
 which Fi-component is easy
 
Index NrEasy
 how many easy Fi-components are there
 
StabFun stab
 tell the stabilization term adopted
 
HpNum t
 proximal parameter term
 
double MaxTime
 max running time for each call to SolveMP
 
double OptEps
 optimality tolerance
 
double FsbEps
 feasibility tolerance
 
double * RhoCol
 entries of the column of rho
 
int * RhoColBse
 ... in sparse format: indices ...
 
Index RhoColBDm
 ... and size
 
Index RhoColSgPos
 position in RhoColBse of the "subgradient"
 
SgRow NewItem
 the new item
 
Index NewItemFi
 the function component relative to new item
 
cIndex_Set NewItemBse
 index vector of the new item
 
Index NewItemBDm
 dimension of new item
 
Index MaxNZ
 max dimension of the item
 
bool NewItemisSG
 subgradient vs constraint
 
double NewItemprice
 price of the item
 
double NewItemScPri
 Gid product.
 
Index_Set NSubG
 number of subgradients (per component)
 
Index_Set NConst
 number of constraints (per component)
 
Index_Set comp_row
 row vocabulary
 
Index_Set dict_item
 item vocabulary
 
Index_Set dict_slack
 slack vocabulary
 
Index_Set dict_stab
 stabilization variables vocabulary
 
Index_Set wcomp
 
Index_Set comp_col
 column vocabulary
 
Index item_maxname
 maximum name among current items + 1
 
HpRow tempHP
 temporary vector of HpNum
 
int * tempI
 temporary vector of ints
 
Index tempHP_size
 size of tempHP
 
Index tempI_size
 size of tempI
 
LMRow Upper
 upper bounds on the primal variable
 
LMRow Lower
 lower bounds on the primal variable
 
Index NNVars
 variables with >= 0 constraints
 
Index BxdVars
 variables with >= 0 or <= UB constraints
 
cIndex_Set Aset
 active set structure
 
Index Asetdim
 size of active set
 
HpRow LwrBnds
 values of lower bounds for each component
 
HpRow GiPerd
 
double * csol
 column (primal) solution
 
int csols
 length of csol
 
double * rsol
 row (dual) solution
 
int rsols
 length of rsol
 
double * rcst
 row reduced cost
 
CoinMessageHandler * derhand
 
Index NRCall
 number of calls to SolveMP procedure
 

Additional Inherited Members

- Protected Attributes inherited from MPSolver
Index MaxBSize
 
Index MaxSGLen
 
Index CrrSGLen
 
Index NrFi
 
std::ostream * MPLog
 
char MPLLvl
 
OPTtimersMPt
 

Detailed Description

The OSIMPSolver class implements a Master Problem Solver for Bundle algorithms, using a generic OSISolverInterface object. This class derives from MPSolver, and therefore conforms to the generic interface defined there for Master Problem Solvers [see MPSolver.h].

Member Enumeration Documentation

◆ OsiAlg

enum OsiAlg

Public enum which is used to define the algorithm used in Solve() [see below] by the OSISolverInterface object for solving the Master problem, which may be either LP or QP problem:

  • kAuto: automatic
  • kPrim: Primal simplex
  • kDual: Dual simplex
  • kNet: Network simplex
  • kBar: Barrier
  • kSif: Sifting
  • kCon: Concurrent (Dual, Barrier, and Primal)

◆ OsiRed

enum OsiRed

Public enum which is used to tell wich reductions (if any) are used by the OSISolverInterface object while solving the Master problem:

  • rNo: no reductions
  • rPrim: primal reductions only
  • rDual: dual reductions only
  • rBoth: both primal and dual reduction

Constructor & Destructor Documentation

◆ ~OSIMPSolver()

virtual ~OSIMPSolver ( )
virtual

Note: the destructor does not delete the OSISolverInterface object, this has to be done by whomever created it in the first place.

Member Function Documentation

◆ AddActvSt()

void AddActvSt ( cIndex_Set Addd,
cIndex AdDm,
cIndex_Set AVrs )
overridevirtual

Adds ‘AdDm’ new variables to the current active set, whose names are contained in the first AdDm positions of the vector Addd (ordered in increasing sense and Inf< Index >()-terminated); those names must not already belong to the active set.

AVrs must point to a vector containing the new active set; a pointer to that vector is kept and used from now on as the new active set, replacing the old vector set with the latest call to one of SetActvSt(), AddActvSt() or RmvActvSt(). Thus, the vector pointed by the "old" pointer can be freely modified after the return of AddActvSt(), while the one pointed by AVrs must not be modified until it is in turn substituted by a new one.

Thus, AVrs must usually be different from the "old" pointer; one exception is permitted, however, in the case when all the variables to be added are to be found after all the variables that were present previously. In this case, the caller is authorized to use the old pointer as the new pointer (and checking this is exactly how the MPSolver detects the fact), overwriting the "Inf< Index >()" which previously terminated the active set with the names of the variables to be added before the call to AddActvSt() (i.e., Addd = AVrs + <n. of variables in the old active set>).

This method can be used only if an active set has been previously defined with SetActvSt( <non-0> ) and not yet removed with SetActvSt( 0 ).

Important note: after a call to AddActvSt(), all the solution information corresponding to the last call to SolveMP() is lost, and the results returned by all the query methods are unpredictable. This is true even for CheckSubG(), CheckCnst() and ChangesMPSol(), so adding items to the bundle after a call to this method is not advised.

Implements MPSolver.

◆ AddVars()

void AddVars ( cIndex NNwVrs)
overridevirtual

This method adds ‘NNwVrs’ new variables to the MP, in response to a call to NDOSolver::AddVariables( NNwVrs , IVs ) [see NDOSlvr.h]. The derived classes are required to use GetUC(), GetUB(), GetGi() and (possibly) GetADesc() [see FiOracle.h] to acquire directly from the FiOracle all the information about the newly created variables they need.

The method is used when NDOSolver::AddVariables() is called, i.e., when new variables Lambda[ i ] of the original NDO problem are created. Note that adding variables Lambda[ i ] changes, in principle, the current point, and therefore the linearization errors of the items. For instance, in the Lagrangian case, it is easy to check that the linearization error alfa_{Lambda}[ i ] of the subgradient named ‘i’ with respect to the current point Lambda is simply

alfa_{Lambda}[ i ] = Fi( Lambda )

  • c( x[ i ] ) + Lambda * ( b - A( x[ i ] ) )

where x[ i ] is the solution of the Lagrangian problem corresponding to the subgradient. Thus, adding a new variable Lambda_j changes the linearization error of Lambda_j * ( b_j - A_j( x[ i ] ) ). Hence, the linearization errors of all items change except if Lambda_j == 0. In order to avoid the recomputation of the linearization errors, it is assumed that all the variables d[ i ] of the Primal Master Problem added with AddVars() correspond to variables Lambda[ i ] that are set to 0; in other words, it is assumed that adding the variables does not change the linearization errors (changing them is always possible with ChangeCurrPoint()).

Note that, if an "active set" has been defined [see SetDim() above], then the newly added variables are not added to the current active set; this must* be done explicitly, if needed, with [Set/Add]ActvSt() [see above]. If the active set option has not been selected, then the variables are automatically added to the MP.

Important note: after a call to AddVars(), all the solution information corresponding to the last call to SolveMP() is lost, and the results returned by all the query methods are unpredictable. This is true even for CheckSubG(), CheckCnst() and ChangesMPSol(), so adding items to the bundle after a call to this method is not advised.

Implements MPSolver.

◆ BCSize()

Index BCSize ( cIndex wFi = InfIndex >())
overridevirtual

returns the current number of constraints in the bundle.

Implements MPSolver.

◆ BSize()

Index BSize ( cIndex wFi = InfIndex >())
overridevirtual

returns the current number of items (of either type) in the bundle.

Implements MPSolver.

◆ ChangeCurrPoint() [1/2]

void ChangeCurrPoint ( cHpNum Tau,
cHpRow DFi )
overridevirtual

Equivalent to ChangeCurrPoint( ( Tau / t ) * d* , DFi ), where d* is the optimal solution of the direction finding subproblem (Tau is a relative step along d, relative w.r.t. t). Here, however, no problem with the active set can arise, since d* is nonzero only for variables of the active set.

Implements MPSolver.

◆ ChangeCurrPoint() [2/2]

void ChangeCurrPoint ( cLMRow DLambda,
cHpRow DFi )
overridevirtual

Updates the data of the Master Problem in response to a change of the "current point" Lambda:

DLambda = newLambda - Lambda,

i.e., DLambda is the direction that is followed (with step 1) to go from the old current point Lambda to the new current point newLambda. DLambda is in "dense" format [see MakeLambda1() above], i.e., ‘DLambda[ i ]’ contains the displacement relative to variable ‘i’ for i = 0 .. CrrSGLen - 1. Note that the set of variables with nonzero displacement may not be a subset of the "active set" [see [Set/Add/Rmv]ActvSt() above], i.e., that the entry of the current point for variable ‘i’ may vary even if i is currently forced to be zero by being out of the active set, In order for the Master problem to become "aware" of the change for these variables, the active set has to be changed.

DFi[ k ] must be the change in the value of the function between the two points, i.e.,

DFi[ k ] = Fi[ k ]( newLambda ) - Fi[ k ]( Lambda )

for 1 <= k <= NrFi, while

DFi[ 0 ] = Fi( newLambda ) - Fi( Lambda )

(i.e., the value of the whole function). Note that the values for "easy" components are provided as well, altough the MPSolver in principle knows them already.

The formulae for updating the data are the following. If z[ i ] is a alfa[ i ]-subgradient to the k-th component of Fi in Lambda, then it is a newalfa[ i ]-subgradient (to the k-th component of Fi) in newLambda, where

newalfa[ i ] = alfa[ i ] + DFi[ k ] - DLambda * z[ i ] .

If z[ i ] is a constraint

z[ i ] * ( Lambda + d ) <= RHSi

its right-hand side in the d-variables is RHSi - z[ i ] * Lambda. Hence, for newLambda = Lambda + DLambda the constraint has to become

z[ i ] * ( newLambda + d ) <= RHSi ,

that is, its right-hand side in the d-variables is

RHSi - z[ i ] * newLambda = ( RHSi - z[ i ] * Lambda ) - z[ i ] * Dlambda .

In particular, for box constraints on the Lambda variables [see GetUC() and GetUB() in FiOracle.h], the Master Problem has to contain box constraints on the d variables (assuming for simplicity every variable is box constrained)

LB - Lambda <= d <= UB - Lambda .

Clearly, when the current point is changed to newLambda = Lambda + DLambda these constraints have to become

LB - newLambda <= d <= UB - newLambda ,

that is, ( LB - Lambda ) - Dlambda <= d <= ( UB - Lambda ) - Dlambda.

Implements MPSolver.

◆ ChangesMPSol()

bool ChangesMPSol ( void )
overridevirtual

After that Check[SubG/Cnst]() [see above] has been called, this method can be used to check whether or not the introduction of the new item changes the solution of the MP. If not, the Bundle algorithm may decide not to insert it in the bundle, and however knows if the newly obtained information is enough for obtaining a different d* at the next iteration. This allows to cope with "non-exact functions", i.e., with a FiOracle that is not capable of - or is instructed not to - computing the Fi-values and/or the subgradients exactly. Note that ChangesMPSol() must not be called for the constant subgradient of the linear 0-th component of Fi.

Implements MPSolver.

◆ CheckCnst()

Index CheckCnst ( HpNum & Ai,
HpNum & ScPri,
cHpRow CrrPnt )
overridevirtual

After that SetItemBse() has been called, CheckCnst() [if the item is a constraint, otherwise see CheckSubG() above] must be called to:

  • verify that the new item is not already in the bundle, and
  • compute some information about the item, i.e., the Right Hand Side of the item and the scalar product with d*.

The return value must be the "name" of a constraint in the bundle that is identical (element-wise) to GS, if such an item exists, or Inf< Index >() otherwise. This allows to check for useless copies (but this check is costly and can be deactivated, see CheckIdentical() above). Note, however, that the component to which the constraint belongs is now irrelevant, and therefore that the checks have to be performed across all the components (comprised the 0-th). The reason is that while a constraint nominally belongs to one component (in Lagrangian terms, the corresponding primal unbounded ray is for one particular subproblem) this has no impact on the master problem: unlike subgradients, two identical constraints (columns) corresponding to two different components actually are the same constraint (column) in the master problem, and therefore one is redundant. Thus one can (and therefore should) safely declare a constraint a copy of another even if they belong to different components. Note that this corresponds, from the primal side, to choose between two different primal rays belonging to different subproblems, and therefore it does have an impact on the generated primal (convexified) solution. However, if one would keep the two constraints then it would be the Master Problem solver which would arbitrarily "split the corresponding primal variables among them", so a decision would still be somewhat taken in this respect that is entirely out of control; better, then, to at least save some running time and memory. If this turns out to be a problem for some weird reason, the check can be deactivated with CheckIdentical() and the problem solved.

As for CheckSubG(), the scalar product between the primal optimal solution d* of the MP and SG must be calculated and returned in ScPri.

Ai must contain the Right Hand Side, i.e., the constraint on the original variables is

SG * Lambda <= Ai.

However, the MP is defined over the d variables, that is the constraint in the MP must be

SG * ( CrrPnt + d ) <= Ai

where CrrPnt is the current point of the Bundle algorithm. The format of CrrPnt depends on whether or not the "active set" strategy has been activated with SetDim() [see above]. If the active set is used, then CrrPnt[] is in "sparse" format, i.e., ‘CrrPnt[ i ]’ is the value for the variable with name ‘AVrs[ i ]’ for i = 0 .. AVDm - 1, all the other values being zero. Otherwise, CrrPnt[] is in "dense" format. Upon return, Ai must be updated to the "scaled" Right Hand Side w.r.t. the current point, that is

Ai = Ai - SG * CurrPnt.

The Ai computed by CheckCnst() are typically used in the MP. Note that CheckCnst() must be called before SetItem() for all constraints, comprised* those of the 0-th component (wFi == 0 in GetItem()).

Implements MPSolver.

◆ CheckIdentical()

void CheckIdentical ( const bool Chk = true)
overridevirtual

When a new item is inserted in the bundle, checks can be done to ensure that it is not identical to other items already present, i.e., an useless duplicate [see CheckItem() below]. Since these checks can be costly, they are done only if CheckIdentical( true ) is called; by default, or if CheckIdentical( false ) is called, they are not.

Reimplemented from MPSolver.

◆ CheckSubG()

Index CheckSubG ( cHpNum DFi,
cHpNum Tau,
HpNum & Ai,
HpNum & ScPri )
overridevirtual

After that SetItemBse() has been called, CheckSubG() (if the item is a subgradient, otherwise see CheckCnst() below] must be called to:

  • verify that the new item is not already in the bundle, and
  • compute some information about the item, i.e., the linearization error of the item and the scalar product with d*.

The return value must be the "name" of an item in the bundle that is identical (element-wise) to GS, if such an item exists, or Inf< Index >() otherwise. This allows to check for useless copies (but this check is costly and can be deactivated, see CheckIdentical() above), which is not an unlikely possibility e.g. when Fi is the Lagrangian function of a combinatorial optimization problem (and therefore the set of possible subgradients is finite). Note that subgradients that are element-wise identical but belong to different components are not identical and should not be substituted for one another; in other words, the check only has to be performed with subgradients of the same component as SG* (which, besides being correct, is also likely to be less costly). Also, subgradients should not be checked against constraints*, as they are two entirely different kinds of items.

Also, the scalar product between the primal optimal solution d* of the MP and SG must be calculated and returned in ScPri.

DFi must contain the valueof Fi[ wFi ]( Lambda1 ) - Fi[ wFi ]( Lambda ), where Lambda1 = Lambda + ( Tau / t ) * d*, and must Ai contain its linearization error w.r.t. Lambda1, i.e., SG is an Ai-subgradient of Fi[ wFi ]() in Lambda1. Ai must be updated to the linearization error of SG w.r.t. the current point Lambda, that is

Ai = Ai - ( DFi - ( Tau / t ) * SG * d* ).

A special case, that has a separate treatment, is that of Tau <= 0. This is taken to mean that SG is an Ai-subgradient in Lambda (typically, DFi then is also 0). This happens when adding subgradients for whom the linearization error in Lambda is already known, an important case being that of aggregation, i.e., adding to the bundle the dual optimal solution z* as provided by ReadZ() and ReadSigma(). In such a case, Ai is not modified, and the original one is directly used.

A further specific difference is whether CheckSubG() with Tau <= 0 is called inside the normal loop "compute function value and subgradient at next iterate, put them in the master problem, re-solve" or outside it; in other words, if it is possible that ChangeCurrPoint() is called before the next master problem solution. If not, then the computation of the scalar product (that anyway is not used to update Ai) can be avoided, saving some cost: this is done if Tau == 0. If, instead, ChangeCurrPoint() may be called, then the scalar product need be computed anyway, even if it is not used to update Ai: this is done if Tau < 0.

The Ai computed by CheckSubG() are typically used in the MP. Note that CheckSubG() must be called before SetItem() for all items except the constant subgradient of the "linear part" of Fi (wFi == 0 in GetItem()), where CheckSubG() must not be called.

Implements MPSolver.

◆ ChgAlfa() [1/3]

void ChgAlfa ( cHpRow DeltaAlfa)
overridevirtual

Changes the vector of the linearization errors of the subgradients.

DeltaAlfa[] is a (NrFi + 1)-vector; the linearization error of each subgradient corresponding to a component wFi is increased by DeltaAlfa[ wFi ]. DeltaAlfa[ 0 ] is ignored, and the right hand sides of the constraints are unchanged.

This is the typical correction that ought to be done when "nonexact" functions have to be handled.

Implements MPSolver.

◆ ChgAlfa() [2/3]

void ChgAlfa ( cHpRow NewAlfa,
cIndex wFi )
overridevirtual

Changes the vector of the linearization errors of the subgradients/right hand sides of the constraints.

For every item i in the bundle, the new value for its linearization error/ right hand side is taken from NewAlfa[ i ]. 0 <= wFi <= NrFi means that only the items corresponding to the wFi-th component of Fi are affected by the change. In particular, if wFi == 0 then what changes is the constant term ‘b0’ in the linear 0-th component of Fi; the new value for b0 is to be found in NewAlfa[ MaxBSize ]. If Inf< Index >() > wFi > NrFi, then all the components except the 0-th have changed, while if wFi == Inf< Index >() then all the components, comprised the 0-th, have changed.

Implements MPSolver.

◆ ChgAlfa() [3/3]

void ChgAlfa ( cIndex i,
cHpNum Ai )
overridevirtual

Changes the linearization error of the item with name i in the bundle to the new value Ai.

Implements MPSolver.

◆ ChgCosts()

void ChgCosts ( Index wFi,
cLMRow Lambda )
overridevirtual

This method allows to change the costs of the "easy" component wFi. It is illegal if it is called with wFi not an "easy" component. The method is supposed to call GetBDesc() (with all parameters 0 save that of costs) to retrieve the new costs and change the master problem accordingly.

IMPORTANT NOTE: the number of columns in the "easy" component must not have changed.

An issue is that the costs that need be stored in the objective of the master problem are not the original ones, but rather the reduced ones w.r.t. the stability centre Lambda. Since the MPSolver typically has no direct access to it, the second parameter is provided; it must point to a "dense" representation of the stability centre.

Since not all MPSolver may implement this advanced feature, a default implementation is given throwing exception.

Reimplemented from MPSolver.

◆ ChgLUBD()

void ChgLUBD ( Index wFi)
overridevirtual

This method allows to change the lower and upper bounds of the variables of the "easy" component wFi. It is illegal if it is called with wFi not an "easy" component. The method calls GetBDesc() (with all parameters 0 save these of lower and upper bounds) to retrieve the new lower and upper bounds and change the master problem accordingly.

IMPORTANT NOTE: the number of columns in the "easy" component must not have changed.

IMPORTANT NOTE: to make the implementation simple, the method only works if the "finite nonzero" status of the bounds is preserved by any change. This (in its simplest terms) means that if the bound was finite and nonzero before the change it must have remained so after the change. That is, bounds that are either 0 or -INF/+INF must remain in one of the two statuses (but they can swap, i.e., 0 can become -INF/+INF and vice/versa), and bounds that are finite and nonzero must remain finite and nonzero.

This is rather restrictive but allowing more flexibility would require a substantial rework of the method that is not possible right now.

Failure to comply with these rules may result in an exception being thrown, or, even worse, wrong results being reported with no warning about it.

Reimplemented from MPSolver.

◆ ChgRLHS()

void ChgRLHS ( Index wFi)
overridevirtual

This method allows to change the LHS/RHS of the constraints of the "easy" component wFi. It is illegal if it is called with wFi not an "easy" component. The method calls GetBDesc() (with all parameters 0 save these of the LHS/RHS) to retrieve the new LHS/RHS and change the master problem accordingly.

IMPORTANT NOTE: the number of rows in the "easy" component must not have changed.

IMPORTANT NOTE: to make the implementation simple, the method only works if the "two sided status" of each constraint is preserved by any change. This (in its simplest terms) means that for each constraint, if a LHS/RHS was 0, finite or -INF/+INF before the change it must have remained 0, finite or -INF/+INF; in other words, only finite nonzero LHS/RHS can change. Also, if the LHS and RHS of a row were equal before the change (an equality constraint) they must be equal (to each other, but possibly different from before) after the change.

This is rather restrictive but allowing more flexibility would require a substantial rework of the method that is not possible right now.

Failure to comply with these rules may result in an exception being thrown, or, even worse, wrong results being reported with no warning about it.

Reimplemented from MPSolver.

◆ ChgSubG()

void ChgSubG ( cIndex strt,
Index stp,
cIndex wFi )
overridevirtual

This method allows to change the entries corresponding to variables with names between strt (included) and stp (excluded) for all the items in the bundle corresponding to the wFi-th component of Fi, in response to a call to NDOSolver::ChgSbG( strt , stp, wFi ) [see NDOSlvr.h]. The method has to use GetGi() and (possibly) GetADesc() [see FiOracle.h] to acquire directly from the FiOracle the new information about the changed variables. If 0 <= wFi <= NrFi, the NDO solver is told that only the wFi-th component of Fi has changed; if Inf< Index >() > wFi > NrFi then all the components except the 0-th have changed, while if wFi == Inf< Index >() then all the components, comprised the 0-th, have changed.

Note that Fi-values, and therefore the Alfas, are assumed not to be involved in the changes.

Implements MPSolver.

◆ EpsilonD()

HpNum EpsilonD ( void )
overridevirtual

Must return the threshold used by the solver to decide when a variable is zero; Lambda[ i ] <= EpsilonD() is considered as == 0.

Implements MPSolver.

◆ GetItem()

SgRow GetItem ( cIndex wFi = InfIndex >())
overridevirtual

GetItem( wFi ) must return a pointer SG to a vector of SgNum's where a new item has to be written. wFi tells to which of the components of Fi() [see SetDim() above] the new item belongs, with the usual format:

  • when wFi == 0, the new item is relative to the linear 0-th component of Fi, i.e., it is either its constant (sub)gradient, or a constraint which defines its domain Feas0;
  • when 1 <= wFi <= NrFi, the new item is relative to the wFi-th component of Fi;
  • when Inf< Index >() > wFi > NrFi, the new item is relative to all the components except the 0-th;
  • when wFi == Inf< Index >(), the new item is relative to all the components, comprised* the 0-th.

This information allows the MPSolver to properly dimension the vector to be returned [see GetMaxNZ() in FiOracle.h].

Implements MPSolver.

◆ IsNN()

bool IsNN ( cIndex i)
overridevirtual

This method shall return true if the variable ‘i’ is constrained to be NonNegative, and false otherwise; the standard implementation is fine for MP Solvers not accepting constrained problems.

Reimplemented from MPSolver.

◆ IsSubG()

bool IsSubG ( cIndex i)
overridevirtual

Returns true if item ‘i’ is in the Bundle and it is a subgradient, and false otherwise; hence, an item ‘i’ for which wFi == WComponent( i ) > 0 [see above] and IsSubG( i ) == false is a constraint relative to the wFi-th component of Fi (again, wFi here can be zero).

Implements MPSolver.

◆ MakeLambda1()

void MakeLambda1 ( cHpRow Lmbd,
HpRow Lmbd1,
cHpNum Tau )
overridevirtual

Writes in (the array pointed by) Lmbd1 the vector Lmbd + ( Tau / t ) * d*, where d* is the optimal primal solution of the MP.

The "format" of Lmbd and Lmbd1 depends on whether or not the "active set" is being used [see SetActvSt() below]; if not, then Lmbd is (and Lmbd1 is expected to be) just a "dense" CrrSGLen-vectors. Otherwise, both Lmbd and Labd1 are intended to be "restricted" to the set of "active" variables that are defined by SetActvSt(); that is, only the first k entries of Lmbd/Lmbd1 are significant, where k is the size of the active set as resulting from all the previous calls to *ActvSt() [see below], and Lmbd/Lmbd1[ i ] contain the value of the i-th "active" variable (all the other variables being fixed to zero).

Note that, if an "active set" is defined, this method can compute only the entries of d* that are strictly necessary for the purpose (like Readd() does, see above), thereby possibly saving some time.

Implements MPSolver.

◆ MaxName()

Index MaxName ( cIndex wFi = InfIndex >())
overridevirtual

Returns 1 + the maximum "name" of an item in the bundle if the bundle is nonempty, and 0 otherwise. wFi "restricts" the count to the items corresponding to that particular component, with the usual convention.

Implements MPSolver.

◆ NumBxdVars()

Index NumBxdVars ( void )
overridevirtual

Return the number of variables (<= CrrSGLen) that have either a lower bound constraint or an upper bound constraint, or both. Hence, NumBxdVars() >= NumNNVars(), and the difference gives the number of variables with u.b. but not l.b.. The standard implementation is fine for MP Solvers not accepting constrained problems.

Reimplemented from MPSolver.

◆ NumNNVars()

Index NumNNVars ( void )
overridevirtual

Return the number of variables that are constrained in sign; this must of course be <= CrrSGLen. The standard implementation is fine for MP Solvers not accepting constrained problems.

Reimplemented from MPSolver.

◆ Readd()

cLMRow Readd ( bool Fulld = false)
overridevirtual

Returns a read-only pointer d to a vector containing the primal optimal solution d*; the entry of d* relative to variable ‘i’, if required (see below), can be found in d[ i ] for i in [ 0 .. CrrSGLen - 1 ] (see [Add/Remove/Subst]Var() below for a discussion on variable "names").

If Fulld == true, the solution d[] returned by Readd() is a "full" one, even if an "active set" of variables has been defined [see SetActvSt() below]; that is, even the entries of d* relative to variables that are not "active" are computed in response to a call to Readd( true ). This may be costly, and it is not always necessary: thus, a call to Readd( false ) computes only the entries relative to "active" variables, the contents of d[ i ] for an "inactive" variable ‘i’ being unpredictable.

Note that the returned pointer may point to a temporary data structure rather than to a "persistent" one: the pointer is no more valid (it may even point to no-more-allocated memory) after any call to any other method of the class.

Implements MPSolver.

◆ ReadDStart()

HpNum ReadDStart ( cHpNum tt = 1)
overridevirtual

Returns the value of the dual stabilizing term D*() in the dual optimal solution z* if tt is taken as the proximal parameter.

Implements MPSolver.

◆ ReadDt()

HpNum ReadDt ( cHpNum tt = 1)
overridevirtual

Returns the value of the primal stabilizing term D() in the primal optimal solution d* if tt is taken as the proximal parameter.

Implements MPSolver.

◆ ReadDualEasy()

cHpRow ReadDualEasy ( cIndex wFi)
overridevirtual

returns the dual variables of the constraints of an easy component

In a Bundle algorithm, easy components will never be evaluated since they "live in the Master Problem". However, at termination (or at any point in time) the user may want to get the optimal value of the corresponding primal and dual solutions. ReadMult [see above] does this for primal variables, this method and ReadReducedCostsEasy [see below] does this for dual variables and reduced costs, respectively. That is, this method returns a read-only pointer to a FiOracle::GetBNR( wFi )-vector so that the i-th entry contains the optimal dual variable of the constraint in the i-th row of the matrix returned by FiOracle::GetBDesc( wFi ).

It is of course an error calling this method if wFi is not the "name" of an easy component.

Note
returned pointer may be that of a temporary data structure; hence, it is no more valid (it may even be dangling) after any call to any other method of the class
the returned pointer may be nullptr if the wFi-th easy component has no rows, hence no dual variables
the method is not pure virtual but it has a default implementation doing nothing (returning nullptr) for derived classes that do not handle easy components and therefore have nothing to implement

Reimplemented from MPSolver.

◆ ReadFiBLambda()

HpNum ReadFiBLambda ( cIndex wFi = InfIndex >())
overridevirtual

Returns the value of the model Fi_{B,Lambda} in d*, the optimal primal solution of the MP; this value is often called "v*". This is the value of the translated model

v* = Fi_{B,Lambda}( d* ) = Fi_B( Lambda + d* ) - Fi( Lambda ) ,

hence it is a negative value that indicates the improvement in the Fi-value that the model Fi_B() predicts for a step to Lambda + d*.

When Fi is a decomposable function [see SetDim() above], each of its components has a separate model Fi[ k ]_{B,Lambda}; ReadFiBx() must return:

  • Fi[ wFi ]_{B,Lambda}( d* ), if 0 <= wFi <= NrFi;
  • \sum_{k = 1}^{NrFi} Fi[ k ]_{B,Lambda}( d* ) if Inf< Index >() > wFi > NrFi, i.e., the value of all the models excluding the linear one;
  • the value of the aggregated model if wFi == Inf< Index >().

Note that if there are no subgradients (of a given component) in B, the problem that is solved is different from the standard one, as there is no such thing as a model for that component. Thus the model is temporarily removed from the master problem (like adding a constraint v[ wFi ] >= 0 for the usual representation where v[ wFi ] is the value of the model). In this case, ReadFiBLambda( wFi ) must return Inf< HpNum >().

The above changes if Fi has "easy" components [see FiOracle::GetBNC() ...], because the model of each "easy" component is not translated by value. Hence, if wFi is the index of an easy component then the method returns

Fi[ wFi ]( Lambda + d* )

rather than

Fi[ wFi ]( Lambda + d* ) - Fi[ wFi ]( Lambda ) .

Note that in both cases this is the actual value of the (corresponding component of the) function, because for easy components the model is the actual function (that is, Fi[ wFi ] == Fi_B[ wFi ]). This information is easily extracted out of the solution of the (Dual) Master Problem: if x[ wFi ]* is the optimal solution of the Master Problem in terms of the "original variables" of the wFi-th component, the required value is just ( c[ wFi ] - (Lambda + d* ) * A[ wFi ] ) x[ wFi ]*).

Note that the 0-th component is an "easy" one but it is dealt with as a "difficult" one, in that the translated model (the only possible one, the true linear function) is used.

Implements MPSolver.

◆ ReadGid()

HpNum ReadGid ( cIndex Nm = InfIndex >())
overridevirtual

Must return the scalar product between the item with name ‘Nm’ [see SetItem() below] and the primal optimal solution d*; this information is typically computed anyway, so it may be made available at little cost. When Nm >= MaxBSize, then the scalar product to be returned is that with the constant subgradient of the "linear part" of Fi, the 0-th component.

Implements MPSolver.

◆ ReadLBMult()

HpNum ReadLBMult ( cIndex wFi = InfIndex >())
overridevirtual

Must return the information required by NDOSolver::ReadLBMult(), see NDOSlver.h.

Implements MPSolver.

◆ ReadLinErr()

cHpRow ReadLinErr ( void )
overridevirtual

Must return a read-only pointer to a vector containing in position i the current value of the linearization error/right hand side of item ‘i’.

Implements MPSolver.

◆ ReadLowerBound()

HpNum ReadLowerBound ( cIndex wFi = InfIndex >())
overridevirtual

Reads back the current value of the Lower Bound; see SetLowerBound() for the meaning of wFi. Note that the value is relative to the function values in the current point, and therefore changes when these (say, the current point itself) do. The updating need be done by the MPSolver, which is why it makes sense to be able to read the current value.

Implements MPSolver.

◆ ReadMult()

cHpRow ReadMult ( cIndex_Set & I,
Index & D,
cIndex wFi = InfIndex >(),
const bool IncldCnst = true )
overridevirtual

The dual optimal solution z* lies in the domain of Fi_{B,Lambda}*; at least with the "classical" Cutting Plane model \hat{Fi}_{B,Lambda}, this is the convex hull of all the subgradients in B. Hence, there is another possible forms in which the dual optimal solution can be represented, i.e., as a set of (convex, for \hat{Fi}) multipliers Theta[] (>= 0 and summing to 1) such that

z* = Sum{ i \in B } Theta[ i ] * z[ i ]

where z[ i ] is the item with "name" ‘i’ [see [Set/Get]Item() below]. ReadMult() returns these multipliers. The format is "dense" if I == 0 and "sparse" otherwise: if I != 0, then only the first D entries of the returned vector contain relevant information, the i-th entry being relative to item with name ‘I[ i ]’, all other entries being zero. If I == 0, a full MaxBSize-vector is returned.

If Fi is separable, then there are as many dual optimal solutions z[ i ]* as components [see ReadSigma() above], the optimal solution z being their sum; thus, ReadMult() must return:

  • the multipliers making z[ wFi ]* (those relative to the subgradients of the i-th component) if 1 <= wFi <= NrFi (wFi = 0 makes no sense since z[ 0 ]* is a known constant vector, i.e., its multiplier is always 1);
  • the union of all the things that would be reported by calls with wFi going from 1 to NrFi otherwise.

If the parameter IncldCnst is left to true, then the returned multipliers are both of subgradients and constraints (corresponding to component wFi). In the constrained case, "extra" dual variables w (w[ i ]) corresponds to constraints, and the optimal solution is actually z* + w* (z[ i ]* + w[ i ]*): the multipliers corresponding to constraints are just nonnegative—they are not constrained to sum to 1—as they represent "extreme rays" in the dual feasible set. If IncldCnst == false, only the multipliers corresponding to subgradients are returned.

Note
The pointers returned by ReadMult() may point to a temporary data structure; the pointers are no more valid (they may even be dangling) after any call to any other method of the class.
No subgradient of "easy" components of Fi() [see FiOracle::GetBNC() ...] are ever produced, so a call to ReadMult( wFi ) with wFi the index of an "easy" component may appear to make no sense. However, something else naturally "takes the place" of the multipliers: the optimal solution x[ wFi ]* of the Master Problem in terms of the "original variables" of the wFi-th component. Thus, calls to ReadMult( wFi ) with an "easy" wFi are indeed allowed, with the following different meaning: the vector returned by the method (that may be "dense" or "sparse" as in the standard case) represents x[ wFi ]*, and it is therefore a FiOracle::GetBNC( wFi )-vector. This information, however, can only be accessed when calling ReadMult( wFi ) for wFi <= NrFi: in "global" calls (wFi > NrFi) only the multipliers corresponding to non-easy components are returned.

Implements MPSolver.

◆ ReadReducedCostsEasy()

cHpRow ReadReducedCostsEasy ( cIndex wFi)
overridevirtual

returns the reduced costs of the variables of an easy component

In a Bundle algorithm, easy components will never be evaluated since they "live in the Master Problem". However, at termination (or at any point in time) the user may want to get the optimal value of the corresponding primal and dual solutions. ReadMult [see above] does this for primal variables, ReadDualEasy [see above] and this method does this for dual variables and reduced costs, respectively. That is, this method returns a read-only pointer to a FiOracle::GetBNC( wFi )-vector so that the i-th entry contains the optimal reduced cost of the variable in the i-th column of the matrix returned by FiOracle::GetBDesc( wFi ).

It is of course an error calling this method if wFi is not the "name" of an easy component.

Note
returned pointer may be that of a temporary data structure; hence, it is no more valid (it may even be dangling) after any call to any other method of the class
the returned pointer may be nullptr if the wFi-th easy component has no columns, hence no reduced costs
the method is not pure virtual but it has a default implementation doing nothing (returning nullptr) for derived classes that do not handle easy components and therefore have nothing to implement

Reimplemented from MPSolver.

◆ ReadSigma()

HpNum ReadSigma ( cIndex wFi = InfIndex >())
overridevirtual

Returns the value of the conjugate Fi_{B,Lambda}* of the model Fi_{B,Lambda} in z*, the optimal dual solution of the MP; this value is often called "sigma*".

When Fi is a decomposable function [see SetDim() above], each component (of the model) has a separate conjugate function Fi[ k ]_{B,Lambda}* and a separate dual solution z[ k ]*: the dual stabilizing term is then

D_t*( - \sum_{k = 0}^{NrFi} z[ k ] ).

If primal constraints Lambda \in L are added to the primal, this is viewed as adding the indicator function I_L of the feasible set to the stabilized Primal Master Problem; this corresponds to having the support function \sigma_L appearing in the dual objective function. Thus, in this case the dual objective function is

\sum_{k = 1}^{NrFi} Fi[ k ]_{B,Lambda}*( z[ k ] ) + \sigma_L( w ) + D_t*( - \sum_{k = 0}^{NrFi} z[ k ] )

i.e., it has a term for each component, the stabilizing term and an extra term for constraints (this term is always zero in the primal because the indicator function I_L is zero in L).

Thus, ReadSigma() must return:

  • Fi[ wFi ]_{B,Lambda}*( z[ i ]* ) if 1 <= wFi <= NrFi: note that the value for the 0-th component of Fi is always 0, since the conjugate of a linear function is 0 iif its argument is equal to the gradient and +Infinity otherwise—in other words, z[ 0 ] is always equal to the gradient of the 0-th component;
  • \sum_{k = 1}^{NrFi} Fi[ k ]_{B,Lambda}*( z[ k ] ) + \sigma_L( w ) if wFi > NrFi, i.e., the dual objective function except the (dual) stabilizing term.
Note
For "easy" components of Fi() [see FiOracle::GetBNC() ...] the value that should be returned by this method is the actual value of the conjugate of the (corresponding component of) Fi(). Although it may be possible to extract that information out of the solution of the MP, a Bundle algorithm typically will never use it, so the MPSolver is allowed to assume that ReadSigma( wFi ) will never be called with wFi the index of an "easy" component. Furthermore, as in the case of ReadFiBLambda() [see above], in an "aggregated" call (that is, with wFi > NrFi) only the difficult components are counted, that is, the value of the conjugate of the model of an "easy" component never has to be computed.

Implements MPSolver.

◆ ReadZ()

void ReadZ ( LMRow tz,
cIndex_Set & I,
Index & D,
cIndex wFi = InfIndex >() )
overridevirtual

Writes in tz[] the dual optimal solution z* in its "natural" format, i.e., a single vector z*, with as many components as d*. The format is "dense" if I == 0 and "sparse" otherwise: if I != 0, then only the first D entries of tz[] contain relevant information, the i-th entry being relative to the variable ‘I[ i ]’, all other entries being zero. If I == 0, a full CrrSGLen-vector is returned. Note that all the entries of z are "significant", even the ones corresponding to variables which are not in the "active set" [see SetActvSt() below].

If Fi is separable, then there are as many dual optimal solutions z[ i ]* as components [see ReadSigma() above], the optimal solution z being their sum. Thus, ReadZ() must return:

  • z[ 0 ]* == the gradient of the 0-th component if wFi == 0;
  • z[ wFi ]* if 1 <= wFi <= NrFi;
  • the sum of all the things that would be reported by calls with wFi going from 1 to NrFi if NrFi < wFi < Inf< Index >();
  • the sum of all the things that would be reported by calls with wFi going from 0 to NrFi if wFi == INF (this is the same as summing what is separately obtained with wFi == 0 and wFi == NrFi + 1).
Note
In the constrained case, "extra" dual variables w (w[ i ]) corresponds to constraints, and the optimal solution is actually z* + w* (z[ i ]* + w[ i ]*); ReadZ() only returns the "z* part" of the dual solution, leaving the "w* part" out.
The vector pointed by I[] may be that of a temporary data structure; hence, it is no more valid (it may even be dangling) after any call to any other method of the class. There is one notable exception: I[] can be passed to SetItemBase(), and be used as set of indices af a sparse subgradient to be inserted in the Master Problem. In fact, this is the typical reason why a z[ i ]* is computed: to perform "aggregation".
For "easy" components of Fi() [see FiOracle::GetBNC() ...] the vector that should be returned by this method is a subgradient of Fi[ wFi ] in Lambda. Although it is possible to extract that information out of the solution of the Master Problem (if x[ wFi ]* is the optimal solution of the Master Problem in terms of the "original variables" of the wFi-th component, this is just - A[ wFi ] x[ wFi ]*), a Bundle algorithm typically will never use it except in its "aggregated" form (that is, in a call with wFi > NrFi), so the MPSolver is allowed to assume that ReadZ( wFi ) will never be called with wFi the index of an "easy" component.

Implements MPSolver.

◆ RmvActvSt()

void RmvActvSt ( cIndex_Set Rmvd,
cIndex RmDm,
cIndex_Set AVrs )
overridevirtual

Removes ‘RmDm’ variables from the current active set, whose names are contained in the first RmDm positions of the vector Rmvd (ordered in increasing sense and Inf< Index >()-terminated); those names must be in the current active set.

AVrs must point to a vector containing the new active set; a pointer to that vector is kept and used from now on as the new active set, replacing the old vector set with the latest call to one of SetActvSt(), AddActvSt() or RmvActvSt(). Thus, the vector pointed by the "old" pointer can be freely modified after the return of RmvActvSt(), while the one pointed by AVrs must not be modified until it is in turn substituted by a new one.

This method can be used only if an active set has been previously defined with SetActvSt( <non-0> ) and not yet removed with SetActvSt( 0 ).

Important note: after a call to RmvActvSt(), all the solution information corresponding to the last call to SolveMP() is lost, and the results returned by all the query methods are unpredictable. This is true even for CheckSubG(), CheckCnst() and ChangesMPSol(), so adding items to the bundle after a call to this method is not advised.

Implements MPSolver.

◆ RmvItem()

void RmvItem ( cIndex i)
overridevirtual

Removes the item with name ‘i’ from the bundle. Items removed are unused, and their names can be used again in GetItem() [see above]. Initially, all the items are unused.

Implements MPSolver.

◆ RmvItems()

void RmvItems ( void )
overridevirtual

Removes all items from the bundle. Note that RmvItems() does not remove the constant subgradient of the linear 0-th component of Fi (RmvItem() [see above] also cannot be used for the purpose as that subgradient has no name). In principle, any Fi has a linear part, possibly all-0; thus, to remove that item one has just to use the [Get/Check/Set]Item() sequence passing an all-0 vector.

Implements MPSolver.

◆ RmvVars()

void RmvVars ( cIndex_Set whch = nullptr,
Index hwmny = 0 )
overridevirtual

This method removes a set of variables from the MP, in response to a call to NDOSolver::RemoveVariables( which ) [see NDOSlvr.h]; the remaining variables in the MS have to be renamed as described in the comments to NDOSolver::RemoveVariables(). As for NDOSolver::RemoveVariables(), whch == 0 means "remove all variables" (hwmny is ignored in this case).

As opposed to RmvActvSt() [see above], the effect of this removal is permanent; the variables are eliminated from the whole MP, rather than just "temporarily set aside" outside of the active set. Indeed, this "permanent" removal can only happen after the removal from the active set: no variable in the active set can be removed, i.e., RmvVars( 0 ) can only be called if the active set is empty (if, of course, the active set has been at all defined).

As for AddVars() [see above], it is assumed that the removal of the variables does not change the linearization errors of the items, i.e., that all the variables d[ i ] of the Primal Master Problem that are removed correspond to variables Lambda[ i ] of the NDO problem that have value zero in the current point; it is responsibility of the caller to ensure this, using ChangeCurrPoint() if required.

Important note: after a call to RmvVars(), all the solution information corresponding to the last call to SolveMP() is lost, and the results returned by all the query methods are unpredictable. This is true even for CheckSubG(), CheckCnst() and ChangesMPSol(), so adding items to the bundle after a call to this method is not advised.

Implements MPSolver.

◆ SensitAnals()

void SensitAnals ( HpNum & lp,
HpNum & cp )
overridevirtual

< Gives a linear lower approximation on the optimal value of v (the predicted decrease from the Cutting Plane model) as a function of the parameter t in the subproblem. That is, returns two numbers lp and cp such that

v( t ) = Fi_{B,x}( x + t * d ) >= t * lp + cp.

Implements MPSolver.

◆ SetActvSt()

void SetActvSt ( cIndex_Set AVrs = nullptr,
cIndex AVDm = 0 )
overridevirtual

Solving the MP can be costly; a possible way of decreasing the cost is to resort to an "active set" strategy over the primal variables d[ i ]. That is, one may guess that a (large) set of the variables will have 0 as the optimal value, and only let the others (the "active set") to vary; if it turns out that the guess has been incorrect, the active set must be revised by inserting or deleting variables [see [Add/Rmv]ActvSt() below]. This way, a smaller MP is solved at each step.

Whether or not an "active set" is used is decided in SetDim() [see above]; if it is, then SetActvSt() sets it to the variables whose names are in the first AVDm entries of the vector AVrs[] (which must be ordered in increasing sense and Inf< Index >()-terminated, i.e. AVrs[ AVDm ] == Inf< Index >()). It is an error to call SetActvSt() if the active set option has not previously been selected with SetDim().

Note that changing the active set has no impact on the linearization errors of the items. Setting the active set dictates which among the primal variables d[ i ] of the MP can be different from zero, hence, which primal variables Lambda[ i ] of the original NDO problem may change their value w.r.t. the value that they have in the current point. The linearization errors of the items depend on the current point only, and therefore they change if and only if the current point changes [see ChangeCurrPoint() and [Add/Rmv]Vars() below]. Since setting the active set affects only the variables d[ i ] of the MP and does not affect the value of the primal variables Lambda[ i ] of the Bundle algorithm, setting (or changing) the active set has no impact on the linearization errors.

Note that the pointer passed to SetActvSt() will be retained and used until it is changed, either by a new call to SetActvSt() - which just substitutes any previous active set with the new one - or by calls to [Add/Rmv]ActvSt() [see below]. Thus, the caller should not change it while it is "in use" by the MPSolver object; the vector can only be changed after* the return of the call to SetActvSt() or [Add/Rmv]ActvSt() which makes it unused. This implies that any call to SetActvSt() which sets an entirely new active set must pass a pointer that is different from the one which contained the old active set [see also [Add/Rmv]ActvSt() below]. Exception is of course AVrs == 0, that sets an empty active set: that is, the active set option is used (only SetDim() can change this), and it currently contains no variables. This is the default if SetActvSt() has not been called yet.

Important note: after a call to SetActvSt(), all the solution information corresponding to the last call to SolveMP() is lost, and the results returned by all the query methods are unpredictable. This is true even for CheckSubG(), CheckCnst() and ChangesMPSol(), so adding items to the bundle after a call to this method is not advised.

Implements MPSolver.

◆ SetAlgo()

void SetAlgo ( const OsiAlg algo = kAuto,
const OsiRed reduc = rNo )

Change "int" algorithmic parameters of the NDO solver. This method is used to set the method for solving either LP or QP.

◆ SetDim()

void SetDim ( cIndex MxBSz = 0,
FiOracle * Oracle = nullptr,
const bool UsAvSt = false )
overridevirtual

Provides the MPSolver with basic information about the "size" of the Master problem:

  • MxBSz is the maximum number of different items (subgradients and constraints) that can be managed by the class, i.e. the maximum size of the bundle: the protected field MaxBSize is offered by the base class to store this information;
  • Oracle is a pointer to the FiOracle object that computes the function Fi to be minimized: the other data describing the size of the problem can be obtained by calls to methods of the FiOracle class, and three protected fields are offered by the base class to store some of this information:

    = MaxSGLen the maximum length of any item, i.e., the maximum number of variables in the Primal Master Problem (PMP);

    = CrrSGLen (<= MaxSGLen) is the current length of any item: the interface supports on-line creation and destruction of primal variables;

    = NrFi the function Fi to be minimized may be "decomposable", i.e., the sum of k functions; in this case, its model Fi_{B,Lambda} is usually decomposable as well (for instance, the Cutting Plane model is), For decomposable functions, it is possible (although not necessary) to keep k "separate bundles", each one describing a model Fi[ k ]_{B,Lambda} of each component Fi[ k ] of Fi, which may be beneficial to improve the "exactness" of the model. Special support is offered for the case when one of the components is a linear function, by considering it as the "0-th component" of Fi [see ReadFiBLambda() and SetItem() and GetItem()].

  • UsAvSt is true if the Active Set technique will be exploited by the solver, and false otherwise. In the first case, the initial "active set" of variables is empty, and a nonempty set of active variables must be set with SetActvSt(). In the second case, SetActvSt() must not be called, and all the variables are considered to be always "active".

This information may be necessary for allocating the data structures of the class, so this method must be called prior to any other method. Other information from the FiOracle is likely to be useful for the MPSolver, such as the max expected density of the subgradients [see FiOracle::GetMaxNZ()], which variables are constrained in sign or has upper bounds [see FiOracle::Get[UC/UB]()], and which of the components of Fi(), if any, is "easy" [see FiOracle::GetBNC() ...].

Note
"easy" components of Fi() have basically to be dealt with by the MPSolver, by inserting their description in the Master Problem; if this is not possible, the Master Problem has to signal it (e.g. by throwing an exception) because the Bundle algorithm relies on this and the FiOracle is not going to provide "ordinary" black-box information (function values, subgradients, ...) for these components.
Important: "easy" components of Fi() in the Master Problem are not translated by value (but they are by argument), meaning that the function that is included in the MP is

Fi_{Lambda}[ k ]( d ) = Fi[ k ]( Lambda + d )

and not

Fi_{Lambda}[ k ]( d ) = Fi[ k ]( Lambda + d ) - Fi[ k ]( Lambda )

as one would expect by analogy with the ordinary "difficult" components. This is done to spare the Bundle algorithm with the need to compute Fi[ k ]( Lambda ) for all "easy" components k and each current point Lambda, which may be problematic especially for the initial current point—consider the case where Lambda is not feasible, so Fi[ k ]( Lambda ) = +INF! However, the value of Fi[ k ] is actually computed by the MPSolver at Lambda + d* (d* being the optimal solution of the Primal MP), so the value of Fi[ k ] is known for the current point at least after every Serious Step with step 1 along d*. Yet, this choice has some impact on the "output" methods [see ReadFiBLambda() and ReadSigma() below].

This method can be called more than once to modify the settings, but expect the implementation to be quite expensive in time and/or memory. There are essentially three different ways for calling SetDim():

  • SetDim( 0 , ... ) makes the MPSolver to deallocate all its memory and to quietly wait for new instructions.
  • SetDim( n , 0 , ASV ) with n != 0 sets the max bundle size to n and activate/deactivate the Active Set Mechanism without changing anything else. The existing items in the bundle (if any) are all kept if n is larger than the previous setting, but a smaller value will force deletion of all the items with "name" [see [Get/Set]Item() and RmvItem() below] greater than or equal to n. Also, if the active set technique is initialized (ASV == true while it was false previously), the set of "active" variables is set to be empty.
  • SetDim( n , Orcl , ASV ) with n != 0 and Orcl != 0 discards all the previous settings and re-allocates everything; all the existing items in the bundle are lost. Note that such a call also resets every parameter of the algorithm, such as the starting point (which is set to 0).

In general, calling this method with a non-empty bundle could be costly, so if the items are to be discarded anyway, this should be done before the call to SetDim() [see RmvItems() below].

Reimplemented from MPSolver.

◆ SetItem()

void SetItem ( cIndex Nm = InfIndex >())
overridevirtual

If – after all the above checks – the Bundle code decides to actually insert the item in the bundle, it has just to call SetItem( Nm ). ‘Nm’ is the "name" of the new item: the names of items are the integers between 0 and MaxBSize - 1. ‘Nm’ must be a currently unused name, that is, either a name that has not yet been used or a name that has been used in RmvItem() [see below] and not yet used in GetItem() ever since. However, if wFi == 0 has been passed to GetItem(), then Nm == Inf< Index >() has to be passed, as the only subgradient corresponding to the (linear) 0-th component of Fi has no name.

Implements MPSolver.

◆ SetItemBse()

void SetItemBse ( cIndex_Set SGBse = nullptr,
cIndex SGBDm = 0 )
overridevirtual

When the correct values have been written in the vector SG returned by GetItem() [see above], SetItemBse() can be used to tell the MPSolver the "format" of the information. If SGBse == 0, then SG is in "dense" format, i.e., SG[ i ] contains the entry of the subgradient relative to the variable with name ‘i’, 0 <= i < CrrSGLen - 1. If SGBse != 0 instead, then SG is in "sparse" format, i.e., SG[ i ] is the entry of the subgradient relative to the variable with name ‘SGBse[ i ]’, 0 <= i < SGBDm - 1. SGBse is assumed to be ordered in increasing sense and Inf< Index >()-terminated, i.e., SGBse[ SGBDm ] == Inf< Index >(); if SGBse == 0, then the value of SGBDm is ignored. If SetItemBse() is not called, then the "dense" format is assumed.

Note that all the CrrSGLen entries of the item must be given (although in the sparse format those that are 0 can be given only implicitly); that is, the entry relative to the variable ‘i’ must be given even if ‘i’ is currently out of the "active set" [see SetActvSt() below].

The MPSolver is allowed not to copy the memory pointed by SGBse (if any) and to keep using it throughout all the sequence of Check[SubG/Const](), ChangesMPSol() and SetItem() [see below] that is needed to decide the "fate" of the item. As soon as SetItem() returns, or GetItem() [see above] is called again, the memory pointed by SGBse must not be accessed again.

Implements MPSolver.

◆ SetLowerBound()

void SetLowerBound ( cHpNum LwBnd = InfHpNum >(),
cIndex wFi = InfIndex >() )
overridevirtual

Sets a lower bound on the value of the (translated) model Fi_{B,Lambda}; since Fi_{B,Lambda}( 0 ) = 0, LwBnd must clearly be a nonpositive value.

Note
The translated model is Fi_{B,Lambda} is only translated w.r.t. the function value in Lambda of the non-easy components, while the easy components retain their original value. This means that LwBnd must be the original value of the lower bound minus the value \sum_{ k is not an easy component } Fi[ k ]( Lambda )
In the above formula, the 0-th linear component is considered to be a non-easy component, i.e., its value has to be subtracted from the original value of the lower bound. There can hardly be an easier component than the 0-th one, but the variable in the master problem is d = Lambda1 - Lambda and the term in the objective function is RHS * d, therefore the 0-th component is naturally translated w.r.t. Fi[ 0 ]( Lambda ) like the non-easy ones.

When Fi is a decomposable function [see SetDim() above], each of its components may have a separate model Fi[ k ]_{B,Lambda}, thus this method has different meanings according to the value of wFi. For 1 <= wFi <= NrFi, the lower bound is on the individual model Fi[ wFi ]_{B,Lambda}, while for wFi > NrFi the lower bound is on the global model Fi_{B,Lambda} = sum_{k = 0 ... NrFi} Fi[ k ]_{B,Lambda}. Note that the 0-th component is taken into account for the global model, while, being a linear function, it clearly cannot have any "individual" lower bound.

Note
Individual lower bounds are meaningless for "easy" components of Fi() [see FiOracle::GetBNC() ...], and therefore they are not allowed.

These lower bounds are automatically updated when the current point Lambda is changed [see ChangeCurrPoint() below], and they are equivalent to inserting in the bundle of the proper component an all-0 subgradient with alfa_{Lambda}[ i ] = - LwBnd (which is nonnegative). The same holds for the "global" lower bound, except that it is an "aggregated" all-zero subgradient.

Passing - Inf< HpNum >() as the lower bound means that no such bound is available; this is the default.

Implements MPSolver.

◆ SetMPLog()

void SetMPLog ( ostream * outs = 0,
const char lvl = 0 )
override

lvl controls the "level of verbosity" of the code. The first four bits of lvl have the following meaning:

0 => no log at all (also assumed if log = 0);

1 => some messages are printed;

2 => as, 1 plus the current Master Problem is also saved as a MPS file just prior to being solved;

3 ==> as 2, plus very verbose information is also printed

4 .. 15 unused, available to derived classes.

◆ SetOsi()

void SetOsi ( OsiSolverInterface * osi = nullptr)

Provides OSIMPSolver with the actual OsiSolverInterface.

Provides OSIMPSolver with an object of any class derived from OsiSolverInterface, which is used to actually solve the Master Problem. Note that the object does not become "property" of the OSIMPSolver, which will not delete it when it is deleted. However, calling SetOsi() does delete the existing solver to replace it with the new one (if any), hence to delete the current OsiSolverInterface one just has to call it with default == nullptr argument. This is a bit awkward, but it'll do for now.

◆ SetPar()

void SetPar ( const int wp,
cHpNum value )
overridevirtual

Change "float" algorithmic parameters of the MPSolver. This method is pure virtual, i.e., it has to be defined in derived classes, which may extend this it to allow setting the algorithmic-specific parameters.

The enum MPParam is used for selecting the parameter to be set, with the following intended meaning:

  • kMaxTme: maximum running time (in seconds) for each call to SolveMP() [see below]. A 0 or negative value means "no time limit", which should be the default. If more time required, SolveMP() should stop returning kStppd. Note that this check can only be performed if timing of the code is activated [see SetMPTime() below].
  • kOptEps: tells the solver to solve the MP only with at least the specified relative precision (>= 0), the exact meaning of this being solver-dependent. The solver can always decide to use higher precision than that, if so it prefers. By default a solver-specific "minimal" precision is used.
  • kFsbEps: tells the solver to tolerate the specified relative violation of in the (primal) constraints; the solver can always decide to use higher precision than that, if so it prefers. By default a solver-specific "minimal" precision is used.

Note that the type of ‘wp’ is int rather than MPParam, which is not a problem as enums can always be used where intsare required (i.e., SetPar( MPSolver::kMaxTime , ... ) is a valid call). This is done in order to allow derived classes to extend the set of parameters they can handle while keeping the same signature for SetPar().

Implements MPSolver.

◆ SetStabType()

void SetStabType ( const StabFun sf = none)

Sets the stabilizing term in the Master Problem.

Sets the stabilizing term in the Master Problem. Possible values are:

  • none: no stabilization;
  • boxstep: a box in the infinity-norm in the primal, a 1-norm penalty in the dual;
  • quadratic: a 2-norm penalty term both in the primal and in the dual.

◆ Sett()

void Sett ( cHpNum tt = 1)
overridevirtual

Sets the value of the proximal parameter, t, to tt. Prior to the first call to Sett(), t is taken to be equal to 1.

Implements MPSolver.

◆ SetThreads()

void SetThreads ( int nthreads)
overridevirtual

Set the maximum number of threads that the MPSolver is allowed to use. This method is given a default empty implementation (doing nothing) for those MPSolver that do no have any parallel processing capability.

Reimplemented from MPSolver.

◆ SolveMP()

MPStatus SolveMP ( void )
overridevirtual

Attempts to solve the current Master Problem.

Note
The bundle for some components of Fi may be empty; if this happens, those components must just be temporarily ignored (like they did not exist) in the solution of the MP. If all the "sub-bundles" for each of the components are empty, what is solved is just a "pure feasibility" MP which returns just any point in the polyhedron defined by the constraints, if any. If there are neither subgradient nor constraints, then the optimal (both primal and dual) solution must be all-0.
The last sentence is not valid if there there are "easy" components of Fi() [see FiOracle::GetBNC() ...], since in this case the primal solution d* of the MP is typically non-0 even if there is no subgradient for any of the "difficult" components; what happens is that the bundle is not empty since "those of the easy components always contain all the possible subgradients".

The method must return:

kOK if the solver found an optimal solution;

kUnbndd if the MP is primal unbounded: this can only happen with some "not-really-stabilizing" terms D_t, and should always be solvable by properly increasing the proximal parameter t;

kUnfsbl if the MP is primal unfeasible (dual unbounded or empty, see ReadMult() below);

kStppd the algorithm has been stopped because the time/resources allotted to solve the MP have been depleted; this is not a "serious" condition, in that it can be "solved" by providing the MPSolver with more time/resources (usually, just calling it again, possibly multiple times);

kError in case of a failure in the MP solver.

If the "active set technique" is being used [see SetActvSt() below], the problem that is solved is actually a restriction of the "real" MP, where (possibly many) primal variables d[ i ] are forced to 0; this "restricted" MP may be empty (so that the solver might have to return kUnfsbl) while the full MP is nonempty. Hence, after a kUnfsbl return a "price in" [see Readd() below] should be done to find if new variables can to be added to the "restricted" MP in order to try to make it nonempty. Thus, in the unfeasible case ReadMult() [see below] must return the multipliers that make a feasible ascent extreme ray in the dual, so that the corresponding direction returned by Readd() [see below] is nonzero on variable ‘i’ (strictly positive if ‘i’ is constrained) if and only if that variable has to be added.

Implements MPSolver.

◆ SubstItem()

void SubstItem ( cIndex Nm)
overridevirtual

Although having multiple copies of the same item in the bundle is legal, it is also useless and it should be avoided. If CheckItem() returned a "valid" name (< Inf< Index >()), then SubstItem( Nm ) can be used to replace the current item in the bundle with name ‘Nm’ with SG. It is assumed that ‘Nm’ is the name of an item element-wise identical to SG (but possibly with different Ai), typically the name returned by CheckItem(). Typically, this operation makes sense when the Ai (Right Hand Side/linearization error) of the new item is smaller than that of the existing copy in the bundle. Thus, SubstItem( Inf< Index >() ) - that would operate on the (only) (sub)gradient of the 0-th component of Fi - cannot be called.

Implements MPSolver.

◆ WComponent()

Index WComponent ( cIndex i)
overridevirtual

Returns the name of the component that the item ‘i’ is relative to; a return value of Inf< Index >() means that no item with name ‘i’ is in the bundle. Note that a return value of 0 is possible only if the item is a constraint, since the 0-th component may have constraints but no subgradients: its only (sub)gradient is dealt with in a special way and it has no name.

Implements MPSolver.

Member Data Documentation

◆ first

bool first
private

true if it's the first master problem to solve

◆ GiPerd

HpRow GiPerd
private

keep the scalar products between each item and the direction: while it is not used at all iterations, its computation uses information that may change (e.g. if Alphas change) so it's better to save it

◆ RhoColSgPos

Index RhoColSgPos
private

position in RhoColBse of the "subgradient"

position in RhoColBse of the first entry corresponding to the "subgradient part", i.e., >= comp_row[ NrFi ]; if there are no easy component then the first part of RhoCol is dense and therefore RhoColSgPos == comp_row[ NrFi ], but if there are easy components some entries of RhoCol before comp_row[ NrFi ] may be 0, and RhoCol[] is stored in sparse format.

◆ wcomp

Index_Set wcomp
private

Has several different uses:

  • if InINF, it means the item is not there

the value masked by the two most significant bits is the component of the function that the item belongs to

  • the value of the most significant bit is 1 if the item is a subgradient, 0 if it is a constraint
  • the value of the second-most significant bit is 1 if the item was not present at the time that the last Master Problem has been solved, 0 otherwise

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