opt_model

class opt_model

Bases: mp_idx_manager

opt_model - MP-Opt-Model optimization model class.

OM = OPT_MODEL
OM = OPT_MODEL(S)

This class implements the optimization model object used to encapsulate
a given optimization problem formulation. It allows for access to
optimization variables, constraints and costs in named blocks, keeping
track of the ordering and indexing of the blocks as variables,
constraints and costs are added to the problem.

Below are the list of available methods for use with the Opt Model class.
Please see the help on each individual method for more details:

Modify the OPF formulation by adding named blocks of costs, constraints,
or variables:
    add_quad_cost
    add_nln_cost
    add_lin_constraint
    add_nln_constraint
    add_var
    init_indexed_name

Modify parameters of existing costs, constraints, or variables:
    set_params

Return the number of linear constraints, nonlinear constraints, or
variables, optionally for a single named block:
    getN

Return index structure for variables, linear and nonlinear constraints
and costs:
    get_idx

Return the intial values, bounds and type for optimization variables:
    params_var

Build and return constraint parameters and evaluate constraints:
    params_lin_constraint
    params_nln_constraint
    eval_nln_constraint
    eval_nln_constraint_hess

Build and return cost parameters and evaluate costs:
    params_nln_cost
    params_quad_cost
    eval_nln_cost
    eval_quad_cost

Determine model type:
    problem_type

Solve the model, access, and display the solution:
    solve
    display_soln
    get_soln
    parse_soln
    is_solved
    has_parsed_soln

Retreive user data in the model object:
    get_userdata

Display the object (called automatically when you omit the semicolon
at the command-line):
    display

Utility methods for dealing with varsets:
    varsets_cell2struct
    varsets_idx
    varsets_len
    varsets_x

Return the value of an individual field:
    get

Indentify variable, constraint or cost row indices:
    describe_idx

Make a shallow copy of the object:
    copy

The following is the structure of the data in the Opt-Model object.
Each field of .idx or .data is a struct whose field names are the names
of the corresponding blocks of vars, constraints or costs (found in
order in the corresponding .order field). The description next to these
fields gives the meaning of the value for each named sub-field.
E.g. om.var.data.v0.Pg contains a vector of initial values for the 'Pg'
block of variables.

om
    .var        - data for optimization variable sets that make up
                  the full optimization variable x
        .idx
            .i1 - starting index within x
            .iN - ending index within x
            .N  - number of elements in this variable set
        .N      - total number of elements in x
        .NS     - number of variable sets or named blocks
        .data   - bounds and initial value data
            .v0 - vector of initial values
            .vl - vector of lower bounds
            .vu - vector of upper bounds
            .vt - scalar or vector of variable types
                    'C' = continuous
                    'I' = integer
                    'B' = binary
        .order  - struct array of names/indices for variable
                  blocks in the order they appear in x
            .name   - name of the block, e.g. Pg
            .idx    - indices for name, {2,3} => Pg(2,3)
        .params - cache for previously assembled aggregate parameters
            .v0  - aggregate vector of variable initial values
            .vl  - aggregate vector of variable lower bounds
            .vu  - aggregate vector of variable upper bounds
            .vt  - aggregate vector of variable types
    .nle        - data for nonlinear equality constraints that make up the
                  full set of nonlinear constraints ghne(x)
        .idx
            .i1 - starting index within ghne(x)
            .iN - ending index within ghne(x)
            .N  - number of elements in this constraint set
        .N      - total number of elements in ghne(x)
        .NS     - number of nonlinear constraint sets or named blocks
        .data   - data for nonlinear constraints
            .fcn - function handle for constraints/gradient evaluation
            .hess - function handle for Hessian evaluation
            .vs - cell array of variable sets that define the xx for
                  this constraint block
        .order  - struct array of names/indices for nonlinear constraint
                  blocks in the order they appear in ghne(x)
            .name   - name of the block, e.g. Pmis
            .idx    - indices for name, {2,3} => Pmis(2,3)
    .nli        - data for nonlinear inequality constraints that make up the
                  full set of nonlinear constraints ghni(x)
        .idx
            .i1 - starting index within ghni(x)
            .iN - ending index within ghni(x)
            .N  - number of elements in this constraint set
        .N      - total number of elements in ghni(x)
        .NS     - number of nonlinear constraint sets or named blocks
        .data   - data for nonlinear constraints
            .fcn - function handle for constraints/gradient evaluation
            .hess - function handle for Hessian evaluation
            .vs - cell array of variable sets that define the xx for
                  this constraint block
        .order  - struct array of names/indices for nonlinear constraint
                  blocks in the order they appear in ghni(x)
            .name   - name of the block, e.g. Pmis
            .idx    - indices for name, {2,3} => Pmis(2,3)
    .lin        - data for linear constraints that make up the
                  full set of linear constraints ghl(x)
        .idx
            .i1 - starting index within ghl(x)
            .iN - ending index within ghl(x)
            .N  - number of elements in this constraint set
        .N      - total number of elements in ghl(x)
        .NS     - number of linear constraint sets or named blocks
        .data   - data for l <= A*xx <= u linear constraints
            .A  - sparse linear constraint matrix
            .l  - left hand side vector, bounding A*x below
            .u  - right hand side vector, bounding A*x above
            .vs - cell array of variable sets that define the xx for
                  this constraint block
        .order  - struct array of names/indices for linear constraint
                  blocks in the order they appear in ghl(x)
            .name   - name of the block, e.g. Pmis
            .idx    - indices for name, {2,3} => Pmis(2,3)
        .params - cache for previously assembled aggregate parameters
            .A  - aggregate sparse linear constraint matrix
            .l  - aggregate left hand side vector, bounding A*x below
            .u  - aggregate right hand side vector, bounding A*x above
    .qdc       - quadratic costs
        .idx
            .i1 - starting row index within quadratic costs
            .iN - ending row index within quadratic costs
            .N  - number of rows in this quadratic cost block
        .N      - total number of rows in quadratic costs
        .NS     - number of quadratic cost blocks
        .data   - data for each quadratic cost block
            .Q  - sparse matrix (or vector) of quadratic cost coefficients
            .c  - column vector of linear cost coefficients
            .k  - constant term
            .vs - cell array of variable sets that define xx for this
                  quadratic cost block, where sizes of Q, c, k conform to xx
        .order  - struct array of names/indices for quadratic cost blocks
                  in the order the were added
            .name   - name of the block, e.g. R
            .idx    - indices for name, {2,3} => R(2,3)
        .params - cache for previously assembled aggregate parameters
            .Q  - aggregate sparse matrix of quadratic cost coefficients
            .c  - aggregate column vector of linear cost coefficients
            .k  - aggregate constant term
    .nlc       - general nonlinear costs
        .idx
            .i1 - starting row index within nonlinear costs
            .iN - ending row index within nonlinear costs
            .N  - number of rows in this nonlinear cost block
                  (always equal to 1 for nonlinear cost blocks)
        .N      - total number of rows in nonlinear costs
        .NS     - number of nonlinear cost blocks
        .data   - data for each nonlinear cost block
            .fcn - function handle for cost, gradient, Hessian evaluation
            .vs - cell array of variable sets that define xx for this
                  nonlinear cost block, where xx is the input to the
                  evaluation function
        .order  - struct array of names/indices for nonlinear cost blocks
                  in the order they were added
            .name   - name of the block, e.g. R
            .idx    - indices for name, {2,3} => R(2,3)
    .prob_type  - used to cache the return value of problem_type()
    .soln       - struct containing the output of the solve() method
                  with the following fields
        .eflag  - exit flag
        .output - output struct with the following fields
            .alg     - algorithm code of solver used
            (others) - solver specific fields
        .x      - solution vector
        .f      - final (objective) function value
        .jac    - final Jacobian matrix (if available, for LEQ/NLEQ probs)
        .lambda - Lagrange and Kuhn-Tucker multipliers on constraints
    .userdata   - any user defined data
        .(user defined fields)

