API Reference
AbstractTree Functions
Defining Trees
JuDGE.narytree
— Functionnarytree(depth::Int, degree::Int)
Given the depth
and degree
, this function returns an N-ary tree. Note that a depth of 0 return a single Leaf
node (which is also the root node of the tree).
Required Arguments
depth
is the maximum number of arcs from the root node any Leaf
node
degree
is the number of children of all nodes, other than the Leaf
nodes
Example
tree = narytree(2,2)
JuDGE.tree_from_file
— Functiontree_from_file(filename::String)
Construct tree from a file, each line in the file is of the form B,A,... representing an arc in the tree, from node "A" to node "B". The total number of columns is arbitrary. The first row of the file should be n,p,... these column headers are converted into symbols used to index the data
. Each column itself converted into a dictionary, indexed by the node.
Required Arguments
string
is the full path of the file containing the tree
Example
tree, data = tree_from_file(joinpath(@__DIR__,"tree.csv"))
JuDGE.tree_from_leaves
— Functiontree_from_leaves(leafnodes::Vector{Vector{Int}}, probs::Vector{Float64})
Construct tree from Array of leaf nodes, and (optionally) the corresponding probabilities
Required Arguments
leafnodes
is an array of arrays defining the set of leaf nodes
Optional Arguments
probs
is an array of probabilities for the leaf nodes
Example
(tree,prob) = tree_from_leaves([[1,1,1],[1,1,2],[1,2,1],[1,2,2]],[0.25,0.25,0.25,0.25])
tree = tree_from_leaves([[1,1,1],[1,1,2],[1,2,1],[1,2,2]])
JuDGE.print_tree
— Methodprint_tree(some_tree::AbstractTree, data::Dict{AbstractTree,Any})
Given some_tree
, this function prints a simple representation of the tree to the REPL.
Required Arguments
some_tree
is the tree we wish to visualise
Optional Arguments
data
is a dictionary indexed by each node in some_tree
Nodes of Trees
Base.collect
— Functioncollect(tree::Tree;order::Symbol=:depth,truncate::Int=-1)
Given tree
, this function returns an array of corresponding nodes. By default this will be in a depth-first order.
Required Arguments
tree
is the tree from which we wish to collect the nodes
Optional Arguments
order
can be set to :depth
or :breadth
to specify the order that the nodes are listed in the array.
truncate
limits the nodes returned to a maximum depth of truncate
.
Examples
nodes = collect(tree) #gets an array of nodes from tree in depth-first order
nodes = collect(tree,order=:breadth) #gets an array of nodes from tree in breadth-first order
JuDGE.get_leafnodes
— Functionget_leafnodes(tree::AbstractTree; truncate::Int = -1)
Given tree
, this function returns an array of corresponding Leaf
nodes.
Required Arguments
tree
is the tree from which we wish to collect leaf nodes
Example
leafnodes = JuDGE.get_leafnodes(tree)
JuDGE.get_node
— Functionget_node(tree::AbstractTree, indices::Vector{Int})
Given a tree
, and an array of indices
, this function returns the corresponding node in the tree.
Required Arguments
tree
is the tree from which we are finding the node
indicies
is an array of integer indices identifying a node within tree
.
Examples
node = get_node(tree,[1]) #get the root node
node = get_node(tree,[1,1]) #get the first child of the root node
node = get_node(tree,[1,2]) #get the second child of the root node
Tree Probabilities
JuDGE.convert_probabilities
— Functionconvert_probabilities(tree::AbstractTree, probabilities::Dict{NodeID,Float64})
Given a dictionary of conditional probabilities for each node in tree
, this function returns a dictionary that maps the NodeID
of each node of tree
to the corresponding unconditional probability.
Required Arguments
tree
is the tree that the probabilities pertain to
probabilities
is a dictionary of condition probabilities for each node in tree
Example
probs = JuDGE.convert_probabilities(tree,probabilities)
JuDGE.ConditionallyUniformProbabilities
— FunctionConditionallyUniformProbabilities(tree::AbstractTree)
Given a tree, this function returns a dictionary which maps nodes of the tree to probabilities, given that there are conditionally uniform probabilities over the children of any node.
Required Arguments
tree
is the tree for which the probability distribution will be generated
Example
probs = ConditionallyUniformProbabilities(tree)
JuDGE.UniformLeafProbabilities
— FunctionUniformLeafProbabilities(tree::AbstractTree)
Given a tree, this function returns a dictionary which maps nodes of the tree to probabilities, given that there is a uniform distribution over the leaf nodes
Required Arguments
tree
is the tree for which the probability distribution will be generated
Example
probs = UniformLeafProbabilities(tree)
Other Tree functions
JuDGE.depth
— Functiondepth(tree::AbstractTree)
Given tree
, this function returns the depth. The root node has a depth of 0.
Required Arguments
tree
is the node we wish to find the depth for.
Example
depth = JuDGE.depth(tree) #returns the depth of a node in a tree
JuDGE.history
— Functionhistory(tree::AbstractTree)
Given tree
, this function returns the history back up to the root node of tree
.
Required Arguments
tree
is the node that we wish to find the history for.
Example
history = JuDGE.history(tree) #get a vector of nodes that precede tree
JuDGE.visualize_tree
— Functionvisualize_tree(some_tree::AbstractTree,
data::Dict{AbstractTree,Dict{Symbol,Any}};
scale_edges = nothing,
scale_nodes::Float64 = 0.0,
max_size::Float64 = 50.0,
custom::Union{Nothing,Dict{Symbol,Tuple{String,String,String}}} = nothing,
truncate::Int = -1)
Given some_tree
, this function generates a html/js visualization of the tree.
Required Arguments
some_tree
is the tree we wish to visualise.
data
is a dictionary of the data we wish to display, each element is a dictionary for a node of the tree, indexed by the Symbols.
Optional Arguments
scale_edges
this is the scale factor for edges as the network gets deeper.
scale_nodes
this is the scale factor for nodes as the network gets deeper.
max_size
this is the size of the root node.
JuDGE.get_groups
— Functionget_groups(tree::AbstractTree; combine=0)
Given a tree
, this function will split it up into an array of subtrees. These can be provided as blocks
for the JuDGE.solve()
function to perform partial pricing.
Required Arguments
tree
is the tree from which we are finding the node
Optional Arguments
combine
this parameter determines the size of the subtrees (higher creates larger subtrees). If set to 0, the subtrees will be the sets of paths to the leaf nodes.
Examples
blocks = get_groups(tree)
JuDGE Functions
JuDGE solving functions
JuDGE.JuDGEModel
— TypeJuDGEModel(tree::AbstractTree,
probabilities,
sub_problem_builder::Function,
solver;
discount_factor=1.0,
risk=RiskNeutral(),
sideconstraints=nothing,
check=true,
perfect_foresight=false)
Define a JuDGE model.
Required arguments
tree
is a reference to a scenario tree
probabilities
is either a function, which returns a dictionary of the probabilities of all nodes in a tree, or simply the dictionary itself
sub_problem_builder
is a function mapping a node to a JuMP model for each subproblems
solver
is a reference to the optimizer used for the master problem (with appropriate settings); this can also be a tuple containing two optimizers (one for solving the relaxation, and one for solving the binary model)
Optional arguments
discount_factor
is a number between 0 and 1 defining a constant discount factor along each arc in the scenario tree
risk
can be either a Risk
object, or a vector of such objects.
sideconstraints
is a function which specifies side constraints in the master problem, see Tutorial 9: Side-constraints for further details
check
is a boolean, which can be set to false
to disable the validation of the JuDGE model.
perfect_foresight
is a boolean; this is an experimental feature, which creates an array of JuDGE models, one for each leaf node. This will enable users to easily compute the EVPI for the stochastic program. Also can be used for regret-based risk implementations.
Examples
judge = JuDGEModel(tree, ConditionallyUniformProbabilities, sub_problems,
Gurobi.Optimizer)
judge = JuDGEModel(tree, probabilities, sub_problems, CPLEX.Optimizer,
discount_factor=0.9, risk=Risk(0.5,0.1)))
JuDGE.solve
— Methodsolve(judge::JuDGEModel;
termination::Termination=Termination(),
max_no_int::Int=typemax(Int),
blocks::Union{Nothing,Vector{Vector{AbstractTree}}}=nothing,
warm_starts::Bool=false,
optimizer_attributes::Union{Nothing,Function}=nothing,
mp_callback::Union{Nothing,Function}=nothing,
prune::Float64=Inf,
heuristic::Union{Nothing,Function}=nothing,
verbose::Int=2)
Solve a JuDGEModel judge
without branch-and-price.
Required Arguments
judge
is the JuDGE model that we wish to solve.
Optional Arguments
termination
is a Termination
object containing all the stopping conditions.
max_no_int
is the maximum number of iterations yielding a fractional solution before a MIP solve is performed on the master. By default, the MIP solves will not occur until the relaxed bound gap is less than the relgap
/ absgap
stopping conditions. To override this, set max_no_int
to the negative of the number desired value.
blocks
specifies the groups of nodes to solve in each iteration (these groups can be generated using JuDGE.get_groups()
, or created manually), after all nodes have been solved, a full pricing iteration is used to compute an updated lower bound. See advanced.jl
for more details.
warm_starts
boolean specifing whether to use warm starts for subproblems and binary solves of master problem.
optimizer_attributes
can be set to a specific function that dynamically changes optimizer attributes for the subproblems; this should only be used by people who have examined the advanced.jl
example.
mp_callback
is a user-defined function that specifies termination conditions for MIP solves of the master problem. See examples/advanced.jl.
prune
is used to stop the algorithm before convergence, if a known upper bound for the problem is specified.
heuristic
is a user-defined function that typically would perform an improvement heuristic when a new incumbent is found.
verbose
if 0, all output from solve will be suppressed, if 1, the subproblem solve process will be suppressed. Default is 2.
Examples
JuDGE.solve(jmodel, termination=Termination(rlx_abstol=10^-6))
JuDGE.solve(jmodel, termination=Termination(rlx_abstol=10^-6), max_no_int=-5)
JuDGE.branch_and_price
— Functionbranch_and_price(models::Union{JuDGEModel,Vector{JuDGEModel}};
branch_method::Function=JuDGE.variable_branch,search::Symbol=:lowestLB,
termination::Termination=Termination(),
max_no_int::Int=typemax(Int),
blocks::Union{Nothing,Vector{Vector{AbstractTree}}}=nothing,
warm_starts::Bool=false,
optimizer_attributes::Union{Nothing,Function}=nothing,
mp_callback::Union{Nothing,Function}=nothing,
bp_callback::Union{Nothing,Function}=nothing,
heuristic::Union{Nothing,Function}=nothing,
verbose::Int=2)
Solve a JuDGEModel judge
without branch and price.
Required Arguments
judge
is the JuDGE model that we wish to solve.
Optional Arguments
branch_method
is a function specifies the way that constraints are added to create new nodes in the branchandprice tree.
search
specifies the order in which nodes are solved in the (branch-and-price) tree. Options are: :lowestLB
, :depth_first_dive
, :depth_first_resurface
, :breadth_first
.
termination
is a Termination
object containing all the stopping conditions.
max_no_int
is the maximum number of iterations yielding a fractional solution before a MIP solve is performed on the master. By default, the MIP solves will not occur until the relaxed bound gap is less than the relgap
/ absgap
stopping conditions. To override this, set max_no_int
to the negative of the number desired value.
blocks
specifies the groups of nodes to solve in each iteration (these groups can be generated using JuDGE.get_groups()
, or created manually), after all nodes have been solved, a full pricing iteration is used to compute an updated lower bound. See advanced.jl
for more details.
warm_starts
boolean specifing whether to use warm starts for subproblems and binary solves of master problem.
optimizer_attributes
can be set to a specific function that dynamically changes optimizer attributes for the subproblems; this should only be used by people who have examined the advanced.jl
example.
mp_callback
is a user-defined function that specifies termination conditions for MIP solves of the master problem. See examples/advanced.jl.
bp_callback
is a user-defined function that allows you to modify the Termination
conditions for JuDGE.solve
and the search
policy during the branch-and-price process.
heuristic
is a user-defined function that typically would perform an improvement heuristic when a new incumbent is found.
verbose
if 0, most output from JuDGE.solve
will be suppressed, if 1, the subproblem solve process will be suppressed. Default is 2.
Examples
JuDGE.branch_and_price(jmodel, termination=Termination(abstol=10^-6))
JuDGE.branch_and_price(jmodel, search=:depth_first_dive, verbose=0)
JuDGE.Termination
— TypeTermination(;abstol::Float64=10^-10, reltol::Float64=10^-10, rlxabstol::Float64=10^-10, rlxreltol::Float64=10^-10, timelimit::Float64=Inf, maxiter::Int=typemax(Int), inttol::Float64=10^-9, allowfrac::Symbol=:binarysolve)
Define the stopping conditions for JuDGE.solve()
/ JuDGE.branch_and_price()
.
Optional Arguments
abstol
is the absolute tolerance for the best integer-feasible objective value and the lower bound.
reltol
is the relative tolerance for the best integer-feasible objective value and the lower bound.
rlx_abstol
is the absolute tolerance for the relaxed master objective value and the lower bound.
rlx_reltol
is the relative tolerance for the relaxed master objective value and the lower bound.
time_limit
is the maximum duration in seconds.
max_iter
is the maximum number of iterations.
inttol
is the maximum deviation from 0 or 1 for any binary/integer variable for integer feasible solutions.
allow_frac
indicates whether a fractional solution will be returned; possible values are: :binary_solve
a binary solve of master will be performed (if needed) prior to the solution being returned; :binary_solve_return_relaxation
a binary solve of master will be performed (if needed), updating the upper bound, but the master problem relation will be returned; :first_fractional
will return the first fractional master solution found; :no_binary_solve
will simply return the solution to the relaxed master problem when terminated.
Examples
Termination(rlx_abstol=10^-6) Termination(abstol=10^-3,inttol=10^-6)
JuDGE.variable_branch
— Functionvariable_branch(jmodel::JuDGEModel, inttol::Float64)
This is an in-built function that is called during branch-and-price to perform a branch. Users can define their own functions that follow this format to create new branching strategies.
Required Arguments
jmodel
is the JuDGE model
inttol
is the maximum permitted deviation from binary/integer for a value to still be considered binary/integer feasible.
JuDGE.resolve_subproblems
— Functionresolve_subproblems(judge::JuDGEModel)
Once a JuDGE model has converged, it is necessary to re-solve the subproblems to find the optimal decisions within each node.
Required Arguments
jmodel
is the JuDGE model that we wish to solve.
Examples
resolve_subproblems(judge)
JuDGE macros for subproblems
JuDGE.@expansion
— Macroexpansion(model, variable, args...)
Defines an expansion variable variable
within a subproblem model
. Note that all subproblems must have the same set of expansion variables.
Required Arguments
model
is the JuDGE subproblem that we are adding the expansion variable to
variable
is the name of the variable being created, this will be continuous by default; follows JuMP syntax if defining a set of variables.
Optional Arguments
This macro has a third, unnamed, argument which can be set to Con, Bin, or Int, similar to the @variable
macro.
lag
is the number of nodes in the scenario between an expansion being decided, and it becoming available.
duration
is the number of consecutive nodes in the scenario over which an expansion is available.
lb
is the lower bound for this variable in the master problem (typically omitted).
ub
is the upper bound for this variable in the master problem (typically omitted).
Examples
@expansion(model, expand[1:5], Bin) #defines an array of 5 binary variables with no lag, and unlimited lifespan
@expansion(model, expand[1:5,1:2]>=0, lag=1) #defines a matrix of 10 continuous variables with a lag of 1, and unlimited duration
@expansion(model, 0<=expand<=10, Int, duration=2) #defines a single integer variable with a lag of 0, and a duration of 2
JuDGE.@shutdown
— Macroshutdown(model, variable, args...)
Defines an shutdown variable variable
within a subproblem model
. Note that all subproblems must have the same set of shutdown variables.
Required Arguments
model
is the JuDGE subproblem that we are adding the shutdown variable to
variable
is the name of the variable being created, this will be continuous by default; follows JuMP syntax if defining a set of variables.
Optional Arguments
This macro has a third, unnamed, argument which can be set to Con, Bin, or Int, similar to the @variable
macro.
lag
is the number of nodes in the scenario between an shutdown being decided, and it becoming unavailable.
duration
is the number of consecutive nodes in the scenario over which the shutdown will last.
lb
is the lower bound for this variable in the master problem (typically omitted).
ub
is the upper bound for this variable in the master problem (typically omitted).
Examples
@shutdown(model, shut[1:5], Bin) #defines an array of 5 binary variables with no lag, and unlimited lifespan
@shutdown(model, shut[1:5,1:2]>=0, lag=1) #defines a matrix of 10 continuous variables with a lag of 1, and unlimited duration
@shutdown(model, 0<=shut<=10, Int, duration=2) #defines a single integer variable with a lag of 0, and a duration of 2
JuDGE.@enforced
— Macroenforced(model, variable, args...)
Defines an enforced variable variable
within a subproblem model
. Note that all subproblems must have the same set of enforced variables. These variables can be used as either expansion or shutdown variables, but since the constraint in the master problem is an equality, convergence is more difficult since there is less flexibility when solving the master problem.
Required Arguments
model
is the JuDGE subproblem that we are adding the expansion variable to
variable
is the name of the variable being created, this will be continuous by default; follows JuMP syntax if defining a set of variables.
Optional Arguments
This macro has a third, unnamed, argument which can be set to Con, Bin, or Int, similar to the @variable
macro.
lag
is the number of nodes in the scenario between an expansion being decided, and it becoming available.
duration
is the number of consecutive nodes in the scenario over which an expansion is available.
lb
is the lower bound for this variable in the master problem (typically omitted).
ub
is the upper bound for this variable in the master problem (typically omitted).
penalty
is a placeholder for a future feature, which may allow the violation of master/subproblem equality constraint, at a cost.
Examples
@expansion(model, forced[1:5], Bin) #defines an array of 5 binary variables with no lag, and unlimited lifespan
@expansion(model, forced[1:5,1:2]>=0, lag=1) #defines a matrix of 10 continuous variables with a lag of 1, and unlimited duration
@expansion(model, 0<=forced<=10, Int, duration=2) #defines a single integer variable with a lag of 0, and a duration of 2
JuDGE.@state
— Macrostate(model, variable, args...)
Defines a state variable variable
within a subproblem model
. Note that all subproblems must have the same set of state variables. These variables can be used to model inventory that is carried forward between the subproblems.
Required Arguments
model
is the JuDGE subproblem that we are adding the expansion variable to
variable
is the name of the variable being created in the subproblem, this will be continuous by default; follows JuMP syntax if defining a set of variables. The subproblem variable corresponds to the change in the state.
Optional Arguments
This macro has a third, unnamed, argument which can be set to Con, Bin, or Int, similar to the @variable
macro.
state_name
is the name for the state variable in the master problem. If omitted, the name of the master problem variable will match the subproblem variable (but this may cause confusion, since only the master problem variable is the state). See the inventory.jl
example to see how this should be implemented.
initial
is the initial value for the master problem's state variable at the root node.
lb
is the lower bound for the variable in the master problem (typically omitted).
ub
is the upper bound for the variable in the master problem (typically omitted).
Examples
@state(sp, -50<=Δstock<=50, state_name=stock, lb=0, ub=200, initial=0) #defines a state variable called stock in the master
#(starting at 0, and able to take values 0 to 200),
#and Δstock in the subproblem (able to change the stock level by ±50).
JuDGE.@capitalcosts
— Macrocapitalcosts(model, expr)
Defines a linear expression specifying the capital cost of expansions and shutdowns at the current node
Required Arguments
model
is the JuDGE subproblem corresponding to the node in the scenario tree that we are adding specifying the costs for
expr
is an AffExpr
which gives the total cost of choosing expansion and shutdown variables at the current node
Example
@capitalcosts(model, sum(expand[i]*cost[node][i] for i in 1:5))
JuDGE.@ongoingcosts
— Macroongoingcosts(model, expr)
Defines a linear expression specifying the ongoing costs of expansions and shutdowns available at the current node
Required Arguments
model
is the JuDGE subproblem corresponding to the node in the scenario tree that we are adding specifying the costs for
expr
is an AffExpr
which gives the ongoing cost of expansions and shutdowns available at the current node
Example
@ongoingcosts(model, sum(expand[i]*ongoingcosts[node][i] for i in 1:5))
JuDGE Output
JuDGE.write_solution_to_file
— Methodwrite_solution_to_file(jmodel::JuDGEModel,filename::String)
Given a JuDGE model and a filename, this function writes the entire solution to a CSV.
Required Arguments
jmodel
is the JuDGE model whose solution we wish to write to a file
filename
is the output filename
JuDGE.print_expansions
— Methodprint_expansions(jmodel::JuDGEModel;
onlynonzero::Bool=true,
inttol=10^-9,
format=nothing)
Given a solved JuDGE model, this function will write the optimal capacity expansion decisions to the REPL.
Required Arguments
jmodel
is the JuDGE model whose solution we wish to write to a file
Optional Arguments
onlynonzero
is a boolean, if set to true
the function will only print expansions with a non-zero value.
inttol
is the integrality tolerance; any expansion variable value less than this will be treated as 0, and any value greater than 1-inttol
will be treated as 1
format
is a function that specifies customised printing of expansion values. See Tutorial 2: Formatting output for more details.
Deterministic Equivalent
Define and solve DetEq model
JuDGE.DetEqModel
— TypeDetEqModel(tree::AbstractTree,
probabilities,
sub_problem_builder::Function,
solver
discount_factor=1.0,
risk=RiskNeutral,
sideconstraints=nothing,
parallel=false,
check=true)
Define a deterministic equivalent model for the stochastic capacity expansion problem.
Required arguments
tree
is a reference to a scenario tree
probabilities
is either a function, which returns a dictionary of the probabilities of all nodes in a tree, or simply the dictionary itself
sub_problem_builder
is a function mapping a node to a JuMP model for each subproblems
solver
is a reference to the optimizer used for this problem (with appropriate settings)
Optional arguments
discount_factor
is a number between 0 and 1 defining a constant discount factor along each arc in the scenario tree
risk
is a tuple with the two CVaR parameters: (λ, α)
sideconstraints
is a function which specifies side constraints in the master problem, see Tutorial 9: Side-constraints for further details.
parallel
is a boolean, setting whether the sub-problems will be formulated in parallel
check
is a boolean, which can be set to false
to disable the validation of the JuDGE model.
Examples
deteq = DetEqModel(tree, ConditionallyUniformProbabilities, sub_problems,
Gurobi.Optimizer)
judge = DetEqModel(tree, probabilities, sub_problems, CPLEX.Optimizer,
discount_factor=0.9, risk=(0.5,0.1)))
JuDGE.solve
— Methodsolve(deteq::DetEqModel)
Solve a determinisitc equivalent model.
Required Arguments
deteq
is the determinisitc equivalent model that we wish to solve.
Example
JuDGE.solve(deteq)
Deterministic Equivalent Output
JuDGE.write_solution_to_file
— Methodwrite_solution_to_file(deteq::DetEqModel,filename::String)
Given a deterministic equivalent model and a filename, this function writes the entire solution to a CSV.
Required Arguments
deteq
is the deterministic equivalent model whose solution we wish to write to a file
filename
is the output filename
JuDGE.print_expansions
— Methodprint_expansions(deteq::DetEqModel;
onlynonzero::Bool=true,
inttol=10^-9,
format=nothing)
Given a solved deterministic equivalent model, this function will write the optimal capacity expansion decisions to the REPL.
Required Arguments
deteq
is the deterministic equivalent model whose solution we wish to write to a file
Optional Arguments
onlynonzero
is a boolean, if set to true
the function will only print expansions with a non-zero value.
inttol
is the integrality tolerance; any expansion variable value less than this will be treated as 0, and any value greater than 1-inttol
will be treated as 1
format
is a function that specifies customised printing of expansion values. See Tutorial 2: Formatting output for more details.
Risk
JuDGE.RiskNeutral
— MethodRiskNeutral()
Create a risk-neutral risk measure.
JuDGE.Risk
— MethodRisk(λ::Float64,
α::Float64;
offset::Union{Dict{Leaf,Float64},Nothing}=nothing,
bound::Union{Float64,Nothing}=nothing,
penalty::Float64=100000.0)
Define the CVaR risk measure to be applied to the accumulated profits at the leaf nodes.
Required Arguments
λ
is weighting applied for the risk measure (max sum of weightings should be 1.0), if sum of weightings is less than 1.0, expected value will make up the rest.
α
is the probability in the tail of the distribution
Optional Arguments
offset
applies a negative offset to each leaf node. This can be used to reorder the outcomes prior to applying the risk measure.
bound
if used, this will create a constraint on CVaR(α) with this as the upper bound.
penalty
if a constraint on CVaR is applied, then the marginal cost of violating the constraint is penalty
.
JuDGE.Risk
— MethodRisk(α::Float64;
offset::Union{Dict{Leaf,Float64},Nothing}=nothing,
bound::Union{Float64,Nothing}=nothing,
penalty::Float64=100000.0)
Define the CVaR risk constraint to be applied to the accumulated profits at the leaf nodes.
Required Arguments
α
is the probability in the tail of the distribution
Optional Arguments
offset
applies a negative offset to each leaf node. This can be used to reorder the outcomes prior to applying the risk measure.
bound
if used, this will create a constraint on CVaR(α) with this as the upper bound.
penalty
if a constraint on CVaR is applied, then the marginal cost of violating the constraint is penalty
.