#include <MatCBSolver.hxx>
Public Member Functions | |
MatrixCBSolver (bool no_bundle=false) | |
if the input parameter no_bundle is set to true then a minimal bundle variant is used consisting of just one aggregate and one new subgradient in each iteration (this is basically "no bundle", try it whenever fast iterations and/or low memory consumption seem advantageous) | |
Initialization | |
void | clear () |
Clears all data structures and problem information but keeps ouptut settings and algorithmic parameter settings. | |
void | set_defaults () |
Sets default values for algorithmic parameters that are not function specific (e.g., relative precision, weight and weight bounds for the augmentedproblem, etc.). | |
int | init_problem (int dim, const CH_Matrix_Classes::Matrix *lbounds=0, const CH_Matrix_Classes::Matrix *ubounds=0, const CH_Matrix_Classes::Matrix *costs=0) |
Initializes the problem by setting up the design space (the dimension and possible box constraints of the variables). | |
int | add_function (FunctionObject &function) |
Adds a function, typically derived from ConicBundle::FunctionOracle; all functions added must have the same argument dimension set in init_problem(). | |
int | set_lower_bound (int i, double lb) |
Sets lower bound for variable i, use ConicBundle::CB_minus_infinity for unbounded from below. | |
int | set_upper_bound (int i, double ub) |
Sets upper bound for variable i, use ConicBundle::CB_plus_infinity for unbounded from below. | |
int | append_variables (int n_append, const CH_Matrix_Classes::Matrix *lbounds=0, const CH_Matrix_Classes::Matrix *ubounds=0, const CH_Matrix_Classes::Matrix *costs=0) |
Append new variables (always in last postions in this order). | |
int | delete_variables (const CH_Matrix_Classes::Indexmatrix &del_indices, CH_Matrix_Classes::Indexmatrix &map_to_old) |
Deletes variables corresponding to the specified indices. | |
int | reassign_variables (const CH_Matrix_Classes::Indexmatrix &assign_new_from_old) |
Reassigns variables to new index positions by mapping to position i the variable that previously had index assign_new_from_old[i]. | |
Basic algorithmic routines and parameters | |
int | do_descent_step (int maxsteps=0) |
Does a descent step for the current center point. | |
int | termination_code () const |
Returns the termination code of the bundle algorithm for the latest descent step. | |
std::ostream & | print_termination_code (std::ostream &out) |
Outputs a text version of termination code, see termination_code(). | |
double | get_objval () const |
Returns the objective value resulting from last descent step (initially undefined). If no problem modification routines were called since then, it is the objective value at the point returned by get_center(). | |
int | get_center (CH_Matrix_Classes::Matrix ¢er) const |
Returns the next center point that was produced by the latest call to do_descent_step (in some problem modification routines the center point may be updated immediately, in others the center point will be corrected automatically directly before starting the next descent step and its values may be infeasible till then). | |
double | get_sgnorm () const |
Returns Euclidean norm of the latest aggregate subgradient. | |
int | get_subgradient (CH_Matrix_Classes::Matrix &subgradient) const |
Returns the latest aggregate subgradient. | |
double | get_cutval () const |
Returns the cutting model value resulting from last call to do_descent_step() (initially undefined). | |
double | get_candidate_value () const |
Returns the objective value computed in the last step of do_descent_step(), independent of whether this was a descent step or a null step (initially undefined). | |
int | get_candidate (CH_Matrix_Classes::Matrix ¢er) const |
Returns the last point, the "candidate", at which the function was evaluated in do_descent_step(). | |
Advanced algorithmic routines and parameters | |
int | set_term_relprec (const double term_relprec) |
Sets the relative precision requirements for successful termination (default 1e-5). | |
int | set_new_center_point (const CH_Matrix_Classes::Matrix ¢er_point) |
Set the starting point/center that will be used in the next call to do_descent_step(). Each call to this routine causes an immediate evaluation of all oracles. | |
int | get_function_status (const FunctionObject &function) const |
Returns the return value of the latest evaluation call to this function. | |
int | get_approximate_slacks (CH_Matrix_Classes::Matrix &) const |
Returns the multipliers for the box constraints on the design variables; in Lagrangean relaxation they may be interpreted as primal slacks for inequality constraints. | |
int | get_approximate_primal (const FunctionObject &function, PrimalData &primal) const |
returns the current approximate primal solution corresponding to the aggregate subgradient of the specified function. | |
int | get_center_primal (const FunctionObject &function, PrimalData &primal) const |
Returns the primal solution corresponding to the best epsilon subgradient returned in the evaluation of the specified function at the current center point. | |
int | get_candidate_primal (const FunctionObject &function, PrimalData &primal) const |
Returns the primal solution returned by the last evaluation of the specified function in the point get_candidate(). | |
int | set_max_bundlesize (const FunctionObject &function, int max_bundlesize) |
Sets the maximum number of subgradients used in forming the cutting model of the specified function. | |
int | set_max_new_subgradients (const FunctionObject &function, int max_new_subgradients) |
Sets the maximum number of new subgradients to be used in the next bundle update of the cutting modle for the specified . | |
int | set_bundle_parameters (const FunctionObject &function, const BundleParameters ¶ms) |
Sets the maximum bundlesize and the maximum number of new subgradients added in a bundle update of the cutting model for the specified function. The meaning of this routine may differ from standard for predefined special functions with special bundle types. | |
int | get_bundle_parameters (const FunctionObject &function, BundleParameters ¶ms) const |
Retrieves current bundle parameters (not the actual size in use!) as set for the cutting model of the specified function. | |
int | get_bundle_values (const FunctionObject &function, BundleParameters ¶ms) const |
Returns the current bundle values: the current bundle_size and the number of subgradients added in the latest update of the cutting model of the specified function. | |
int | reinit_function_model (const FunctionObject &function) |
Clears cutting model, subgradients and stored function values for the specified function. | |
int | clear_aggregates (const FunctionObject &function) |
Clears the aggregate parts of the cutting model of this function. | |
int | call_primal_extender (const FunctionObject &function, PrimalExtender &primal_extender) |
Asks function to call primal_extender for each of its primal objects (see also FunctionOracle::evaluate() ). | |
double | get_last_weight () const |
Returns the current weight for the quadratic term in the augmented subproblem (may be interpreted as 1./step_size or 1./trustregion-radius). | |
int | set_next_weight (const double weight) |
Sets the weight (>0) to be used in the quadratic term of the next augmented subproblem (may be interpreted as 1./step_size or 1./trustregion-radius). | |
int | set_min_weight (const double min_weight) |
Sets a lower bound on the weight for the quadratic term of the augmented subproblem. | |
int | set_max_weight (const double max_weight) |
Sets an upper bound on the weight for the quadratic term of the augmented subproblem. | |
int | adjust_multiplier (void) |
Adjusts on all conic functions the penalty parameter for conic violations to twice the trace of the primal approximation. | |
int | set_scaling (bool do_scaling) |
Use a scaling heuristic or switch off scaling alltogether. | |
int | set_scaling (const CH_Matrix_Classes::Matrix &scale) |
user defined diagonal scaling, values greater than 1 allow more movement for this variable, values smaller than 1 allow less movement. | |
void | set_active_bounds_fixing (bool allow_fixing) |
If set to true (the default is false), some variables will be fixed automatically to the center value if their bounds are strongly active (i.e., the corresponding multipliers are big). | |
void | clear_fail_counts (void) |
clears all fail counts on numerical function oder model failures, may be useful if this caused premature termination. | |
void | set_eval_limit (int eval_limit) |
Sets an upper bound on the number of calls to the oracle (use negative numbers for no limit). | |
void | set_inner_update_limit (int update_limit) |
Set an upper bound on the number of inner updates for the cutting model with primal slacks within one null step (use negative numbers for no limit). | |
Look up basic paramaters (dimension, number of functions, ...) | |
int | get_dim () const |
Returns the current dimension of the design space/argument or -1 if no dimension is set. | |
int | get_n_functions () const |
Returns the current number of functions in the problem. | |
int | get_n_oracle_calls () const |
Returns the number of function evaluations. | |
int | get_n_descent_steps () const |
Returns the number of function descent setps. | |
int | get_n_inner_iterations () const |
Returns the number of inner iterations of the bundle method. | |
int | get_n_inner_updates () const |
Returns the number of inner multiplier updates for the box constraints. | |
const CH_Matrix_Classes::Matrix & | get_lbounds () const |
Returns the vector of lower bounds. | |
const CH_Matrix_Classes::Matrix & | get_ubounds () const |
Returns the vector of upper bounds. | |
const CH_Matrix_Classes::Indexmatrix & | get_active_bounds_indicator () const |
Returns the indicator vector of variables temporarily fixed to the center value due to significantly positive multipliers for the box constraints. | |
Output | |
void | set_out (std::ostream *out=0, int print_level=1) |
Specifies the output level (out==NULL: no output at all, out!=NULL and level=0: errors and warnings, level>0 increasingly detailed information). | |
std::ostream & | print_line_summary (std::ostream &out) const |
print a one line summary of important evaluation data | |
Private Member Functions | |
MatrixCBSolver (const CBSolver &) | |
not available, blocked deliberately | |
MatrixCBSolver & | operator= (const CBSolver &) |
not available, blocked deliberately | |
Private Attributes | |
MatrixBSolver * | solver |
pointer to internal solver data |
Minimizes the sum of convex functions that are given via ConicBundle::MatrixFunctionOracle interfaces, see the text explaining the C++ interface for Matrix Classes for a quick overview.
It provides special support for Lagrangean relaxation by generating primal approximate solutions if such information is provided in the function oracles.
Based on these primal approximations it is also possible to implement cutting plane schemes. Routines for adding and deleting corresponding dual variables as well as a framework for extending subgradients in order not to loose the cutting model are available.
int ConicBundle::MatrixCBSolver::init_problem | ( | int | dim, | |
const CH_Matrix_Classes::Matrix * | lbounds = 0 , |
|||
const CH_Matrix_Classes::Matrix * | ubounds = 0 , |
|||
const CH_Matrix_Classes::Matrix * | costs = 0 | |||
) |
Initializes the problem by setting up the design space (the dimension and possible box constraints of the variables).
Clears all data structures and sets the dimension @ m for a new problem. for solving min_{y in R^m} f_0(y) + f_1(y) + ... Box constraints may be specified for y. (The functions f_i must be added by add_function()).
Lower and/or upper bounds must be speicified for all variables or for none of them. To specify no bounds at all, give Null pointers. Otherwise use ConicBundle::CB_minus_infinity for unbounded below and ConicBundle::CB_plus_infinity for unbounded above. For NULL pointers, unbounded will be used as default for all variables. Specifying bounds selectively is also possible by set_lower_bound() or set_upper_bound().
[in] | dimm | (int) the dimension of the argument/design space/the number of Lagrange multipliers |
[in] | lbounds | (const Matrix*) If NULL, all variables are considered unbounded below, otherwise lowerb[i] gives the minimum feasible value for variable y[i], use ConicBundle::CB_minus_infinity for unbounded below. |
[in] | ubounds | (const Matrix*) If NULL, all variables are considered unbounded above, otherwise upperb[i] gives the maximum feasible value for variable y[i], use ConicBundle::CB_plus_infinity for unbounded above. |
[in] | costs | (const Matrix*) Use this in order to specify linear costs on the variables in addition to the functions (may be convenient in Lagrangean relaxation for the right hand side of coupling contsraints); NULL is equivalent to costs zero. |
int ConicBundle::MatrixCBSolver::add_function | ( | FunctionObject & | function | ) |
Adds a function, typically derived from ConicBundle::FunctionOracle; all functions added must have the same argument dimension set in init_problem().
Besides the standard ConicBundle::MatrixFunctionOracle the interface only accepts a few other prespecified derivations of the class FunctionObject that come along with the CH_Matrix_Classes interface (e.g. for semidefinite and second order cones). Functions not derived from these will fail to be added and return a value !=0.
int ConicBundle::MatrixCBSolver::set_lower_bound | ( | int | i, | |
double | lb | |||
) |
Sets lower bound for variable i, use ConicBundle::CB_minus_infinity for unbounded from below.
The algorithm may have to adapt the center point aftwards. In this case the old function values will be marked as outdated and will be recomputed at the next call to e.g. do_descent_step().
int ConicBundle::MatrixCBSolver::set_upper_bound | ( | int | i, | |
double | ub | |||
) |
Sets upper bound for variable i, use ConicBundle::CB_plus_infinity for unbounded from below.
The algorithm may have to adapt the center point aftwards. In this case the old function values will be marked as outdated and will be recomputed at the next call to e.g. do_descent_step().
int ConicBundle::MatrixCBSolver::append_variables | ( | int | n_append, | |
const CH_Matrix_Classes::Matrix * | lbounds = 0 , |
|||
const CH_Matrix_Classes::Matrix * | ubounds = 0 , |
|||
const CH_Matrix_Classes::Matrix * | costs = 0 | |||
) |
Append new variables (always in last postions in this order).
If 0 is feasible for the new coordinates then this is selected as starting value for the new coordinates; otherwise, the number closest to zero is used. If all new coordinates can be set to zero then it is assumed that for an existing center point the function values need not be recomputed (this is e.g. the case in Lagrangean relaxation; if this is not correct call reinit_function_model() below). Otherwise the old function values will be marked as outdated and will be recomputed at the next call to e.g. do_descent_step().
[in] | n_append | (int) number of variables to append (always in last position in the same order) |
[in] | lbounds | (const Matrix*) If NULL, all appended variables are considered unbounded below, otherwise lowerb[i] gives the minimum feasible value for variable y[i], use ConicBundle::CB_minus_infinity for unbounded below. |
[in] | ubounds | (const Matrix*) If NULL, all appended variables are considered unbounded above, otherwise upperb[i] gives the maximum feasible value for variable y[i], use ConicBundle::CB_plus_infinity for unbounded above. |
int ConicBundle::MatrixCBSolver::delete_variables | ( | const CH_Matrix_Classes::Indexmatrix & | del_indices, | |
CH_Matrix_Classes::Indexmatrix & | map_to_old | |||
) |
Deletes variables corresponding to the specified indices.
The indices of the remaining variables are reassigned so that they are consecutive again, the routine returns in map_to_old a vector giving for each new index of these remaining variables the old coordinate.
If all of the deleted variables are zero, function values are assumed to remain correct (if this is not so, call reinit_function_model() below) Otherwise the old function values will be marked as outdated and will be recomputed at the next call to e.g. do_descent_step().
[in] | delete_indices | (const Indexmatrix&) the entries delete_indices[i] specify the indices of the variables to be deleted |
[out] | map_to_old | (Indexmatrix&) after the call, element map_to_old[i] gives the old index (before the call) of the variable that now has index position i. |
int ConicBundle::MatrixCBSolver::reassign_variables | ( | const CH_Matrix_Classes::Indexmatrix & | assign_new_from_old | ) |
Reassigns variables to new index positions by mapping to position i the variable that previously had index assign_new_from_old[i].
Old variables, that are not mapped to any position will be deleted. It is allowed to generate several copies of old variables.
If all of the deleted variables as well as new multiple copies are zero, function values are assumed to remain correct (if this is not so, call reinit_function_model() below). Otherwise the old function values will be marked as outdated and will be recomputed at the next call to e.g. do_descent_step().
[in] | assign_new_from_old | (const IVector&) entry assign_new_from_old[i] specifies the old index of the variable, that has to be copied to index position i. |
int ConicBundle::MatrixCBSolver::do_descent_step | ( | int | maxsteps = 0 |
) |
Does a descent step for the current center point.
A descent step may consist of several function evaluations (null steps), that lead to no immediate progress but serve for building a cutting model of the objective function close to the current center point. A minimizer to the model is accepted as descent step if the function value at this point satisfies a sufficient decrease criterion in comparison to the decrease predicted by the model. Having found a descent step, the next center is automatically shifted to this successful candidate. Termination criteria may stop the process of seeking for a descent step, in which case the current center is kept and the routine termination_code() returns the termination code.
Restarting, after each descent step, the bundle method from scratch with the new center as starting point does not endanger convergence. Therefore, a descent step is the smallest unit, after which user interaction can take place safely and this is the default choice.
If you know what your are doing, you may also use the input parameter maxsteps to force the algorithm to return after at most maxsteps null steps. Calling do_descent_step again without any intermediate problem configurations will then simply continue the process where it stopped and convergence is save. During null steps one may not decrease the weight or delete nonzero variables of the center or the current candidate!
In a Lagrangean relaxation cutting plane approach one may want to separate and enlarge the dimension after a certain number of null steps. In this case the code will try to preserve the model, given appropriate subgradient extension routines have been provided. If the model cannot be extended, it has to be discarded (if subgradient extension is not successful this is done automatically), and the algorithm will be restarted from the current center point.
[in] | maxsteps | (int) if maxsteps>0 the code returns after at most so many null steps |
int ConicBundle::MatrixCBSolver::termination_code | ( | ) | const |
Returns the termination code of the bundle algorithm for the latest descent step.
For resetting all counters relevant for termination see clear_fail_counts() .
std::ostream& ConicBundle::MatrixCBSolver::print_termination_code | ( | std::ostream & | out | ) |
Outputs a text version of termination code, see termination_code().
int ConicBundle::MatrixCBSolver::get_center | ( | CH_Matrix_Classes::Matrix & | center | ) | const |
Returns the next center point that was produced by the latest call to do_descent_step (in some problem modification routines the center point may be updated immediately, in others the center point will be corrected automatically directly before starting the next descent step and its values may be infeasible till then).
int ConicBundle::MatrixCBSolver::get_subgradient | ( | CH_Matrix_Classes::Matrix & | subgradient | ) | const |
Returns the latest aggregate subgradient.
double ConicBundle::MatrixCBSolver::get_candidate_value | ( | ) | const |
Returns the objective value computed in the last step of do_descent_step(), independent of whether this was a descent step or a null step (initially undefined).
If no problem modification routines were called since then, it is the objective value at the point returned by get_candidate(). If this last evaluation led to a descent step, then it is the same value as in get_objval().
int ConicBundle::MatrixCBSolver::get_candidate | ( | CH_Matrix_Classes::Matrix & | center | ) | const |
Returns the last point, the "candidate", at which the function was evaluated in do_descent_step().
If this evaluation lead to a descent step, it is the same point as in get_center().
int ConicBundle::MatrixCBSolver::set_term_relprec | ( | const double | term_relprec | ) |
Sets the relative precision requirements for successful termination (default 1e-5).
[in] | term_relprec | (double) The algorithm stops with termination code 1, if predicted progress for the next step is less than term_relprec times absolute function value plus one. |
int ConicBundle::MatrixCBSolver::set_new_center_point | ( | const CH_Matrix_Classes::Matrix & | center_point | ) |
Set the starting point/center that will be used in the next call to do_descent_step(). Each call to this routine causes an immediate evaluation of all oracles.
int ConicBundle::MatrixCBSolver::get_approximate_slacks | ( | CH_Matrix_Classes::Matrix & | ) | const |
Returns the multipliers for the box constraints on the design variables; in Lagrangean relaxation they may be interpreted as primal slacks for inequality constraints.
int ConicBundle::MatrixCBSolver::get_approximate_primal | ( | const FunctionObject & | function, | |
PrimalData & | primal | |||
) | const |
returns the current approximate primal solution corresponding to the aggregate subgradient of the specified function.
PrimalData solutions must have been supplied in all previous calls to evaluate; In this case it returns the current approximate primal solution aggregated alongside with the aggregate subgradient. A primal solution may not be available after addition of constraints, if extension of the aggregate subgradient to the new coordinates failed.
int ConicBundle::MatrixCBSolver::get_center_primal | ( | const FunctionObject & | function, | |
PrimalData & | primal | |||
) | const |
Returns the primal solution corresponding to the best epsilon subgradient returned in the evaluation of the specified function at the current center point.
PrimalData solutions must have been supplied in all previous calls to evaluate; It may not be available or may correspond to an aggregate primal after addition or deletion of design variables/primal constraints.
int ConicBundle::MatrixCBSolver::get_candidate_primal | ( | const FunctionObject & | function, | |
PrimalData & | primal | |||
) | const |
Returns the primal solution returned by the last evaluation of the specified function in the point get_candidate().
It will only be available if also supplied by the function
int ConicBundle::MatrixCBSolver::set_max_bundlesize | ( | const FunctionObject & | function, | |
int | max_bundlesize | |||
) |
Sets the maximum number of subgradients used in forming the cutting model of the specified function.
Quite often a very small model, e.g., 2, yields very fast iterations and good progress in time (sometimes at the cost of more evaluations). By limited numerical experience, a significant reduction in the number of evaluations can only be expected if the bundle is large enough to wrap the function rather tightly. Quite frequently, unfortunately, this entails that solving the quadratic subproblems is more expensive than function evaluation.
The meaning of this routine may differ from standard for predefined special functions with special bundle types.
[in] | function | (const FunctionObject&) the function added in add_function() |
[in] | bundlesize | (int) maximum number of subgradients to be used in forming the cutting model |
int ConicBundle::MatrixCBSolver::set_max_new_subgradients | ( | const FunctionObject & | function, | |
int | max_new_subgradients | |||
) |
Sets the maximum number of new subgradients to be used in the next bundle update of the cutting modle for the specified .
The meaning of this routine may differ from standard for predefined special functions with special bundle types.
[in] | function | (const FunctionObject&) the function added in add_function() |
[in] | max_new_subgradients | (int) maximum number of new epsilon subgradients to be used in bundle updates |
int ConicBundle::MatrixCBSolver::set_bundle_parameters | ( | const FunctionObject & | function, | |
const BundleParameters & | params | |||
) |
Sets the maximum bundlesize and the maximum number of new subgradients added in a bundle update of the cutting model for the specified function. The meaning of this routine may differ from standard for predefined special functions with special bundle types.
[in] | function | (const FunctionObject&) the function added in add_function() |
[in] | params | (const BundleParameters&) some update parameters for the cutting model, see e.g. ConicBundle::BundleParameters |
int ConicBundle::MatrixCBSolver::get_bundle_parameters | ( | const FunctionObject & | function, | |
BundleParameters & | params | |||
) | const |
Retrieves current bundle parameters (not the actual size in use!) as set for the cutting model of the specified function.
This may differ for predefined special functions with derived BundleParameter classes.
[in] | function | (const FunctionObject&) the function added in add_function() |
[out] | params | (BundleParameters&) |
int ConicBundle::MatrixCBSolver::get_bundle_values | ( | const FunctionObject & | function, | |
BundleParameters & | params | |||
) | const |
Returns the current bundle values: the current bundle_size and the number of subgradients added in the latest update of the cutting model of the specified function.
This may differ for predefined special functions with derived BundleParameter classes.
[in] | function | (const FunctionObject&) the function added in add_function() |
[out] | params | (BundleParameters&) |
int ConicBundle::MatrixCBSolver::reinit_function_model | ( | const FunctionObject & | function | ) |
Clears cutting model, subgradients and stored function values for the specified function.
This has to be called whenever the specified function was modified so that the old subgradients and/or primal generators are no longer valid.
[in] | function | (const FunctionObject&) the function added in add_function() |
int ConicBundle::MatrixCBSolver::clear_aggregates | ( | const FunctionObject & | function | ) |
Clears the aggregate parts of the cutting model of this function.
This has to be called whenever the specified function was modified so that the old aggregate subgradients and/or primal generators are no longer valid.
[in] | function | (const FunctionObject&) the function added in add_function() |
int ConicBundle::MatrixCBSolver::call_primal_extender | ( | const FunctionObject & | function, | |
PrimalExtender & | primal_extender | |||
) |
Asks function to call primal_extender for each of its primal objects (see also FunctionOracle::evaluate() ).
If the function is the Lagrangian dual of a primal problem and primal_data returned previous calls to the oracle has now to be updated due to changes in the primal problem -- e.g., this may happen in column generation -- the call causes updates of all internally stored primal_data objects by calling PrimalExtender::extend on each of these.
[in] | function | (const FunctionObject&) the function added in add_function() |
[in] | primal_extender | (PrimalExtender&) the object holding the extension function for primal_data |
int ConicBundle::MatrixCBSolver::set_next_weight | ( | const double | weight | ) |
Sets the weight (>0) to be used in the quadratic term of the next augmented subproblem (may be interpreted as 1./step_size or 1./trustregion-radius).
Independent of whether the weight violates current min- and max-bounds set in set_min_weight() and set_max_weight(), the next model will be computed for this value. Thereafter, however, it will be updated as usual; in particular, it may be truncated by min and max bounds immediately after the first subproblem.
In order to guarantee a constant weight (e.g. 1 is frequently a reasonable choice if the automatic default heuristic performs poorly), set the min and max bounds to the same value, too.
[in] | weight | (double) |
int ConicBundle::MatrixCBSolver::set_min_weight | ( | const double | min_weight | ) |
Sets a lower bound on the weight for the quadratic term of the augmented subproblem.
Nonpositive values indicate no bound. The new value shows its effect only at first dynamic change of the weight.
[in] | min_weight | (double) |
int ConicBundle::MatrixCBSolver::set_max_weight | ( | const double | max_weight | ) |
Sets an upper bound on the weight for the quadratic term of the augmented subproblem.
Nonpositive values indicate no bound. The new value shows its effect only at first dynamic change of the weight.
[in] | max_weight | (double) |
int ConicBundle::MatrixCBSolver::adjust_multiplier | ( | void | ) |
Adjusts on all conic functions the penalty parameter for conic violations to twice the trace of the primal approximation.
This routine is only needed for conic function objects such as the nonnegative cone, the second order cone and the semidefinite cone if no good upper bound on the trace of feasible points is known and has to be determined automatically.
If after some time, the trace values settle, the upper bounds on the trace may be way to high and can then be reset with this call.
int ConicBundle::MatrixCBSolver::set_scaling | ( | bool | do_scaling | ) |
Use a scaling heuristic or switch off scaling alltogether.
int ConicBundle::MatrixCBSolver::set_scaling | ( | const CH_Matrix_Classes::Matrix & | scale | ) |
user defined diagonal scaling, values greater than 1 allow more movement for this variable, values smaller than 1 allow less movement.
It is the users responsibility to guarantee that the scaling vector fits in dimension to the runnig problem data, in particular if routines such as append_variables(), delete_variables(), and reassign_variables() are used.
[in] | scale | (const Matrix&) |
void ConicBundle::MatrixCBSolver::set_active_bounds_fixing | ( | bool | allow_fixing | ) |
If set to true (the default is false), some variables will be fixed automatically to the center value if their bounds are strongly active (i.e., the corresponding multipliers are big).
The coordinates to be fixed are redetermined in each call following a descent step or a change of the function. An indicator vector of the variables fixed in the last call can be obtained via the routine get_active_bounds_indicator().
Setting this value to true might improve the performance of the algorithm in some instances but there is no convergence theory. It might be particularly helpful within Lagrangian relaxation if a primal cutting plane approach is used and non-tight inequalities should be eliminated quickly (fixing then indicates large primal slack values).
[in] | allow_fixing | (bool) |
void ConicBundle::MatrixCBSolver::clear_fail_counts | ( | void | ) |
clears all fail counts on numerical function oder model failures, may be useful if this caused premature termination.
void ConicBundle::MatrixCBSolver::set_eval_limit | ( | int | eval_limit | ) |
Sets an upper bound on the number of calls to the oracle (use negative numbers for no limit).
If this number is reached, the algorithm will terminate independently of whether the last step was a descent or a null step. A negative number will be interepreted as no limit.
[in] | eval_limit | (Integer) |
void ConicBundle::MatrixCBSolver::set_inner_update_limit | ( | int | update_limit | ) |
Set an upper bound on the number of inner updates for the cutting model with primal slacks within one null step (use negative numbers for no limit).
A negative number will be interepreted as no limit, i.e., the updates will be done till a certain precision of the cutting model is achieved.
[in] | update_limit | (Integer) |
const CH_Matrix_Classes::Indexmatrix& ConicBundle::MatrixCBSolver::get_active_bounds_indicator | ( | ) | const |
Returns the indicator vector of variables temporarily fixed to the center value due to significantly positive multipliers for the box constraints.
Such a fixing indicates that the corresponding variables would like to stay at their bounds. If no variables were fixed, the dimension of the vector is zero.
void ConicBundle::MatrixCBSolver::set_out | ( | std::ostream * | out = 0 , |
|
int | print_level = 1 | |||
) |
Specifies the output level (out==NULL: no output at all, out!=NULL and level=0: errors and warnings, level>0 increasingly detailed information).
[in] | out | (ostream*) direct all output to (*out). If out==NULL, there will be no output at all. |
[in] | print_level | (int) |
Example for level 1:
00:00:00.00 endit 1 1 1 563. 563. 39041.188 39043.162 00:00:00.00 endit 2 2 2 563. 559. 38488.165 38490.200 00:00:00.00 endit 3 3 3 56.3 555. 33014.533 33211.856 00:00:00.00 endit 4 4 4 5.63 517. -14306.459 2738.0343 00:00:00.00 endit 5 5 5 4.04 148. -2692.1131 2.2150883 00:00:00.00 endit 6 6 6 4.01 1.29 1.7908952 2.0000581 00:00:00.00 endit 7 7 7 3.95 0.0213 1.9999387 2.0000000 00:00:00.00 _endit 8 8 8 3.95 2.94e-05 2.0000000 2.0000000 Column 1 2 3 4 5 6 7 8 9