See also mp_idx_manager.

Constructor Summary
opt_model(varargin)

Constructor.

om = opt_model()
om = opt_model(a_struct)
om = opt_model(an_om)
Property Summary
var

variables

lin

linear constraints

nle

nonlinear equality constraints

nli

nonlinear inequality constraints

qdc

quadratic costs

nlc

general nonlinear costs

prob_type = ''

problem type

soln = struct('eflag',[], ... 'lambda',[])

results of solve()

Method Summary
def_set_types(om)

Define set types var, lin, nle, nli, qdc, nlc.

init_set_types(om)

Initialize data structures for each set type.

is_solved(om)

Return true if model has been solved.

has_parsed_soln(om)

Return true if model has a parsed solution.

params_nln_cost(om, name, idx)

params_nln_cost() - Returns cost parameters for general nonlinear costs.

[N, FCN] = OM.PARAMS_NLN_COST(NAME)
[N, FCN] = OM.PARAMS_NLN_COST(NAME, IDX_LIST)
[N, FCN, VS] = OM.PARAMS_NLN_COST(...)

Returns the parameters N and FCN provided when the corresponding
named general nonlinear cost set was added to the model. Likewise
for indexed named sets specified by NAME and IDX_LIST.

An optional 3rd output argument VS indicates the variable sets used by
this cost set.

See also opt_model, add_nln_cost(), eval_nln_cost().

params_quad_cost(om, name, idx)

params_quad_cost() - Returns the cost parameters for quadratic costs.

[Q, C] = OM.PARAMS_QUAD_COST()
[Q, C] = OM.PARAMS_QUAD_COST(NAME)
[Q, C] = OM.PARAMS_QUAD_COST(NAME, IDX_LIST)
[Q, C, K] = OM.PARAMS_QUAD_COST(...)
[Q, C, K, VS] = OM.PARAMS_QUAD_COST(...)

With no input parameters, it assembles and returns the parameters
for the aggregate quadratic cost from all quadratic cost sets added
using ADD_QUAD_COST. The values of these parameters are cached
for subsequent calls. The parameters are Q, C, and optionally K,
where the quadratic cost is of the form
    F(X) = 1/2 * X'*Q*X + C'*X + K

If a NAME is provided then it simply returns the parameters for the
corresponding named set. Likewise for indexed named sets specified
by NAME and IDX_LIST. In this case, Q and K may be vectors, corresponding
to a cost function of the form
    F(X) = 1/2 * Q .* X.^2 + C .* X + K

An optional 4th output argument VS indicates the variable sets used by
this cost set. The size of Q and C will be consistent with VS.

See also opt_model, add_quad_cost().

display(om, varargin)

display() - Displays the object.

Called when semicolon is omitted at the command-line. Displays the details of the variables, constraints, costs included in the model.

See also opt_model.

is_mixed_integer(om)

is_mixed_integer() - Return true if model is mixed integer, false otherwise.

TorF = OM.IS_MIXED_INTEGER()

Outputs:
    TorF : 1 or 0, indicating whether any of the variables are
           binary or integer

See also opt_model.

parse_soln(om, stash)

parse_soln() - Parse solution vector and shadow prices by all named sets.

PS = OM.PARSE_SOLN()
OM.PARSE_SOLN(STASH)

For a solved model, OM.PARSE_SOLN() returns a struct of parsed
solution vector and shadow price values for each named set of
variables and constraints. The returned PS (parsed solution) struct
has the following format, where each of the terminal elements is a
struct with fields corresponding to the respective named sets:

Output:
    PS
        .var
            .val
            .mu_l
            .mu_u
        .lin
            .mu_l
            .mu_u
        .nle
            .lam
        .nli
            .mu

The value of each element in the returned struct can be obtained
via the GET_SOLN method as well, but using PARSE_SOLN is generally
more efficient if a complete set of values is needed.

If the optional STASH input argument is present and true, the fields
of the return struct are copied to OM.SOLN.

See also get_soln().

eval_nln_constraint_hess(om, x, lam, iseq)

eval_nln_constraint_hess() - Builds and returns Hessian of nonlinear constraints.

D2G = OM.EVAL_NLN_CONSTRAINT_HESS(X, LAM, ISEQ)
Builds the Hessian of the full set of nonlinear equality or inequality
constraints (ISEQ equal to 1 or 0, respectively) for given values of
the optimization vector X and dual variables LAM, based on constraints
added by ADD_NLN_CONSTRAINT.

    g(X) = 0
    h(X) <= 0

Example:
    d2G = om.eval_nln_constraint_hess(x, lam, 1)
    d2H = om.eval_nln_constraint_hess(x, lam, 0)

See also opt_model, add_nln_constraint(), eval_nln_constraint().

eval_nln_cost(om, x, name, idx)

eval_nln_cost() - Evaluates individual or full set of general nonlinear costs.

[F, DF, D2F] = OM.EVAL_NLN_COST(X)
[F, DF, D2F] = OM.EVAL_NLN_COST(X, NAME)
[F, DF, D2F] = OM.EVAL_NLN_COST(X, NAME, IDX_LIST)
Evaluates the cost function and its derivatives for an individual named
set or the full set of general nonlinear costs for a given value of the
optimization vector X, based on costs added by ADD_NLN_COST.

