Functions
EAGO._add_constraint_store_ci! — Method
_add_constraint_store_ci!(
d,
m::InputProblem,
_::Type{F},
_::Type{S}
)
Add a constraint to the local problem, storing the new constraint index and the associated index in the input problem.
EAGO._add_constraint_store_ci_linear! — Method
_add_constraint_store_ci_linear!(d, ip::InputProblem)
Add linear constraints to the local problem, storing the new constraint indices and the associated indices in the input problem.
EAGO._add_constraint_store_ci_quadratic! — Method
_add_constraint_store_ci_quadratic!(d, ip::InputProblem)
Add quadratic constraints to the local problem, storing the new constraint indices and the associated indices in the input problem.
EAGO._constraint_primal — Method
_constraint_primalHelper function which simplifies finding constraints of different types. See also: _constraints
EAGO._constraints — Method
_constraintsHelper function which simplifies finding constraints of different types. See also: _constraint_primal
EAGO._diam — Method
_diam(::BranchVar, ::GlobalOptimizer, ::Int)
_diam(::FullVar, ::GlobalOptimizer, ::Int)Return the diameter of a variable (upper bound - lower bound).
EAGO._extract_primal! — Method
_extract_primal!(d, m::InputProblem, _::Type{F}, _::Type{S})
Extract primal constraint value from the local problem and save the result to the appropriate field of the optimizer.
EAGO._extract_primal_linear! — Method
_extract_primal_linear!(d, ip::InputProblem)
Extract linear primal constraint values from the local problem and save the result to the appropriate field of the optimizer.
EAGO._extract_primal_quadratic! — Method
_extract_primal_quadratic!(d, ip::InputProblem)
Extract quadratic primal constraint values from the local problem and save the result to the appropriate field of the optimizer.
EAGO._is_branch_var — Method
_is_branch_var(m::GlobalOptimizer, i) -> Any
Check if m._branch_variables[i] is true or false.
EAGO._is_input_min — Method
_is_input_min(m::GlobalOptimizer) -> Bool
Check to see if the sense of the optimization problem is Min (true) or Max (false).
EAGO._relaxed_optimizer — Method
_relaxed_optimizer(::GlobalOptimizer)
_relaxed_optimizer(::SubSolvers)Return the relaxed optimizer (<: MOI.AbstractOptimizer) from the SubSolvers field or object.
EAGO._upper_optimizer — Method
_upper_optimizer(::GlobalOptimizer)
_upper_optimizer(::SubSolvers)Return the upper optimizer (<: MOI.AbstractOptimizer) from the SubSolvers field or object.
EAGO.add_nonlinear! — Method
add_nonlinear!(m::GlobalOptimizer)
Add an Evaluator and nonlinear functions and populate each appropriately.
EAGO.add_nonlinear_evaluator! — Method
add_nonlinear_evaluator!Add an Evaluator structure if nonlinear terms are attached.
EAGO.add_objective! — Method
add_objective!Add objective function (if any) to the parsed problem.
EAGO.affine_relax_quadratic! — Method
affine_relax_quadratic!Default routine for relaxing quadratic constraint func < 0.0 on node n. Takes affine bounds of convex part at point x0 and secant line bounds on concave parts. This approach comes from Section 2.4.3.2 of: Vigerske, S. and Gleixner, A. "SCIP: global optimization of mixed-integer nonlinear programs in a branch-and-cut framework". Optimization Methods and Software, 33:3, 563-593 (2018). DOI: 10.1080/10556788.2017.1335312.
EAGO.aggressive_filtering! — Method
aggressive_filtering!(
m::GlobalOptimizer{R, S, Q<:ExtensionType},
n::NodeBB
) -> Bool
Excludes OBBT on variable indices after a search in a filtering direction.
EAGO.bool_indx_diff! — Method
bool_indx_diff!(
z::Vector{Bool},
x::Vector{Bool},
y::Vector{Bool}
) -> Vector{Bool}
Utility function used to set vector of booleans z to x & ~y. Avoids the generation of conversion of the BitArray created by broadcasting logical operators.
EAGO.bound_objective — Method
bound_objective(::GlobalOptimizer, ::T)
bound_objective(::ExtensionType, ::GlobalOptimizer)Compute a tuple representing the lower and upper bounds for an objective function. Note: bound_objective(::GlobalOptimizer) dispatches to bound_objective(::ExtensionType, ::GlobalOptimizer).
Options for T:
- AffineFunctionIneq
- BufferedNonlinearFunction
- BufferedQuadraticIneq
- MOI.VariableIndex
EAGO.branch_node! — Method
branch_node!(t::ExtensionType, m::GlobalOptimizer)
Create two nodes from current_node and store them on the stack. Call select_branch_variable(t, m) and select_branch_point(t, m, k) to determine the variable that should be branched on and the point at which branching should occur, respectively.
EAGO.build_model — Method
build_modelCreate the model and variables used with extension t::EAGO.ExtensionType in algorithm a::AbstractSIPAlgo in subproblem s::AbstractSubproblemType via the command build_model(t::ExtensionType, a::AbstractSIPAlgo, s::AbstractSubproblemType, p::SIPProblem).
EAGO.check_set_affine_nl! — Method
check_set_affine_nl!(
m::GlobalOptimizer{R, S, Q<:ExtensionType},
f::EAGO.BufferedNonlinearFunction{V, N, T<:RelaxTag},
finite_cut::Bool,
check_safe::Bool
) -> Bool
EAGO.convergence_check — Method
convergence_check(
t::ExtensionType,
m::GlobalOptimizer
) -> Bool
Check for problem convergence.
By default, check if the lower and upper bounds have converged to within absolute and/or relative tolerances.
EAGO.create_buffer_dict — Method
create_buffer_dict(
func::MathOptInterface.ScalarQuadraticFunction{Float64}
) -> Dict{Int64, Float64}
Create a buffer dictionary from a ScalarQuadraticFunction.
EAGO.cut — Method
cut
Intersects the new set valued operator with the prior and performs affine bound tightening
- First forward pass:
postshould be set by user option,is_intersectshould be false so that the tape overwrites existing values, and theinterval_intersectflag could be set to either value. - Forward CP pass (assumes same reference point):
postshould be set by user option,is_intersectshould be true so that the tape intersects with existing values, and theinterval_intersectflag should be false. - Forward CP pass (assumes same reference point):
postshould be set by user option,is_intersectshould be true so that the tape intersects with existing values, and theinterval_intersectflag should be false. - Subsequent forward passes at new points: post
should be set by user option,isintersectshould be true so that the tape intersects with existing values, and theintervalintersectflag should betrue` as predetermined interval bounds are valid but the prior values may correspond to different points of evaluation.
EAGO.cut_condition — Method
cut_condition(t::ExtensionType, m::GlobalOptimizer) -> Bool
Returns true if a cut should be added and computes a new reference point to add the cut at. By default, checks that cut_max_iterations are not exceeded and that the improvement in the objective value associated with the previous cut is greater than both an absolute tolerance cut_ϵ_abs and a relative tolerance cut_ϵ_rel. Returns false otherwise.
EAGO.default_upper_heuristic — Method
default_upper_heuristic(m::GlobalOptimizer) -> Bool
Default check to see if the upper bounding problem should be run. By default, The upper bounding problem is run on every node up to depth upper_bounding_depth and is triggered with a probability of 0.5^(depth - upper_bounding_depth) afterwards for continuous problems. For integral problems, the upper_bounding_depth approach is used as well as running on every node up to depth upper_bounding_depth + cont_depth with another trigger of probability 0.5^(depth - upper_bounding_depth - cont_depth).
EAGO.eliminate_fixed_variables! — Function
eliminate_fixed_variables!(::T, ::Vector{VariableInfo})Eliminate fixed variables by rearrangment or restructuring of the AbstractEAGOConstraint.
Options for T (all are subtypes of AbstractEAGOConstraint):
- AffineFunctionEq
- AffineFunctionIneq
- BufferedQuadraticIneq
- BufferedNonlinearFunction{N,T} where {N, T<:RelaxTag}
- NonlinearExpression{V,N,T} where {V, N, T<:RelaxTag}
EAGO.f_init! — Method
Initializes information in cache c for each node in g that may be used in a forward-pass of attribute t.
EAGO.fathom! — Method
fathom!(t::ExtensionType, m::GlobalOptimizer)
Remove nodes from the stack. By default, delete nodes from the stack if their lower bounds are greater than the current global upper bound.
EAGO.fbbt! — Function
fbbt!(m::GlobalOptimizer, f::T)Performs feasibility-based bound tightening on a back-end constraint and returns true if it is feasible or false if it is infeasible.
Options for T (all are subtypes of AbstractEAGOConstraint):
- AffineFunctionIneq
- AffineFunctionEq
EAGO.fprop! — Method
Populates information associated with attribute t for a constant v at index k in cache c associated with graph g using information at index k.
EAGO.fprop! — Method
Populates information associated with attribute t for a expression v at index k in cache c associated with graph g using information taken from the children of k.
EAGO.fprop! — Method
Populates information associated with attribute t for a parameter v at index k in cache c associated with graph g using information at index k.
EAGO.fprop! — Method
Populates information associated with attribute t for a subexpression v at index k in cache c associated with graph g using information taken from the children of k.
EAGO.fprop! — Method
Populates information associated with attribute t for a variable v at index k in cache c associated with graph g using information at index k.
EAGO.get_sip_optimizer — Method
get_sip_optimizerSpecifices the optimizer to be used in extension t::EAGO.ExtensionType with algorithm alg::AbstractSIPAlgo in subproblem s::AbstractSubproblemType via the command get_sip_optimizer(t::ExtensionType, alg::AbstractSIPAlgo, s::AbstractSubproblemType).
EAGO.global_solve! — Method
global_solve!(m::GlobalOptimizer)
Solves the branch-and-bound problem with the input EAGO.GlobalOptimizer object.
Pseudocode description of the algorithm, as implemented here:
-I) Prepare optimizers and stack for branch-and-bound
-II) While no reason to terminate the algorithm has occurred:
–-II.A) Fathom nodes from the stack
–-II.B) Select the new "current node" from the stack
–-II.C) Perform preprocessing on current node
–-II.D) If preprocessing result is feasible:
––-II.D.1) Solve lower problem for current node
––-II.D.2) If lower problem result is feasible and lower/upper bounds have not converged:
–––-II.D.2.a) Solve upper problem for current node
–––-II.D.2.b) Update the global upper bound if necessary
–––-II.D.2.c) Perform postprocessing
–––-II.D.2.d) If postprocessing result is feasible:
––––-II.D.2.d.α) Branch and add the nodes back to the stack
–-II.E) Update the global lower bound if necessary
–-II.F) Update log information
-III) Set termination and result statuses
-IV) Print solution
EAGO.initial_parse! — Method
initial_parse!(m::Optimizer{R, S, T})
Translate the input problem to the working problem. Any checks or optional manipulation are left to the presolve stage.
EAGO.initialize! — Method
initialize!(_::AbstractCache, _::AbstractDirectedGraph)
Function used to initialize the storage cache d::AbstractCache for a given type of directed acyclic graph g::AbstractDirectedGraph.
EAGO.initialize_stack! — Method
initialize_stack!(t::ExtensionType, m::GlobalOptimizer)
Prepare the stack for the branch-and-bound routine. By default, create an initial node with the variable bounds as box constraints and add it to the stack.
EAGO.interval_bound — Function
interval_bound(::GlobalOptimizer, ::T)Compute a tuple representing the lower and upper interval bounds for an AbstractEAGOConstraint representing an equality constraint.
Options for T (all are subtypes of AbstractEAGOConstraint):
- AffineFunctionEq
- AffineFunctionIneq
- BufferedQuadraticEq
- BufferedQuadraticIneq
- BufferedNonlinearFunction{V,N,T} where {V,N,T}
EAGO.interval_objective_bound! — Function
EAGO.is_feasible — Method
is_feasible(::GlobalOptimizer, ::T)Check if a given bound is feasible.
Options for T (all are subtypes of AbstractEAGOConstraint):
- AffineFunctionIneq
- AffineFunctionEq
- BufferedQuadraticIneq
- BufferedQuadraticEq
- BufferedNonlinearFunction{V,N,T}
EAGO.is_integer_feasible_local — Method
is_integer_feasible_local(m::GlobalOptimizer, d) -> Bool
Checks that the solution of a local solve is integer feasible to within the tolerances specified by integer_abs_tol and integer_rel_tol.
EAGO.is_integer_feasible_relaxed — Method
is_integer_feasible_relaxed(m::GlobalOptimizer) -> Bool
Check that the solution of the lower (relaxed problem) is integer feasible to within tolerances specified by the parameters: integer_abs_tol (absolute tolerance) and integer_rel_tol (relative tolerance).
EAGO.is_integer_subproblem — Method
is_integer_subproblem(m)
Returns true that the subproblem at the current node n has participating integer variables that have not been fixed to constant valued as the branch-and-bound algorithm progresses. Returns false otherwise.
EAGO.is_safe_cut! — Method
is_safe_cut!(
m::GlobalOptimizer,
f::MathOptInterface.ScalarAffineFunction{Float64}
) -> Bool
Applies the safe cut checks detailed in Khajavirad, 2018 [Khajavirad, Aida, and Nikolaos V. Sahinidis. "A hybrid LP/NLP paradigm for global optimization relaxations." Mathematical Programming Computation 10.3 (2018): 383-421] to ensure that only numerically safe affine relaxations are added. Checks that:
|b| <= safe b,safe_l <= abs(ai) <= safe u, andsafe_l <= abs(ai/aj) <= safe_u.
EAGO.label_branch_variables! — Method
label_branch_variables!(m::GlobalOptimizer)
Detect any variables participating in nonconvex terms and populate the _branch_variables storage array.
EAGO.label_fixed_variables! — Method
label_fixed_variables!(m::GlobalOptimizer) -> Vector{Bool}
Detect any variables set to a fixed value by equality or inequality constraints and populate the _fixed_variable storage array.
EAGO.load_fbbt_buffer! — Method
load_fbbt_buffer!(m::GlobalOptimizer)
EAGO.load_relaxed_problem! — Method
load_relaxed_problem!(
m::GlobalOptimizer{R, S, Q<:ExtensionType}
)
Load variables, linear constraints, and empty storage space for the first NLP and quadratic cut into the relaxed optimizer.
EAGO.local_problem_status — Method
local_problem_status(t, r)
Takes an MOI.TerminationStatusCode and a MOI.ResultStatusCode and returns true if this corresponds to a solution that is proven to be feasible. Returns false otherwise.
EAGO.log_iteration! — Method
log_iteration!(m::GlobalOptimizer)
If log_on is true, the global_lower_bound, global_upper_bound, run_time, and node_count are stored every log_interval. If log_subproblem_info then the lower bound, feasibility and run times of the subproblems are logged every log_interval.
EAGO.lower_interval_bound — Function
lower_interval_bound(::GlobalOptimizer, ::T)Compute the lower interval bound for an AbstractEAGOConstraint representing an inequality constraint.
Options for T (all are subtypes of AbstractEAGOConstraint):
- AffineFunctionIneq
- BufferedQuadraticIneq
- BufferedSOC
- BufferedNonlinearFunction{V,N,T} where {V,N,T}
EAGO.lower_problem! — Method
lower_problem!(
t::ExtensionType,
m::GlobalOptimizer{R, S, Q<:ExtensionType}
)
Constructs a relaxation of the MINLP on node y and solves it using the default EAGO relaxation scheme. By default, EAGO applies Kelley's algorithm (from Kelley Jr., J.E.: The cutting-plane method for solving convex programs. J. Soc. Ind. Appl. Math. 8(4), 703 to 712 (1960)) while cut_condition(m) returns true then activates the integrality constraints of the relaxed problems and solves the resulting MILP relaxation. results are stored to the _lower_solution, _lower_termination_status, _lower_primal_status, _lower_dual_status, _lower_objective_value, and _lower_feasibility. Further, lower and upper variable duals are stored _lower_lvd and _lower_uvd, respectively, for use in duality based bound tightening. If relaxation-based bounds are weaker or cutting-planes are numerically poorly ill-posed, then interval bounds are used instead. If the problem is dual feasible but the primal status is ambiguous the dual objective value is used for the lower bound to avoid numerical issues.
EAGO.mc_type — Method
mc_type(::RelaxCache{V, N, T} where {V, N, T<:RelaxTag}
mc_type(::NonlinearExpression{V, N, T} where {V, N, T<:RelaxTag}
mc_type(::BufferedNonlinearFunction{V, N, T} where {V, N, T<:RelaxTag}Returns a McCormick structure of type MC{N, T<:RelaxTag}.
EAGO.node_selection! — Method
node_selection!(t::ExtensionType, m::GlobalOptimizer)
Select the next node in the stack to evaluate. By default, perform best-first node selection (select the node with the lowest lower bound in the stack).
EAGO.obbt! — Method
obbt!(m::GlobalOptimizer{R, S, Q<:ExtensionType}) -> Bool
Performs OBBT with filtering and greedy ordering as detailed in: Gleixner, A.M., Berthold, T., Müller, B. et al. J Glob Optim (2017) 67: 731. https://doi.org/10.1007/s10898-016-0450-4
EAGO.objective_cut! — Method
objective_cut!
Add linear objective cut constraint to the m._subsolvers.relaxed_optimizer.
EAGO.optimize_hook! — Method
optimize_hook!(t::ExtensionType, m::Optimizer)Provide a hook for extensions to EAGO.
The user-defined extension of optimize_hook! is used in EAGO's overloading of MOI.optimize! (see EAGO.jl/src/eago_optimizer/optimize/optimize.jl). Without the optimize_hook! specified, EAGO will run initial_parse!, parse_classify_problem!, and then optimize! using the parsed problem type. The user-specified optimize_hook! should thus take the new extension and Optimizer as inputs and will execute when the user writes optimize!(model).
Example
Here, optimize_hook! is used to bypass EAGO's problem parsing and treat every problem using its branch-and-bound routine. This is done in this example by telling EAGO to treat the problem as a mixed integer nonconvex problem, which normally dispatches to branch-and-bound.
struct MyNewExtension <: EAGO.ExtensionType end
import EAGO: optimize_hook!
function EAGO.optimize_hook!(t::MyNewExtension, m::Optimizer)
initial_parse!(m)
optimize!(EAGO.MINCVX(), m)
endThe same functionality could be accomplished by setting the EAGOParameter field force_global_solve to be true.
EAGO.parse_classify_problem! — Method
parse_classify_problem!(m::GlobalOptimizer)Interprets the type/number of constraints and the type of objective function to infer a problem type. Current possible types include:
LP: Linear program; sent tooptimize_lp.jl- 'MILP
: Mixed integer linear program; sent tooptimize_lp.jl` SOCP: Second-order cone program; sent tooptimize_conic.jlMINCVX: Mixed-integer nonconvex; sent tooptimize_nonconvex.jl
If the force_global_solve parameter is set to true, parse_classify_problem! will set the problem type to MINCVX to pass the problem to optimize_nonconvex.jl.
EAGO.parse_global! — Method
Basic parsing for global solutions (no extensive manipulation). By default, does nothing.
EAGO.postprocess! — Method
postprocess!(t::ExtensionType, m::GlobalOptimizer)
Default postprocess perfoms duality-based bound tightening (Tawarmalani, M., Sahinidis, N.V.: Global optimization of mixed-integer nonlinear programs: a theoretical and computational study. Math. Progr. 99, 563–591 (2004).) up to an iteration limit set by dbbt_depth.
EAGO.preprocess! — Method
preprocess!(
t::ExtensionType,
m::GlobalOptimizer{R, S, Q<:ExtensionType}
)
Runs contractor methods prior to solving lower bounding problem. By default linear and quadratic contractor methods followed by interval constraint propagation then optimization-based bound tightening for a specified number of iterations while the subproblem at current node n has not been proven infeasible.
EAGO.presolve_global! — Method
presolve_global!(t::ExtensionType, m::GlobalOptimizer)
Perform any necessary work prior to running branch-and-bound.
- Set subsolver configs using values in
EAGOParameters; - Load variables, linear constraints, and empty storage space into the relaxed optimizer;
- Prepare the stack for the start of branch-and-bound;
- Fill fields of the
GlobalOptimizerwith zeros of the proper dimensions; - Pass necessary flags from the
GlobalOptimizerto the working problem.
EAGO.print_iteration! — Method
print_iteration!
Print status information based on iteration count. The header print frequency is based on the header_iterations setting, and the data print frequency is based on the output_iterations setting.
EAGO.print_node! — Method
print_node!
Print information about the current node. Includes node ID, lower bound, upper bound, and interval box.
EAGO.print_preamble! — Method
print_preamble!
Print noteworthy information prior to running branch-and-bound. Currently prints a note about flipping max(f) to -min(-f) internally, if the input is a maximization problem and verbosity >= 3.
EAGO.print_results! — Method
print_results!(m::GlobalOptimizer, lower_flag::Bool)
Print the results of a single (lower or upper) bounding problem. lower_flag=true prints information for the lower problem, lower_flag=false prints information for the upper problem.
EAGO.print_solution! — Method
print_solution!
Print solution information for the B&B problem. Display node with the best solution, solution value, solution, and time spent solving subproblems. This print occurs following termination of the B&B algorithm.
EAGO.r_init! — Method
Initializes information in cache c for each node in g that may be used in a reverse-pass of attribute t.
EAGO.reform_epigraph_min! — Method
reform_epigraph_min!(m::GlobalOptimizer)
Perform an epigraph reformulation assuming the working_problem is a minimization problem.
EAGO.rel_diam — Method
rel_diam(m::GlobalOptimizer, i::Int64) -> Float64
Return the relative diameter of a variable. In the case of diam(X)=Inf, rel_diam returns 0.0 to prevent "branching" on this variable endlessly.
EAGO.relax! — Method
relax!(
m::GlobalOptimizer,
f::EAGO.BufferedQuadraticEq,
i::Int64,
check_safe::Bool
) -> Bool
EAGO.relax! — Method
relax!(
m::GlobalOptimizer,
f::EAGO.BufferedQuadraticIneq,
k::Int64,
check_safe::Bool
) -> Bool
EAGO.relax! — Method
relax!(
m::GlobalOptimizer{R, S, Q<:ExtensionType},
f::EAGO.BufferedNonlinearFunction{V, N, T<:RelaxTag},
k::Int64,
check_safe::Bool
) -> Tuple{Bool, Bool}
EAGO.relax! — Function
relax!
Relax the constraint by adding an affine constraint to the model.
EAGO.relax_all_constraints! — Method
relax_all_constraints!(
t::ExtensionType,
m::GlobalOptimizer,
k::Int64
) -> Tuple{Bool, Bool}
A routine that adds relaxations for all nonlinear constraints and quadratic constraints corresponding to the current node to the relaxed problem. This adds an objective cut (if specified by objective_cut_on) and then sets the _new_eval_constraint flag to false indicating that an initial evaluation of the constraints has occurred. If the objective_cut_on flag is true then the _new_eval_objective flag is also set to false indicating that the objective expression was evaluated.
EAGO.relaxed_problem_status — Method
relaxed_problem_status(t, p, d)
Take an MOI.TerminationStatusCode and two MOI.ResultStatusCodes (one each for the primal and dual status) and return a RelaxResultStatus. Returns RRS_OPTIMAL if the codes prove that the subproblem solution was solved to global optimality. Returns RRS_INFEASIBLE if the codes prove that the subproblem solution is infeasible. Returns RRS_DUAL_FEASIBLE if subproblem solution is not optimal and not proven infeasible, but the dual status is MOI.FEASIBLE_POINT. Returns RRS_INVALID otherwise.
EAGO.repeat_check — Method
repeat_check(t::ExtensionType, m::GlobalOptimizer) -> Bool
Check to see if current node should be reprocessed. Without any custom extension, return false by default.
EAGO.rprop! — Method
Populates information associated with attribute t for a constant v at index k in cache c associated with graph g using information at index k taken from the parents of k.
EAGO.rprop! — Method
Populates information associated with attribute t for a expressions v at index k in cache c associated with graph g using information at index k taken from the parents of k.
EAGO.rprop! — Method
Populates information associated with attribute t for a parameters v at index k in cache c associated with graph g using information at index k taken from the parents of k.
EAGO.rprop! — Method
Populates information associated with attribute t for a subexpressions v at index k in cache c associated with graph g using information at index k taken from the parents of k.
EAGO.rprop! — Method
Populates information associated with attribute t for a variable v at index k in cache c associated with graph g using information at index k taken from the parents of k.
EAGO.rprop! — Method
rprop!
Updates storage tapes with reverse evalution of node representing n = *(x,y,z...) which updates x, y, z and so on.
EAGO.rprop! — Method
rprop!
Updates storage tapes with reverse evalution of node representing n = +(x,y,z...) which updates x, y, z and so on.
EAGO.rprop_2! — Method
rprop_2!
Updates storage tapes with reverse evalution of node representing n = x * y which updates x and y.
EAGO.rprop_2! — Method
rprop_2!
Updates storage tapes with reverse evalution of node representing n = x + y which updates x and y.
EAGO.rprop_n! — Method
rprop_n!
Updates storage tapes with reverse evalution of node representing n = *(x,y,z...) which updates x, y, z and so on.
EAGO.rprop_n! — Method
rprop_n!
Updates storage tapes with reverse evalution of node representing n = +(x,y,z...) which updates x, y, z and so on.
EAGO.same_box — Method
same_box(x::NodeBB, y::NodeBB, r::Float64) -> Bool
Check that node x and y have equal domains within an absolute tolerance of r.
EAGO.select_branch_point — Method
select_branch_point(
t::ExtensionType,
m::GlobalOptimizer,
i
) -> Float64
Select a point xb within the domain of the ith branching variable. By default, this point is a convex combination of the solution to the relaxation and the midpoint of the node (branch_cvx_factor*xmid + (1-branch_cvx_factor)*xsol). If the solution lies within branch_offset of a bound, then the branch point is moved to a distance of branch_offset from that bound.
EAGO.select_branch_variable — Method
select_branch_variable(
t::ExtensionType,
m::GlobalOptimizer
) -> Any
Choose a variable to branch on. A maximum relative width branching rule is used by default.
EAGO.set_constraint_propagation_fbbt! — Method
set_constraint_propagation_fbbt!(
m::GlobalOptimizer{R, S, Q<:ExtensionType}
) -> Bool
Performs bound tightening based on forward/reverse interval and/or McCormick passes. This routine resets the current node with new interval bounds.
EAGO.set_default_config! — Method
set_default_config!Configures subsolver tolerances based on tolerance parameters provided to EAGO (provided that a specialized subsolver configuration routine has been provided and m.user_solver_config = false).
EAGO.set_dual! — Method
set_dual!(m)
Retrieves the lower and upper duals for variable bounds from the relaxed_optimizer and sets the appropriate values in the _lower_lvd and _lower_uvd storage fields.
EAGO.set_first_relax_point! — Method
set_first_relax_point!(m::GlobalOptimizer)
EAGO.set_global_lower_bound! — Method
set_global_lower_bound!(m::GlobalOptimizer)
If the previous best-known global lower bound is lower than the lowest lower bound in the stack, set the global lower bound equal to the lowest lower bound in the stack.
EAGO.set_node! — Method
set_node!
Sets the current node in the Evaluator structure.
EAGO.set_reference_point! — Method
set_reference_point!(m::GlobalOptimizer)
EAGO.set_result_status! — Method
set_result_status!(m::GlobalOptimizer)
Convert EAGO's ending status code into an MOI.ResultStatusCode.
EAGO.set_termination_status! — Method
set_termination_status!(m::GlobalOptimizer)
Convert EAGO's ending status code into an MOI.TerminationStatusCode.
EAGO.single_storage! — Method
single_storage!(t::ExtensionType, m::GlobalOptimizer)
Store the current node to the stack, without branching, after updating lower/upper bounds.
EAGO.sip_bnd! — Method
sip_bnd!Solves the bounding problem for the ith-SIP used with extension t::EAGO.ExtensionType in algorithm a::AbstractSIPAlgo in subproblem s::AbstractSubproblemType via the command sip_bnd!(t::ExtensionType, a::AbstractSIPAlgo, s::AbstractSubproblemType, ..., i, tol).
EAGO.sip_llp! — Method
sip_llp!Solves the lower level problem for the ith-SIP used with extension t::EAGO.ExtensionType in algorithm a::AbstractSIPAlgo in subproblem s::AbstractSubproblemType via the command sip_llp!(t::ExtensionType, a::AbstractSIPAlgo, s::AbstractSubproblemType, ..., i, tol).
EAGO.sip_res! — Method
sip_res!Solves the restriction problem for extension t::EAGO.ExtensionType in algorithm a::AbstractSIPAlgo in subproblem s::AbstractSubproblemType via the command sip_res!(t::ExtensionType, a::AbstractSIPAlgo, ...).
EAGO.sip_solve — Method
sip_solveSolve an SIP with decision variable bounds x_l to x_u, uncertain variable bounds p_l to p_u, an objective function of f, and gSIP seminfiniite constraint(s).
EAGO.solve_local_nlp! — Method
Constructs and solves the problem locally on node y updated the upper solution informaton in the optimizer.
EAGO.store_candidate_solution! — Method
store_candidate_solution!(m::GlobalOptimizer)
If the most recent upper problem returned a feasible result, and the upper objective value is less than the previous best-known global upper bound, set the most recent upper problem result to be the new global upper bound. Update the _feasible_solution_found, _solution_node, _global_upper_bound, and _continuous_solution fields of the GlobalOptimizer accordingly.
EAGO.stored_adjusted_upper_bound! — Method
stored_adjusted_upper_bound!(d, v)
Shifts the resulting local nlp objective value f* by (1.0 + relative_tolerance/100.0)*f* + absolute_tolerance/100.0. This assumes that the local solvers relative tolerance and absolute tolerance is significantly lower than the global tolerance (local problem is minimum).
EAGO.termination_check — Method
termination_check(m::GlobalOptimizer)
termination_check(t::ExtensionType, m::GlobalOptimizer) -> BoolCheck for termination of the branch-and-bound algorithm.
If only the GlobalOptimizer is given as an argument, termination_check dispatches to the other form using the ExtensionType given in the SubSolvers. If there is no user-defined extension, then by default, this will check for satisfaction of absolute or relative tolerances, solution infeasibility, and other specified limits. Returns true if any conditions are met and branch-and-bound should end, and false otherwise.
EAGO.trivial_filtering! — Method
trivial_filtering!(
m::GlobalOptimizer{R, S, Q<:ExtensionType},
n::NodeBB
)
Excludes OBBT on variable indices that are tight for the solution of the relaxation.
EAGO.unbounded_check! — Method
unbounded_check!(m::GlobalOptimizer) -> Union{Nothing, Bool}
Check the optimization problem for unbounded branching variables, which would interfere with EAGO's branch-and-bound routine since there are no well-defined branching rules for cases where the interval bounds contain -Inf or Inf. If any branching variables are missing bounds, add the missing bound at +/- 1E6 and warn the user.
EAGO.unpack_fbbt_buffer! — Method
unpack_fbbt_buffer!(m::GlobalOptimizer)
EAGO.unpack_global_solution! — Method
unpack_global_solution!(
m::Optimizer{R, S, Q<:ExtensionType}
)
If global optimization was performed, much of the work happened within the _global_optimizer::GlobalOptimizer. The unpack_global_solution! function extracts results from the GlobalOptimizer and puts them in the correct fields of the Optimizer.
EAGO.unsafe_check_fill! — Method
unsafe_check_fill!(f, y::Array{T, 1}, x, n::Int64)
Performs map!(f, y, x) in an unsafe manner if y[i] is true, else no-op. Assumes n == length(x) == length(y). About 2x faster for small arrays (n < 1000).
[Unused]
EAGO.update_relaxed_problem_box! — Method
update_relaxed_problem_box!(m)
Update the relaxed constraint by setting the constraint set of v == x* , xL_i <= x_i, and x_i <= xU_i for each such constraint added to the relaxed optimizer. Resets integral valued constraints to either EqualTo or Interval constraints.
EAGO.upper_problem! — Method
upper_problem!(t::ExtensionType, m::GlobalOptimizer)
Default upper bounding problem which simply calls solve_local_nlp! to solve the NLP locally.
EAGO.variable_dbbt! — Method
variable_dbbt!(
n::NodeBB,
mult_lo::Vector{Float64},
mult_hi::Vector{Float64},
LBD::Float64,
UBD::Float64,
nx::Int64
)
Tighten the bounds of the _current_node using the current global upper bound and the duality information obtained from the relaxation.
EAGO.variable_load_parse! — Method
variable_load_parse!Parse constraint information of a given constraint type to fill in related variable information.