ConicBundle
|
abstract interface extending BundleModel so that any such model can be used alone or within SumModel and so that it supports AffineFunctionTransformation as well as switching to a SumBundle with some SumBundleHandler More...
#include <SumBlockModel.hxx>
Public Member Functions | |
virtual | ~SumBlockModel () |
destructor | |
virtual void | clear () |
resets all data to the initial status of this class, also aftmodel and vm_selection are deleted | |
int | initialize_aft (AffineFunctionTransformation *aft=0) |
first it discards an old affine function transformation if there is one, then it sets the new one, if there ins one. | |
virtual int | set_variable_metric_selection (VariableMetricSelection *vms=0) |
delete old selector and set a new one (0 is allowed resulting in no local selector) | |
virtual VariableMetricSelection * | get_variable_metric_selection () const |
delete old selector and set a new one (0 is allowed resulting in no local selector) | |
bool | call_aftmodel_first (const FunObjModMap &funmdfmap) |
in apply_modification this routine is needed to check whether the aftmodel is modified already More... | |
SumBlockModel (CBout *cb=0, int cbinc=-1) | |
calls clear | |
virtual ModifiableOracleObject * | get_oracle_object ()=0 |
return the function oracle this provides a model for, or some dummy oracle | |
virtual const SumBlockModel * | model (const FunctionObject *) const |
returns the submodel for FunctionObject fo if it in this model, otherwise 0 | |
virtual CH_Matrix_Classes::Integer | nsubmodels () const |
returns the number of submodels in this model (direct ones, not all descendants; most have none, but see e.g. SumModel) | |
virtual int | add_model (SumBlockModel *) |
adds the model as submodel to this model (if this model may have submodels) | |
virtual SumBlockModel * | remove_model (const FunctionObject *) |
remove the submodel identified by fo from this model, this does NOT destruct the model. It returns the pointer to the model if there is one, otherwise 0 | |
virtual SumBlockModel * | remove_model (SumBlockModel *model) |
remove the submodel identified by the given models FunctionObject from this model, this does NOT destruct the model. It returns the pointer to the model if there is one, otherwise 0 | |
some default implementations of abstract class BundleModel | |
virtual int | start_augmodel (QPModelDataPointer &blockp, CH_Matrix_Classes::Integer cand_id, const CH_Matrix_Classes::Matrix &cand_y, const CH_Matrix_Classes::Indexmatrix *indices=0) |
see BundleModel::start_augmodel, here it just moves on to start_sumaugmodel | |
virtual int | get_model_aggregate (CH_Matrix_Classes::Integer &model_aggregate_id, MinorantPointer &model_aggregate) |
see BundleModel::get_model_aggregate | |
virtual int | update_model (ModelUpdate model_update, CH_Matrix_Classes::Integer center_id, const CH_Matrix_Classes::Matrix ¢er_y, CH_Matrix_Classes::Integer cand_id, const CH_Matrix_Classes::Matrix &cand_y, CH_Matrix_Classes::Real model_maxviol, BundleProxObject &H)=0 |
generate the next cutting model and store the center information in the case of a descent step More... | |
virtual int | synchronize_ids (CH_Matrix_Classes::Integer &new_center_ub_fid, CH_Matrix_Classes::Integer new_center_id, CH_Matrix_Classes::Integer old_center_id, CH_Matrix_Classes::Integer &new_cand_ub_fid, CH_Matrix_Classes::Integer new_cand_id, CH_Matrix_Classes::Integer old_cand_id, CH_Matrix_Classes::Integer &new_aggregate_id) |
see BundleModel::synchronize_ids | |
virtual bool | center_modified (CH_Matrix_Classes::Integer ¢er_fid, CH_Matrix_Classes::Integer center_id) |
see BundleModel::center_modified | |
virtual bool | model_aggregate_modified (CH_Matrix_Classes::Integer old_model_aggregate_id) |
see BundleModel::model_aggregate_modified | |
virtual int | apply_modification (bool &no_changes, const GroundsetModification &gsmdf, const FunObjModMap &funmdfmap, CH_Matrix_Classes::Integer new_center_id, const CH_Matrix_Classes::Matrix &new_center, CH_Matrix_Classes::Integer old_center_id, const CH_Matrix_Classes::Matrix &old_center) |
passes on modification information about function and ground set changes More... | |
int | add_variable_metric (VariableMetric &H, CH_Matrix_Classes::Integer y_id, const CH_Matrix_Classes::Matrix &y, bool descent_step, CH_Matrix_Classes::Real weightu, CH_Matrix_Classes::Real model_maxviol, const CH_Matrix_Classes::Indexmatrix *indices=0) |
see VariableMetricModel::add_variable_metric; this implements the non local option | |
int | get_function_minorant (CH_Matrix_Classes::Integer &function_modification_id, MinorantPointer &minorant) |
see BundleModel::get_function_minorant() | |
int | get_center_minorant (CH_Matrix_Classes::Integer &function_modification_id, MinorantPointer &minorant) |
see BundleModel::get_center_minorant() | |
BundleModel * | transform () |
if an affine function transformation is defined for this model, return it, otherwise use this | |
new virtual functions needed for the classes AFTModel and SumModel | |
SumBlockModel * | sbm_transform () |
like BundleModel::transform() this allows SumBlockModel routines to use a given transformation automatically | |
virtual int | get_model_aggregate (CH_Matrix_Classes::Integer &model_aggregate_id, MinorantPointer &model_aggregate, bool all_parts, const AffineFunctionTransformation *aft=0) |
returns the model aggregate if available. More... | |
virtual CH_Matrix_Classes::Real | lb_function (CH_Matrix_Classes::Integer y_id, const CH_Matrix_Classes::Matrix &y) |
returns a quick lower bound for the function value at y (eg by a previous subgradient) | |
virtual int | get_function_minorant (MinorantPointer &function_minorant, const AffineFunctionTransformation *aft=0)=0 |
returns the minorant corresponding to the subgradient inequality returned by the last function evaluation (if the function has not changed since then). If the input minorant is not empty, the local minorant is added. | |
virtual int | get_center_minorant (MinorantPointer ¢er_minorant, const AffineFunctionTransformation *aft=0)=0 |
returns the minorant corresponding to the subgradient inequality returned by the function evaluation for the current center (if the function has not changed since then). If the input minorant is not empty, the local minorant is added. | |
virtual int | adjust_multiplier (bool &values_may_have_changed)=0 |
for conic subproblems with adjustable multipliers, reset the multiplier to twice the current trace and set values_may_have_changed to true in this case. Otherwise leave values_may_have_changed unmodified. More... | |
virtual int | sumbundle_mode (SumBundle::Mode &mode, SumBundleHandler *bh=0, AffineFunctionTransformation *aft=0)=0 |
called by start_sumaugmodel in order to suggest (or to force) *this to opt in or out of a common SumBundle handled by the SumBundleHandler More... | |
virtual int | start_sumaugmodel (QPModelDataPointer &blockp, CH_Matrix_Classes::Integer cand_id, const CH_Matrix_Classes::Matrix &cand_y, const CH_Matrix_Classes::Indexmatrix *indices=0, SumBundleHandler *bh=0, SumBundle::Mode mode=SumBundle::inactive, AffineFunctionTransformation *aft=0)=0 |
see BundleModel::start_augmodel() for the first four parameters; for the others see sumbundle_mode() | |
virtual int | update_model (ModelUpdate model_update, CH_Matrix_Classes::Integer center_id, const CH_Matrix_Classes::Matrix ¢er_y, CH_Matrix_Classes::Integer cand_id, const CH_Matrix_Classes::Matrix &cand_y, CH_Matrix_Classes::Real model_maxviol, BundleProxObject &H, CH_Matrix_Classes::Real &model_deviation, CH_Matrix_Classes::Real &model_curvature)=0 |
generate the next cutting model and store the center information in the case of a descent step. This version also allows for updating a common model owned by parent SumModel. More... | |
messages for direct get/set requestss | |
AFTModel * | get_aftmodel () |
returns aftmodel | |
const AFTModel * | get_aftmodel () const |
returns aftmodel (as a const variant) | |
virtual BundleData * | get_data ()=0 |
returns the essential information on the bundle date for continuing from the same point lateron, see BundleData | |
virtual const BundleData * | get_data () const =0 |
const version returning the essential information on the bundle date for continuing from the same point lateron, see BundleData | |
virtual int | set_data (BundleData *) |
reinstalls the bundle date from get_data() for continuing from this previous point, see BundleData | |
virtual const PrimalData * | get_approximate_primal () const |
if primal data is provided by the oracle then the primal corresponding to the current aggregate is formed in primal (this may differ if primal is a recognized derived class) | |
virtual const PrimalData * | get_center_primal () const |
if primal data is provided by the oracle then the primal corresponding to the best eps-subgradient of the evaluation in the current center is returned (this may differ if primal is a recognized derived class) | |
virtual const PrimalData * | get_candidate_primal () const |
if primal data is provided by the oracle then the primal corresponding to the best eps-subgradient of the evaluation in the current center is returned (this may differ if primal is a recognized derived class) | |
virtual int | call_primal_extender (PrimalExtender &pext) |
if the function is the Lagrangian dual and primal_data of previous calls has now to be updated due to changes in the primal problem – e.g., this may happen in column generation – the problem updates all its internally stored primal_data objects by calling PrimalExtender::extend on each of these. returns 0 on success, 1 if not applicable to this function, 2 if it would be applicable but there is no primal data. | |
virtual int | set_bundle_parameters (const BundleParameters &) |
set max_bundle_size and max_model_size (this may differ if the parameter is a recognized derived class); model blocks without bundle return 1. | |
virtual const BundleParameters * | get_bundle_parameters () const |
returns the current parameter settings; model blocks without bundle parameters return NULL. | |
virtual int | set_sumbundle_parameters (const BundleParameters &) |
set max_bundle_size and max_model_size (this may differ if the parameter is a derived class, in particular a SumBunldeParametersObject); model blocks without bundle return 1. | |
virtual const SumBundleParametersObject * | get_sumbundle_parameters () const |
returns the current parameter settings; model blocks without bundle parameters return NULL. | |
virtual void | clear_model (bool discard_minorants_only=false) |
modifications of this specific problem were such that old subgradient data and function values have to be removed completely; do so. In the special case of some particular support functions it may be possile to regenerate some of the minorants and keep the core of the model of the support set; in this case set discard_minorants_only=true. | |
virtual void | clear_aggregates () |
remove all current aggregate cutting planes, but not the other parts of the model; returns 0 on success, 1 on failure | |
virtual int | get_ret_code () const |
for functions given by an oracle: return the return value that was produced by the last call to the function evaluation routine; model blocks without oracle return 0. | |
output and statistics on solution times | |
CH_Tools::Microseconds | get_evalmodel_time () const |
return time spent in total for evaluating the model in eval_model() | |
CH_Tools::Microseconds | get_updatemodel_time () const |
return time spent in total for updating the model in update_model() | |
virtual CH_Tools::Microseconds | get_preeval_time () const |
return time spent in total for the oracle in eval_function() | |
virtual CH_Tools::Microseconds | get_eval_time () const |
return time spent in total for the oracle in eval_function() | |
virtual CH_Tools::Microseconds | get_posteval_time () const |
return time spent in total for the oracle in eval_function() | |
CH_Tools::Microseconds | get_metric_time () const |
return time spent in total for the add_variable_metric routine | |
std::ostream & | print_statistics (std::ostream &out) const |
output the timing statistics | |
void | set_out (std::ostream *o=0, int pril=1) |
set output and outputlevel of warnings and errors recursively, see CBout | |
void | set_cbout (const CBout *cb, int incr) |
set output and outputlevel of warnings and errors recursively with CBout | |
Public Member Functions inherited from ConicBundle::BundleModel | |
BundleModel (CBout *cb=0, int cbinc=-1) | |
constructor (cb allows to set output options) | |
virtual | ~BundleModel () |
virtual destructor | |
virtual int | eval_function (CH_Matrix_Classes::Integer &ub_fid, CH_Matrix_Classes::Real &ub, CH_Matrix_Classes::Integer y_id, const CH_Matrix_Classes::Matrix &y, CH_Matrix_Classes::Real nullstep_bound, CH_Matrix_Classes::Real relprec)=0 |
evaluates the objective function in y and returns an upper bound in ub within relative precision relprec. More... | |
virtual int | eval_model (CH_Matrix_Classes::Real &lb, CH_Matrix_Classes::Integer y_id, const CH_Matrix_Classes::Matrix &y, CH_Matrix_Classes::Real relprec)=0 |
evaluate the current cutting model in the given point More... | |
virtual int | make_model_aggregate (bool &penalty_parameter_increased, bool keep_penalty_fixed)=0 |
after the common QP is solved, this call asks to form the new aggregate from the solution. If keep_penalty_fixed==false the model may decide to increase some internal penalty parameter, has then to report this in penalty_parameter_increased but need not form the aggregate. More... | |
virtual int | recompute_center (CH_Matrix_Classes::Integer &new_center_ub_fid, CH_Matrix_Classes::Real &new_center_ub, CH_Matrix_Classes::Integer center_id, const CH_Matrix_Classes::Matrix ¢er_y, bool accept_only_higher_values=false, CH_Matrix_Classes::Real relprec=-1.)=0 |
after modifications of the problem the center information may have to be recomputed partially or completely More... | |
virtual int | provide_model_aggregate (CH_Matrix_Classes::Integer y_id, const CH_Matrix_Classes::Matrix &y)=0 |
makes sure that the model_aggregate returned by add_model_aggregate is actually a minorant contained in the next cutting model. More... | |
virtual int | check_center_validity_by_candidate (bool &cand_minorant_is_below, CH_Matrix_Classes::Integer center_id, const CH_Matrix_Classes::Matrix ¢er_y)=0 |
consistency check for oracle computations: test if the subgradient inequality arising out of the last eval_function holds for center_y. More... | |
VariableMetricModel * | variable_metric_transform () |
replaces variable_metric_transform by transform | |
Public Member Functions inherited from ConicBundle::VariableMetricModel | |
VariableMetricModel (CBout *cb=0, int cbincr=-1) | |
constructor for passing on ouptut information | |
virtual | ~VariableMetricModel () |
virtual destructor | |
Public Member Functions inherited from ConicBundle::CBout | |
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 Attributes | |
CH_Tools::Clock | clock |
for collecting time statistics | |
CH_Tools::Microseconds | evalmodel_time |
total time spent in evaluating the model in eval_model() | |
CH_Tools::Microseconds | updatemodel_time |
total time spent in updating the bundle in update_model() | |
CH_Tools::Microseconds | eval_time |
total time spent in the oracle in eval_function() or calls to eval_function of children | |
CH_Tools::Microseconds | preeval_time |
total time spent in eval_function() before the oracle call | |
CH_Tools::Microseconds | posteval_time |
total time spent in eval_function() after the oracle call | |
CH_Tools::Microseconds | metric_time |
total time spent in add_variable_metric() | |
AFTModel * | aftmodel |
if not NULL this points to an AFTModel describing an AffineFunctionTransformation | |
VariableMetricSelection * | vm_selection |
if not NULL this points to the selection routine for computing the local contribution to the metric | |
SumBundleHandler * | bundlehandler |
if the sumbundle is used, this points to the handler that operates it | |
SumBundleParametersObject * | sumbundle_parameters |
this holds the default settings whether and how to start or contribute to a SumBundle | |
MinorantPointer | old_model_aggregate |
for testing purposes | |
Additional Inherited Members | |
Public Types inherited from ConicBundle::BundleModel | |
enum | ModelUpdate { new_subgradient, descent_step, null_step } |
for informing update_model() at what stage it is called to update the bundle so that the amount of information available is clear More... | |
abstract interface extending BundleModel so that any such model can be used alone or within SumModel and so that it supports AffineFunctionTransformation as well as switching to a SumBundle with some SumBundleHandler
The BundleSolver will alway assume that there is one cutting model and it does not need to know how this looks like. A SumBlockModel like the SumModel is used to present this single model view to the BundleSolver, yet it allows e.g. in the case of SumModel for a sum of convex functions to use a separate cutting model for each of the functions. Of course it is also possible to set up a common single model via an appropriate oracle or to recursively build a sum of models of sums of models, the latter with some loss of efficiency. The joint quadratic bundle subproblem will consist of a quadratic cost matrix coupling all models (start_augmodel() serves to collect the bundle model, BundleProxObject forms the quadratic data from this), but the generating descriptions for the linear minorants of each model will each be a separate block of constraints (obtained in start_augmodel()). After the quadratic subproblem is solved the local aggregates are formed via make_model_aggregate(), afterwards the joint aggregate of the entire sum will need to be collected by recursive calls to get_model_aggregate().
If an objective is a sum of convex functions, most functions will just work on part of the coordinates or on an affine transformation of the coordinates. Here we implement this by supplying each model with the possibility to be combined with an AFTModel, which first performs an AffineFunctionTransformation of the input arguments.
Depending on its ConicBundle::FunctionTask each function may be used as a simple cost function or it may be of penalty type for its level set to the value zero. Functions of penalty type may have a prespecified constant penalty parameter or an adaptive penalty parameter that may be increased or decreased upon need (dynamic updates of the penalty parameter take place in make_model_aggregate() and adjust_multiplier()). Intuitively speeking, the penalty parameter allows to increase the length of subgradients, thereby increasing their slope so that the iterates get forced to lie within the level set in the limit. Whenever the feasible set is compact and the functions are finite valued, a finite penalty value is sufficient (exact penalty) if on the boundary of the level set of the penalty function the slope is bounded away from zero.
If SumModel contains many separate functions, having a separate model for each of the functions may be unnecessary and inefficient. Making use of additional information provided in the extended routine update_model(ModelUpdate,CH_Matrix_Classes::Integer,const CH_Matrix_Classes::Matrix&,CH_Matrix_Classes::Real&,CH_Matrix_Classes::Real&,CH_Matrix_Classes::Real&,CH_Matrix_Classes::Real&,CH_Matrix_Classes::Real&) the class SumModel offers a heuristic that suggests which of the functions should keep a separate model and which should only contribute its subgradient information to a common polyhedral model for the remaining functions. The heuristic gives a new suggestion on the fly whenever the next model is formed. Thus each function needs to keep track of the subgradients and the aggregate in the global polyhedral model should it want to participate in the common global model, so that switching back and forth does not endanger convergence. Note that no model is forced to contribute and that models may even decide to contribute one part but keep another part as a local model. Due to the different types of functions (fixed, bounded, dynamic) SumModel has to use a different common model for each type appearing. For keeping consistent with the recursive structure it also needs to assume that each of its functions has these three different types of models. This makes the whole business a bit clumsy. For uniform treatment, each Model is equipped with a SumBundle, and this is employed only if the Model is given a SumBundleHandler and informed to use it via the routine sumbundle_contribution(). If a model/function wants to stop contributing to a SumBundle of some parent model without destroying the bundle information of the others e.g. because it is modified, it can do so by calling sumbundle_remove_contributions().
The class also offers a number of information collection routines that are not needed for running BundleSolver but are requested in many applications, e.g., for providing the common subgradient in the center, information about the PrimalData generating the aggregate, or some computation time measurements.
|
pure virtual |
for conic subproblems with adjustable multipliers, reset the multiplier to twice the current trace and set values_may_have_changed to true in this case. Otherwise leave values_may_have_changed unmodified.
This message is recursively passed on. Afterwards a recompute_center has to be called if values_may_have_changed is true.
[in,out] | values_may_have_changed |
|
Implemented in ConicBundle::SumModel, ConicBundle::AFTModel, and ConicBundle::ConeModel.
Referenced by remove_model().
|
virtual |
passes on modification information about function and ground set changes
In addition to the explanation given in BundleModel::apply_modification() the variable no_changes also acts as an input variable.
If on input no_change==false and the model contributes as a child to sumbundle, it has to remove these contributions before excuting any changes, even if there are no changes within.
Implements ConicBundle::BundleModel.
Reimplemented in ConicBundle::SumModel, and ConicBundle::AFTModel.
Referenced by remove_model().
bool ConicBundle::SumBlockModel::call_aftmodel_first | ( | const FunObjModMap & | funmdfmap | ) |
in apply_modification this routine is needed to check whether the aftmodel is modified already
If an aftmodel is present or fundmdfmap contains no data for the aftmodel of this, the routine returns false. If, however, no aftmodel is present so far and in spite of this funmdfmap contains modifications for it, these were no yet carried out but need to be carried out before the data is passed on to *this. So in this case the routine initializes an aftmodel and returns true (otherwise it always returns false). Whenever this routine returns true, the apply_modification routine has to pass the call on to the new aftmodell first so that it is called itself again with the correct parameters.
Referenced by get_variable_metric_selection().
|
virtual |
returns the model aggregate if available.
The current routine differs from BundelModle::get_model_aggregate in the two additional parameters all_parts, and aft. If add==false, the aggregate has to be initiliazed. If == true, the aggregate has to be added. If all_parts==true, the entire aggregate is used. If all_parts==false, then the local parts are used that correspnd to model parts handled by *this rather than a sumbundle of a parent. If aft!=0, the AffineFunctionTransformation pointed to has first to be applied and the result has to be stored/added in/to the aggregate.
[out] | model_aggregate_id | nonnegative counter for quickly checking validity of the aggregate |
[in,out] | model_aggregate | the aggregate linear model minorant |
[in] | all_parts |
|
[in] | aft |
|
Reimplemented in ConicBundle::AFTModel.
|
pure virtual |
called by start_sumaugmodel in order to suggest (or to force) *this to opt in or out of a common SumBundle handled by the SumBundleHandler
The routine decides on contributing to the parents sumbundle or not (unless force==true). Its outcome may be, that this does not contribute or that it contributes just some part of its subgradient information. As the parent can only discern between contributing and not contributing, it is the task of *this to feed the parent's Bundlehandler *bh with the correct information. In particular, if the model is willing to contribute at all, it should always update its own SumBundle by the same rules in order to allow for switching back and forth between common and local model without siginificant loss of information.
The parent's bundlehandler and its corresponding aft is only passed here and has to be memorized by the local handler if it intends to contribute or communicate with the parent about the sumbundle at some point in time.
The basic assumption is that when the calling parent passes the next (the same or another) handler bh (with corresponding aft) while *this is already contributing, this next handler and aft is consistent with the previous information in that the contribution of *this is not changed and *this will be able to continue by only updating its contribution. If this is not the case, the caller has to first force this to not contribute to the old handler (i.e. to remove its contribution) and then call the routine again with its new handler.
If the SumBundleHandler bh==NULL, the parent does not maintain a sumbundle and the model will be strictly local.
On input:
On output:
It returns 0 if things worked out as they should, 1 otherwise.
Implemented in ConicBundle::SumModel, ConicBundle::AFTModel, and ConicBundle::ConeModel.
Referenced by remove_model().
|
pure virtual |
generate the next cutting model and store the center information in the case of a descent step
If model_update==null_step, the next model has to contain at least all convex combinations of the current subgradient of eval_function() and the aggregate minorant of (re)eval_augmodel(). update_model may only be called if at least one of both is available.
The point passed in the argument is the candidate where the last function evaluation took place, its id is passed on so that some consistency check is possible without having to store the entire vector. The coordinates of the point may help in setting up a model of good quality close to this point
If model_upate is descent_step and the last function evaluation took place for this candidate, the function evaluation data and point identifier (not nec. the point itself) needs to be stored as the new function value of the center.
If model_update is new_subgradient the model may wish to use new subgradient information available in the candidate for inclusion in the model; if the modle is empty it has to include it so that the modle is initialized. The center is not moved in this case.
[in] | model_update |
|
[in] | center_id | the point id of the center, so that the model can check consistency of the evaluation data in moving the center for a descent step |
[in] | center_y | this may help to estimate the model changes relative to this point |
[in] | cand_id | the point id of the candidate, so that the model can check consistency of the evaluation data in moving the center for a descent step |
[in] | cand_y | the next model should be good close to this point |
[in] | model_maxviol | a minorant violated by this ammount would have resulted in a null step |
[in] | H | the weight intended for the proximal term of the next model. it should not be changed here, but may use computation routines of H |
Implements ConicBundle::BundleModel.
Implemented in ConicBundle::SumModel, ConicBundle::AFTModel, and ConicBundle::ConeModel.
Referenced by remove_model().
|
pure virtual |
generate the next cutting model and store the center information in the case of a descent step. This version also allows for updating a common model owned by parent SumModel.
If model_update==null_step , the next model has to contain at least all convex combinations of the current subgradient of eval_function() and the aggregate minorant of (re)eval_augmodel(). update_model may only be called if at least one of both is available.
The point passed in the argument is the candidate where the last function evaluation took place, its id is passed on so that some consistency check is possible without having to store the entire vector. The coordinates of the point may help in setting up a model of good quality close to this point
If model_upate is descent_step and the last function evaluation took place for this candidate, the function evaluation data and point identifier (not nec. the point itself) has to be stored as the function value of the center.
If model_update is new_subgradient the model may wish to use new subgradient information available in the candidate for inclusion in the model; if the modle is empty it has to include it so that the modle is initialized. The center is not moved in this case.
Whenever the model is already initialized before the call to this routine, the output variable model_deviation returns the difference between the value of the local aggregate and the true function in cand_y, otherwise it is set to zero.
If there is a parent SumBundle, it has already been updated except that – if *this is indeed contributing – it needs the information about the new subgradient from *this and the bundlehandler passed in sumbundle_contribution() knows where to add this. The update of the common and local sumbundle now requires the following steps:
If model_update is new_subgradient the model may wish to use new subgradient information available in the candidate for inclusion in the model; if the modle is empty it has to include it so that the modle is initialized. The center is not moved in this case.
The default implementation sets default choices for model_deviation and model_curvature and then calls the routine update_model(ModelUpdate,CH_Matrix_Classes::Integer,const CH_Matrix_Classes::Matrix&,CH_Matrix_Classes::Real,CH_Matrix_Classes::Real).
[in] | model_update |
|
[in] | center_id | the point id of the current center |
[in] | center_y | the coordinates of the current center, its availability may help to compute the relevance of minorants |
[in] | cand_id | the point id of the candidate, so that the model can check consistency of the evaluation data in moving the center for a descent step |
[in] | cand_y | the next model should be good close to this point |
[in] | model_maxviol | a minorant violated by this would have resulted in a null step |
[in] | H | the weight intended for the proximal term of the next model |
[out] | model_deviation | the difference of function value and aggregate in cand_y |
[out] | model_curvature | an estimate of the curvature when violating minorants is allowed by model_maxviol |
Implemented in ConicBundle::SumModel, ConicBundle::AFTModel, and ConicBundle::ConeModel.