Example:
    [f, df, d2f] = om.eval_nln_cost(x)
    [f, df, d2f] = om.eval_nln_cost(x, name)
    [f, df, d2f] = om.eval_nln_cost(x, name, idx_list)

See also opt_model, add_nln_cost(), params_nln_cost().

params_lin_constraint(om, name, idx)

params_lin_constraint() - Builds and returns linear constraint parameters.

[A, L, U] = OM.PARAMS_LIN_CONSTRAINT()
[A, L, U] = OM.PARAMS_LIN_CONSTRAINT(NAME)
[A, L, U] = OM.PARAMS_LIN_CONSTRAINT(NAME, IDX_LIST)
[A, L, U, VS] = OM.PARAMS_LIN_CONSTRAINT(...)
[A, L, U, VS, I1, IN] = OM.PARAMS_LIN_CONSTRAINT(...)
[A, L, U, VS, I1, IN, TR] = OM.PARAMS_LIN_CONSTRAINT(NAME ...)

With no input parameters, it assembles and returns the parameters
for the aggregate linear constraints from all linear constraint sets
added using ADD_LIN_CONSTRAINT. The values of these parameters are
cached for subsequent calls. The parameters are A, L and U where the
linear constraint is of the form
    L <= A * x <= U

If a NAME is provided then it simply returns the parameters for the
corresponding named set. Likewise for indexed named sets specified
by NAME and IDX_LIST. An optional 4th output argument VS indicates the
variable sets used by this constraint set. The size of A will be
consistent with VS. Optional 5th and 6th output arguments I1 and IN
indicate the starting and ending row indices of the corresponding
constraint set in the full aggregate constraint matrix. Finally, TR
will be true if it was the transpose of the A matrix that was
supplied/stored for this set (NAME or NAME/IDX_LIST must be supplied).

Examples:
    [A, l, u] = om.params_lin_constraint();
    [A, l, u, vs, i1, iN] = om.params_lin_constraint('Pmis');

See also opt_model, add_lin_constraint().

add_quad_cost(om, name, idx, Q, c, k, varsets)

add_quad_cost() - Adds a set of user costs to the model.

OM.ADD_QUAD_COST(NAME, Q, C);
OM.ADD_QUAD_COST(NAME, Q, C, K);
OM.ADD_QUAD_COST(NAME, Q, C, K, VARSETS);
OM.ADD_QUAD_COST(NAME, IDX_LIST, Q, C);
OM.ADD_QUAD_COST(NAME, IDX_LIST, Q, C, K);
OM.ADD_QUAD_COST(NAME, IDX_LIST, Q, C, K, VARSETS);

Adds a named block of quadratic costs to the model. Costs are of the
form
    F(X) = 1/2 * X'*Q*X + C'*X + K
where Q is an NX x NX matrix (possibly sparse), C is an NX x 1 vector,
K is a scalar and NX is the number of elements in X. Here X is the vector
formed by combining the specified VARSETS (the full optimization vector
by default). Alternatively, if Q is an NX x 1 vector or empty, then F(X)
is also NX x 1, and K can be either NX + 1 or scalar.
    F(X) = 1/2 * Q .* X.^2 + C .* X + K

Indexed Named Sets
    A cost set can be identified by a single NAME, as described
    above, such as 'PgCost', or by a name that is indexed by one
    or more indices, such as 'PgCost(3,4)'. For an indexed named
    set, before adding the cost sets themselves, the dimensions
    of the indexed set must be set by calling INIT_INDEXED_NAME.

    The constraints are then added using the following, where
    all arguments are as described above, except IDX_LIST is a cell
    array of the indices for the particular cost set being added.

    OM.ADD_QUAD_COST(NAME, IDX_LIST, Q, C, K);
    OM.ADD_QUAD_COST(NAME, IDX_LIST, Q, C, K, VARSETS);

Examples:
    om.add_quad_cost('quad_cost1', Q1, c1, 0);
    om.add_quad_cost('lin_cost2',  [], c2, k2, {'Vm', 'Pg', 'z'});

    om.init_indexed_name('c', {2, 3});
    for i = 1:2
      for j = 1:3
        om.add_quad_cost('c', {i, j}, Q{i,j}, ...);
      end
    end

See also opt_model, params_quad_cost(), eval_quad_cost().

display_soln(om, set_type, name, idx)

display_soln() - Display solution values.

OM.DISPLAY_SOLN()
OM.DISPLAY_SOLN(SET_TYPE)
OM.DISPLAY_SOLN(SET_TYPE, NAME)
OM.DISPLAY_SOLN(SET_TYPE, NAME, IDX)

Displays the model solution, including values, bounds and shadow
prices for variables and linear constraints, values and shadow
prices for nonlinear constraints, and individual cost components.

Results are displayed for each SET_TYPE or specified SET_TYPE and
for each named/indexed set or a specified NAME/IDX.

Inputs:
    SET_TYPE - one of the following, specifying the type of set:
        'var' - variables
        'lin' - linear constraints
        'nle' - nonlinear equality constraints
        'nli' - nonlinear inequality constraints
        'nlc' - nonlinear costs
        'qdc' - quadratic costs
      or
        a cell array of one or more of the above
      or
        '' or 'all' - indicating to display all
    NAME - (optional) char array specifying the name of the set
    IDX  - (optional) cell array specifying the indices of the set

Examples:
    om.display_soln('var');
    om.display_soln({'nle', 'nli'});
    om.display_soln('var', 'P');
    om.display_soln('lin', 'lin_con_1');
    om.display_soln('nle', 'nle_con_b', {2,3});

See also get_soln(), parse_soln().

varsets_x(om, x, vs, want_vector)

varsets_x() - Returns a cell array of sub-vectors of X specified by VARSETS

X = OM.VARSETS_X(X, VARSETS)
X = OM.VARSETS_X(X, VARSETS, 'vector')

Returns a cell array of sub-vectors of X specified by VARSETS, or
the full optimization vector X, if VARSETS is empty.

If a 3rd argument is present (value is ignored) the returned value is
a single numeric vector with the individual components stacked vertically.

See also varsets_len().

eval_lin_constraint(om, x, name, idx)

eval_lin_constraint() - Builds and returns full set of linear constraints.

AX_U = OM.EVAL_LIN_CONSTRAINT(X)
AX_U = OM.EVAL_LIN_CONSTRAINT(X, NAME)
AX_U = OM.EVAL_LIN_CONSTRAINT(X, NAME, IDX_LIST)
[AX_U, L_AX] = OM.EVAL_LIN_CONSTRAINT(...)
[AX_U, L_AX, A] = OM.EVAL_LIN_CONSTRAINT(...)
Builds the linear constraints for the full set of constraints or an
individual named subset for a given value of the vector X, based on
constraints added by ADD_LIN_CONSTRAINT.

    l <= A * x <= u

