Public reference
Building a problem
RandomizedProgressiveHedging.AbstractScenario
— TypeAbstractScenario
Abstract type that user scenario concrete types should descend from.
Example
struct MyScenario <: AbstractScenario
value::Vector{Float64}
end
RandomizedProgressiveHedging.build_fs!
— Functionbuild_fs!(model::JuMP.Model, s::S, id_scen::ScenarioId) where S<:AbstractScenario
Define variables, build the objective function, build and attach constraints relative to the scenario s
, of identifier id_scen
, into the given JuMP model
.
Return value: a tuple of
- array of subproblem variables
- expression of objective function
- array of references to constraints that define the objective, as opposed to constraints that define the feasible space.
Example
function build_fs!(model::JuMP.Model, s::MyScenario, id_scen::ScenarioId)
nstages = length(s.value)
Y = @variable(model, [1:nstages], base_name="y_s$id_scen")
m = @variable(model)
# feasibility constraints
@constraint(model, Y[:] .>= 0)
# objective constraints, enforcing m=maximum(Y)
max_ctrs = @constraint(model, [i in 1:stages], m .>= Y[i])
objexpr = sum( (Y[i]-s.value[i])^2 for i in 1:nstages) + m
return Y, objexpr, max_ctrs
end
RandomizedProgressiveHedging.ScenarioTree
— TypeScenarioTree
A tree structure to hold the non-anticipativity structure of a Problem
.
All nodes are stored in the vecnodes
vector, and referenced by their index.
Note: The tree nodes should be indexed in such a way that all sets of equivalent scenarios be made of adjacent indices (n1:n2
).
RandomizedProgressiveHedging.Problem
— TypeProblem{T}
Describe the problem exactly. Attributes are:
scenarios::Vector{T}
: vector of all scenario objects of typeT <: AbstractScenario
build_subpb::Function
: function specializingbuild_fs!
probas::Vector{Float64}
: vector of probability of scenarios, defining objective expectation.nscenarios::Int
: number of scenarios.nstages::Int
: number of stages.stage_to_dim::Vector{UnitRange{Int}}
: map stage to set of corresponding indices of a scenario vector variable.scenariotree::ScenarioTree
: hold the non-anticipatory structure as aScenarioTree
object.
Remarks:
- indexing of scenarios and stages should take values in
1:nscenarios
and1:nstages
respectively. - for any stage, the set of indices of scenarios indistinguishable up to this point should be packed, that is look like
n1:n2
.
Solving a problem
RandomizedProgressiveHedging.solve_direct
— Methodx = solve_direct(pb::Problem; optimizer = GLPK.Optimizer, printlev=1)
Build the progressive hedging problem by explicitly laying out non-anticipatory constraints, and solve globally.
Keyword arguments:
optimizer
: optimizer used for solve. Default isGLPK.Optimizer
.optimizer_params
: aDict{Symbol, Any}
storing parameters for the optimizer.printlev
: if 0, mutes output from the function (not solver). Default value is 1.
RandomizedProgressiveHedging.solve_progressivehedging
— Methodx = solve_progressivehedging(pb::Problem)
Run the classical Progressive Hedging scheme on problem pb
.
Stopping criterion is based on primal dual residual, maximum iterations or time can also be set. Return a feasible point x
.
Keyword arguments:
ϵ_primal
: Tolerance on primal residual.ϵ_dual
: Tolerance on dual residual.μ
: Regularization parameter.maxtime
: Limit on time spent insolve_progressivehedging
.maxcomputingtime
: Limit time spent in computations, excluding computation of initial feasible point and computations required by plottting / logging.maxiter
: Maximum iterations.printlev
: if 0, mutes output.printstep
: number of iterations skipped between each print and logging.hist
: if notnothing
, will record::functionalvalue
: array of functional value indexed by iteration,:time
: array of time at the end of each iteration, indexed by iteration,:dist_opt
: if dict has entry:approxsol
, array of distance between iterate andhist[:approxsol]
, indexed by iteration.:logstep
: number of iteration between each log.
optimizer
: an optimizer for subproblem solve.optimizer_params
: aDict{Symbol, Any}
storing parameters for the optimizer.callback
: either nothing or a functioncallback(pb, x, hist)::nothing
called at each log phase.x
is the current feasible global iterate.
RandomizedProgressiveHedging.solve_randomized_sync
— Methodx = solve_randomized_sync(pb::Problem)
Run the Randomized Progressive Hedging scheme on problem pb
.
Stopping criterion is maximum iterations or time. Return a feasible point x
.
Keyword arguments:
μ
: Regularization parameter.qdistr
: strategy of scenario sampling for subproblem selection,:pdistr
: samples according to objective scenario probability distribution,:unifdistr
: samples according to uniform distribution,- otherwise, an
Array
specifying directly the probablility distribution.
maxtime
: Limit on time spent insolve_progressivehedging
.maxcomputingtime
: Limit time spent in computations, excluding computation of initial feasible point and computations required by plottting / logging.maxiter
: Maximum iterations.printlev
: if 0, mutes output.printstep
: number of iterations skipped between each print and logging.seed
: if notnothing
, specifies the seed used for scenario sampling.hist
: if notnothing
, will record::functionalvalue
: array of functional value indexed by iteration,:time
: array of time at the end of each iteration, indexed by iteration,:dist_opt
: if dict has entry:approxsol
, array of distance between iterate and
hist[:approxsol]
, indexed by iteration.:logstep
: number of iteration between each log.
optimizer
: an optimizer for subproblem solve.optimizer_params
: aDict{Symbol, Any}
storing parameters for the optimizer.callback
: eithernothing
or a functioncallback(pb, x, hist)::nothing
called at each
log phase. x
is the current feasible global iterate.
RandomizedProgressiveHedging.solve_randomized_par
— Functionsolve_randomized_par(pb::Problem{T}) where T<:AbstractScenario
Run the Randomized Progressive Hedging scheme on problem pb
. All workers should be available.
Stopping criterion is maximum iterations or time. Return a feasible point x
.
Keyword arguments:
μ
: Regularization parameter.c
: parameter for step length.qdistr
: strategy of scenario sampling for subproblem selection,:pdistr
: samples according to objective scenario probability distribution,:unifdistr
: samples according to uniform distribution,- otherwise, an
Array
specifying directly the probablility distribution.
maxtime
: Limit on time spent insolve_progressivehedging
.maxcomputingtime
: Limit time spent in computations, excluding computation of initial feasible point and computations required by plottting / logging.maxiter
: Maximum iterations.printlev
: if 0, mutes output.printstep
: number of iterations skipped between each print and logging.seed
: if notnothing
, specifies the seed used for scenario sampling.hist
: if notnothing
, will record::functionalvalue
: array of functional value indexed by iteration,:time
: array of time at the end of each iteration, indexed by iteration,:number_waitingworkers
: array of number of wainting workers, indexed by iteration,:maxdelay
: array of maximal delay among done workers, indexed by iteration,:dist_opt
: if dict has entry:approxsol
, array of distance between iterate andhist[:approxsol]
, indexed by iteration.:logstep
: number of iteration between each log.
optimizer
: an optimizer for subproblem solve.optimizer_params
: aDict{Symbol, Any}
storing parameters for the optimizer.callback
: eithernothing
or a functioncallback(pb, x, hist)::nothing
called at each
log phase. x
is the current feasible global iterate.
RandomizedProgressiveHedging.solve_randomized_async
— Methodsolve_randomized_async(pb::Problem{T}) where T<:AbstractScenario
Run the Randomized Progressive Hedging scheme on problem pb
. All workers should be available.
Stopping criterion is maximum iterations or time. Return a feasible point x
.
Keyword arguments:
μ
: Regularization parameter.c
: parameter for step length.stepsize
: ifnothing
uses theoretical formula for stepsize, otherwise uses constant numerical value.qdistr
: strategy of scenario sampling for subproblem selection,:pdistr
: samples according to objective scenario probability distribution,:unifdistr
: samples according to uniform distribution,- otherwise, an
Array
specifying directly the probablility distribution.
maxtime
: Limit on time spent insolve_progressivehedging
.maxcomputingtime
: Limit time spent in computations, excluding computation of initial feasible point and computations required by plottting / logging.maxiter
: Maximum iterations.printlev
: if 0, mutes output.printstep
: number of iterations skipped between each print and logging.seed
: if notnothing
, specifies the seed used for scenario sampling.hist
: if notnothing
, will record::functionalvalue
: array of functional value indexed by iteration,:time
: array of time at the end of each iteration, indexed by iteration,:number_waitingworkers
: array of number of wainting workers, indexed by iteration,:maxdelay
: array of maximal delay among done workers, indexed by iteration,:dist_opt
: if dict has entry:approxsol
, array of distance between iterate andhist[:approxsol]
, indexed by iteration.:logstep
: number of iteration between each log.
optimizer
: an optimizer for subproblem solve.optimizer_params
: aDict{Symbol, Any}
storing parameters for the optimizer.callback
: eithernothing
or a functioncallback(pb, x, hist)::nothing
called at each
log phase. x
is the current feasible global iterate.
Other
JuMP.objective_value
— Functionobjective_value(pb, x)
Evaluate the objective of problem pb
at point x
.
Note: This function discards all subproblem constraints not explicitly returned by the build_fs!
function.
objective_value(pb, x, id_scen)
Evaluate the objective of the subproblem associated with scenario id_scen
of problem pb
at point x
.
Note: This function discards all subproblem constraints not explicitly returned by the build_fs!
function.
RandomizedProgressiveHedging.cvar_problem
— Functionpb_cvar = cvar_problem(pb::Problem, cvarlevel::Real)
Build the problem with risk measure corresponding to cvar
.