Measure Operators
A technical manual for measures in InfiniteOpt. See the respective guide for more information.
Measure Data
InfiniteOpt.AbstractMeasureData — TypeAbstractMeasureDataAn abstract type to define data for measures to define the behavior of Measure.
InfiniteOpt.DiscreteMeasureData — MethodDiscreteMeasureData(pref::GeneralVariableRef,
coefficients::Vector{<:Real},
supports::Vector{<:Real};
[label::Type{<:AbstractSupportLabel} = generate_unique_label(),
weight_function::Function = [`default_weight`](@ref),
lower_bound::Real = NaN,
upper_bound::Real = NaN,
is_expect::Bool = false]
)::DiscreteMeasureDataReturns a 1-dimensional DiscreteMeasureData object that can be utilized to define measures using measure. This accepts input for a scalar (single) infinite parameter. A description of the other arguments is provided in the documentation for DiscreteMeasureData. Errors if supports are out bounds or an unequal number of supports and coefficients are given. Note that by default a unique label is generated via generate_unique_label to ensure the supports can be located in the infinite parameter support storage. Advanced implementations, may choose a different behavior but should do so with caution.
Example
julia> data = DiscreteMeasureData(pref, [0.5, 0.5], [1, 2]);InfiniteOpt.DiscreteMeasureData — MethodDiscreteMeasureData(
prefs::Array{GeneralVariableRef},
coefficients::Vector{<:Real},
supports::Vector{<:Array{<:Real}};
label::Type{<:AbstractSupportLabel} = generate_unique_label(),
weight_function::Function = [`default_weight`](@ref),
lower_bounds::Array{<:Real} = [NaN...],
upper_bounds::Array{<:Real} = [NaN...],
is_expect::Bool = false
)::DiscreteMeasureDataReturns a DiscreteMeasureData object that can be utilized to define measures using measure. This accepts input for an array (multi) parameter. The inner arrays in the supports vector need to match the formatting of the array used for parameter_refs. A description of the other arguments is provided in the documentation for DiscreteMeasureData. Errors if supports are out bounds, an unequal number of supports and coefficients are given, the array formats do not match, or if mixed infinite parameter types are given. Note that by default a unique label is generated via generate_unique_label to ensure the supports can be located in the infinite parameter support storage. Advanced implementations, may choose a different behavior but should do so with caution.
Example
julia> data = DiscreteMeasureData(prefs, [0.5, 0.5], [[1, 1], [2, 2]]);InfiniteOpt.DiscreteMeasureData — TypeDiscreteMeasureData{P <: Union{JuMP.AbstractVariableRef,
Vector{<:JuMP.AbstractVariableRef}},
N, B <: Union{Float64, Vector{Float64}},
F <: Function
} <: AbstractMeasureDataA DataType for immutable measure abstraction data where the abstraction is of the form: $measure = \int_{\tau \in T} f(\tau) w(\tau) d\tau \approx \sum_{i = 1}^N \alpha_i f(\tau_i) w(\tau_i)$. The supports and coefficients are immutable (i.e., they will not change even if supports are changed for the underlying infinite parameter.) This type can be used for both 1-dimensional and multi-dimensional measures.
Fields
parameter_refs::P: The infinite parameter(s) over which the integration occurs. These can be comprised of multiple independent parameters, but dependent parameters cannot be mixed with other types.coefficients::Vector{Float64}: Coefficients $\alpha_i$ for the above measure abstraction.supports::Array{Float64, N}: Supports points $\tau_i$. This is aVectorif only one parameter is given, otherwise it is aMatrixwhere the supports are stored column-wise.label::DataType: Label for the support points $\tau_i$ when stored in the infinite parameter(s), stemming fromAbstractSupportLabel.weight_function::F: Weighting function $w$ must map an individual support value to aRealscalar value.lower_bounds::B: Lower bound in accordance with $T$, this denotes the intended interval of the measure and should beNaNif ignoredupper_bounds::B: Same as above but the upper bound.is_expect::Bool: Is this data associated with an expectation call?
InfiniteOpt.FunctionalDiscreteMeasureData — MethodFunctionalDiscreteMeasureData(pref::GeneralVariableRef,
coeff_func::Function,
min_num_supports::Int,
label::Type{<:AbstractSupportLabel};
[weight_function::Function = [`default_weight`](@ref),
lower_bound::Real = NaN,
upper_bound::Real = NaN,
is_expect::Bool = false,
generative_support_info::AbstractGenerativeInfo = NoGenerativeSupports()]
)::FunctionalDiscreteMeasureDataReturns a 1-dimensional FunctionalDiscreteMeasureData object that can be utilized to define measures using measure. This accepts input for a scalar (single) infinite parameter. A description of the other arguments is provided in the documentation for FunctionalDiscreteMeasureData. Errors if pref is not an infinite parameter. Built-in choices for label include:
All: Use all of the supports stored inprefMCSample: Use Monte Carlo samples associated withprefWeightedSample: Use weighted Monte Carlo samples associated withprefUniformGrid: Use uniform grid points associated withpref.
Example
julia> data = FunctionalDiscreteMeasureData(pref, my_func, 20, UniformGrid);InfiniteOpt.FunctionalDiscreteMeasureData — MethodFunctionalDiscreteMeasureData(prefs::Array{GeneralVariableRef},
coeff_func::Function,
min_num_supports::Int,
label::Type{<:AbstractSupportLabel};
[weight_function::Function = [`default_weight`](@ref),
lower_bounds::Array{<:Real} = [NaN...],
upper_bounds::Array{<:Real} = [NaN...],
is_expect::Bool = false]
)::FunctionalDiscreteMeasureDataReturns a multi-dimensional FunctionalDiscreteMeasureData object that can be utilized to define measures using measure. This accepts input for an array of infinite parameters. A description of the other arguments is provided in the documentation for FunctionalDiscreteMeasureData. Errors if prefs are not infinite parameters or if the mixed parameter types are provided. Built-in choices for label include:
All: Use all of the supports stored inprefsMCSample: Use Monte Carlo samples associated withprefsWeightedSample: Use weighted Monte Carlo samples associated withprefsUniformGrid: Use uniform grid points associated withprefs.
Example
julia> data = FunctionalDiscreteMeasureData(prefs, my_func, 20, MCSample);InfiniteOpt.FunctionalDiscreteMeasureData — TypeFunctionalDiscreteMeasureData{P <: Union{JuMP.AbstractVariableRef,
Vector{<:JuMP.AbstractVariableRef}},
B <: Union{Float64, Vector{Float64}},
I <: AbstractGenerativeInfo,
F1 <: Function,
F2 <: Function
} <: AbstractMeasureDataA DataType for mutable measure abstraction data where the abstraction is of the form: $measure = \int_{\tau \in T} f(\tau) w(\tau) d\tau \approx \sum_{i = 1}^N \alpha_i f(\tau_i) w(\tau_i)$. This abstraction is equivalent to that of DiscreteMeasureData, but the difference is that the supports are not fully known at the time of measure creation. Thus, functions are stored that will be used to generate the concrete support points $\tau_i$ and their coefficients $\alpha_i$ when the measure is evaluated (expanded). These supports are identified/generated in accordance with the label with a gaurantee that at least num_supports are generated. For example, if label = MCSample and num_supports = 100 then the measure will use all of the supports stored in the parameter_refs with the label MCSample and will ensure there are at least 100 are generated. This type can be used for both 1-dimensional and multi-dimensional measures.
For 1-dimensional measures over independent infinite parameters, the generative_supp_info specifies the info needed to make generative supports based on those with that exist with label. Note that only 1 kind of generative supports are allowed for each infinite parameter.
Fields
parameter_refs::P: The infinite parameter(s) over which the integration occurs. These can be comprised of multiple independent parameters, but dependent parameters cannot be mixed with other types.coeff_function::F1: Coefficient generation function making $\alpha_i$ for the above measure abstraction. It should take all the supports as input (formatted as an Array) and return the corresponding vector of coefficients.min_num_supports::Int: Specifies the minimum number of supports $\tau_i$ desired in association withparameter_refsandlabel.label::DataType: Label for the support points $\tau_i$ which are/will be stored in the infinite parameter(s), stemming fromAbstractSupportLabel.generative_supp_info::I: Information needed to generate supports based on other existing ones.weight_function::F2: Weighting function $w$ must map an individual support value to aRealscalar value.lower_bounds::B: Lower bounds in accordance with $T$, this denotes the intended interval of the measure and should beNaNif ignoredupper_bounds::B: Same as above but the upper bounds.is_expect::Bool: Is this data associated with an expectation call?
InfiniteOpt.parameter_refs — Methodparameter_refs(
data::AbstractMeasureData
)::Union{GeneralVariableRef, Array{GeneralVariableRef}}Return the infinite parameter reference(s) in data. This is intended as an internal function to be used with measure addition. User-defined measure data types will need to extend this function otherwise an error is thrown.
InfiniteOpt.support_label — Methodsupport_label(data::AbstractMeasureData)::Type{<:AbstractSupportLabel}Return the label stored in data associated with its supports. This is intended as en internal method for measure creation and ensures any new supports are added to parameters with such a label. User-defined measure data types should extend this functionif supports are used, otherwise an error is thrown.
InfiniteOpt.generative_support_info — Methodgenerative_support_info(data::AbstractMeasureData)::AbstractGenerativeInfoReturn the generative support creation info that corresponds to data. This is intended as an internal method and only needs to be extended for user-defined measure data types that use generative supports.
JuMP.lower_bound — MethodJuMP.lower_bound(data::AbstractMeasureData)::Union{Float64, Vector{Float64}}Return the lower bound associated with data that defines its domain. This is intended as an internal method, but may be useful for extensions. User-defined measure data types should extend this function if desired, otherwise NaN is returned
JuMP.upper_bound — MethodJuMP.upper_bound(data::AbstractMeasureData)::Union{Float64, Vector{Float64}}Return the lower bound associated with data that defines its domain. This is intended as an internal method, but may be useful for extensions. User-defined measure data types should extend this function if desired, otherwise NaN is returned.
InfiniteOpt.supports — Methodsupports(data::AbstractMeasureData)::Array{Float64}Return the supports associated with data and its infinite parameters. This is intended as en internal method for measure creation and ensures any new supports are added to parameters. User-defined measure data types should extend this function if appropriate, otherwise an empty vector is returned.
InfiniteOpt.num_supports — Methodnum_supports(data::AbstractMeasureData)::IntReturn the number supports associated with data and its infinite parameters. This is intended as an internal method for measure creation. User-defined measure data types should extend this function if appropriate, otherwise 0 is returned.
InfiniteOpt.min_num_supports — Methodmin_num_supports(data::AbstractMeasureData)::IntReturn the minimum number of supports associated with data. By fallback, this will just return num_supports(data). This is primarily intended for internal queries of FunctionalDiscreteMeasureData, but can be extended for other measure data types if needed.
InfiniteOpt.coefficient_function — Methodcoefficient_function(data::AbstractMeasureData)::FunctionReturn the coefficient function stored in data associated with its expansion abstraction is there is such a function. This is intended as an internal method for measure creation. User-defined measure data types should extend this function if appropriate, otherwise an error is thrown for unsupported types.
InfiniteOpt.coefficients — Methodcoefficients(data::AbstractMeasureData)::Vector{<:Real}Return the coefficients associated with data associated with its expansion abstraction. This is intended as en internal method for measure creation. User-defined measure data types should extend this function if appropriate, otherwise an empty vector is returned.
InfiniteOpt.weight_function — Methodweight_function(data::AbstractMeasureData)::FunctionReturn the weight function stored in data associated with its expansion abstraction. This is intended as en internal method for measure creation. User-defined measure data types should extend this function if appropriate, otherwise an error is thrown.
InfiniteOpt.default_weight — Functiondefault_weight(t) = 1Default weight function for DiscreteMeasureData and FunctionalDiscreteMeasureData. Returns 1 regardless of the input value.
Definition
General
InfiniteOpt.measure — Functionmeasure(expr::JuMP.AbstractJuMPScalar,
data::AbstractMeasureData;
[name::String = "measure"])::GeneralVariableRefReturn a measure reference that evaluates expr using according to data. The measure data data determines how the measure is to be evaluated. Typically, the DiscreteMeasureData and the FunctionalDiscreteMeasureData constructors can be used to for data. The variable expression expr can contain InfiniteOpt variables, infinite parameters, other measure references (meaning measures can be nested), and constants. Typically, this is called inside of JuMP.@expression, JuMP.@objective, and JuMP.@constraint in a manner similar to sum. Note measures are not explicitly evaluated until build_transformation_backend! is called or unless they are expanded via expand or expand_all_measures!.
Example
julia> tdata = DiscreteMeasureData(t, [0.5, 0.5], [1, 2]);
julia> xdata = DiscreteMeasureData(xs, [0.5, 0.5], [[-1, -1], [1, 1]]);
julia> constr_RHS = @expression(model, measure(g - s + 2, tdata) + s^2)
measure{t}[g(t) - s + 2] + s²
julia> @objective(model, Min, measure(g - 1 + measure(T, xdata), tdata))
measure{xs}[g(t) - 1 + measure{xs}[T(t, x)]]InfiniteOpt.@measure — Macro@measure(expr::JuMP.AbstractJuMPScalar,
data::AbstractMeasureData;
[name::String = "measure"])::GeneralVariableRefAn efficient wrapper for measure, please see its doc string for more information.
InfiniteOpt.build_measure — Functionbuild_measure(expr::JuMP.AbstractJuMPScalar,
data::AbstractMeasureData)::MeasureBuild and return a Measure given the expression to be measured expr using measure data data. This principally serves as an internal method for measure definition. Errors if the supports associated with data violate an finite variable parameter bounds of finite variables that are included in the measure.
InfiniteOpt.Measure — TypeMeasure{T <: JuMP.AbstractJuMPScalar, V <: AbstractMeasureData}A DataType for measure abstractions. The abstraction is determined by data and is enacted on func when the measure is evaluated (expended).
Fields
func::TTheInfiniteOptexpression to be measured.data::VData of the abstraction as described in aAbstractMeasureDataconcrete subtype.group_int_idxs::Vector{Int}: The parameter group integer indices of the evaluated measure expression (i.e., the group integer indices offuncexcluding those that belong todata).parameter_nums::Vector{Int}: The parameter numbers that parameterize the evaluated measure expression. (i.e., the parameter numbers offuncexcluding those that belong todata).constant_func::Bool: Indicates iffuncis not parameterized by the infinite parameters indata. (i.e., do the group integer indices offuncanddatahave no intersection?) This is useful to enable analytic evaluations if possible.
InfiniteOpt.add_measure — Functionadd_measure(model::InfiniteModel, meas::Measure,
name::String = "measure")::GeneralVariableRefAdd a measure to model and return the corresponding measure reference. This operates in a manner similar to JuMP.add_variable. Note this intended as an internal method.
InfiniteOpt.add_supports_to_parameters — Methodadd_supports_to_parameters(data::AbstractMeasureData)::NothingAdd supports as appropriate with data to the underlying infinite parameters. This is an internal method with by add_measure and should be defined for user-defined measure data types.
InfiniteOpt.MeasureIndex — TypeMeasureIndex <: ObjectIndexA DataType for storing the index of a Measure.
Fields
value::Int64: The index value.
InfiniteOpt.MeasureData — TypeMeasureData{M <: Measure} <: AbstractDataObjectA mutable DataType for storing Measures and their data.
Fields
measure::M: The measure structure.name::String: The base name used for printingname(meas_expr d(par)).measure_indices::Vector{MeasureIndex}: Indices of dependent measures.constraint_indices::Vector{InfOptConstraintIndex}: Indices of dependent constraints.derivative_indices::Vector{DerivativeIndex}: Indices of dependent derivatives.in_objective::Bool: Is this used in objective?
InfiniteOpt.MeasureRef — TypeMeasureRef <: DispatchVariableRefA DataType for referring to measure abstractions.
Fields
model::InfiniteModel: Infinite model.index::MeasureIndex: Index of the measure in model.
Integrals
InfiniteOpt.MeasureToolbox.integral — Methodintegral(expr::JuMP.AbstractJuMPScalar,
pref::GeneralVariableRef,
[lower_bound::Real = lower_bound(pref),
upper_bound::Real = upper_bound(pref);
kwargs...])::GeneralVariableRefReturns a measure reference that evaluates the integral of expr with respect to infinite parameter pref from lower_bound to upper_bound. This thus considers integrals of the form: $\int_{p \in P} expr(p) w(p) dp$ where $p$ is an infinite parameter and $w$ is the weight function is 1 by default. This function provides a high-level interface that ultimately constructs an appropriate concrete form of AbstractMeasureData via generate_integral_data in accordance with the keyword arugment eval_method that is then used with measure. Note that it is preferred to call @integral when expr is not just a single variable reference. Errors for bad bound input.
The keyword arguments are as follows:
eval_method::AbstractUnivariateMethod: Used to determine the numerical evaluation scheme. Possible choices include:num_supports: The minimum number of supports to be generated (if used byeval_method)weight_func: $w(p)$ above with parameter value inputs and scalar output
See set_uni_integral_defaults to update the default keyword argument values for all one-dimensional integral calls.
Example
julia> @infinite_parameter(model, x in [0, 1])
x
julia> @variable(model, f, Infinite(x))
f(x)
julia> int = integral(f, x)
∫{x ∈ [0, 1]}[f(x)]
julia> expand(int)
0.2 f(0.8236475079774124) + 0.2 f(0.9103565379264364) + 0.2 f(0.16456579813368521) + 0.2 f(0.17732884646626457) + 0.2 f(0.278880109331201)InfiniteOpt.MeasureToolbox.∫ — Method∫(expr::JuMP.AbstractJuMPScalar,
pref::GeneralVariableRef,
[lower_bound::Real = NaN,
upper_bound::Real = NaN;
kwargs...])::GeneralVariableRefA convenient wrapper for integral. The ∫ unicode symbol is produced via \int.
InfiniteOpt.MeasureToolbox.uni_integral_defaults — Functionuni_integral_defaults()::Dict{Symbol, Any}Get the default keyword argument values for defining one-dimensional integrals.
julia> uni_integral_defaults()
Dict{Symbol,Any} with 1 entry:
:eval_method => Automatic()InfiniteOpt.MeasureToolbox.set_uni_integral_defaults — Functionset_uni_integral_defaults(; kwargs...)::NothingSet the default keyword argument settings for one-dimensional integrals. The keyword arguments of this function will be recorded in the default keyword argument dictionary. These will determine the default keyword argument values when calling integral with a single infinite parameter.
Example
julia> uni_integral_defaults()
Dict{Symbol,Any} with 1 entry:
:eval_method => Automatic()
julia> set_uni_integral_defaults(num_supports = 5, eval_method = Quadrature(),
new_kwarg = true)
julia> uni_integral_defaults()
Dict{Symbol,Any} with 3 entries:
:new_kwarg => true
:num_supports => 5
:eval_method => Quadrature()InfiniteOpt.MeasureToolbox.clear_uni_integral_defaults — Functionclear_uni_integral_defaults()::NothingClears and resets the keyword argument defaults for univariate integrals to their default state.
Example
julia> uni_integral_defaults()
Dict{Symbol,Any} with 3 entries:
:new_kwarg => true
:num_supports => 5
:eval_method => Quadrature()
julia> clear_uni_integral_defaults()
julia> uni_integral_defaults()
Dict{Symbol,Any} with 1 entry:
:eval_method => Automatic()InfiniteOpt.MeasureToolbox.integral — Methodintegral(expr::JuMP.AbstractJuMPScalar,
prefs::Array{GeneralVariableRef},
[lower_bounds::Union{Real, Array{<:Real}} = [lower_bound(pref)...],
upper_bounds::Union{Real, Array{<:Real}} = [upper_bound(pref)...];
kwargs...])::GeneralVariableRefReturns a measure reference that evaluates the integral of expr with respect to infinite parameters prefs from lower_bounds to upper_bounds. This thus considers integrals of the form: $\int_{p \in P} expr(p) w(p) dp$ where $p$ is an infinite parameter and $w$ is the weight function is 1 by default. This function provides a high-level interface that ultimately constructs an appropriate concrete form of AbstractMeasureData via generate_integral_data in accordance with the keyword arugment eval_method that is then used with measure. Note that it is preferred to call @integral when expr is not just a single variable reference. Errors when the container types and dimensions do not match or the bounds are invalid.
The keyword arguments are as follows:
eval_method::AbstractMultivariateMethod: Used to determine the numerical evaluation scheme. Possible choices include:num_supports: The minimum number of supports to be generated (if used byeval_method)weight_func: $w(p)$ above with parameter value inputs and scalar output
See set_multi_integral_defaults to update the default keyword argument values for all multi-dimensional integral calls.
Example
julia> @infinite_parameter(model, x[1:2] in [0, 1], independent = true);
julia> @variable(model, f, Infinite(x));
julia> int = integral(f, x)
∫{x ∈ [0, 1]^2}[f(x)]InfiniteOpt.MeasureToolbox.∫ — Method∫(expr::JuMP.AbstractJuMPScalar,
prefs::Array{GeneralVariableRef},
[lower_bounds::Union{Real, Array{<:Real}} = NaN,
upper_bounds::Union{Real, Array{<:Real}} = NaN;
kwargs...])::GeneralVariableRefA convenient wrapper for integral. The unicode symbol ∫ is produced via \int.
InfiniteOpt.MeasureToolbox.multi_integral_defaults — Functionmulti_integral_defaults()::Dict{Symbol, Any}Get the default keyword argument values for defining multi-dimensional integrals.
julia> multi_integral_defaults()
Dict{Symbol,Any} with 1 entry:
:eval_method => Automatic()InfiniteOpt.MeasureToolbox.set_multi_integral_defaults — Functionset_multi_integral_defaults(; kwargs...)::NothingSet the default keyword argument settings for multi-dimesnional integrals. The keyword arguments of this function will be recorded in the default keyword argument dictionary. These will determine the default keyword argument values when calling integral with an array of infinite parameters.
Example
julia> multi_integral_defaults()
Dict{Symbol,Any} with 1 entry:
:eval_method => Automatic()
julia> set_multi_integral_defaults(num_supports = 5, new_kwarg = true)
julia> multi_integral_defaults()
Dict{Symbol,Any} with 3 entries:
:new_kwarg => true
:num_supports => 5
:eval_method => Automatic()InfiniteOpt.MeasureToolbox.clear_multi_integral_defaults — Functionclear_multi_integral_defaults()::NothingClears and resets the keyword argument defaults for multivariate integrals to their default state.
Example
julia> multi_integral_defaults()
Dict{Symbol,Any} with 3 entries:
:new_kwarg => true
:num_supports => 5
:eval_method => Automatic()
julia> clear_multi_integral_defaults()
julia> multi_integral_defaults()
Dict{Symbol,Any} with 1 entry:
:eval_method => Automatic()InfiniteOpt.MeasureToolbox.@integral — Macro@integral(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, Array{GeneralVariableRef}},
[lower_bounds::Union{Real, Array{<:Real}} = default_bounds,
upper_bounds::Union{Real, Array{<:Real}} = default_bounds;
kwargs...])::GeneralVariableRefAn efficient wrapper for integral and integral. Please see the above doc strings for more information.
InfiniteOpt.MeasureToolbox.@∫ — Macro@∫(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, Array{GeneralVariableRef}},
[lower_bounds::Union{Real, Array{<:Real}} = default_bounds,
upper_bounds::Union{Real, Array{<:Real}} = default_bounds;
kwargs...])::GeneralVariableRefA convenient wrapper for @integral. The unicode symbol ∫ is produced via \int.
InfiniteOpt.MeasureToolbox.AbstractIntegralMethod — TypeAbstractIntegralMethodAn abstract type for integral evaluation methods use in combination with integral and generate_integral_data.
InfiniteOpt.MeasureToolbox.Automatic — TypeAutomatic <: AbstractIntegralMethodAn integral evaluation type for automically selecting an appropriate integral evaluation method. Contains no fields.
InfiniteOpt.MeasureToolbox.AbstractUnivariateMethod — TypeAbstractUnivariateMethod <: AbstractIntegralMethodAn abstract type for integral evaluation methods for 1-dimensional integrals.
InfiniteOpt.MeasureToolbox.UniTrapezoid — TypeUniTrapezoid <: AbstractUnivariateMethodAn integral evalution method that uses the trapezoid rule to in combination with all parameter supports available when the integral is expanded and/or when the infinite model is optimized, whichever comes first. Note this method will ignore the num_supports keyword argument. The upper and lower bounds of the integral will automatically be added as supports. Note this is valid only for finite integral domains. Contains no fields.
InfiniteOpt.MeasureToolbox.UniMCSampling — TypeUniMCSampling <: AbstractUnivariateMethodAn integral evaluation method that uses uniform Monte Carlo sampling to approximate the integral. This variant will add more supports to the model as needed to satisfy num_supports and it will include all supports with the MCSample label up till the integral is expanded and/or when the infinite model is optimized, whichever comes first. Note this is valid only for finite integral domains. Contains no fields.
InfiniteOpt.MeasureToolbox.UniIndepMCSampling — TypeUniIndepMCSampling <: AbstractUnivariateMethodAn integral evaluation method that uses uniform Monte Carlo sampling to approximate the integral similar to UniMCSampling. However, this variant will generate its own set of supports and ignore all other supports with the MCSample label. Note this is valid only for finite integral domains. This is not compatible with individual dependent parameters. Contains no fields.
InfiniteOpt.MeasureToolbox.Quadrature — TypeQuadrature <: AbstractUnivariateMethodA general integral evaluation method that will automatically select the appropriate quadrature method to approximate the integral. Please note that this will generate a unique set of parameter supports and will ignore existing supports when the integral is evaluated and thus should be used with caution. However, this method is able to handle infinite and semi-infinite integral domains. This is not compatible with individual dependent parameters. Contains no fields.
InfiniteOpt.MeasureToolbox.GaussHermite — TypeGaussHermite <: AbstractUnivariateMethodAn integral evaulation method that uses Gauss-Hermite quadrature to evaluate integrals. This is valid for infinite integral domains. It will take this form:
$\int_{-∞}^{∞} f(x) e^{-x^2} \approx \sum_{i=1}^{n} \alpha_i f(x_i)$
Using the weight function: $w(x)$ = $e^{-x^2}$
Note this will generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters.
InfiniteOpt.MeasureToolbox.GaussLegendre — TypeGaussLegendre <: FiniteGaussQuadAn integral evaulation method that uses Gauss-Legendre quadrature to evaluate integrals. This is valid for finite integral domains. Note this will generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters. Contains no fields.
InfiniteOpt.MeasureToolbox.GaussRadau — TypeGaussRadau <: FiniteGaussQuadAn integral evaulation method that uses Gauss-Radau quadrature to evaluate integrals. This is valid for finite integral domains. Note this will generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters. Contains no fields.
InfiniteOpt.MeasureToolbox.GaussLobatto — TypeGaussLobatto <: FiniteGaussQuadAn integral evaulation method that uses Gauss-Lobatto quadrature to evaluate integrals. This is valid for finite integral domains. Note this will generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters. Contains no fields.
InfiniteOpt.MeasureToolbox.GaussJacobi — TypeGaussJacobi <: FiniteGaussQuadAn integral evaulation method that uses Gauss-Jacobi quadrature to evaluate integrals. It will take this form:
$\int_{-1}^{1} f(x) (1-x)^\alpha (1+x)^\beta dx \approx \sum_{i=1}^{n} \alpha_i f(x_i)$
Where,
$(1-x)^\alpha (1+x)^\beta$
is the weight function. This is valid for finite integral domains. This requires the user to input the alpha and beta shape parameters for their function. This will then generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters. If α or β is < -1, an error will be returned.
Fields
α::Float64: Shape parameter that must be > -1β::Float64: Shape parameter that must be > -1
InfiniteOpt.MeasureToolbox.FEGaussLobatto — TypeFEGaussLobatto <: AbstractUnivariateMethodIntegral evaluation method that allows for the user to specify supports to be included in quadrature evaluation. The upper and lower bounds of the integral will automatically be added as supports. This method uses Gauss Lobatto quadrature to decompose the overall Integral into smaller integrals that span the user defined supports as follows:
$\int_{x_1}^{x_3} f(x) dx = \int_{x_1}^{x_2} f(x) dx + \int_{x_2}^{x_3} f(x) dx$
where the integrals are evaluated using Gauss Lobatto quadrature:
$\int f(x) dx \approx \sum_{i=1}^{n} \alpha_i f(x_i)$
InfiniteOpt.MeasureToolbox.GaussChebyshev — TypeGaussChebyshev <: FiniteGaussQuadAn integral evaulation method that uses Gauss-Chebyshev quadrature to evaluate integrals. This is valid for finite integral domains. This requires the user to input the order of Guass-Chebyshev Quadrature they want to use. If the order is not between 1 and 4 an error will be returned. The integral evaluated is as follows:
$\int_{-1}^{1} f(x) w(x) \approx \sum_{i=1}^{n} \alpha_i f(x_i)$
The weight functions are as follows:
- 1st order: $w(x) = \frac{1}{\sqrt{1-x^2}}$
- 2nd order: $w(x) = {\sqrt{1-x^2}}$
- 3rd order: $w(x) = \sqrt{(1+x)/(1-x)}$
- 4th order: $w(x) = \sqrt{(1-x)/(1+x)}$
This will then generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters.
Fields
order::Int: Specifies the order of Gauss-Chebyshev Quadrature. Must be between 1 and 4.
InfiniteOpt.MeasureToolbox.GaussLaguerre — TypeGaussLaguerre <: AbstractUnivariateMethodAn integral evaulation method that uses Gauss-Laguerre quadrature to evaluate integrals. This is valid for semi-infinite integral domains.
This method evaluates the following integral:
$\int_{0}^{+∞} f(x) e^{-x} \approx \sum_{i=1}^{n} \alpha_i f(x_i)$
Using the weight function:
$w(x) = e^{-x}$
Note this will generate its own set of supports and will ignore other parameter supports. This is not compatible with individual dependent parameters.
InfiniteOpt.MeasureToolbox.AbstractMultivariateMethod — TypeAbstractMultivariateMethod <: AbstractIntegralMethodAn abstract type for integral evaluation methods for multi-dimensional integrals.
InfiniteOpt.MeasureToolbox.MultiMCSampling — TypeMultiMCSampling <: AbstractMultivariateMethodAn integral evaluation method that uses uniform Monte Carlo sampling to approximate the integral. This variant will add more supports to the model as needed to satisfy num_supports and it will include all supports with the MCSample label up till the integral is expanded and/or when the infinite model is optimized, whichever comes first. Note this is valid only for finite integral domains. If an array of independent infinite parameters is specified, they must use the same amount of supports. Contains no fields.
InfiniteOpt.MeasureToolbox.MultiIndepMCSampling — TypeMultiIndepMCSampling <: AbstractMultivariateMethodAn integral evaluation method that uses uniform Monte Carlo sampling to approximate the integral similar to MultiMCSampling. However, this variant will generate its own set of supports and ignore all other supports with the MCSample label. Note this is valid only for finite integral domains. Contains no fields.
InfiniteOpt.MeasureToolbox.generate_integral_data — Functiongenerate_integral_data(
prefs::Union{InfiniteOpt.GeneralVariableRef, Vector{InfiniteOpt.GeneralVariableRef}},
lower_bounds::Union{Real, Vector{<:Real}},
upper_bounds::Union{Real, Vector{<:Real}},
method::V; [num_supports::Int = InfiniteOpt.DefaultNumSupports,
weight_func::Function = InfiniteOpt.default_weight,
extra_kwargs...]
)::InfiniteOpt.AbstractMeasureData where {V <: AbstractIntegralMethod}Generate the appropriate concrete realization of AbstractMeasureData using method. Here prefs, lower_bounds, and upper_bounds will always have a 1 to 1 correspondence when this is called from integral. Please refer to the method docstrings for an explanation of each one.
User-defined method extensions should first define a concrete method type inheriting from AbstractUnivariateMethod or AbstractMultivariateMethod as appropriate and then implement extend this method using that type for method.
Expectations
InfiniteOpt.MeasureToolbox.expect — Functionexpect(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, Array{GeneralVariableRef};
[num_supports::Int = DefaultNumSupports])::GeneralVariableRefMakes a measure for expr based on its expectation with respect to prefs. For prefs with distribution domains this is essentially equivalent to
1/total_num_supports * support_sum(expr, prefs, label = WeightedSample)Thus, for these domain types it only considers supports that are added to prefs via generation on creation (i.e., specifying the num_supports keyword when creating prefs). For incorporating other supports consider calling integral and using the weight_func argument to specify the probability density function.
For a single infinite parameter defined over a bounded interval domain the syntax becomes:
expect(expr::JuMP.AbstractJuMPScalar,
prefs::GeneralVariableRef;
[num_supports::Int = DefaultNumSupports,
pdf::Function = (supp) -> 1 / (ub - lb)])::GeneralVariableRefThe behavior with the default pdf is equivalent to evaluating the mean value theorem for integrals for expr with respect to pref using UniTrapezoid. Other density functions can be given via pdf. Errors if the interval domain is not bounded.
Note that num_supports should be 0 if a single dependent parameter is given. Also, note that it is preferred to call @expect when expr is not just a single variable reference.
Example
julia> @infinite_parameter(model, x ~ Normal(), num_supports = 2)
x
julia> @variable(model, f, Infinite(x))
f(x)
julia> meas = expect(f, x)
𝔼{x}[f(x)]
julia> expand(meas)
0.5 f(0.6791074260357777) + 0.5 f(0.8284134829000359)InfiniteOpt.MeasureToolbox.@expect — Macro@expect(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, AbstractArray{GeneralVariableRef};
[num_supports::Int = DefaultNumSupports,
kwargs...]
)::GeneralVariableRefAn efficient wrapper for expect. Please see its doc string more information.
InfiniteOpt.MeasureToolbox.𝔼 — Function𝔼(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, AbstractArray{GeneralVariableRef}};
[num_supports::Int = DefaultNumSupports,
kwargs...]
)::GeneralVariableRef)A convenient wrapper for expect. The unicode symbol 𝔼 is produced by \bbE.
InfiniteOpt.MeasureToolbox.@𝔼 — Macro@𝔼(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, AbstractArray{GeneralVariableRef};
[num_supports::Int = DefaultNumSupports],
kwargs...)::GeneralVariableRefA convenient wrapper for @expect. The unicode symbol 𝔼 is produced by \bbE.
InfiniteOpt.MeasureToolbox.generate_expect_data — Functiongenerate_expect_data(domain::AbstractInfiniteDomain,
prefs::Union{GeneralVariableRef, Vector{GeneralVariableRef}},
num_supports::Int;
[kwargs...])::AbstractMeasureDataGenerate a concrete instance of AbstractMeasureData in accordance with the domain and infinite parameter(s) prefs given for computing the expectation. This is intended as an internal method, but should be extended for user defined infinite domain types.
Support Sum
InfiniteOpt.MeasureToolbox.@support_sum — Macro@support_sum(expr::JuMP.AbstractJuMPScalar,
prefs::Union{GeneralVariableRef, Array{GeneralVariableRef}};
label = All
)::GeneralVariableRefAn efficient wrapper for support_sum please see its doc string for more information.
InfiniteOpt.MeasureToolbox.support_sum — Functionsupport_sum(expr::JuMP.AbstractJuMPScalar,
params::Union{GeneralVariableRef, Array{GeneralVariableRef}};
label = All
)::GeneralVariableRefCreates a measure that represents the sum of the expression over a parameter(s) using all of its supports corresponding to label. Also, note that it is preferred to call @support_sum when expr is not just a single variable reference.
Example
julia> @infinite_parameter(model, x in [0, 1], supports = [0.3, 0.7])
x
julia> @variable(model, f, Infinite(x))
f(x)
julia> meas = support_sum(f, x)
support_sum{x}[f(x)]
julia> expand(meas)
f(0.3) + f(0.7)Queries
JuMP.name — MethodJuMP.name(mref::MeasureRef)::StringExtend JuMP.name to return the name associated with a measure reference.
InfiniteOpt.num_measures — Functionnum_measures(model::InfiniteModel)::IntReturn the number of measures defined in model.
Example
julia> num_measures(model)
2InfiniteOpt.all_measures — Functionall_measures(model::InfiniteModel)::Vector{GeneralVariableRef}Return the list of all measures added to model.
Examples
julia> all_measures(model)
2-element Array{GeneralVariableRef,1}:
∫{t ∈ [0, 6]}[w(t, x)]
𝔼{x}[w(t, x)]InfiniteOpt.measure_function — Functionmeasure_function(mref::MeasureRef)::JuMP.AbstractJuMPScalarReturn the function associated with mref.
Example
julia> measure_function(meas)
y(x, t) + 2measure_function(mref::GeneralVariableRef)Define measure_function for general variable references. Errors if mref does not correspond to a MeasureRef. See the underlying docstrings for more information.
InfiniteOpt.measure_data — Functionmeasure_data(mref::MeasureRef)::AbstractMeasureDataReturn the measure data associated with mref.
Example
julia> data = measure_data(meas);
julia> typeof(data)
FunctionalDiscreteMeasureData{Vector{GeneralVariableRef},Vector{Float64}}measure_data(mref::GeneralVariableRef)Define measure_data for general variable references. Errors if mref does not correspond to a MeasureRef. See the underlying docstrings for more information.
InfiniteOpt.is_analytic — Functionis_analytic(mref::MeasureRef)::BoolReturn if mref is evaluated analytically.
Example
julia> is_analytic(meas)
falseis_analytic(mref::GeneralVariableRef)Define is_analytic for general variable references. Errors if mref does not correspond to a MeasureRef. See the underlying docstrings for more information.
InfiniteOpt.parameter_refs — Methodparameter_refs(mref::MeasureRef)::TupleReturn the tuple of infinite parameters that the measured expression associated mref depends on once the measure has been evaluated. Note that this will correspond to the parameter dependencies of the measure function excluding those included in the measure data.
Example
julia> parameter_refs(meas)
(t,)InfiniteOpt.is_used — Methodis_used(mref::MeasureRef)::BoolReturn a Bool indicating if mref is used in the model.
Example
julia> is_used(mref)
trueInfiniteOpt.used_by_derivative — Methodused_by_derivative(mref::MeasureRef)::BoolReturn a Bool indicating if mref is used by a derivative.
Example
julia> used_by_derivative(mref)
trueInfiniteOpt.used_by_constraint — Methodused_by_constraint(mref::MeasureRef)::BoolReturn a Bool indicating if mref is used by a constraint.
Example
julia> used_by_constraint(mref)
falseInfiniteOpt.used_by_measure — Methodused_by_measure(mref::MeasureRef)::BoolReturn a Bool indicating if mref is used by a measure.
Example
julia> used_by_measure(mref)
trueInfiniteOpt.used_by_objective — Methodused_by_objective(mref::MeasureRef)::BoolReturn a Bool indicating if mref is used by the objective.
Example
julia> used_by_objective(mref)
trueInfiniteOpt.core_object — Methodcore_object(mref::MeasureRef)::MeasureRetrieve the underlying core Measure object for vref. This is intended as an advanced method for developers.
InfiniteOpt.parameter_group_int_indices — Methodparameter_group_int_indices(mref::MeasureRef)::Vector{Int}Return the list of infinite parameter group integer indices used by mref.
Modification
JuMP.set_name — MethodJuMP.set_name(mref::MeasureRef, name::String)::NothingExtend JuMP.set_name to specify the name of a measure reference.
JuMP.delete — MethodJuMP.delete(model::InfiniteModel, mref::MeasureRef)::NothingExtend JuMP.delete to delete measures. Errors if measure is invalid, meaning it does not belong to the model or it has already been deleted.
Example
julia> print(model)
Min ∫{t ∈ [0, 6]}[g(t)] + z
Subject to
z ≥ 0.0
∫{t ∈ [0, 6]}[g(t)] = 0
g(t) + z ≥ 42.0, ∀ t ∈ [0, 6]
g(0.5) = 0
julia> delete(model, meas)
julia> print(model)
Min z
Subject to
z ≥ 0.0
0 = 0
g(t) + z ≥ 42.0, ∀ t ∈ [0, 6]
g(0.5) = 0Expansion
InfiniteOpt.expand — Functionexpand(mref::MeasureRef)::JuMP.AbstractJuMPScalarReturn a JuMP scalar function containing the explicit expansion of the measure mref. This expansion is done according to the measure data. Note that variables are added to the model as necessary to accomodate the expansion (i.e., point variables and semi-infinite variables are made as needed). Errors if expansion is undefined for the measure data and/or the measure expression. If desired this can be used in combination with measure to expand measures on the fly.
This is useful for extensions that employ a custom transformation backend since it can be used evaluate measures before expressions are translated to the new model. This method can also be extended to handle custom measure data types by extending expand_measure. Optionally, analytic_expansion can also be extended which is triggered by is_analytic for such types if analytic expansion is possible in certain cases.
Example
julia> tdata = DiscreteMeasureData(t, [0.5, 0.5], [0, 1])
julia> expr = expand(measure(g + z + T - h - 2, tdata))
0.5 g(0) + 0.5 g(1) + z + 0.5 T(0, x) + 0.5 T(1, x) - h(x) - 2expand(mref::GeneralVariableRef)Define expand for general variable references. Errors if mref does not correspond to a MeasureRef. See the underlying docstrings for more information.
InfiniteOpt.expand_all_measures! — Functionexpand_all_measures!(model::InfiniteModel)::NothingExpand all of the measures used in the objective and/or constraints of model. The objective and constraints are updated accordingly. Note that variables are added to the model as necessary to accomodate the expansion (i.e., point variables and semi-infinite variables are made as needed). Errors if expansion is undefined for the measure data and/or the measure expression.
This is useful for extensions that employ a custom transformation backend since it can be used evaluate measures before model is translated into the new model. This method can also be extended to handle custom measure data types by extending expand_measure. Note that this method leverages expand_measure via expand_measures. Optionally, analytic_expansion can also be extended which is triggered by is_analytic for such types if analytic expansion is possible in certain cases.
Example
julia> print(model)
Min integral{t ∈ [0, 6]}[g(t)*t] + z
Subject to
T(t, x) ≥ 0.0, ∀ t ∈ [0, 6], xi ∈ [-1, 1]
z ≥ 0.0
g(t) + z ≥ 42.0, ∀ t ∈ [0, 6]
integral{t ∈ [0, 6]}[T(t, x)] ≥ 0.0, ∀ x ∈ [-1, 1]
julia> expand_all_measures!(model)
julia> print(model)
Min 3 g(6) + z
Subject to
T(t, x) ≥ 0.0, ∀ t ∈ [0, 6], xi ∈ [-1, 1]
z ≥ 0.0
g(t) + z ≥ 42.0, ∀ t ∈ [0, 6]
0.5 T(0, x) + 0.5 T(6, xi) ≥ 0.0, ∀ x ∈ [-1, 1]InfiniteOpt.expand_measure — Functionexpand_measure(
expr,
data::AbstractMeasureData,
write_model::Union{InfiniteModel, AbstractTransformationBackend}
)::JuMP.AbstractJuMPScalarReturn the finite reformulation of a measure containing a variable/parameter expression expr with measure data data. Here write_model is the target model where this expanded expression will be used. Thus, any variables that need to be created will be added to write_model. The methods make_point_variable_ref and make_semi_infinite_variable_ref should be used as appropriate to create these variables. Note this is intended as an internal function, but will need to be extended for unsupported expr types and for user-defined measure data types. Principally, this is leveraged to enable the user methods expand and expand_all_measures!.
InfiniteOpt.analytic_expansion — Functionanalytic_expansion(
expr,
data::AbstractMeasureData,
write_model::Union{InfiniteModel, AbstractTransformationBackend}
)::JuMP.AbstractJuMPScalarAnalytically, evaluate measure in the simple case where the measure expression expr doesn't depend on data and thus expr can be treated as a constant in conjunction with an analytic result of the data. This is intended as an internal method that is used by expand and expand_measures. For unrecognized data types, expand_measure is called instead. User defined measure data type may choose to extend this method if desired. This is triggered when is_analytic(mref) = true.
InfiniteOpt.expand_measures — Functionexpand_measures(
expr,
write_model::Union{InfiniteModel, AbstractTransformationBackend}
)Expand all MeasureRefs in expr in-place via expand_measure and return the expanded expression. This is an internal method used by expand_all_measures! and TranscriptionOpt but can be useful for user-defined transformation backend extensions that add implement add_point_variable/add_semi_infinite_variable in combination with expand_measure. write_model is the model that the measure variables are added to as described in expand_measure.
InfiniteOpt.make_point_variable_ref — Functionmake_point_variable_ref(
write_model::Union{InfiniteModel, AbstractTransformationBackend},
ivref::GeneralVariableRef,
support::Vector{Float64}
)::GeneralVariableRefMake a point variable for infinite variable/derivative ivref at support, add it to the write_model, and return the GeneralVariableRef. This is an internal method for point variables produced by expanding measures via expand_measure. This is also useful for those writing extension transformation backends and wish to expand measures without modifiying the InfiniteModel. In such cases, write_model should be the transformation backend and add_point_variable should be extended appropriately for point variables. Errors if write_model is an transformation backend and add_point_variable is not properly extended.
Note this is also accomodates infinite parameter functions, in which case the infinite parameter function is called with the support as input.
InfiniteOpt.make_semi_infinite_variable_ref — Functionmake_semi_infinite_variable_ref(
write_model::Union{InfiniteModel, AbstractTransformationBackend},
ivref::GeneralVariableRef,
eval_support::Vector{Float64}
)::GeneralVariableRefMake a semi-infinite variable for infinite variable/derivative/parameter function ivref at eval_support (see SemiInfiniteVariable), add it to the write_model, and return the GeneralVariableRef. This is an internal method for semi-infinite variables produced by expanding measures via expand_measure. This is also useful for those writing extension transformation backends and wish to expand measures without modifiying the InfiniteModel. In such cases, write_model should be the transformation backend and add_semi_infinite_variable should be extended appropriately for semi-infinite variables. Errors if write_model is an transformation backend and add_semi_infinite_variable is not properly extended. Note this is only intended for transformation backends that are currently stored in InfiniteModel.backend.
InfiniteOpt.add_point_variable — Methodadd_point_variable(
backend::AbstractTransformationBackend,
var::PointVariable
)::GeneralVariableRefAdd a point variable var to the tranformation backend backend and return the correct InfiniteOpt variable reference. This is an internal method used by make_point_variable_ref to make point variables when the write_model is an transformation backend. This is useful for extensions that wish to expand measures, but without changing the original InfiniteModel. An error is thrown for unextended backend types.
InfiniteOpt.add_semi_infinite_variable — Methodadd_semi_infinite_variable(
backend::AbstractTransformationBackend,
var::SemiInfiniteVariable,
)::GeneralVariableRefAdd a semi-infinite variable var to the transformation backend backend and return the correct InfiniteOpt variable reference. This is an internal method used by make_semi_infinite_variable_ref to make semi-infinite variables when the write_model is a transformation backend. This is useful for extensions that wish to expand measures, but without changing the original InfiniteModel. An error is thrown for new transformation backend types. Note if this is extended, than internal_semi_infinite_variable should also be extended in order to direct semi-infinite variables references to the underlying SemiInfiniteVariable.
InfiniteOpt.internal_semi_infinite_variable — Functioninternal_semi_infinite_variable(
vref::SemiInfiniteVariableRef,
backend::AbstractTransformationBackend
)::SemiInfiniteVariableReturn the semi-infinite variable object of vref assuming it is an internal variable made during measure expansion within a transformation backend. This will apply to transformation backend extensions that utilize add_measure_variable in combination with expand_measure.