Returns A*X - U, and optionally L - A*X and A.

Example:
    [Ax_u, l_Ax, A] = om.eval_lin_constraint(x)

See also opt_model, add_lin_constraint(), params_lin_constraint().

get_soln(om, set_type, tags, name, idx)

get_soln() - Fetch solution values for specific named/indexed sets.

VALS = OM.GET_SOLN(SET_TYPE, NAME)
VALS = OM.GET_SOLN(SET_TYPE, NAME, IDX)
VALS = OM.GET_SOLN(SET_TYPE, TAGS, NAME)
VALS = OM.GET_SOLN(SET_TYPE, TAGS, NAME, IDX)

Returns named/indexed results for a solved model, evaluated at
the solution found.

Inputs:
    SET_TYPE - one of the following, specifying the type of set:
        'var' - variables
        'lin' - linear constraints
        'nle' - nonlinear equality constraints
        'nli' - nonlinear inequality constraints
        'nlc' - nonlinear costs
        'qdc' - quadratic costs
    TAGS - char array or cell array of char arrays specifying the
        desired output(s). Valid tags vary by SET_TYPE as follows:
        'var' - default is {'x', 'mu_l', 'mu_u'}
            'x' - value of solution variable
            'mu_l' - shadow price on variable lower bound
            'mu_u' - shadow price on variable upper bound
        'lin' - default is {'g', 'mu_l', 'mu_u'}
            'g' - 1 x 2 cell array of upper and lower constraint
                values, {A*x - u, l - A*x}
            'Ax_u' - upper constraint value, A*x - u
            'l_Ax' - lower constraint value, l - A*x
            'mu_l' - shadow price on constraint lower bound
            'mu_u' - shadow price on constraint upper bound
        'nle' - default is {'g', 'lam', 'dg'}
            'g' - constraint value g(x)
            'lam' - shadow price on constraint
            'dg' - Jacobian of constraint
        'nli' - default is {'h', 'mu', 'dh'}
            'h' - constraint value h(x)
            'mu' - shadow price on constraint
            'dh' - Jacobian of constraint
        'nlc' and 'qdc' - default is {'f', 'df', 'd2f'}
            'f' - cost function value f(x) (for 'qdc' can return a vector)
            'df' - gradient of cost function
            'd2f' - Hessian of cost function
    NAME - char array specifying the name of the set
    IDX  - cell array specifying the indices of the set

Outputs:
    Variable number of outputs corresponding to TAGS input. If TAGS
    is empty or not specified, the calling context will define the
    number of outputs, returned in order of default tags for the
    specified SET_TYPE.

Examples:
    [P, muPmin, muPmax] = om.get_soln('var', 'P');
    [mu_u, mu_l] = om.get_soln('lin', {'mu_u', 'mu_l'}, 'lin_con_1');
    dg_b_2_3 = om.get_soln('nle', 'dg', 'nle_con_b', {2,3});

For a complete set of solution vector values and shadow prices, using
the PARSE_SOLN method may be more efficient.

See also parse_soln().

add_lin_constraint(om, name, idx, A, l, u, varsets, tr)

add_lin_constraint() - Adds a set of linear constraints to the model.

OM.ADD_LIN_CONSTRAINT(NAME, A, L, U);
OM.ADD_LIN_CONSTRAINT(NAME, A, L, U, VARSETS);
OM.ADD_LIN_CONSTRAINT(NAME, A, L, U, VARSETS, TR);

OM.ADD_LIN_CONSTRAINT(NAME, IDX_LIST, A, L, U);
OM.ADD_LIN_CONSTRAINT(NAME, IDX_LIST, A, L, U, VARSETS);
OM.ADD_LIN_CONSTRAINT(NAME, IDX_LIST, A, L, U, VARSETS, TR);

Linear constraints are of the form L <= A * x <= U, where x is a
vector made of the vars specified in VARSETS (in the order given).
This allows the A matrix to be defined only in terms of the relevant
variables without the need to manually create a lot of zero columns.
If VARSETS is empty, x is taken to be the full vector of all
optimization variables. If L or U are empty, they are assumed to be
appropriately sized vectors of -Inf and Inf, respectively. If TR is
present and true, it means that A' is supplied/stored rather than A.
In some contexts this can be used to save significant memory.

Indexed Named Sets
    A constraint set can be identified by a single NAME, as described
    above, such as 'Pmismatch', or by a name that is indexed by one
    or more indices, such as 'Pmismatch(3,4)'. For an indexed named
    set, before adding the constraint sets themselves, the dimensions
    of the indexed set must be set by calling INIT_INDEXED_NAME.

    The constraints are then added using the following, where
    all arguments are as described above, except IDX_LIST is a cell
    array of the indices for the particular constraint set being added.

    OM.ADD_LIN_CONSTRAINT(NAME, IDX_LIST, A, L, U);
    OM.ADD_LIN_CONSTRAINT(NAME, IDX_LIST, A, L, U, VARSETS);

Examples:
    %% linear constraint
    om.add_lin_constraint('vl', Avl, lvl, uvl, {'Pg', 'Qg'});

    %% linear constraints with indexed named set 'R(i,j)'
    om.init_indexed_name('lin', 'R', {2, 3});
    for i = 1:2
      for j = 1:3
        om.add_lin_constraint('R', {i, j}, A{i,j}, ...);
      end
    end

See also opt_model, params_lin_constraint().

varsets_idx(om, vs)

varsets_idx() - Returns a vector of indices into X specified by VARSETS

K = OM.VARSETS_IDX(VARSETS)

Returns a vector of indices into X corresponding to the variable
sets specified by VARSETS.

See also varsets_x().

varsets_len(om, vs)

varsets_len() - Returns the total number of variables in VARSETS

NV = OM.VARSETS_LEN(VARSETS)

Returns the total number of elements in the optimization sub-vector
specified by VARSETS.

See also varsets_cell2struct().

varsets_cell2struct(om, vs)

varsets_cell2struct() - Converts VARSETS from cell array to struct array.

VARSETS = OM.VARSETS_CELL2STRUCT(VARSETS)

Converts VARSETS from a cell array to a struct array, if necessary.

See also varsets_len().

add_named_set(om, set_type, name, idx, N, varargin)

add_named_set() - Adds a named set of variables/constraints/costs to the model.

-----  PRIVATE METHOD  -----

This method is intended to be a private method, used internally by
the public methods ADD_VAR, ADD_LIN_CONSTRAINT, ADD_NLN_CONSTRAINT
ADD_QUAD_COST and ADD_NLN_COST.

Variable Set
    OM.ADD_NAMED_SET('var', NAME, IDX_LIST, N, V0, VL, VU, VT);

Linear Constraint Set
    OM.ADD_NAMED_SET('lin', NAME, IDX_LIST, N, A, L, U, VARSETS, TR);

Nonlinear Inequality Constraint Set
    OM.ADD_NAMED_SET('nle', NAME, IDX_LIST, N, FCN, HESS, COMPUTED_BY, VARSETS);

Nonlinear Inequality Constraint Set
    OM.ADD_NAMED_SET('nli', NAME, IDX_LIST, N, FCN, HESS, COMPUTED_BY, VARSETS);

Quadratic Cost Set
    OM.ADD_NAMED_SET('qdc', NAME, IDX_LIST, N, CP, VARSETS);

See also opt_model, add_var(), add_lin_constraint(), add_nln_constraint() add_quad_cost(), add_nln_cost().

params_var(om, name, idx)

params_var() - Returns initial value, lower bound and upper bound for opt variables.

[V0, VL, VU] = OM.PARAMS_VAR()
[V0, VL, VU] = OM.PARAMS_VAR(NAME)
[V0, VL, VU] = OM.PARAMS_VAR(NAME, IDX_LIST)
[V0, VL, VU, VT] = PARAMS_VAR(...)

Returns the initial value V0, lower bound VL and upper bound VU for
the full set of optimization variables, or for a specific named or named
and indexed variable set. Values for the full set are cached for
subsequent calls. Optionally also returns a corresponding char
vector VT of variable types, where 'C', 'I' and 'B' represent continuous
integer and binary variables, respectively.

Examples:
    [x0, xmin, xmax] = om.params_var();
    [Pg0, Pmin, Pmax] = om.params_var('Pg');
    [zij0, zijmin, zijmax, ztype] = om.params_var('z', {i, j});

See also opt_model, add_var().

add_var(om, name, idx, varargin)

add_var() - Adds a set of variables to the model.

OM.ADD_VAR(NAME, N, V0, VL, VU, VT)
OM.ADD_VAR(NAME, N, V0, VL, VU)
OM.ADD_VAR(NAME, N, V0, VL)
OM.ADD_VAR(NAME, N, V0)
OM.ADD_VAR(NAME, N)
OM.ADD_VAR(NAME, DIM_LIST) (deprecated, use INIT_INDEXED_NAME instead)
OM.ADD_VAR(NAME, IDX_LIST, N, V0, VL, VU, VT)
OM.ADD_VAR(NAME, IDX_LIST, N, V0, VL, VU)
OM.ADD_VAR(NAME, IDX_LIST, N, V0, VL)
OM.ADD_VAR(NAME, IDX_LIST, N, V0)
OM.ADD_VAR(NAME, IDX_LIST, N)

Adds a set of variables to the model, where N is the number of
variables in the set, V0 is the initial value of those variables,
VL and VU are the lower and upper bounds on the variables and VT
is the variable type. The accepted values for elements of VT are:
    'C' - continuous
    'I' - integer
    'B' - binary
V0, VL and VU are N x 1 column vectors, VT is a scalar or a 1 x N row
vector. The defaults for the last four arguments, which are all optional,
are for all values to be initialized to zero (V0 = 0), unbounded
(VL = -Inf, VU = Inf), and continuous (VT = 'C').

Examples:
    om.add_var('V', nb, V0, Vmin, Vmax, 'C');

    om.init_indexed_name('x', {2, 3});
    for i = 1:2
      for j = 1:3
        om.add_var('x', {i, j}, nx(i,j), ...);
      end
    end

See also opt_model, params_var().

solve(om, opt)

solve() - Solve the optimization model.

X = OM.SOLVE()
[X, F] = OM.SOLVE()
[X, F, EXITFLAG] = OM.SOLVE()
[X, F, EXITFLAG, OUTPUT] = OM.SOLVE()
[X, F, EXITFLAG, OUTPUT, JAC] = OM.SOLVE()      (LEQ/NLEQ problems)
[X, F, EXITFLAG, OUTPUT, LAMBDA] = OM.SOLVE()   (other problem types)
[X ...] = OM.SOLVE(OPT)

Solves the optimization model using one of the following, depending
on the problem type: QPS_MASTER, MIQPS_MASTER, NLPS_MASTER, NLEQS_MASTER.

Inputs:
    OPT : optional options structure with the following fields,
        all of which are also optional (default values shown in
        parentheses)
        alg ('DEFAULT') : determines which solver to use, list of relevant
                problem types are listed in parens next to each
            'DEFAULT' : automatic, depending on problem type, uses the
                    the first available of:
                LP - Gurobi, CPLEX, MOSEK, linprog (if MATLAB), GLPK,
                        BPMPD, MIPS
                QP - Gurobi, CPLEX, MOSEK, quadprog (if MATLAB), BPMPD,
                        MIPS
                MILP - Gurobi, CPLEX, MOSEK, Opt Tbx (intlingprog), GLPK
                MIQP - Gurobi, CPLEX, MOSEK
                NLP - MIPS
                MINLP - Artelys Knitro (not yet implemented)
                LEQ - built-in backslash operator
                NLEQ - Newton's method
            'BPMPD'   : (LP, QP) BPMPD_MEX
            'CLP'     : (LP, QP) CLP
            'CPLEX'   : (LP, QP, MILP, MIQP) CPLEX
            'FD'      : (NLEQ) fast-decoupled Newon's method
            'FMINCON' : (NLP) FMINCON, MATLAB Optimization Toolbox
            'FSOLVE'  : (NLEQ) FSOLVE, MATLAB Optimization Toolbox
            'GLPK'    : (LP, MILP) GLPK
            'GS'      : (NLEQ) Gauss-Seidel
            'GUROBI'  : (LP, QP, MILP, MIQP) Gurobi
            'IPOPT'   : (LP, QP, NLP) IPOPT, requires MEX interface to IPOPT solver
                        https://github.com/coin-or/Ipopt
            'KNITRO'  : (NLP, MINLP) Artelys Knitro, requires Artelys Knitro solver
                        https://www.artelys.com/solvers/knitro/
            'MIPS'    : (LP, QP, NLP) MIPS, MATPOWER Interior Point Solver
                     pure MATLAB implementation of a primal-dual
                     interior point method, if mips_opt.step_control = 1
                     it uses MIPS-sc, a step controlled variant of MIPS
            'MOSEK'   : (LP, QP, MILP, MIQP) MOSEK
            'NEWTON'  : (NLEQ) Newton's method
            'OSQP'    : (LP, QP) OSQP, https://osqp.org
            'OT'      : (LP, QP, MILP) MATLAB Optimization Toolbox,
                        LINPROG, QUADPROG or INTLINPROG
        verbose (0) - controls level of progress output displayed
            0 = no progress output
            1 = some progress output
            2 = verbose progress output
        bp_opt      - options vector for BP (BPMPD)
        clp_opt     - options vector for CLP
        cplex_opt   - options struct for CPLEX
        fd_opt      - options struct for fast-decoupled Newton's method,
                        nleqs_fd_newton()
        fmincon_opt - options struct for FMINCON
        fsolve_opt  - options struct for FSOLVE
        glpk_opt    - options struct for GLPK
        grb_opt     - options struct for GUROBI
        gs_opt      - options struct for Gauss-Seidel method,
                        nleqs_gauss_seidel()
        intlinprog_opt - options struct for INTLINPROG
        ipopt_opt   - options struct for IPOPT
        knitro_opt  - options struct for Artelys Knitro
        leq_opt     - options struct for MPLINSOLVE, with optional fields
            'solver' and 'opt' corresponding to respective MPLINSOLVE args,
            and 'thresh' specifying a threshold on the absolute value of
            any element X, above which EXITFLAG will be set to 0
        linprog_opt - options struct for LINPROG
        mips_opt    - options struct for MIPS
        mosek_opt   - options struct for MOSEK
        newton_opt  - options struct for Newton method, NLEQS_NEWTON
        osqp_opt    - options struct for OSQP
        quadprog_opt - options struct for QUADPROG
        parse_soln (0) - flag that specifies whether or not to call
            the PARSE_SOLN method and place the return values in OM.soln.
        price_stage_warn_tol (1e-7) - tolerance on the objective fcn
            value and primal variable relative match required to avoid
            mis-match warning message if mixed integer price computation
            stage is not skipped
        skip_prices (0) - flag that specifies whether or not to skip the
            price computation stage for mixed integer problems, in which
            the problem is re-solved for only the continuous variables,
            with all others being constrained to their solved values
        x0 (empty)  - optional initial value of x, overrides value
            stored in model (ignored by some solvers)

Outputs:
    X : solution vector
    F : final (objective) function value
    EXITFLAG : exit flag
        1 = converged
        0 or negative values = solver specific failure codes
    OUTPUT : output struct with the following fields:
        alg - algorithm code of solver used
        et  - elapsed time (sec)
        (others) - solver specific fields
    JAC : final Jacobian matrix (if available, for LEQ/NLEQ problems)
    LAMBDA : (for all non-NLEQ problem types) struct containing the
        Langrange and Kuhn-Tucker multipliers on the constraints, with
        fields:
        eqnonlin - nonlinear equality constraints
        ineqnonlin - nonlinear inequality constraints
        mu_l - lower (left-hand) limit on linear constraints
        mu_u - upper (right-hand) limit on linear constraints
        lower - lower bound on optimization variables
        upper - upper bound on optimization variables

See also opt_model, qps_master(), miqps_master(), nlps_master(), nleqs_master(), pnes_master(), mp_linsolve.

add_nln_cost(om, name, idx, N, fcn, varsets)

add_nln_cost() - Adds a set of general nonlinear costs to the model.

OM.ADD_NLN_COST(NAME, N, FCN);
OM.ADD_NLN_COST(NAME, N, FCN, VARSETS);
OM.ADD_NLN_COST(NAME, IDX_LIST, N, FCN);
OM.ADD_NLN_COST(NAME, IDX_LIST, N, FCN, VARSETS);

Adds a named block of general nonlinear costs to the model. FCN is
a handle to function that evaluates the cost, its gradient and Hessian
as described below.

The N parameter specifies the dimension for vector valued cost
functions, which are not yet implemented. Currently N must equal 1
or it will throw an error.

For a cost function f(x), FCN should point to a function with the
following interface:
    F = FCN(X)
    [F, DF] = FCN(X)
    [F, DF, D2F] = FCN(X)
where F is a scalar with the value of the function, DF is the NX x 1
gradient, and D2F is the NX x NX Hessian and NX is the number of
elements in X.

The first input argument X can take two forms. If the constraint set
is added with VARSETS empty or missing, then X will be the full
optimization vector. Otherwise it will be a cell array of vectors
corresponding to the variable sets specified in VARSETS.

Indexed Named Sets
    A cost set can be identified by a single NAME, as described
    above, such as 'PgCost', or by a name that is indexed by one
    or more indices, such as 'PgCost(3,4)'. For an indexed named
    set, before adding the cost sets themselves, the dimensions
    of the indexed set must be set by calling INIT_INDEXED_NAME.

    The constraints are then added using the following, where
    all arguments are as described above, except IDX_LIST is a cell
    array of the indices for the particular cost set being added.

    OM.ADD_NLN_COST(NAME, IDX_LIST, FCN);
    OM.ADD_NLN_COST(NAME, IDX_LIST, FCN, VARSETS);

Examples:

    fcn1 = @(x)my_cost_function1(x, other_args)
    fcn2 = @(x)my_cost_function2(x, other_args)
    om.add_nln_cost('mycost1', 1, fcn1);
    om.add_nln_cost('mycost2', 1, fcn2, {'Vm', 'Pg', 'Qg', 'z'});

    om.init_indexed_name('c', {2, 3});
    for i = 1:2
      for j = 1:3
        om.add_nln_cost('c', {i, j}, 1, fcn(i,j), ...);
      end
    end

See also opt_model, eval_nln_cost().

params_nln_constraint(om, iseq, name, idx)

params_nln_constraint() - Returns parameters for general nonlinear constraints.

N = OM.PARAMS_NLN_CONSTRAINT(ISEQ, NAME)
N = OM.PARAMS_NLN_CONSTRAINT(ISEQ, NAME, IDX_LIST)
[N, FCN] = OM.PARAMS_NLN_CONSTRAINT(...)
[N, FCN, HESS] = OM.PARAMS_NLN_CONSTRAINT(...)
[N, FCN, HESS, VS] = OM.PARAMS_NLN_CONSTRAINT(...)
[N, FCN, HESS, VS, INCLUDE] = OM.PARAMS_NLN_CONSTRAINT(...)

Returns the parameters N, and optionally FCN, and HESS provided when
the corresponding named nonlinear constraint set was added to the
model. Likewise for indexed named sets specified by NAME and IDX_LIST.
The ISEQ input should be set to 1 for equality constrainst and to 0
for inequality constraints.

An optional 4th output argument VS indicates the variable sets used by
this constraint set.
problem_type(om, recheck)

problem_type() - Return a string identifying the type of mathematical program

PROB_TYPE = OM.PROBLEM_TYPE()
PROB_TYPE = OM.PROBLEM_TYPE(RECHECK)

Returns a string identifying the type of mathematical program
represented by the current model, based on the variables, costs,
and constraints that have been added to the model. Used to
automatically select an appropriate solver.

Linear and nonlinear equations are models with no costs, no inequality
constraints, and an equal number of continuous variables and equality
constraints. If the number of variables in a nonlinear equation model
is one more than the number of constraints, it is a parameterized
nonlinear equation.

Outputs:
    PROB_TYPE : problem type, one of the following strings:
        'LEQ'   - linear equations
        'NLEQ'  - nonlinear equations
        'PNE'   - parameterized nonlinear equations
        'LP'    - linear program
        'QP'    - quadratic program
        'NLP'   - nonlinear program
        'MILP'  - mixed-integer linear program
        'MIQP'  - mixed-integer quadratic program
        'MINLP' - mixed-integer nonlinear program

The output value is cached for future calls, but calling with a true
value for the optional RECHECK argument will force it to recheck in
case the problem type has changed due to modifying the variables,
constraints or costs in the model.

See also opt_model.

eval_nln_constraint(om, x, iseq, name, idx)

eval_nln_constraint() - Builds and returns full set of nonlinear constraints.

G = OM.EVAL_NLN_CONSTRAINT(X, ISEQ)
G = OM.EVAL_NLN_CONSTRAINT(X, ISEQ, NAME)
G = OM.EVAL_NLN_CONSTRAINT(X, ISEQ, NAME, IDX_LIST)
[G, DG] = OM.EVAL_NLN_CONSTRAINT(...)
Builds the nonlinear equality or inequality constraints (ISEQ equal to
1 or 0, respectively) and optionally their gradients for the full set
of constraints or an individual named subset for a given value of the
vector X, based on constraints added by ADD_NLN_CONSTRAINT.

    g(X) = 0
    h(X) <= 0

Example:
    [g, dg] = om.eval_nln_constraint(x, 1)
    [h, dh] = om.eval_nln_constraint(x, 0)

See also opt_model, add_nln_constraint(), eval_nln_constraint_hess().

init_indexed_name(om, set_type, name, dim_list)

init_indexed_name() - Initializes the dimensions for an indexed named set.

OM.INIT_INDEXED_NAME(SET_TYPE, NAME, DIM_LIST)

Initializes the dimensions for an indexed named variable, constraint
or cost set.

Variables, constraints and costs are referenced in OPT_MODEL in terms
of named sets. The specific type of named set being referenced is
given by SET_TYPE, with the following valid options:
    SET_TYPE = 'var'   => variable set
    SET_TYPE = 'lin'   => linear constraint set
    SET_TYPE = 'nle'   => nonlinear equality constraint set
    SET_TYPE = 'nli'   => nonlinear inequality constraint set
    SET_TYPE = 'qdc'   => quadratic cost set
    SET_TYPE = 'nlc'   => nonlinear cost set

Indexed Named Sets

A variable, constraint or cost set can be identified by a single NAME,
such as 'Pmismatch', or by a name that is indexed by one or more indices,
such as 'Pmismatch(3,4)'. For an indexed named set, before adding the
indexed variable, constraint or cost sets themselves, the dimensions of
the indexed set must be set by calling INIT_INDEXED_NAME, where
DIM_LIST is a cell array of the dimensions.

Examples:
    %% linear constraints with indexed named set 'R(i,j)'
    om.init_indexed_name('lin', 'R', {2, 3});
    for i = 1:2
      for j = 1:3
        om.add_lin_constraint('R', {i, j}, A{i,j}, ...);
      end
    end

See also opt_model, add_var(), add_lin_constraint(), add_nln_constraint(), add_quad_cost(), add_nln_cost().

add_nln_constraint(om, name, idx, N, iseq, fcn, hess, varsets)

add_nln_constraint() - Adds a set of nonlinear constraints to the model.

OM.ADD_NLN_CONSTRAINT(NAME, N, ISEQ, FCN, HESS);
OM.ADD_NLN_CONSTRAINT(NAME, N, ISEQ, FCN, HESS, VARSETS);

OM.ADD_NLN_CONSTRAINT(NAME, IDX_LIST, N, ISEQ, FCN, HESS);
OM.ADD_NLN_CONSTRAINT(NAME, IDX_LIST, N, ISEQ, FCN, HESS, VARSETS);

N specifies the number of constraints in the set, ISEQ is a 1 for
an equality constraint set, 0 for inequality, FCN is the handle of a
function that evaluates the constraint and its gradients, and HESS is
the handle of a function that evaluates the Hessian of the constraints.

For a constraint G(x) = 0, FCN should point to a function with the
following interface:
    G = FCN(X)
    [G, DG] = FCN(X)
where G is an N x 1 vector and DG is the N x NX Jacobian, where
    DG(i, j) = dG(i)/dX(j)
and NX is the number of elements in X. Note: DG is the transpose of
what is expected from an input function for NLPS_MASTER and friends.

HESS should point to a function that returns an NX x NX matrix of
derivatives of DG * LAMBDA, with the following interface:
    D2G = HESS(X, LAMBDA)

For both functions, the first input argument X can take two forms. If
the constraint set is added with VARSETS empty or missing, then X will
be the full optimization vector. Otherwise it will be a cell array of
vectors corresponding to the variable sets specified in VARSETS.

For simple (not indexed) named sets, NAME can be a cell array of
constraint set names, in which case N is a vector, specifying the number
of constraints in each corresponding set. FCN and HESS are each still
a single function handle, but the values computed by each correspond
to the entire stacked collection of constraint sets together, as if
they were a single set.

Likewise, if FCN or HESS are empty, it also indicates a placeholder in
the indexing for a constraint set whose implementation is included in
another constraint set. This functionality is only intended to be used
internally to handle constraint/gradient and Hessian functions that
compute the values for more than one constraint set simultaneously.

Indexed Named Sets
    A constraint set can be identified by a single NAME, as described
    above, such as 'Pmismatch', or by a name that is indexed by one
    or more indices, such as 'Pmismatch(3,4)'. For an indexed named
    set, before adding the constraint sets themselves, the dimensions
    of the indexed set must be set by calling INIT_INDEXED_NAME.

    The constraints are then added using the following, where
    all arguments are as described above, except IDX_LIST is a cell
    array of the indices for the particular constraint set being added.

    OM.ADD_NLN_CONSTRAINT(NAME, IDX_LIST, N, ISEQ, FCN, HESS);
    OM.ADD_NLN_CONSTRAINT(NAME, IDX_LIST, N, ISEQ, FCN, HESS, VARSETS);

Examples:
    %% nonlinear equality constraint with constraint/gradient and Hessian
    %% evaluation functions provided
    om.add_nln_constraint('Qmis', nb, 1, fcn, hess);

    %% nonlinear inequality constraints with indexed named set 'S(i,j)'
    om.init_indexed_name('nli', 'S', {2, 3});
    for i = 1:2
      for j = 1:3
        om.add_nln_constraint('S', {i, j}, N{i,j}, ...);
      end
    end

See also opt_model, eval_nln_constraint().

get_idx(om, varargin)

get_idx() - Returns the idx struct for vars, lin/nonlin constraints, costs.

VV = OM.GET_IDX()
[VV, LL] = OM.GET_IDX()
[VV, LL, NNE] = OM.GET_IDX()
[VV, LL, NNE, NNI] = OM.GET_IDX()
[VV, LL, NNE, NNI, QQ] = OM.GET_IDX()
[VV, LL, NNE, NNI, QQ, NNC] = OM.GET_IDX()

Returns a structure for each with the beginning and ending
index value and the number of elements for each named block.
The 'i1' field (that's a one) is a struct with all of the
starting indices, 'iN' contains all the ending indices and
'N' contains all the sizes. Each is a struct whose fields are
the named blocks.

Alternatively, you can specify the type of named set(s) directly
as inputs ...

[IDX1, IDX2, ...] = OM.GET_IDX(SET_TYPE1, SET_TYPE2, ...);
VV = OM.GET_IDX('var');
[LL, NNE, NNI] = OM.GET_IDX('lin', 'nle', 'nli');

The specific type of named set being referenced is
given by the SET_TYPE inputs, with the following valid options:
    SET_TYPE = 'var'   => variable set
    SET_TYPE = 'lin'   => linear constraint set
    SET_TYPE = 'nle'   => nonlinear equality constraint set
    SET_TYPE = 'nli'   => nonlinear inequality constraint set
    SET_TYPE = 'qdc'   => quadratic cost set
    SET_TYPE = 'nnc'   => nonlinear cost set

Examples:
    [vv, ll, nne] = om.get_idx();
    [vv, ll, qq] = om.get_idx('var', 'lin', 'qdc');

    For a variable block named 'z' we have ...
        vv.i1.z - starting index for 'z' in optimization vector x
        vv.iN.z - ending index for 'z' in optimization vector x
        vv.N.z  - number of elements in 'z'

    To extract a 'z' variable from x:
        z = x(vv.i1.z:vv.iN.z);

    To extract the multipliers on a linear constraint set
    named 'foo', where mu_l and mu_u are the full set of
    linear constraint multipliers:
        mu_l_foo = mu_l(ll.i1.foo:ll.iN.foo);
        mu_u_foo = mu_u(ll.i1.foo:ll.iN.foo);

    The number of nonlinear equality constraints in a set named 'bar':
        nbar = nne.N.bar;
      (note: the following is preferable ...
        nbar = om.getN('nle', 'bar');
      ... if you haven't already called get_idx to get nne.)

    If 'z', 'foo' and 'bar' are indexed sets, then you can
    replace them with something like 'z(i,j)', 'foo(i,j,k)'
    or 'bar(i)' in the examples above.

See also opt_model, add_var(), add_lin_constraint(), add_nln_constraint(), add_quad_cost(), add_nln_cost().

eval_quad_cost(om, x, name, idx)

eval_quad_cost() - Evaluates individual or full set of quadratic costs.

F = OM.EVAL_QUAD_COST(X ...)
[F, DF] = OM.EVAL_QUAD_COST(X ...)
[F, DF, D2F] = OM.EVAL_QUAD_COST(X ...)
[F, DF, D2F] = OM.EVAL_QUAD_COST(X, NAME)
[F, DF, D2F] = OM.EVAL_QUAD_COST(X, NAME, IDX_LIST)
Evaluates the cost function and its derivatives for an individual named
set or the full set of quadratic costs for a given value of the
optimization vector X, based on costs added by ADD_QUAD_COST.

Example:
    [f, df, d2f] = om.eval_quad_cost(x)
    [f, df, d2f] = om.eval_quad_cost(x, name)
    [f, df, d2f] = om.eval_quad_cost(x, name, idx_list)

See also opt_model, add_quad_cost(), params_quad_cost().

set_params(om, st, name, idx, params, vals)

set_params() - Modifies parameters for variable, cost or constraint in model

This method can be used to modify parameters for an existing variable,
constraint or cost in the model.

OM.SET_PARAMS(SET_TYPE, NAME, PARAMS, VALS)
OM.SET_PARAMS(SET_TYPE, NAME, IDX, PARAMS, VALS)

Inputs:
    SET_TYPE : one of 'var', 'lin', 'nle', 'nli', 'nlc', 'qdc' for
        variables, linear constraints, nonlinear equality constraints,
        nonlinear inequality constraints, general nonlinear costs,
        and quadratic costs, respectively
    NAME : name of set
    IDX : index of named set (for an indexed set)
    PARAMS : can be one of three options:
        1 - 'all', indicating that VALS is a cell array whose elements
            correspond to the input parameters of the respective
            add_*() method
        2 - the name of a PARAM, VAL is the value of that parameter
        3 - a cell array of PARAM names, VALS is a cell array of
            corresponding values
        Note: Changing the dimension of a 'var' is not allowed
            and changing the #1 ('all') is the only option for 'nle', 'nli', and 'nlc'
    VALS : new value or cell array of new values for PARAMS

Valid PARAM names:
    var - N, v0, vl, vu, vt
    lin - A, l, u, vs
    nle - N, fcn, hess, include, vs
    nli - N, fcn, hess, include, vs
    nlc - N, fcn, vs
    qdc - Q, c, k, vs

Examples:
    om.set_params('var', 'Pg', 'v0', Pg0);
    om.set_params('lin', 'y', {2,3}, {'l', 'u'}, {l, u});
    om.set_params('nle', 'Pmis', 'all', {N, @fcn, @hess, vs});

See also opt_model, add_var(), add_lin_constraint(), add_nln_constraint(), add_quad_cost(), add_nln_cost().