Home > matpower7.0 > lib > mpoption.m

mpoption

PURPOSE ^

MPOPTION Used to set and retrieve a MATPOWER options struct.

SYNOPSIS ^

function opt = mpoption(varargin)

DESCRIPTION ^

MPOPTION  Used to set and retrieve a MATPOWER options struct.

   OPT = MPOPTION
       Returns the default options struct.

   OPT = MPOPTION(OVERRIDES)
       Returns the default options struct, with some fields overridden
       by values from OVERRIDES, which can be a struct or the name of
       a function that returns a struct.

   OPT = MPOPTION(NAME1, VALUE1, NAME2, VALUE2, ...)
       Same as previous, except override options are specified by NAME,
       VALUE pairs. This can be used to set any part of the options
       struct. The names can be individual fields or multi-level field
       names with embedded periods. The values can be scalars or structs.

       For backward compatibility, the NAMES and VALUES may correspond
       to old-style MATPOWER option names (elements in the old-style
       options vector) as well.

   OPT = MPOPTION(OPT0)
       Converts an old-style options vector OPT0 into the corresponding
       options struct. If OPT0 is an options struct it does nothing.

   OPT = MPOPTION(OPT0, OVERRIDES)
       Applies overrides to an existing set of options, OPT0, which
       can be an old-style options vector or an options struct.

   OPT = MPOPTION(OPT0, NAME1, VALUE1, NAME2, VALUE2, ...)
       Same as above except it uses the old-style options vector OPT0
       as a base instead of the old default options vector.

   OPT_VECTOR = MPOPTION(OPT, [])
       Creates and returns an old-style options vector from an
       options struct OPT.

   Note: The use of old-style MATPOWER options vectors and their
         names and values has been deprecated and will be removed
         in a future version of MATPOWER. Until then, all uppercase
         option names are not permitted for new top-level options.

   Examples:
       mpopt = mpoption('pf.alg', 'FDXB', 'pf.tol', 1e-4);
       mpopt = mpoption(mpopt, 'opf.dc.solver', 'CPLEX', 'verbose', 2);

The currently defined options are as follows:

   name                    default     description [options]
----------------------    ---------   ----------------------------------
Model options:
   model                   'AC'        AC vs. DC power flow model
       [ 'AC' - use nonlinear AC model & corresponding algorithms/options  ]
       [ 'DC' - use linear DC model & corresponding algorithms/options     ]

Power Flow options:
   pf.alg                  'NR'        AC power flow algorithm
       [ 'NR'    - Newton's method (formulation depends on values of       ]
       [           pf.current_balance and pf.v_cartesian options)          ]
       [ 'NR-SC' - Newton's method (power mismatch, cartesian)             ]
       [ 'NR-IP' - Newton's method (current mismatch, polar)               ]
       [ 'NR-IC' - Newton's method (current mismatch, cartesian)           ]
       [ 'FDXB'  - Fast-Decoupled (XB version)                             ]
       [ 'FDBX'  - Fast-Decoupled (BX version)                             ]
       [ 'GS'    - Gauss-Seidel                                            ]
       [ 'PQSUM' - Power Summation method (radial networks only)           ]
       [ 'ISUM'  - Current Summation method (radial networks only)         ]
       [ 'YSUM'  - Admittance Summation method (radial networks only)      ]
   pf.current_balance     0           type of nodal balance equation
       [ 0 - use complex power balance equations                           ]
       [ 1 - use complex current balance equations                         ]
   pf.v_cartesian         0           voltage representation
       [ 0 - bus voltage variables represented in polar coordinates        ]
       [ 1 - bus voltage variables represented in Cartesian coordinates    ]
   pf.tol                  1e-8        termination tolerance on per unit
                                       P & Q mismatch
   pf.nr.max_it            10          maximum number of iterations for
                                       Newton's method
   pf.nr.lin_solver        ''          linear solver passed to MPLINSOLVE to
                                       solve Newton update step
       [ ''      - default to '\' for small systems, 'LU3' for larger ones ]
       [ '\'     - built-in backslash operator                             ]
       [ 'LU'    - explicit default LU decomposition and back substitution ]
       [ 'LU3'   - 3 output arg form of LU, Gilbert-Peierls algorithm with ]
       [           approximate minimum degree (AMD) reordering             ]
       [ 'LU4'   - 4 output arg form of LU, UMFPACK solver (same as 'LU')  ]
       [ 'LU5'   - 5 output arg form of LU, UMFPACK solver, w/row scaling  ]
       [       (see MPLINSOLVE for complete list of all options)           ]
   pf.fd.max_it            30          maximum number of iterations for
                                       fast decoupled method
   pf.gs.max_it            1000        maximum number of iterations for
                                       Gauss-Seidel method
   pf.radial.max_it        20          maximum number of iterations for
                                       radial power flow methods
   pf.radial.vcorr         0           perform voltage correction procedure
                                       in distribution power flow
       [  0 - do NOT perform voltage correction                            ]
       [  1 - perform voltage correction                                   ]
   pf.enforce_q_lims       0           enforce gen reactive power limits at
                                       expense of |V|
       [  0 - do NOT enforce limits                                        ]
       [  1 - enforce limits, simultaneous bus type conversion             ]
       [  2 - enforce limits, one-at-a-time bus type conversion            ]

Continuation Power Flow options:
   cpf.parameterization    3           choice of parameterization
       [  1 - natural                                                      ]
       [  2 - arc length                                                   ]
       [  3 - pseudo arc length                                            ]
   cpf.stop_at             'NOSE'      determins stopping criterion
       [ 'NOSE'     - stop when nose point is reached                      ]
       [ 'FULL'     - trace full nose curve                                ]
       [ <lam_stop> - stop upon reaching specified target lambda value     ]
   cpf.enforce_p_lims      0           enforce gen active power limits
       [  0 - do NOT enforce limits                                        ]
       [  1 - enforce limits, simultaneous bus type conversion             ]
   cpf.enforce_q_lims      0           enforce gen reactive power limits at
                                       expense of |V|
       [  0 - do NOT enforce limits                                        ]
       [  1 - enforce limits, simultaneous bus type conversion             ]
   cpf.enforce_v_lims      0           enforce bus voltage magnitude limits
       [  0 - do NOT enforce limits                                        ]
       [  1 - enforce limits, termination on detection                     ]
   cpf.enforce_flow_lims   0           enforce branch flow MVA limits
       [  0 - do NOT enforce limits                                        ]
       [  1 - enforce limits, termination on detection                     ]
   cpf.step                0.05        continuation power flow step size
   cpf.adapt_step          0           toggle adaptive step size feature
       [  0 - adaptive step size disabled                                  ]
       [  1 - adaptive step size enabled                                   ]
   cpf.step_min            1e-4        minimum allowed step size
   cpf.step_max            0.2         maximum allowed step size
   cpf.adapt_step_damping  0.7         damping factor for adaptive step
                                       sizing
   cpf.adapt_step_tol      1e-3        tolerance for adaptive step sizing
   cpf.target_lam_tol      1e-5        tolerance for target lambda detection
   cpf.nose_tol            1e-5        tolerance for nose point detection (pu)
   cpf.p_lims_tol          0.01        tolerance for generator active
                                       power limit enforcement (MW)
   cpf.q_lims_tol          0.01        tolerance for generator reactive
                                       power limit enforcement (MVAR)
   cpf.v_lims_tol          1e-4        tolerance for bus voltage
                                       magnitude enforcement (p.u)
   cpf.flow_lims_tol       0.01        tolerance for line MVA flow
                                       enforcement (MVA)
   cpf.plot.level          0           control plotting of noze curve
       [  0 - do not plot nose curve                                       ]
       [  1 - plot when completed                                          ]
       [  2 - plot incrementally at each iteration                         ]
       [  3 - same as 2, with 'pause' at each iteration                    ]
   cpf.plot.bus            <empty>     index of bus whose voltage is to be
                                       plotted
   cpf.user_callback       <empty>     string containing the name of a user
                                       callback function, or struct with
                                       function name, and optional priority
                                       and/or args, or cell array of such
                                       strings and/or structs, see
                                       'help cpf_default_callback' for details

Optimal Power Flow options:
   name                    default     description [options]
----------------------    ---------   ----------------------------------
   opf.ac.solver           'DEFAULT'   AC optimal power flow solver
       [ 'DEFAULT' - choose default solver, i.e. 'MIPS'                    ]
       [ 'MIPS'    - MIPS, MATPOWER Interior Point Solver, primal/dual     ]
       [             interior point method (pure MATLAB/Octave)            ]
       [ 'FMINCON' - MATLAB Optimization Toolbox, FMINCON                  ]
       [ 'IPOPT'   - IPOPT, requires MEX interface to IPOPT solver         ]
       [             available from:                                       ]
       [                   https://github.com/coin-or/Ipopt                ]
       [ 'KNITRO'  - Artelys Knitro, requires Artelys Knitro solver,       ]
       [             available from:https://www.artelys.com/solvers/knitro/]
       [ 'MINOPF'  - MINOPF, MINOS-based solver, requires optional         ]
       [             MEX-based MINOPF package, available from:             ]
       [                   http://www.pserc.cornell.edu/minopf/            ]
       [ 'PDIPM'   - PDIPM, primal/dual interior point method, requires    ]
       [             optional MEX-based TSPOPF package, available from:    ]
       [                   http://www.pserc.cornell.edu/tspopf/            ]
       [ 'SDPOPF'  - SDPOPF, solver based on semidefinite relaxation of    ]
       [             OPF problem, requires optional packages:              ]
       [               SDP_PF, available in extras/sdp_pf                  ]
       [               YALMIP, available from:                             ]
       [                   https://yalmip.github.io                        ]
       [               SDP solver such as SeDuMi, available from:          ]
       [                   http://sedumi.ie.lehigh.edu/                    ]
       [ 'TRALM'   - TRALM, trust region based augmented Langrangian       ]
       [             method, requires TSPOPF (see 'PDIPM')                 ]
   opf.dc.solver           'DEFAULT'   DC optimal power flow solver
       [ 'DEFAULT' - choose solver based on availability in the following  ]
       [             order: 'GUROBI', 'CPLEX', 'MOSEK', 'OT',              ]
       [             'GLPK' (linear costs only), 'BPMPD', 'MIPS'           ]
       [ 'MIPS'    - MIPS, MATPOWER Interior Point Solver, primal/dual     ]
       [             interior point method (pure MATLAB/Octave)            ]
       [ 'BPMPD'   - BPMPD, requires optional MEX-based BPMPD_MEX package  ]
       [             available from: http://www.pserc.cornell.edu/bpmpd/   ]
       [ 'CLP'     - CLP, requires interface to COIN-OP LP solver          ]
       [             available from:https://github.com/coin-or/Clp         ]
       [ 'CPLEX'   - CPLEX, requires CPLEX solver available from:          ]
       [             https://www.ibm.com/analytics/cplex-optimizer         ]
       [ 'GLPK'    - GLPK, requires interface to GLPK solver               ]
       [             available from: https://www.gnu.org/software/glpk/    ]
       [             (GLPK does not work with quadratic cost functions)    ]
       [ 'GUROBI'  - GUROBI, requires Gurobi optimizer (v. 5+)             ]
       [             available from: https://www.gurobi.com/               ]
       [ 'IPOPT'   - IPOPT, requires MEX interface to IPOPT solver         ]
       [             available from:                                       ]
       [                 https://github.com/coin-or/Ipopt                  ]
       [ 'MOSEK'   - MOSEK, requires MATLAB interface to MOSEK solver      ]
       [             available from: https://www.mosek.com/                ]
       [ 'OT'      - MATLAB Optimization Toolbox, QUADPROG, LINPROG        ]
   opf.current_balance     0           type of nodal balance equation
       [ 0 - use complex power balance equations                           ]
       [ 1 - use complex current balance equations                         ]
   opf.v_cartesian         0           voltage representation
       [ 0 - bus voltage variables represented in polar coordinates        ]
       [ 1 - bus voltage variables represented in Cartesian coordinates    ]
   opf.violation           5e-6        constraint violation tolerance
   opf.use_vg              0           respect gen voltage setpt     [ 0-1 ]
       [ 0 - use specified bus Vmin & Vmax, and ignore gen Vg              ]
       [ 1 - replace specified bus Vmin & Vmax by corresponding gen Vg     ]
       [ between 0 and 1 - use a weighted average of the 2 options         ]
   opf.flow_lim            'S'         quantity limited by branch flow
                                       constraints
       [ 'S' - apparent power flow (limit in MVA)                          ]
       [ 'P' - active power flow, implemented using P   (limit in MW)      ]
       [ '2' - active power flow, implemented using P^2 (limit in MW)      ]
       [ 'I' - current magnitude (limit in MVA at 1 p.u. voltage)          ]
   opf.ignore_angle_lim    0           angle diff limits for branches
       [ 0 - include angle difference limits, if specified                 ]
       [ 1 - ignore angle difference limits even if specified              ]
   opf.softlims.default    1           behavior of OPF soft limits for
                                       which parameters are not explicitly
                                       provided
       [ 0 - do not include softlims if not explicitly specified           ]
       [ 1 - include softlims w/default values if not explicitly specified ]
   opf.init_from_mpc       -1          (DEPRECATED: use opf.start instead)
                                       specify whether to use current state
                                       in MATPOWER case to initialize OPF
                                       (currently only supported by fmincon,
                                        Ipopt, Knitro and MIPS solvers,
                                        others always use mpc)
       [  -1 - MATPOWER decides, based on solver/algorithm                 ]
       [   0 - ignore current state in MATPOWER case (only applies to      ]
       [       fmincon, Ipopt, Knitro and MIPS, which use an interior pt   ]
       [       estimate; others use current state as with opf.start = 2)   ]
       [   1 - use current state in MATPOWER case                          ]
   opf.start               0           strategy for initializing OPF start pt
       [   0 - default, MATPOWER decides based on solver                   ]
       [       (currently identical to 1)                                  ]
       [   1 - ignore current state in MATPOWER case (only applies to      ]
       [       fmincon, Ipopt, Knitro and MIPS, which use an interior pt   ]
       [       estimate; others use current state as with opf.start = 2)   ]
       [   2 - use current state in MATPOWER case                          ]
       [   3 - solve power flow and use resulting state                    ]
   opf.return_raw_der      0           for AC OPF, return constraint and
                                       derivative info in results.raw
                                       (in fields g, dg, df, d2f) [ 0 or 1 ]

Output options:
   name                    default     description [options]
----------------------    ---------   ----------------------------------
   verbose                 1           amount of progress info printed
       [   0 - print no progress info                                      ]
       [   1 - print a little progress info                                ]
       [   2 - print a lot of progress info                                ]
       [   3 - print all progress info                                     ]
   out.all                 -1          controls pretty-printing of results
       [  -1 - individual flags control what prints                        ]
       [   0 - do not print anything (overrides individual flags, ignored  ]
       [       for files specified as FNAME arg to runpf(), runopf(), etc.)]
       [   1 - print everything (overrides individual flags)               ]
   out.sys_sum             1           print system summary       [ 0 or 1 ]
   out.area_sum            0           print area summaries       [ 0 or 1 ]
   out.bus                 1           print bus detail           [ 0 or 1 ]
   out.branch              1           print branch detail        [ 0 or 1 ]
   out.gen                 0           print generator detail     [ 0 or 1 ]
   out.lim.all             -1          controls constraint info output
       [  -1 - individual flags control what constraint info prints        ]
       [   0 - no constraint info (overrides individual flags)             ]
       [   1 - binding constraint info (overrides individual flags)        ]
       [   2 - all constraint info (overrides individual flags)            ]
   out.lim.v               1           control voltage limit info
       [   0 - do not print                                                ]
       [   1 - print binding constraints only                              ]
       [   2 - print all constraints                                       ]
       [   (same options for OUT_LINE_LIM, OUT_PG_LIM, OUT_QG_LIM)         ]
   out.lim.line            1           control line flow limit info
   out.lim.pg              1           control gen active power limit info
   out.lim.qg              1           control gen reactive pwr limit info
   out.force               0           print results even if success
                                       flag = 0                   [ 0 or 1 ]
   out.suppress_detail     -1          suppress all output but system summary
       [  -1 - suppress details for large systems (> 500 buses)            ]
       [   0 - do not suppress any output specified by other flags         ]
       [   1 - suppress all output except system summary section           ]
       [       (overrides individual flags, but not out.all = 1)           ]

Solver specific options:
       name                    default     description [options]
   -----------------------    ---------   ----------------------------------
   MIPS:
       mips.linsolver          ''          linear system solver
           [   '' or '\'   build-in backslash \ operator (e.g. x = A \ b)  ]
           [   'PARDISO'   PARDISO solver (if available)                   ]
       mips.feastol            0           feasibility (equality) tolerance
                                           (set to opf.violation by default)
       mips.gradtol            1e-6        gradient tolerance
       mips.comptol            1e-6        complementary condition
                                           (inequality) tolerance
       mips.costtol            1e-6        optimality tolerance
       mips.max_it             150         maximum number of iterations
       mips.step_control       0           enable step-size cntrl [ 0 or 1 ]
       mips.sc.red_it          20          maximum number of reductions per
                                           iteration with step control
       mips.xi                 0.99995     constant used in alpha updates*
       mips.sigma              0.1         centering parameter*
       mips.z0                 1           used to initialize slack variables*
       mips.alpha_min          1e-8        returns "Numerically Failed" if
                                           either alpha parameter becomes
                                           smaller than this value*
       mips.rho_min            0.95        lower bound on rho_t*
       mips.rho_max            1.05        upper bound on rho_t*
       mips.mu_threshold       1e-5        KT multipliers smaller than this
                                           value for non-binding constraints
                                           are forced to zero
       mips.max_stepsize       1e10        returns "Numerically Failed" if the
                                           2-norm of the reduced Newton step
                                           exceeds this value*
           * See the corresponding Appendix in the manual for details.

   CPLEX:
       cplex.lpmethod          0           solution algorithm for LP problems
           [   0 - automatic: let CPLEX choose                             ]
           [   1 - primal simplex                                          ]
           [   2 - dual simplex                                            ]
           [   3 - network simplex                                         ]
           [   4 - barrier                                                 ]
           [   5 - sifting                                                 ]
           [   6 - concurrent (dual, barrier, and primal)                  ]
       cplex.qpmethod          0           solution algorithm for QP problems
           [   0 - automatic: let CPLEX choose                             ]
           [   1 - primal simplex optimizer                                ]
           [   2 - dual simplex optimizer                                  ]
           [   3 - network optimizer                                       ]
           [   4 - barrier optimizer                                       ]
       cplex.opts              <empty>     see CPLEX_OPTIONS for details
       cplex.opt_fname         <empty>     see CPLEX_OPTIONS for details
       cplex.opt               0           see CPLEX_OPTIONS for details

   FMINCON:
       fmincon.alg             4           algorithm used by fmincon() for OPF
                                           for Opt Toolbox 4 and later
            [  1 - active-set (not suitable for large problems)            ]
            [  2 - interior-point, w/default 'bfgs' Hessian approx         ]
            [  3 - interior-point, w/ 'lbfgs' Hessian approx               ]
            [  4 - interior-point, w/exact user-supplied Hessian           ]
            [  5 - interior-point, w/Hessian via finite differences        ]
            [  6 - sqp (not suitable for large problems)                   ]
       fmincon.tol_x           1e-4        termination tol on x
       fmincon.tol_f           1e-4        termination tol on f
       fmincon.max_it          0           maximum number of iterations
                                                           [  0 => default ]

   GUROBI:
       gurobi.method           0           solution algorithm (Method)
           [  -1 - automatic, let Gurobi decide                            ]
           [   0 - primal simplex                                          ]
           [   1 - dual simplex                                            ]
           [   2 - barrier                                                 ]
           [   3 - concurrent (LP only)                                    ]
           [   4 - deterministic concurrent (LP only)                      ]
       gurobi.timelimit        Inf         maximum time allowed (TimeLimit)
       gurobi.threads          0           max number of threads (Threads)
       gurobi.opts             <empty>     see GUROBI_OPTIONS for details
       gurobi.opt_fname        <empty>     see GUROBI_OPTIONS for details
       gurobi.opt              0           see GUROBI_OPTIONS for details

   IPOPT:
       ipopt.opts              <empty>     see IPOPT_OPTIONS for details
       ipopt.opt_fname         <empty>     see IPOPT_OPTIONS for details
       ipopt.opt               0           see IPOPT_OPTIONS for details

   KNITRO:
       knitro.tol_x            1e-4        termination tol on x
       knitro.tol_f            1e-4        termination tol on f
       knitro.maxit            0           maximum number of iterations
                                                           [  0 => default ]
       knitro.opt_fname        <empty>     name of user-supplied native
                                           Knitro options file that overrides
                                           all other options
       knitro.opt              0           if knitro.opt_fname is empty and
                                           knitro.opt is a non-zero integer N
                                           then knitro.opt_fname is auto-
                                           generated as:
                                           'knitro_user_options_N.txt'

   LINPROG:
       linprog                 <empty>     LINPROG options passed to
                                           OPTIMOPTIONS or OPTIMSET.
                                           see LINPROG in the Optimization
                                           Toolbox for details

   MINOPF:
       minopf.feastol          0 (1e-3)    primal feasibility tolerance
                                           (set to opf.violation by default)
       minopf.rowtol           0 (1e-3)    row tolerance
       minopf.xtol             0 (1e-4)    x tolerance
       minopf.majdamp          0 (0.5)     major damping parameter
       minopf.mindamp          0 (2.0)     minor damping parameter
       minopf.penalty          0 (1.0)     penalty parameter
       minopf.major_it         0 (200)     major iterations
       minopf.minor_it         0 (2500)    minor iterations
       minopf.max_it           0 (2500)    iterations limit
       minopf.verbosity        -1          amount of progress info printed
           [  -1 - controlled by 'verbose' option                          ]
           [   0 - print nothing                                           ]
           [   1 - print only termination status message                   ]
           [   2 - print termination status and screen progress            ]
           [   3 - print screen progress, report file (usually fort.9)     ]
       minopf.core             0 (1200*nb + 2*(nb+ng)^2) memory allocation
       minopf.supbasic_lim     0 (2*nb + 2*ng) superbasics limit
       minopf.mult_price       0 (30)      multiple price

   MOSEK:
       mosek.lp_alg            0           solution algorithm
                                               (MSK_IPAR_OPTIMIZER)
           for MOSEK 8.x ...        (see MOSEK_SYMBCON for a "better way")
           [   0 - automatic: let MOSEK choose                             ]
           [   1 - dual simplex                                            ]
           [   2 - automatic: let MOSEK choose                             ]
           [   3 - automatic simplex (MOSEK chooses which simplex method)  ]
           [   4 - interior point                                          ]
           [   6 - primal simplex                                          ]
       mosek.max_it            0 (400)     interior point max iterations
                                               (MSK_IPAR_INTPNT_MAX_ITERATIONS)
       mosek.gap_tol           0 (1e-8)    interior point relative gap tol
                                               (MSK_DPAR_INTPNT_TOL_REL_GAP)
       mosek.max_time          0 (-1)      maximum time allowed
                                               (MSK_DPAR_OPTIMIZER_MAX_TIME)
       mosek.num_threads       0 (1)       max number of threads
                                               (MSK_IPAR_INTPNT_NUM_THREADS)
       mosek.opts              <empty>     see MOSEK_OPTIONS for details
       mosek.opt_fname         <empty>     see MOSEK_OPTIONS for details
       mosek.opt               0           see MOSEK_OPTIONS for details

   QUADPROG:
       quadprog                <empty>     QUADPROG options passed to
                                           OPTIMOPTIONS or OPTIMSET.
                                           see QUADPROG in the Optimization
                                           Toolbox for details

   TSPOPF:
       pdipm.feastol           0           feasibility (equality) tolerance
                                           (set to opf.violation by default)
       pdipm.gradtol           1e-6        gradient tolerance
       pdipm.comptol           1e-6        complementary condition
                                           (inequality) tolerance
       pdipm.costtol           1e-6        optimality tolerance
       pdipm.max_it            150         maximum number of iterations
       pdipm.step_control      0           enable step-size cntrl [ 0 or 1 ]
       pdipm.sc.red_it         20          maximum number of reductions per
                                           iteration with step control
       pdipm.sc.smooth_ratio   0.04        piecewise linear curve smoothing
                                           ratio

       tralm.feastol           0           feasibility tolerance
                                           (set to opf.violation by default)
       tralm.primaltol         5e-4        primal variable tolerance
       tralm.dualtol           5e-4        dual variable tolerance
       tralm.costtol           1e-5        optimality tolerance
       tralm.major_it          40          maximum number of major iterations
       tralm.minor_it          40          maximum number of minor iterations
       tralm.smooth_ratio      0.04        piecewise linear curve smoothing
                                           ratio

Experimental Options:
   exp.sys_wide_zip_loads.pw   <empty>     1 x 3 vector of active load fraction
                                           to be modeled as constant power,
                                           constant current and constant
                                           impedance, respectively, where
                                           <empty> means use [1 0 0]
   exp.sys_wide_zip_loads.qw   <empty>     same for reactive power, where
                                           <empty> means use same value as
                                           for 'pw'

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SUBFUNCTIONS ^

SOURCE CODE ^

0001 function opt = mpoption(varargin)
0002 %MPOPTION  Used to set and retrieve a MATPOWER options struct.
0003 %
0004 %   OPT = MPOPTION
0005 %       Returns the default options struct.
0006 %
0007 %   OPT = MPOPTION(OVERRIDES)
0008 %       Returns the default options struct, with some fields overridden
0009 %       by values from OVERRIDES, which can be a struct or the name of
0010 %       a function that returns a struct.
0011 %
0012 %   OPT = MPOPTION(NAME1, VALUE1, NAME2, VALUE2, ...)
0013 %       Same as previous, except override options are specified by NAME,
0014 %       VALUE pairs. This can be used to set any part of the options
0015 %       struct. The names can be individual fields or multi-level field
0016 %       names with embedded periods. The values can be scalars or structs.
0017 %
0018 %       For backward compatibility, the NAMES and VALUES may correspond
0019 %       to old-style MATPOWER option names (elements in the old-style
0020 %       options vector) as well.
0021 %
0022 %   OPT = MPOPTION(OPT0)
0023 %       Converts an old-style options vector OPT0 into the corresponding
0024 %       options struct. If OPT0 is an options struct it does nothing.
0025 %
0026 %   OPT = MPOPTION(OPT0, OVERRIDES)
0027 %       Applies overrides to an existing set of options, OPT0, which
0028 %       can be an old-style options vector or an options struct.
0029 %
0030 %   OPT = MPOPTION(OPT0, NAME1, VALUE1, NAME2, VALUE2, ...)
0031 %       Same as above except it uses the old-style options vector OPT0
0032 %       as a base instead of the old default options vector.
0033 %
0034 %   OPT_VECTOR = MPOPTION(OPT, [])
0035 %       Creates and returns an old-style options vector from an
0036 %       options struct OPT.
0037 %
0038 %   Note: The use of old-style MATPOWER options vectors and their
0039 %         names and values has been deprecated and will be removed
0040 %         in a future version of MATPOWER. Until then, all uppercase
0041 %         option names are not permitted for new top-level options.
0042 %
0043 %   Examples:
0044 %       mpopt = mpoption('pf.alg', 'FDXB', 'pf.tol', 1e-4);
0045 %       mpopt = mpoption(mpopt, 'opf.dc.solver', 'CPLEX', 'verbose', 2);
0046 %
0047 %The currently defined options are as follows:
0048 %
0049 %   name                    default     description [options]
0050 %----------------------    ---------   ----------------------------------
0051 %Model options:
0052 %   model                   'AC'        AC vs. DC power flow model
0053 %       [ 'AC' - use nonlinear AC model & corresponding algorithms/options  ]
0054 %       [ 'DC' - use linear DC model & corresponding algorithms/options     ]
0055 %
0056 %Power Flow options:
0057 %   pf.alg                  'NR'        AC power flow algorithm
0058 %       [ 'NR'    - Newton's method (formulation depends on values of       ]
0059 %       [           pf.current_balance and pf.v_cartesian options)          ]
0060 %       [ 'NR-SC' - Newton's method (power mismatch, cartesian)             ]
0061 %       [ 'NR-IP' - Newton's method (current mismatch, polar)               ]
0062 %       [ 'NR-IC' - Newton's method (current mismatch, cartesian)           ]
0063 %       [ 'FDXB'  - Fast-Decoupled (XB version)                             ]
0064 %       [ 'FDBX'  - Fast-Decoupled (BX version)                             ]
0065 %       [ 'GS'    - Gauss-Seidel                                            ]
0066 %       [ 'PQSUM' - Power Summation method (radial networks only)           ]
0067 %       [ 'ISUM'  - Current Summation method (radial networks only)         ]
0068 %       [ 'YSUM'  - Admittance Summation method (radial networks only)      ]
0069 %   pf.current_balance     0           type of nodal balance equation
0070 %       [ 0 - use complex power balance equations                           ]
0071 %       [ 1 - use complex current balance equations                         ]
0072 %   pf.v_cartesian         0           voltage representation
0073 %       [ 0 - bus voltage variables represented in polar coordinates        ]
0074 %       [ 1 - bus voltage variables represented in Cartesian coordinates    ]
0075 %   pf.tol                  1e-8        termination tolerance on per unit
0076 %                                       P & Q mismatch
0077 %   pf.nr.max_it            10          maximum number of iterations for
0078 %                                       Newton's method
0079 %   pf.nr.lin_solver        ''          linear solver passed to MPLINSOLVE to
0080 %                                       solve Newton update step
0081 %       [ ''      - default to '\' for small systems, 'LU3' for larger ones ]
0082 %       [ '\'     - built-in backslash operator                             ]
0083 %       [ 'LU'    - explicit default LU decomposition and back substitution ]
0084 %       [ 'LU3'   - 3 output arg form of LU, Gilbert-Peierls algorithm with ]
0085 %       [           approximate minimum degree (AMD) reordering             ]
0086 %       [ 'LU4'   - 4 output arg form of LU, UMFPACK solver (same as 'LU')  ]
0087 %       [ 'LU5'   - 5 output arg form of LU, UMFPACK solver, w/row scaling  ]
0088 %       [       (see MPLINSOLVE for complete list of all options)           ]
0089 %   pf.fd.max_it            30          maximum number of iterations for
0090 %                                       fast decoupled method
0091 %   pf.gs.max_it            1000        maximum number of iterations for
0092 %                                       Gauss-Seidel method
0093 %   pf.radial.max_it        20          maximum number of iterations for
0094 %                                       radial power flow methods
0095 %   pf.radial.vcorr         0           perform voltage correction procedure
0096 %                                       in distribution power flow
0097 %       [  0 - do NOT perform voltage correction                            ]
0098 %       [  1 - perform voltage correction                                   ]
0099 %   pf.enforce_q_lims       0           enforce gen reactive power limits at
0100 %                                       expense of |V|
0101 %       [  0 - do NOT enforce limits                                        ]
0102 %       [  1 - enforce limits, simultaneous bus type conversion             ]
0103 %       [  2 - enforce limits, one-at-a-time bus type conversion            ]
0104 %
0105 %Continuation Power Flow options:
0106 %   cpf.parameterization    3           choice of parameterization
0107 %       [  1 - natural                                                      ]
0108 %       [  2 - arc length                                                   ]
0109 %       [  3 - pseudo arc length                                            ]
0110 %   cpf.stop_at             'NOSE'      determins stopping criterion
0111 %       [ 'NOSE'     - stop when nose point is reached                      ]
0112 %       [ 'FULL'     - trace full nose curve                                ]
0113 %       [ <lam_stop> - stop upon reaching specified target lambda value     ]
0114 %   cpf.enforce_p_lims      0           enforce gen active power limits
0115 %       [  0 - do NOT enforce limits                                        ]
0116 %       [  1 - enforce limits, simultaneous bus type conversion             ]
0117 %   cpf.enforce_q_lims      0           enforce gen reactive power limits at
0118 %                                       expense of |V|
0119 %       [  0 - do NOT enforce limits                                        ]
0120 %       [  1 - enforce limits, simultaneous bus type conversion             ]
0121 %   cpf.enforce_v_lims      0           enforce bus voltage magnitude limits
0122 %       [  0 - do NOT enforce limits                                        ]
0123 %       [  1 - enforce limits, termination on detection                     ]
0124 %   cpf.enforce_flow_lims   0           enforce branch flow MVA limits
0125 %       [  0 - do NOT enforce limits                                        ]
0126 %       [  1 - enforce limits, termination on detection                     ]
0127 %   cpf.step                0.05        continuation power flow step size
0128 %   cpf.adapt_step          0           toggle adaptive step size feature
0129 %       [  0 - adaptive step size disabled                                  ]
0130 %       [  1 - adaptive step size enabled                                   ]
0131 %   cpf.step_min            1e-4        minimum allowed step size
0132 %   cpf.step_max            0.2         maximum allowed step size
0133 %   cpf.adapt_step_damping  0.7         damping factor for adaptive step
0134 %                                       sizing
0135 %   cpf.adapt_step_tol      1e-3        tolerance for adaptive step sizing
0136 %   cpf.target_lam_tol      1e-5        tolerance for target lambda detection
0137 %   cpf.nose_tol            1e-5        tolerance for nose point detection (pu)
0138 %   cpf.p_lims_tol          0.01        tolerance for generator active
0139 %                                       power limit enforcement (MW)
0140 %   cpf.q_lims_tol          0.01        tolerance for generator reactive
0141 %                                       power limit enforcement (MVAR)
0142 %   cpf.v_lims_tol          1e-4        tolerance for bus voltage
0143 %                                       magnitude enforcement (p.u)
0144 %   cpf.flow_lims_tol       0.01        tolerance for line MVA flow
0145 %                                       enforcement (MVA)
0146 %   cpf.plot.level          0           control plotting of noze curve
0147 %       [  0 - do not plot nose curve                                       ]
0148 %       [  1 - plot when completed                                          ]
0149 %       [  2 - plot incrementally at each iteration                         ]
0150 %       [  3 - same as 2, with 'pause' at each iteration                    ]
0151 %   cpf.plot.bus            <empty>     index of bus whose voltage is to be
0152 %                                       plotted
0153 %   cpf.user_callback       <empty>     string containing the name of a user
0154 %                                       callback function, or struct with
0155 %                                       function name, and optional priority
0156 %                                       and/or args, or cell array of such
0157 %                                       strings and/or structs, see
0158 %                                       'help cpf_default_callback' for details
0159 %
0160 %Optimal Power Flow options:
0161 %   name                    default     description [options]
0162 %----------------------    ---------   ----------------------------------
0163 %   opf.ac.solver           'DEFAULT'   AC optimal power flow solver
0164 %       [ 'DEFAULT' - choose default solver, i.e. 'MIPS'                    ]
0165 %       [ 'MIPS'    - MIPS, MATPOWER Interior Point Solver, primal/dual     ]
0166 %       [             interior point method (pure MATLAB/Octave)            ]
0167 %       [ 'FMINCON' - MATLAB Optimization Toolbox, FMINCON                  ]
0168 %       [ 'IPOPT'   - IPOPT, requires MEX interface to IPOPT solver         ]
0169 %       [             available from:                                       ]
0170 %       [                   https://github.com/coin-or/Ipopt                ]
0171 %       [ 'KNITRO'  - Artelys Knitro, requires Artelys Knitro solver,       ]
0172 %       [             available from:https://www.artelys.com/solvers/knitro/]
0173 %       [ 'MINOPF'  - MINOPF, MINOS-based solver, requires optional         ]
0174 %       [             MEX-based MINOPF package, available from:             ]
0175 %       [                   http://www.pserc.cornell.edu/minopf/            ]
0176 %       [ 'PDIPM'   - PDIPM, primal/dual interior point method, requires    ]
0177 %       [             optional MEX-based TSPOPF package, available from:    ]
0178 %       [                   http://www.pserc.cornell.edu/tspopf/            ]
0179 %       [ 'SDPOPF'  - SDPOPF, solver based on semidefinite relaxation of    ]
0180 %       [             OPF problem, requires optional packages:              ]
0181 %       [               SDP_PF, available in extras/sdp_pf                  ]
0182 %       [               YALMIP, available from:                             ]
0183 %       [                   https://yalmip.github.io                        ]
0184 %       [               SDP solver such as SeDuMi, available from:          ]
0185 %       [                   http://sedumi.ie.lehigh.edu/                    ]
0186 %       [ 'TRALM'   - TRALM, trust region based augmented Langrangian       ]
0187 %       [             method, requires TSPOPF (see 'PDIPM')                 ]
0188 %   opf.dc.solver           'DEFAULT'   DC optimal power flow solver
0189 %       [ 'DEFAULT' - choose solver based on availability in the following  ]
0190 %       [             order: 'GUROBI', 'CPLEX', 'MOSEK', 'OT',              ]
0191 %       [             'GLPK' (linear costs only), 'BPMPD', 'MIPS'           ]
0192 %       [ 'MIPS'    - MIPS, MATPOWER Interior Point Solver, primal/dual     ]
0193 %       [             interior point method (pure MATLAB/Octave)            ]
0194 %       [ 'BPMPD'   - BPMPD, requires optional MEX-based BPMPD_MEX package  ]
0195 %       [             available from: http://www.pserc.cornell.edu/bpmpd/   ]
0196 %       [ 'CLP'     - CLP, requires interface to COIN-OP LP solver          ]
0197 %       [             available from:https://github.com/coin-or/Clp         ]
0198 %       [ 'CPLEX'   - CPLEX, requires CPLEX solver available from:          ]
0199 %       [             https://www.ibm.com/analytics/cplex-optimizer         ]
0200 %       [ 'GLPK'    - GLPK, requires interface to GLPK solver               ]
0201 %       [             available from: https://www.gnu.org/software/glpk/    ]
0202 %       [             (GLPK does not work with quadratic cost functions)    ]
0203 %       [ 'GUROBI'  - GUROBI, requires Gurobi optimizer (v. 5+)             ]
0204 %       [             available from: https://www.gurobi.com/               ]
0205 %       [ 'IPOPT'   - IPOPT, requires MEX interface to IPOPT solver         ]
0206 %       [             available from:                                       ]
0207 %       [                 https://github.com/coin-or/Ipopt                  ]
0208 %       [ 'MOSEK'   - MOSEK, requires MATLAB interface to MOSEK solver      ]
0209 %       [             available from: https://www.mosek.com/                ]
0210 %       [ 'OT'      - MATLAB Optimization Toolbox, QUADPROG, LINPROG        ]
0211 %   opf.current_balance     0           type of nodal balance equation
0212 %       [ 0 - use complex power balance equations                           ]
0213 %       [ 1 - use complex current balance equations                         ]
0214 %   opf.v_cartesian         0           voltage representation
0215 %       [ 0 - bus voltage variables represented in polar coordinates        ]
0216 %       [ 1 - bus voltage variables represented in Cartesian coordinates    ]
0217 %   opf.violation           5e-6        constraint violation tolerance
0218 %   opf.use_vg              0           respect gen voltage setpt     [ 0-1 ]
0219 %       [ 0 - use specified bus Vmin & Vmax, and ignore gen Vg              ]
0220 %       [ 1 - replace specified bus Vmin & Vmax by corresponding gen Vg     ]
0221 %       [ between 0 and 1 - use a weighted average of the 2 options         ]
0222 %   opf.flow_lim            'S'         quantity limited by branch flow
0223 %                                       constraints
0224 %       [ 'S' - apparent power flow (limit in MVA)                          ]
0225 %       [ 'P' - active power flow, implemented using P   (limit in MW)      ]
0226 %       [ '2' - active power flow, implemented using P^2 (limit in MW)      ]
0227 %       [ 'I' - current magnitude (limit in MVA at 1 p.u. voltage)          ]
0228 %   opf.ignore_angle_lim    0           angle diff limits for branches
0229 %       [ 0 - include angle difference limits, if specified                 ]
0230 %       [ 1 - ignore angle difference limits even if specified              ]
0231 %   opf.softlims.default    1           behavior of OPF soft limits for
0232 %                                       which parameters are not explicitly
0233 %                                       provided
0234 %       [ 0 - do not include softlims if not explicitly specified           ]
0235 %       [ 1 - include softlims w/default values if not explicitly specified ]
0236 %   opf.init_from_mpc       -1          (DEPRECATED: use opf.start instead)
0237 %                                       specify whether to use current state
0238 %                                       in MATPOWER case to initialize OPF
0239 %                                       (currently only supported by fmincon,
0240 %                                        Ipopt, Knitro and MIPS solvers,
0241 %                                        others always use mpc)
0242 %       [  -1 - MATPOWER decides, based on solver/algorithm                 ]
0243 %       [   0 - ignore current state in MATPOWER case (only applies to      ]
0244 %       [       fmincon, Ipopt, Knitro and MIPS, which use an interior pt   ]
0245 %       [       estimate; others use current state as with opf.start = 2)   ]
0246 %       [   1 - use current state in MATPOWER case                          ]
0247 %   opf.start               0           strategy for initializing OPF start pt
0248 %       [   0 - default, MATPOWER decides based on solver                   ]
0249 %       [       (currently identical to 1)                                  ]
0250 %       [   1 - ignore current state in MATPOWER case (only applies to      ]
0251 %       [       fmincon, Ipopt, Knitro and MIPS, which use an interior pt   ]
0252 %       [       estimate; others use current state as with opf.start = 2)   ]
0253 %       [   2 - use current state in MATPOWER case                          ]
0254 %       [   3 - solve power flow and use resulting state                    ]
0255 %   opf.return_raw_der      0           for AC OPF, return constraint and
0256 %                                       derivative info in results.raw
0257 %                                       (in fields g, dg, df, d2f) [ 0 or 1 ]
0258 %
0259 %Output options:
0260 %   name                    default     description [options]
0261 %----------------------    ---------   ----------------------------------
0262 %   verbose                 1           amount of progress info printed
0263 %       [   0 - print no progress info                                      ]
0264 %       [   1 - print a little progress info                                ]
0265 %       [   2 - print a lot of progress info                                ]
0266 %       [   3 - print all progress info                                     ]
0267 %   out.all                 -1          controls pretty-printing of results
0268 %       [  -1 - individual flags control what prints                        ]
0269 %       [   0 - do not print anything (overrides individual flags, ignored  ]
0270 %       [       for files specified as FNAME arg to runpf(), runopf(), etc.)]
0271 %       [   1 - print everything (overrides individual flags)               ]
0272 %   out.sys_sum             1           print system summary       [ 0 or 1 ]
0273 %   out.area_sum            0           print area summaries       [ 0 or 1 ]
0274 %   out.bus                 1           print bus detail           [ 0 or 1 ]
0275 %   out.branch              1           print branch detail        [ 0 or 1 ]
0276 %   out.gen                 0           print generator detail     [ 0 or 1 ]
0277 %   out.lim.all             -1          controls constraint info output
0278 %       [  -1 - individual flags control what constraint info prints        ]
0279 %       [   0 - no constraint info (overrides individual flags)             ]
0280 %       [   1 - binding constraint info (overrides individual flags)        ]
0281 %       [   2 - all constraint info (overrides individual flags)            ]
0282 %   out.lim.v               1           control voltage limit info
0283 %       [   0 - do not print                                                ]
0284 %       [   1 - print binding constraints only                              ]
0285 %       [   2 - print all constraints                                       ]
0286 %       [   (same options for OUT_LINE_LIM, OUT_PG_LIM, OUT_QG_LIM)         ]
0287 %   out.lim.line            1           control line flow limit info
0288 %   out.lim.pg              1           control gen active power limit info
0289 %   out.lim.qg              1           control gen reactive pwr limit info
0290 %   out.force               0           print results even if success
0291 %                                       flag = 0                   [ 0 or 1 ]
0292 %   out.suppress_detail     -1          suppress all output but system summary
0293 %       [  -1 - suppress details for large systems (> 500 buses)            ]
0294 %       [   0 - do not suppress any output specified by other flags         ]
0295 %       [   1 - suppress all output except system summary section           ]
0296 %       [       (overrides individual flags, but not out.all = 1)           ]
0297 %
0298 %Solver specific options:
0299 %       name                    default     description [options]
0300 %   -----------------------    ---------   ----------------------------------
0301 %   MIPS:
0302 %       mips.linsolver          ''          linear system solver
0303 %           [   '' or '\'   build-in backslash \ operator (e.g. x = A \ b)  ]
0304 %           [   'PARDISO'   PARDISO solver (if available)                   ]
0305 %       mips.feastol            0           feasibility (equality) tolerance
0306 %                                           (set to opf.violation by default)
0307 %       mips.gradtol            1e-6        gradient tolerance
0308 %       mips.comptol            1e-6        complementary condition
0309 %                                           (inequality) tolerance
0310 %       mips.costtol            1e-6        optimality tolerance
0311 %       mips.max_it             150         maximum number of iterations
0312 %       mips.step_control       0           enable step-size cntrl [ 0 or 1 ]
0313 %       mips.sc.red_it          20          maximum number of reductions per
0314 %                                           iteration with step control
0315 %       mips.xi                 0.99995     constant used in alpha updates*
0316 %       mips.sigma              0.1         centering parameter*
0317 %       mips.z0                 1           used to initialize slack variables*
0318 %       mips.alpha_min          1e-8        returns "Numerically Failed" if
0319 %                                           either alpha parameter becomes
0320 %                                           smaller than this value*
0321 %       mips.rho_min            0.95        lower bound on rho_t*
0322 %       mips.rho_max            1.05        upper bound on rho_t*
0323 %       mips.mu_threshold       1e-5        KT multipliers smaller than this
0324 %                                           value for non-binding constraints
0325 %                                           are forced to zero
0326 %       mips.max_stepsize       1e10        returns "Numerically Failed" if the
0327 %                                           2-norm of the reduced Newton step
0328 %                                           exceeds this value*
0329 %           * See the corresponding Appendix in the manual for details.
0330 %
0331 %   CPLEX:
0332 %       cplex.lpmethod          0           solution algorithm for LP problems
0333 %           [   0 - automatic: let CPLEX choose                             ]
0334 %           [   1 - primal simplex                                          ]
0335 %           [   2 - dual simplex                                            ]
0336 %           [   3 - network simplex                                         ]
0337 %           [   4 - barrier                                                 ]
0338 %           [   5 - sifting                                                 ]
0339 %           [   6 - concurrent (dual, barrier, and primal)                  ]
0340 %       cplex.qpmethod          0           solution algorithm for QP problems
0341 %           [   0 - automatic: let CPLEX choose                             ]
0342 %           [   1 - primal simplex optimizer                                ]
0343 %           [   2 - dual simplex optimizer                                  ]
0344 %           [   3 - network optimizer                                       ]
0345 %           [   4 - barrier optimizer                                       ]
0346 %       cplex.opts              <empty>     see CPLEX_OPTIONS for details
0347 %       cplex.opt_fname         <empty>     see CPLEX_OPTIONS for details
0348 %       cplex.opt               0           see CPLEX_OPTIONS for details
0349 %
0350 %   FMINCON:
0351 %       fmincon.alg             4           algorithm used by fmincon() for OPF
0352 %                                           for Opt Toolbox 4 and later
0353 %            [  1 - active-set (not suitable for large problems)            ]
0354 %            [  2 - interior-point, w/default 'bfgs' Hessian approx         ]
0355 %            [  3 - interior-point, w/ 'lbfgs' Hessian approx               ]
0356 %            [  4 - interior-point, w/exact user-supplied Hessian           ]
0357 %            [  5 - interior-point, w/Hessian via finite differences        ]
0358 %            [  6 - sqp (not suitable for large problems)                   ]
0359 %       fmincon.tol_x           1e-4        termination tol on x
0360 %       fmincon.tol_f           1e-4        termination tol on f
0361 %       fmincon.max_it          0           maximum number of iterations
0362 %                                                           [  0 => default ]
0363 %
0364 %   GUROBI:
0365 %       gurobi.method           0           solution algorithm (Method)
0366 %           [  -1 - automatic, let Gurobi decide                            ]
0367 %           [   0 - primal simplex                                          ]
0368 %           [   1 - dual simplex                                            ]
0369 %           [   2 - barrier                                                 ]
0370 %           [   3 - concurrent (LP only)                                    ]
0371 %           [   4 - deterministic concurrent (LP only)                      ]
0372 %       gurobi.timelimit        Inf         maximum time allowed (TimeLimit)
0373 %       gurobi.threads          0           max number of threads (Threads)
0374 %       gurobi.opts             <empty>     see GUROBI_OPTIONS for details
0375 %       gurobi.opt_fname        <empty>     see GUROBI_OPTIONS for details
0376 %       gurobi.opt              0           see GUROBI_OPTIONS for details
0377 %
0378 %   IPOPT:
0379 %       ipopt.opts              <empty>     see IPOPT_OPTIONS for details
0380 %       ipopt.opt_fname         <empty>     see IPOPT_OPTIONS for details
0381 %       ipopt.opt               0           see IPOPT_OPTIONS for details
0382 %
0383 %   KNITRO:
0384 %       knitro.tol_x            1e-4        termination tol on x
0385 %       knitro.tol_f            1e-4        termination tol on f
0386 %       knitro.maxit            0           maximum number of iterations
0387 %                                                           [  0 => default ]
0388 %       knitro.opt_fname        <empty>     name of user-supplied native
0389 %                                           Knitro options file that overrides
0390 %                                           all other options
0391 %       knitro.opt              0           if knitro.opt_fname is empty and
0392 %                                           knitro.opt is a non-zero integer N
0393 %                                           then knitro.opt_fname is auto-
0394 %                                           generated as:
0395 %                                           'knitro_user_options_N.txt'
0396 %
0397 %   LINPROG:
0398 %       linprog                 <empty>     LINPROG options passed to
0399 %                                           OPTIMOPTIONS or OPTIMSET.
0400 %                                           see LINPROG in the Optimization
0401 %                                           Toolbox for details
0402 %
0403 %   MINOPF:
0404 %       minopf.feastol          0 (1e-3)    primal feasibility tolerance
0405 %                                           (set to opf.violation by default)
0406 %       minopf.rowtol           0 (1e-3)    row tolerance
0407 %       minopf.xtol             0 (1e-4)    x tolerance
0408 %       minopf.majdamp          0 (0.5)     major damping parameter
0409 %       minopf.mindamp          0 (2.0)     minor damping parameter
0410 %       minopf.penalty          0 (1.0)     penalty parameter
0411 %       minopf.major_it         0 (200)     major iterations
0412 %       minopf.minor_it         0 (2500)    minor iterations
0413 %       minopf.max_it           0 (2500)    iterations limit
0414 %       minopf.verbosity        -1          amount of progress info printed
0415 %           [  -1 - controlled by 'verbose' option                          ]
0416 %           [   0 - print nothing                                           ]
0417 %           [   1 - print only termination status message                   ]
0418 %           [   2 - print termination status and screen progress            ]
0419 %           [   3 - print screen progress, report file (usually fort.9)     ]
0420 %       minopf.core             0 (1200*nb + 2*(nb+ng)^2) memory allocation
0421 %       minopf.supbasic_lim     0 (2*nb + 2*ng) superbasics limit
0422 %       minopf.mult_price       0 (30)      multiple price
0423 %
0424 %   MOSEK:
0425 %       mosek.lp_alg            0           solution algorithm
0426 %                                               (MSK_IPAR_OPTIMIZER)
0427 %           for MOSEK 8.x ...        (see MOSEK_SYMBCON for a "better way")
0428 %           [   0 - automatic: let MOSEK choose                             ]
0429 %           [   1 - dual simplex                                            ]
0430 %           [   2 - automatic: let MOSEK choose                             ]
0431 %           [   3 - automatic simplex (MOSEK chooses which simplex method)  ]
0432 %           [   4 - interior point                                          ]
0433 %           [   6 - primal simplex                                          ]
0434 %       mosek.max_it            0 (400)     interior point max iterations
0435 %                                               (MSK_IPAR_INTPNT_MAX_ITERATIONS)
0436 %       mosek.gap_tol           0 (1e-8)    interior point relative gap tol
0437 %                                               (MSK_DPAR_INTPNT_TOL_REL_GAP)
0438 %       mosek.max_time          0 (-1)      maximum time allowed
0439 %                                               (MSK_DPAR_OPTIMIZER_MAX_TIME)
0440 %       mosek.num_threads       0 (1)       max number of threads
0441 %                                               (MSK_IPAR_INTPNT_NUM_THREADS)
0442 %       mosek.opts              <empty>     see MOSEK_OPTIONS for details
0443 %       mosek.opt_fname         <empty>     see MOSEK_OPTIONS for details
0444 %       mosek.opt               0           see MOSEK_OPTIONS for details
0445 %
0446 %   QUADPROG:
0447 %       quadprog                <empty>     QUADPROG options passed to
0448 %                                           OPTIMOPTIONS or OPTIMSET.
0449 %                                           see QUADPROG in the Optimization
0450 %                                           Toolbox for details
0451 %
0452 %   TSPOPF:
0453 %       pdipm.feastol           0           feasibility (equality) tolerance
0454 %                                           (set to opf.violation by default)
0455 %       pdipm.gradtol           1e-6        gradient tolerance
0456 %       pdipm.comptol           1e-6        complementary condition
0457 %                                           (inequality) tolerance
0458 %       pdipm.costtol           1e-6        optimality tolerance
0459 %       pdipm.max_it            150         maximum number of iterations
0460 %       pdipm.step_control      0           enable step-size cntrl [ 0 or 1 ]
0461 %       pdipm.sc.red_it         20          maximum number of reductions per
0462 %                                           iteration with step control
0463 %       pdipm.sc.smooth_ratio   0.04        piecewise linear curve smoothing
0464 %                                           ratio
0465 %
0466 %       tralm.feastol           0           feasibility tolerance
0467 %                                           (set to opf.violation by default)
0468 %       tralm.primaltol         5e-4        primal variable tolerance
0469 %       tralm.dualtol           5e-4        dual variable tolerance
0470 %       tralm.costtol           1e-5        optimality tolerance
0471 %       tralm.major_it          40          maximum number of major iterations
0472 %       tralm.minor_it          40          maximum number of minor iterations
0473 %       tralm.smooth_ratio      0.04        piecewise linear curve smoothing
0474 %                                           ratio
0475 %
0476 %Experimental Options:
0477 %   exp.sys_wide_zip_loads.pw   <empty>     1 x 3 vector of active load fraction
0478 %                                           to be modeled as constant power,
0479 %                                           constant current and constant
0480 %                                           impedance, respectively, where
0481 %                                           <empty> means use [1 0 0]
0482 %   exp.sys_wide_zip_loads.qw   <empty>     same for reactive power, where
0483 %                                           <empty> means use same value as
0484 %                                           for 'pw'
0485 
0486 %   MATPOWER
0487 %   Copyright (c) 2013-2017, Power Systems Engineering Research Center (PSERC)
0488 %   by Ray Zimmerman, PSERC Cornell
0489 %   and Shrirang Abhyankar, Argonne National Laboratory
0490 %
0491 %   This file is part of MATPOWER.
0492 %   Covered by the 3-clause BSD License (see LICENSE file for details).
0493 %   See https://matpower.org for more info.
0494 
0495 %% some constants
0496 N = 124;                %% number of options in old-style vector (MATPOWER 4.1)
0497 N40 = 116;              %% dimension of MATPOWER 4.0 options vector
0498 N32 = 93;               %% dimension of MATPOWER 3.2 options vector
0499 v = mpoption_version;   %% version number of MATPOWER options struct
0500 
0501 %% initialize flags and arg counter
0502 have_opt0 = 0;          %% existing options struct or vector provided?
0503 have_old_style_ov = 0;  %% override options using old-style names?
0504 return_old_style = 0;   %% return value as old-style vector?
0505 k = 1;
0506 if nargin > 0
0507     opt0 = varargin{k};
0508     if isstruct(opt0) && isfield(opt0, 'v')         %% options struct
0509         have_opt0 = 1;
0510         k = k + 1;
0511     elseif isnumeric(opt0) && size(opt0, 2) == 1    %% options vector
0512         nn = size(opt0, 1);
0513         if ismember(nn, [N N40 N32])                %% of valid size
0514             %% expand smaller option vectors (from before MATPOWER 4.1)
0515             if nn < N
0516                 optv = mpoption_old();
0517                 opt0(nn+1:N) = optv(nn+1:N);
0518             end
0519             have_opt0 = 1;
0520             k = k + 1;
0521         end
0522     end
0523 end
0524 
0525 %% create base options vector to which overrides are made
0526 if have_opt0
0527     if isstruct(opt0)               %% it's already a valid options struct
0528         if DEBUG, fprintf('OPT0 is a valid options struct\n'); end
0529         if opt0.v < v
0530             %% convert older version to current version
0531             opt_d = mpoption_default();
0532             if opt0.v == 1          %% convert version 1 to 2
0533                 if isfield(opt_d, 'linprog')
0534                     opt0.lingprog = opt_d.linprog;
0535                 end
0536                 if isfield(opt_d, 'quadprog')
0537                     opt0.quadprog = opt_d.quadprog;
0538                 end
0539             end
0540             if opt0.v <= 2          %% convert version 2 to 3
0541                 opt0.out.suppress_detail = opt_d.out.suppress_detail;
0542             end
0543             %if opt0.v <= 3          %% convert version 3 to 4
0544                 %% new mips options were all optional, no conversion needed
0545             %end
0546             if opt0.v <= 4          %% convert version 4 to 5
0547                 opt0.opf.init_from_mpc = opt_d.opf.init_from_mpc;
0548             end
0549             if opt0.v <= 5          %% convert version 5 to 6
0550                 if isfield(opt_d, 'clp')
0551                     opt0.clp = opt_d.clp;
0552                 end
0553             end
0554             if opt0.v <= 6          %% convert version 6 to 7
0555                 if isfield(opt_d, 'intlinprog')
0556                     opt0.intlinprog = opt_d.intlinprog;
0557                 end
0558             end
0559             if opt0.v <= 7          %% convert version 7 to 8
0560                 opt0.mips.linsolver = opt_d.mips.linsolver;
0561             end
0562             if opt0.v <= 8          %% convert version 8 to 9
0563                 opt0.exp.sys_wide_zip_loads = opt_d.exp.sys_wide_zip_loads;
0564             end
0565             if opt0.v <= 9          %% convert version 9 to 10
0566                 opt0.most = opt_d.most;
0567             end
0568             if opt0.v <= 10         %% convert version 10 to 11
0569                 opt0.cpf.enforce_p_lims     = opt_d.cpf.enforce_p_lims;
0570                 opt0.cpf.enforce_q_lims     = opt_d.cpf.enforce_q_lims;
0571                 opt0.cpf.adapt_step_damping = opt_d.cpf.adapt_step_damping;
0572                 opt0.cpf.target_lam_tol     = opt_d.cpf.target_lam_tol;
0573                 opt0.cpf.nose_tol           = opt_d.cpf.nose_tol;
0574                 opt0.cpf.p_lims_tol         = opt_d.cpf.p_lims_tol;
0575                 opt0.cpf.q_lims_tol         = opt_d.cpf.q_lims_tol;
0576                 if (~isempty(opt0.cpf.user_callback_args) && ...
0577                         ~isstruct(opt0.cpf.user_callback_args)) || ...
0578                         (isstruct(opt0.cpf.user_callback_args) && ...
0579                          ~isempty(fields(opt0.cpf.user_callback_args)))
0580                     warning('The ''cpf.user_callback_args'' option has been removed. Please include the args in a struct in ''cpf.user_callback'' instead.')
0581                 end
0582                 opt0.cpf = rmfield(opt0.cpf, 'user_callback_args');
0583             end
0584             if opt0.v <= 11         %% convert version 11 to 12
0585                 opt0.opf.use_vg = opt_d.opf.use_vg;
0586             end
0587             if opt0.v <= 12         %% convert version 12 to 13
0588                 opt0.pf.radial.max_it   = opt_d.pf.radial.max_it;
0589                 opt0.pf.radial.vcorr    = opt_d.pf.radial.vcorr;
0590             end
0591             if opt0.v <= 13         %% convert version 13 to 14
0592                 opt0.pf.nr.lin_solver   = opt_d.pf.nr.lin_solver;
0593             end
0594             if opt0.v <= 14         %% convert version 14 to 15
0595                 opt0.cpf.enforce_v_lims     = opt_d.cpf.enforce_v_lims;
0596                 opt0.cpf.enforce_flow_lims  = opt_d.cpf.enforce_flow_lims;
0597                 opt0.cpf.v_lims_tol         = opt_d.cpf.v_lims_tol;
0598                 opt0.cpf.flow_lims_tol      = opt_d.cpf.flow_lims_tol;
0599             end
0600             if opt0.v <= 15          %% convert version 15 to 16
0601                 opt0.opf.start  = opt_d.opf.start;
0602             end
0603             if opt0.v <= 16          %% convert version 16 to 17
0604                 if isfield(opt_d, 'knitro')
0605                     opt0.knitro.maxit  = opt_d.knitro.maxit;
0606                 end
0607             end
0608             if opt0.v <= 17         %% convert version 17 to 18
0609                 opt0.opf.current_balance = opt_d.opf.current_balance;
0610                 opt0.opf.v_cartesian     = opt_d.opf.v_cartesian;
0611             end
0612             if opt0.v <= 18         %% convert version 18 to 19
0613                 opt0.opf.softlims.default = opt_d.opf.softlims.default;
0614             end
0615             if opt0.v <= 19         %% convert version 19 to 20
0616                 opt0.pf.current_balance = opt_d.pf.current_balance;
0617                 opt0.pf.v_cartesian     = opt_d.pf.v_cartesian;
0618             end
0619             opt0.v = v;
0620         end
0621         opt = opt0;
0622     else                            %% convert from old-style options vector
0623         if DEBUG, fprintf('OPT0 is a old-style options vector\n'); end
0624         opt = mpoption_v2s(opt0);
0625     end
0626 else                                %% use default options struct as base
0627     if DEBUG, fprintf('no OPT0, starting with default options struct\n'); end
0628     opt = mpoption_default();
0629 end
0630 
0631 
0632 %% do we have OVERRIDES or NAME/VALUE pairs
0633 ov = [];
0634 if nargin - k == 0          %% looking at last arg, must be OVERRIDES
0635     if isstruct(varargin{k})        %% OVERRIDES provided as struct
0636         if DEBUG, fprintf('OVERRIDES struct\n'); end
0637         ov = varargin{k};
0638     elseif ischar(varargin{k})      %% OVERRIDES provided as file/function name
0639         if DEBUG, fprintf('OVERRIDES file/function name\n'); end
0640         try
0641             ov = feval(varargin{k});
0642         catch
0643             error('mpoption: Unable to load MATPOWER options from ''%s''', varargin{k});
0644         end
0645         if ~isstruct(ov)
0646             error('mpoption: calling ''%s'' did not return a struct', varargin{k});
0647         end
0648     elseif isempty(varargin{k})
0649         return_old_style = 1;
0650     else
0651         error('mpoption: OVERRIDES must be a struct or the name of a function that returns a struct');
0652     end
0653 elseif nargin - k > 0 && mod(nargin-k, 2)   %% even number of remaining args
0654     if DEBUG, fprintf('NAME/VALUE pairs override defaults\n'); end
0655     %% process NAME/VALUE pairs
0656     if strcmp(varargin{k}, upper(varargin{k}))  %% old-style, all UPPERCASE option pairs
0657         %% NOTE: new top-level option fields cannot be all uppercase
0658         if ~have_opt0
0659             opt_v = mpoption_old(varargin{:});  %% create modified vector ...
0660             opt = mpoption_v2s(opt_v);          %% ... then convert
0661         else
0662             have_old_style_ov = 1;
0663             %% convert pairs to struct
0664             while k < nargin
0665                 name = varargin{k};
0666                 val  = varargin{k+1};
0667                 k = k + 2;
0668                 ov.(name) = val;
0669             end
0670         end
0671     else                                        %% new option pairs
0672         %% convert pairs to struct
0673         while k < nargin
0674             name = varargin{k};
0675             val  = varargin{k+1};
0676             k = k + 2;
0677             c = regexp(name, '([^\.]*)', 'tokens');
0678             s = struct();
0679             for i = 1:length(c)
0680                 s(i).type = '.';
0681                 s(i).subs = c{i}{1};
0682             end
0683             ov = subsasgn(ov, s, val);
0684         end
0685     end
0686 elseif nargin == 0 || nargin == 1
0687     if DEBUG, fprintf('no OVERRIDES, return default options struct or converted OPT0 vector\n'); end
0688 else
0689     error('mpoption: invalid calling syntax, see ''help mpoption'' to double-check the valid options');
0690 end
0691 
0692 %% apply overrides
0693 if ~isempty(ov)
0694     if have_old_style_ov
0695         opt = apply_old_mpoption_overrides(opt, ov);
0696     else
0697         persistent nsc_opt;     %% cache this to speed things up
0698         if ~isstruct(nsc_opt)
0699             vf = nested_struct_copy(mpoption_default(), mpoption_info_mips('V'));
0700             vf = nested_struct_copy(vf, mpoption_optional_fields());
0701             ex = struct(...
0702                 'name', {}, ...
0703                 'check', {}, ...
0704                 'copy_mode', {} ...
0705             );
0706             %% add exceptions for optional packages
0707             opt_pkgs = mpoption_optional_pkgs();
0708             n = length(ex);
0709             for k = 1:length(opt_pkgs)
0710                 fname = ['mpoption_info_' opt_pkgs{k}];
0711                 if exist(fname, 'file') == 2
0712                     opt_ex = feval(fname, 'E');
0713                     nex = length(opt_ex);
0714                     if ~isempty(opt_ex)
0715                         for j = 1:nex
0716                             ex(n+j).name = opt_ex(j).name;
0717                         end
0718                         if isfield(opt_ex, 'check')
0719                             for j = 1:nex
0720                                 ex(n+j).check = opt_ex(j).check;
0721                             end
0722                         end
0723                         if isfield(opt_ex, 'copy_mode')
0724                             for j = 1:nex
0725                                 ex(n+j).copy_mode = opt_ex(j).copy_mode;
0726                             end
0727                         end
0728                         if isfield(opt_ex, 'valid_fields')
0729                             for j = 1:nex
0730                                 ex(n+j).valid_fields = opt_ex(j).valid_fields;
0731                             end
0732                         end
0733                         n = n + nex;
0734                     end
0735                 end
0736             end
0737             nsc_opt = struct('check', 1, 'valid_fields', vf, 'exceptions', ex);
0738         end
0739 %         if have_fcn('catchme')
0740 %             try
0741 %                 opt = nested_struct_copy(opt, ov, nsc_opt);
0742 %             catch me
0743 %                 str = strrep(me.message, 'field', 'option');
0744 %                 str = strrep(str, 'nested_struct_copy', 'mpoption');
0745 %                 error(str);
0746 %             end
0747 %         else
0748             try
0749                 opt = nested_struct_copy(opt, ov, nsc_opt);
0750             catch
0751                 me = lasterr;
0752                 str = strrep(me, 'field', 'option');
0753                 str = strrep(str, 'nested_struct_copy', 'mpoption');
0754                 error(str);
0755             end
0756 %         end
0757     end
0758 end
0759 if return_old_style
0760     opt = mpoption_s2v(opt);
0761 end
0762 
0763 
0764 %%-------------------------------------------------------------------
0765 function opt = apply_old_mpoption_overrides(opt0, ov)
0766 %
0767 %   OPT0 is assumed to already have all of the fields and sub-fields found
0768 %   in the default options struct.
0769 
0770 %% initialize output
0771 opt = opt0;
0772 
0773 errstr = 'mpoption: %g is not a valid value for the old-style ''%s'' option';
0774 fields = fieldnames(ov);
0775 for f = 1:length(fields)
0776     ff = fields{f};
0777     switch ff
0778         case 'PF_ALG'
0779             switch ov.(ff)
0780                 case 1
0781                     opt.pf.alg = 'NR';      %% Newton's method
0782                 case 2
0783                     opt.pf.alg = 'FDXB';    %% fast-decoupled (XB version)
0784                 case 3
0785                     opt.pf.alg = 'FDBX';    %% fast-decoupled (BX version)
0786                 case 4
0787                     opt.pf.alg = 'GS';      %% Gauss-Seidel
0788                 otherwise
0789                     error(errstr, ov.(ff), ff);
0790             end
0791         case 'PF_TOL'
0792             opt.pf.tol = ov.(ff);
0793         case 'PF_MAX_IT'
0794             opt.pf.nr.max_it = ov.(ff);
0795         case 'PF_MAX_IT_FD'
0796             opt.pf.fd.max_it = ov.(ff);
0797         case 'PF_MAX_IT_GS'
0798             opt.pf.gs.max_it = ov.(ff);
0799         case 'ENFORCE_Q_LIMS'
0800             opt.pf.enforce_q_lims = ov.(ff);
0801         case 'PF_DC'
0802             switch ov.(ff)
0803                 case 0
0804                     opt.model = 'AC';
0805                 case 1
0806                     opt.model = 'DC';
0807                 otherwise
0808                     error(errstr, ov.(ff), ff);
0809             end
0810         case 'OPF_ALG'
0811             switch ov.(ff)
0812                 case 0
0813                     opt.opf.ac.solver = 'DEFAULT';
0814                 case 500
0815                     opt.opf.ac.solver = 'MINOPF';
0816                 case 520
0817                     opt.opf.ac.solver = 'FMINCON';
0818                 case {540, 545}
0819                     opt.opf.ac.solver = 'PDIPM';
0820                     if ov.(ff) == 545
0821                         opt.pdipm.step_control = 1;
0822                     else
0823                         opt.pdipm.step_control = 0;
0824                     end
0825                 case 550
0826                     opt.opf.ac.solver = 'TRALM';
0827                 case {560, 565}
0828                     opt.opf.ac.solver = 'MIPS';
0829                     if ov.(ff) == 565
0830                         opt.mips.step_control = 1;
0831                     else
0832                         opt.mips.step_control = 0;
0833                     end
0834                 case 580
0835                     opt.opf.ac.solver = 'IPOPT';
0836                 case 600
0837                     opt.opf.ac.solver = 'KNITRO';
0838                 otherwise
0839                     error(errstr, ov.(ff), ff);
0840             end
0841         case 'OPF_VIOLATION'
0842             opt.opf.violation = ov.(ff);
0843         case 'CONSTR_TOL_X'
0844             opt.fmincon.tol_x = ov.(ff);
0845             opt.knitro.tol_x = ov.(ff);
0846         case 'CONSTR_TOL_F'
0847             opt.fmincon.tol_f = ov.(ff);
0848             opt.knitro.tol_f = ov.(ff);
0849         case 'CONSTR_MAX_IT'
0850             opt.fmincon.max_it = ov.(ff);
0851         case 'OPF_FLOW_LIM'
0852             switch ov.(ff)
0853                 case 0
0854                     opt.opf.flow_lim = 'S';   %% apparent power (MVA)
0855                 case 1
0856                     opt.opf.flow_lim = 'P';   %% real power (MW)
0857                 case 2
0858                     opt.opf.flow_lim = 'I';   %% current magnitude (MVA @ 1 p.u. voltage)
0859                 otherwise
0860                     error(errstr, ov.(ff), ff);
0861             end
0862         case 'OPF_IGNORE_ANG_LIM'
0863             opt.opf.ignore_angle_lim = ov.(ff);
0864         case 'OPF_ALG_DC'
0865             switch ov.(ff)
0866                 case 0
0867                     opt.opf.dc.solver = 'DEFAULT';
0868                 case 100
0869                     opt.opf.dc.solver = 'BPMPD';
0870                 case {200, 250}
0871                     opt.opf.dc.solver = 'MIPS';
0872                     if ov.(ff) == 250
0873                         opt.mips.step_control = 1;
0874                     else
0875                         opt.mips.step_control = 0;
0876                     end
0877                 case 300
0878                     opt.opf.dc.solver = 'OT';     %% QUADPROG, LINPROG
0879                 case 400
0880                     opt.opf.dc.solver = 'IPOPT';
0881                 case 500
0882                     opt.opf.dc.solver = 'CPLEX';
0883                 case 600
0884                     opt.opf.dc.solver = 'MOSEK';
0885                 case 700
0886                     opt.opf.dc.solver = 'GUROBI';
0887                 otherwise
0888                     error(errstr, ov.(ff), ff);
0889             end
0890         case 'VERBOSE'
0891             opt.verbose = ov.(ff);
0892         case 'OUT_ALL'
0893             opt.out.all = ov.(ff);
0894         case 'OUT_SYS_SUM'
0895             opt.out.sys_sum = ov.(ff);
0896         case 'OUT_AREA_SUM'
0897             opt.out.area_sum = ov.(ff);
0898         case 'OUT_BUS'
0899             opt.out.bus = ov.(ff);
0900         case 'OUT_BRANCH'
0901             opt.out.branch = ov.(ff);
0902         case 'OUT_GEN'
0903             opt.out.gen = ov.(ff);
0904         case 'OUT_ALL_LIM'
0905             opt.out.lim.all = ov.(ff);
0906         case 'OUT_V_LIM'
0907             opt.out.lim.v = ov.(ff);
0908         case 'OUT_LINE_LIM'
0909             opt.out.lim.line = ov.(ff);
0910         case 'OUT_PG_LIM'
0911             opt.out.lim.pg = ov.(ff);
0912         case 'OUT_QG_LIM'
0913             opt.out.lim.qg = ov.(ff);
0914         case 'OUT_FORCE'
0915             opt.out.force = ov.(ff);
0916         case 'RETURN_RAW_DER'
0917             opt.opf.return_raw_der = ov.(ff);
0918         case 'FMC_ALG'
0919             opt.fmincon.alg = ov.(ff);
0920         case 'KNITRO_OPT'
0921             opt.knitro.opt = ov.(ff);
0922         case 'IPOPT_OPT'
0923             opt.ipopt.opt = ov.(ff);
0924         case 'MNS_FEASTOL'
0925             opt.minopf.feastol = ov.(ff);
0926         case 'MNS_ROWTOL'
0927             opt.minopf.rowtol = ov.(ff);
0928         case 'MNS_XTOL'
0929             opt.minopf.xtol = ov.(ff);
0930         case 'MNS_MAJDAMP'
0931             opt.minopf.majdamp = ov.(ff);
0932         case 'MNS_MINDAMP'
0933             opt.minopf.mindamp = ov.(ff);
0934         case 'MNS_PENALTY_PARM'
0935             opt.minopf.penalty = ov.(ff);
0936         case 'MNS_MAJOR_IT'
0937             opt.minopf.major_it = ov.(ff);
0938         case 'MNS_MINOR_IT'
0939             opt.minopf.minor_it = ov.(ff);
0940         case 'MNS_MAX_IT'
0941             opt.minopf.max_it = ov.(ff);
0942         case 'MNS_VERBOSITY'
0943             opt.minopf.verbosity = ov.(ff);
0944         case 'MNS_CORE'
0945             opt.minopf.core = ov.(ff);
0946         case 'MNS_SUPBASIC_LIM'
0947             opt.minopf.supbasic_lim = ov.(ff);
0948         case 'MNS_MULT_PRICE'
0949             opt.minopf.mult_price = ov.(ff);
0950         case 'FORCE_PC_EQ_P0'
0951             opt.sopf.force_Pc_eq_P0 = ov.(ff);
0952         case 'PDIPM_FEASTOL'
0953             opt.mips.feastol = ov.(ff);
0954             opt.pdipm.feastol = ov.(ff);
0955         case 'PDIPM_GRADTOL'
0956             opt.mips.gradtol = ov.(ff);
0957             opt.pdipm.gradtol = ov.(ff);
0958         case 'PDIPM_COMPTOL'
0959             opt.mips.comptol = ov.(ff);
0960             opt.pdipm.comptol = ov.(ff);
0961         case 'PDIPM_COSTTOL'
0962             opt.mips.costtol = ov.(ff);
0963             opt.pdipm.costtol = ov.(ff);
0964         case 'PDIPM_MAX_IT'
0965             opt.mips.max_it = ov.(ff);
0966             opt.pdipm.max_it = ov.(ff);
0967         case 'SCPDIPM_RED_IT'
0968             opt.mips.sc.red_it = ov.(ff);
0969             opt.pdipm.sc.red_it = ov.(ff);
0970         case 'TRALM_FEASTOL'
0971             opt.tralm.feastol = ov.(ff);
0972         case 'TRALM_PRIMETOL'
0973             opt.tralm.primaltol = ov.(ff);
0974         case 'TRALM_DUALTOL'
0975             opt.tralm.dualtol = ov.(ff);
0976         case 'TRALM_COSTTOL'
0977             opt.tralm.costtol = ov.(ff);
0978         case 'TRALM_MAJOR_IT'
0979             opt.tralm.major_it = ov.(ff);
0980         case 'TRALM_MINOR_IT'
0981             opt.tralm.minor_it = ov.(ff);
0982         case 'SMOOTHING_RATIO'
0983             opt.pdipm.sc.smooth_ratio = ov.(ff);
0984             opt.tralm.smooth_ratio = ov.(ff);
0985         case 'CPLEX_LPMETHOD'
0986             opt.cplex.lpmethod = ov.(ff);
0987         case 'CPLEX_QPMETHOD'
0988             opt.cplex.qpmethod = ov.(ff);
0989         case 'CPLEX_OPT'
0990             opt.cplex.opt = ov.(ff);
0991         case 'MOSEK_LP_ALG'
0992             opt.mosek.lp_alg = ov.(ff);
0993         case 'MOSEK_MAX_IT'
0994             opt.mosek.max_it = ov.(ff);
0995         case 'MOSEK_GAP_TOL'
0996             opt.mosek.gap_tol = ov.(ff);
0997         case 'MOSEK_MAX_TIME'
0998             opt.mosek.max_time = ov.(ff);
0999         case 'MOSEK_NUM_THREADS'
1000             opt.mosek.num_threads = ov.(ff);
1001         case 'MOSEK_OPT'
1002             opt.mosek.opt = ov.(ff);
1003         case 'GRB_METHOD'
1004             opt.gurobi.method = ov.(ff);
1005         case 'GRB_TIMELIMIT'
1006             opt.gurobi.timelimit = ov.(ff);
1007         case 'GRB_THREADS'
1008             opt.gurobi.threads = ov.(ff);
1009         case 'GRB_OPT'
1010             opt.gurobi.opt = ov.(ff);
1011         otherwise
1012             error('mpoption: ''%s'' is not a valid old-style option name', ff);
1013     end
1014 end
1015 % ov
1016 
1017 
1018 %%-------------------------------------------------------------------
1019 function opt_s = mpoption_v2s(opt_v)
1020 if DEBUG, fprintf('mpoption_v2s()\n'); end
1021 opt_s = mpoption_default();
1022 errstr = 'mpoption: %g is not a valid value for the old-style ''%s'' option';
1023 switch opt_v(1)                                 %% PF_ALG
1024     case 1
1025         opt_s.pf.alg = 'NR';        %% Newton's method
1026     case 2
1027         opt_s.pf.alg = 'FDXB';      %% fast-decoupled (XB version)
1028     case 3
1029         opt_s.pf.alg = 'FDBX';      %% fast-decoupled (BX version)
1030     case 4
1031         opt_s.pf.alg = 'GS';        %% Gauss-Seidel
1032     otherwise
1033         error(errstr, opt_v(1), 'PF_ALG');
1034 end
1035 opt_s.pf.tol                = opt_v(2);         %% PF_TOL
1036 opt_s.pf.nr.max_it          = opt_v(3);         %% PF_MAX_IT
1037 opt_s.pf.fd.max_it          = opt_v(4);         %% PF_MAX_IT_FD
1038 opt_s.pf.gs.max_it          = opt_v(5);         %% PF_MAX_IT_GS
1039 opt_s.pf.enforce_q_lims     = opt_v(6);         %% ENFORCE_Q_LIMS
1040 switch opt_v(10)                                %% PF_DC
1041     case 0
1042         opt_s.model = 'AC';
1043     case 1
1044         opt_s.model = 'DC';
1045     otherwise
1046         error(errstr, opt_v(10), 'PF_DC');
1047 end
1048 switch opt_v(11)                                %% OPF_ALG
1049     case 0
1050         opt_s.opf.ac.solver = 'DEFAULT';
1051     case 500
1052         opt_s.opf.ac.solver = 'MINOPF';
1053     case 520
1054         opt_s.opf.ac.solver = 'FMINCON';
1055     case {540, 545}
1056         opt_s.opf.ac.solver = 'PDIPM';
1057     case 550
1058         opt_s.opf.ac.solver = 'TRALM';
1059     case {560, 565}
1060         opt_s.opf.ac.solver = 'MIPS';
1061     case 580
1062         opt_s.opf.ac.solver = 'IPOPT';
1063     case 600
1064         opt_s.opf.ac.solver = 'KNITRO';
1065     otherwise
1066         error(errstr, opt_v(11), 'OPF_ALG');
1067 end
1068 opt_s.opf.violation         = opt_v(16);        %% OPF_VIOLATION
1069 
1070 opt_s.fmincon.tol_x         = opt_v(17);        %% CONSTR_TOL_X
1071 opt_s.fmincon.tol_f         = opt_v(18);        %% CONSTR_TOL_F
1072 opt_s.fmincon.max_it        = opt_v(19);        %% CONSTR_MAX_IT
1073 
1074 opt_s.knitro.tol_x          = opt_v(17);        %% CONSTR_TOL_X
1075 opt_s.knitro.tol_f          = opt_v(18);        %% CONSTR_TOL_F
1076 
1077 switch opt_v(24)                                %% OPF_FLOW_LIM
1078     case 0
1079         opt_s.opf.flow_lim = 'S';   %% apparent power (MVA)
1080     case 1
1081         opt_s.opf.flow_lim = 'P';   %% real power (MW)
1082     case 2
1083         opt_s.opf.flow_lim = 'I';   %% current magnitude (MVA @ 1 p.u. voltage)
1084     otherwise
1085         error(errstr, opt_v(10), 'PF_DC');
1086 end
1087 
1088 opt_s.opf.ignore_angle_lim  = opt_v(25);        %% OPF_IGNORE_ANG_LIM
1089 
1090 switch opt_v(26)                                %% OPF_ALG_DC
1091     case 0
1092         opt_s.opf.dc.solver = 'DEFAULT';
1093     case 100
1094         opt_s.opf.dc.solver = 'BPMPD';
1095     case {200, 250}
1096         opt_s.opf.dc.solver = 'MIPS';
1097     case 300
1098         opt_s.opf.dc.solver = 'OT';     %% QUADPROG, LINPROG
1099     case 400
1100         opt_s.opf.dc.solver = 'IPOPT';
1101     case 500
1102         opt_s.opf.dc.solver = 'CPLEX';
1103     case 600
1104         opt_s.opf.dc.solver = 'MOSEK';
1105     case 700
1106         opt_s.opf.dc.solver = 'GUROBI';
1107     otherwise
1108         error(errstr, opt_v(26), 'OPF_ALG_DC');
1109 end
1110 
1111 opt_s.verbose               = opt_v(31);        %% VERBOSE
1112 opt_s.out.all               = opt_v(32);        %% OUT_ALL
1113 opt_s.out.sys_sum           = opt_v(33);        %% OUT_SYS_SUM
1114 opt_s.out.area_sum          = opt_v(34);        %% OUT_AREA_SUM
1115 opt_s.out.bus               = opt_v(35);        %% OUT_BUS
1116 opt_s.out.branch            = opt_v(36);        %% OUT_BRANCH
1117 opt_s.out.gen               = opt_v(37);        %% OUT_GEN
1118 opt_s.out.lim.all           = opt_v(38);        %% OUT_ALL_LIM
1119 opt_s.out.lim.v             = opt_v(39);        %% OUT_V_LIM
1120 opt_s.out.lim.line          = opt_v(40);        %% OUT_LINE_LIM
1121 opt_s.out.lim.pg            = opt_v(41);        %% OUT_PG_LIM
1122 opt_s.out.lim.qg            = opt_v(42);        %% OUT_QG_LIM
1123 opt_s.out.force             = opt_v(44);        %% OUT_FORCE
1124 
1125 opt_s.opf.return_raw_der    = opt_v(52);        %% RETURN_RAW_DER
1126 
1127 opt_s.fmincon.alg           = opt_v(55);        %% FMC_ALG
1128 opt_s.knitro.opt            = opt_v(58);        %% KNITRO_OPT
1129 opt_s.ipopt.opt             = opt_v(60);        %% IPOPT_OPT
1130 
1131 opt_s.minopf.feastol        = opt_v(61);        %% MNS_FEASTOL
1132 opt_s.minopf.rowtol         = opt_v(62);        %% MNS_ROWTOL
1133 opt_s.minopf.xtol           = opt_v(63);        %% MNS_XTOL
1134 opt_s.minopf.majdamp        = opt_v(64);        %% MNS_MAJDAMP
1135 opt_s.minopf.mindamp        = opt_v(65);        %% MNS_MINDAMP
1136 opt_s.minopf.penalty        = opt_v(66);        %% MNS_PENALTY_PARM
1137 opt_s.minopf.major_it       = opt_v(67);        %% MNS_MAJOR_IT
1138 opt_s.minopf.minor_it       = opt_v(68);        %% MNS_MINOR_IT
1139 opt_s.minopf.max_it         = opt_v(69);        %% MNS_MAX_IT
1140 opt_s.minopf.verbosity      = opt_v(70);        %% MNS_VERBOSITY
1141 opt_s.minopf.core           = opt_v(71);        %% MNS_CORE
1142 opt_s.minopf.supbasic_lim   = opt_v(72);        %% MNS_SUPBASIC_LIM
1143 opt_s.minopf.mult_price     = opt_v(73);        %% MNS_MULT_PRICE
1144 
1145 opt_s.sopf.force_Pc_eq_P0   = opt_v(80);        %% FORCE_PC_EQ_P0, for c3sopf
1146 
1147 if (opt_v(11) == 565 && opt_v(10) == 0) || (opt_v(26) == 250 && opt_v(10) == 1)
1148     opt_s.mips.step_control = 1;
1149 end
1150 opt_s.mips.feastol          = opt_v(81);        %% PDIPM_FEASTOL
1151 opt_s.mips.gradtol          = opt_v(82);        %% PDIPM_GRADTOL
1152 opt_s.mips.comptol          = opt_v(83);        %% PDIPM_COMPTOL
1153 opt_s.mips.costtol          = opt_v(84);        %% PDIPM_COSTTOL
1154 opt_s.mips.max_it           = opt_v(85);        %% PDIPM_MAX_IT
1155 opt_s.mips.sc.red_it        = opt_v(86);        %% SCPDIPM_RED_IT
1156 
1157 opt_s.pdipm.feastol         = opt_v(81);        %% PDIPM_FEASTOL
1158 opt_s.pdipm.gradtol         = opt_v(82);        %% PDIPM_GRADTOL
1159 opt_s.pdipm.comptol         = opt_v(83);        %% PDIPM_COMPTOL
1160 opt_s.pdipm.costtol         = opt_v(84);        %% PDIPM_COSTTOL
1161 opt_s.pdipm.max_it          = opt_v(85);        %% PDIPM_MAX_IT
1162 opt_s.pdipm.sc.red_it       = opt_v(86);        %% SCPDIPM_RED_IT
1163 opt_s.pdipm.sc.smooth_ratio = opt_v(93);        %% SMOOTHING_RATIO
1164 if opt_v(11) == 545 && opt_v(10) == 0
1165     opt_s.pdipm.step_control = 1;
1166 end
1167 
1168 opt_s.tralm.feastol         = opt_v(87);        %% TRALM_FEASTOL
1169 opt_s.tralm.primaltol       = opt_v(88);        %% TRALM_PRIMETOL
1170 opt_s.tralm.dualtol         = opt_v(89);        %% TRALM_DUALTOL
1171 opt_s.tralm.costtol         = opt_v(90);        %% TRALM_COSTTOL
1172 opt_s.tralm.major_it        = opt_v(91);        %% TRALM_MAJOR_IT
1173 opt_s.tralm.minor_it        = opt_v(92);        %% TRALM_MINOR_IT
1174 opt_s.tralm.smooth_ratio    = opt_v(93);        %% SMOOTHING_RATIO
1175 
1176 opt_s.cplex.lpmethod        = opt_v(95);        %% CPLEX_LPMETHOD
1177 opt_s.cplex.qpmethod        = opt_v(96);        %% CPLEX_QPMETHOD
1178 opt_s.cplex.opt             = opt_v(97);        %% CPLEX_OPT
1179 
1180 opt_s.mosek.lp_alg          = opt_v(111);       %% MOSEK_LP_ALG
1181 opt_s.mosek.max_it          = opt_v(112);       %% MOSEK_MAX_IT
1182 opt_s.mosek.gap_tol         = opt_v(113);       %% MOSEK_GAP_TOL
1183 opt_s.mosek.max_time        = opt_v(114);       %% MOSEK_MAX_TIME
1184 opt_s.mosek.num_threads     = opt_v(115);       %% MOSEK_NUM_THREADS
1185 opt_s.mosek.opt             = opt_v(116);       %% MOSEK_OPT
1186 
1187 opt_s.gurobi.method         = opt_v(121);       %% GRB_METHOD
1188 opt_s.gurobi.timelimit      = opt_v(122);       %% GRB_TIMELIMIT
1189 opt_s.gurobi.threads        = opt_v(123);       %% GRB_THREADS
1190 opt_s.gurobi.opt            = opt_v(124);       %% GRB_OPT
1191 
1192 
1193 %%-------------------------------------------------------------------
1194 function opt_v = mpoption_s2v(opt_s)
1195 if DEBUG, fprintf('mpoption_s2v()\n'); end
1196 %% PF_ALG
1197 old = mpoption_old;
1198 switch upper(opt_s.pf.alg)
1199     case 'NR'
1200         PF_ALG = 1;
1201     case 'FDXB'
1202         PF_ALG = 2;
1203     case 'FDBX'
1204         PF_ALG = 3;
1205     case 'GS'
1206         PF_ALG = 4;
1207 end
1208 
1209 %% PF_DC
1210 if strcmp(upper(opt_s.model), 'DC')
1211     PF_DC = 1;
1212 else
1213     PF_DC = 0;
1214 end
1215 
1216 %% OPF_ALG
1217 switch upper(opt_s.opf.ac.solver)
1218     case 'DEFAULT'
1219         OPF_ALG = 0;
1220     case 'MINOPF'
1221         OPF_ALG = 500;
1222     case 'FMINCON'
1223         OPF_ALG = 520;
1224     case 'PDIPM'
1225         if isfield(opt_s, 'pdipm') && opt_s.pdipm.step_control
1226             OPF_ALG = 545;
1227         else
1228             OPF_ALG = 540;
1229         end
1230     case 'TRALM'
1231         OPF_ALG = 550;
1232     case 'MIPS'
1233         if opt_s.mips.step_control
1234             OPF_ALG = 565;
1235         else
1236             OPF_ALG = 560;
1237         end
1238     case 'IPOPT'
1239         OPF_ALG = 580;
1240     case 'KNITRO'
1241         OPF_ALG = 600;
1242 end
1243 
1244 %% FMINCON, Knitro tol_x, tol_f, max_it
1245 if strcmp(upper(opt_s.opf.ac.solver), 'KNITRO') && isfield(opt_s, 'knitro')
1246     CONSTR_TOL_X = opt_s.knitro.tol_x;
1247     CONSTR_TOL_F = opt_s.knitro.tol_f;
1248 elseif isfield(opt_s, 'fmincon')
1249     CONSTR_TOL_X  = opt_s.fmincon.tol_x;
1250     CONSTR_TOL_F  = opt_s.fmincon.tol_f;
1251 else
1252     CONSTR_TOL_X = old(17);
1253     CONSTR_TOL_F = old(18);
1254 end
1255 if isfield(opt_s, 'fmincon')
1256     CONSTR_MAX_IT   = opt_s.fmincon.max_it;
1257     FMC_ALG         = opt_s.fmincon.alg;
1258 else
1259     CONSTR_MAX_IT   = old(19);
1260     FMC_ALG         = old(55);
1261 end
1262 
1263 %% OPF_FLOW_LIM
1264 switch upper(opt_s.opf.flow_lim)
1265     case 'S'
1266         OPF_FLOW_LIM = 0;
1267     case {'P', '2'}
1268         OPF_FLOW_LIM = 1;
1269     case 'I'
1270         OPF_FLOW_LIM = 2;
1271 end
1272 
1273 %% OPF_ALG_DC
1274 switch upper(opt_s.opf.dc.solver)
1275     case 'DEFAULT'
1276         OPF_ALG_DC = 0;
1277     case 'BPMPD'
1278         OPF_ALG_DC = 100;
1279     case 'MIPS'
1280         if opt_s.mips.step_control
1281             OPF_ALG_DC = 250;
1282         else
1283             OPF_ALG_DC = 200;
1284         end
1285     case 'OT'
1286         OPF_ALG_DC = 300;
1287     case 'IPOPT'
1288         OPF_ALG_DC = 400;
1289     case 'CPLEX'
1290         OPF_ALG_DC = 500;
1291     case 'MOSEK'
1292         OPF_ALG_DC = 600;
1293     case 'GUROBI'
1294         OPF_ALG_DC = 700;
1295 end
1296 
1297 %% KNITRO_OPT
1298 if isfield(opt_s, 'knitro')
1299     KNITRO_OPT  = opt_s.knitro.opt;
1300 else
1301     KNITRO_OPT  = old(58);
1302 end
1303 
1304 %% IPOPT_OPT
1305 if isfield(opt_s, 'ipopt')
1306     IPOPT_OPT  = opt_s.ipopt.opt;
1307 else
1308     IPOPT_OPT  = old(58);
1309 end
1310 
1311 %% MINOPF options
1312 if isfield(opt_s, 'minopf')
1313     MINOPF_OPTS = [
1314         opt_s.minopf.feastol;   %% 61 - MNS_FEASTOL
1315         opt_s.minopf.rowtol;    %% 62 - MNS_ROWTOL
1316         opt_s.minopf.xtol;      %% 63 - MNS_XTOL
1317         opt_s.minopf.majdamp;   %% 64 - MNS_MAJDAMP
1318         opt_s.minopf.mindamp;   %% 65 - MNS_MINDAMP
1319         opt_s.minopf.penalty;   %% 66 - MNS_PENALTY_PARM
1320         opt_s.minopf.major_it;  %% 67 - MNS_MAJOR_IT
1321         opt_s.minopf.minor_it;  %% 68 - MNS_MINOR_IT
1322         opt_s.minopf.max_it;    %% 69 - MNS_MAX_IT
1323         opt_s.minopf.verbosity; %% 70 - MNS_VERBOSITY
1324         opt_s.minopf.core;      %% 71 - MNS_CORE
1325         opt_s.minopf.supbasic_lim;  %% 72 - MNS_SUPBASIC_LIM
1326         opt_s.minopf.mult_price;%% 73 - MNS_MULT_PRICE
1327     ];
1328 else
1329     MINOPF_OPTS = old(61:73);
1330 end
1331 
1332 %% FORCE_PC_EQ_P0
1333 if isfield(opt_s, 'sopf') && isfield(opt_s.sopf, 'force_Pc_eq_P0')
1334     FORCE_PC_EQ_P0 = opt_s.sopf.force_Pc_eq_P0;
1335 else
1336     FORCE_PC_EQ_P0 = 0;
1337 end
1338 
1339 %% PDIPM options
1340 if isfield(opt_s, 'pdipm')
1341     PDIPM_OPTS = [
1342         opt_s.pdipm.feastol;    %% 81 - PDIPM_FEASTOL
1343         opt_s.pdipm.gradtol;    %% 82 - PDIPM_GRADTOL
1344         opt_s.pdipm.comptol;    %% 83 - PDIPM_COMPTOL
1345         opt_s.pdipm.costtol;    %% 84 - PDIPM_COSTTOL
1346         opt_s.pdipm.max_it;     %% 85 - PDIPM_MAX_IT
1347         opt_s.pdipm.sc.red_it;  %% 86 - SCPDIPM_RED_IT
1348     ];
1349 else
1350     PDIPM_OPTS = old(81:86);
1351 end
1352 
1353 %% TRALM options
1354 if isfield(opt_s, 'tralm')
1355     TRALM_OPTS = [
1356         opt_s.tralm.feastol;    %% 87 - TRALM_FEASTOL
1357         opt_s.tralm.primaltol;  %% 88 - TRALM_PRIMETOL
1358         opt_s.tralm.dualtol;    %% 89 - TRALM_DUALTOL
1359         opt_s.tralm.costtol;    %% 90 - TRALM_COSTTOL
1360         opt_s.tralm.major_it;   %% 91 - TRALM_MAJOR_IT
1361         opt_s.tralm.minor_it;   %% 92 - TRALM_MINOR_IT
1362     ];
1363 else
1364     TRALM_OPTS = old(87:92);
1365 end
1366 
1367 %% SMOOTHING_RATIO
1368 if strcmp(upper(opt_s.opf.ac.solver), 'TRALM') && isfield(opt_s, 'tralm')
1369     SMOOTHING_RATIO = opt_s.tralm.smooth_ratio;
1370 elseif isfield(opt_s, 'pdipm')
1371     SMOOTHING_RATIO = opt_s.pdipm.sc.smooth_ratio;
1372 else
1373     SMOOTHING_RATIO = old(93);
1374 end
1375 
1376 %% CPLEX options
1377 if isfield(opt_s, 'cplex')
1378     CPLEX_OPTS = [
1379         opt_s.cplex.lpmethod;   %% 95 - CPLEX_LPMETHOD
1380         opt_s.cplex.qpmethod;   %% 96 - CPLEX_QPMETHOD
1381         opt_s.cplex.opt;        %% 97 - CPLEX_OPT
1382     ];
1383 else
1384     CPLEX_OPTS = old(95:97);
1385 end
1386 
1387 %% MOSEK options
1388 if isfield(opt_s, 'mosek')
1389     MOSEK_OPTS = [
1390         opt_s.mosek.lp_alg;     %% 111 - MOSEK_LP_ALG
1391         opt_s.mosek.max_it;     %% 112 - MOSEK_MAX_IT
1392         opt_s.mosek.gap_tol;    %% 113 - MOSEK_GAP_TOL
1393         opt_s.mosek.max_time;   %% 114 - MOSEK_MAX_TIME
1394         opt_s.mosek.num_threads;%% 115 - MOSEK_NUM_THREADS
1395         opt_s.mosek.opt;        %% 116 - MOSEK_OPT
1396     ];
1397 else
1398     MOSEK_OPTS = old(111:116);
1399 end
1400 
1401 %% Gurobi options
1402 if isfield(opt_s, 'gurobi')
1403     GUROBI_OPTS = [
1404         opt_s.gurobi.method;    %% 121 - GRB_METHOD
1405         opt_s.gurobi.timelimit; %% 122 - GRB_TIMELIMIT
1406         opt_s.gurobi.threads;   %% 123 - GRB_THREADS
1407         opt_s.gurobi.opt;       %% 124 - GRB_OPT
1408     ];
1409 else
1410     GUROBI_OPTS = old(121:124);
1411 end
1412 
1413 opt_v = [
1414         %% power flow options
1415         PF_ALG;                 %% 1  - PF_ALG
1416         opt_s.pf.tol;           %% 2  - PF_TOL
1417         opt_s.pf.nr.max_it;     %% 3  - PF_MAX_IT
1418         opt_s.pf.fd.max_it;     %% 4  - PF_MAX_IT_FD
1419         opt_s.pf.gs.max_it;     %% 5  - PF_MAX_IT_GS
1420         opt_s.pf.enforce_q_lims;%% 6  - ENFORCE_Q_LIMS
1421         0;                      %% 7  - RESERVED7
1422         0;                      %% 8  - RESERVED8
1423         0;                      %% 9  - RESERVED9
1424         PF_DC;                  %% 10 - PF_DC
1425         
1426         %% OPF options
1427         OPF_ALG;                %% 11 - OPF_ALG
1428         0;                      %% 12 - RESERVED12 (was OPF_ALG_POLY = 100)
1429         0;                      %% 13 - RESERVED13 (was OPF_ALG_PWL = 200)
1430         0;                      %% 14 - RESERVED14 (was OPF_POLY2PWL_PTS = 10)
1431         0;                      %% 15 - OPF_NEQ (removed)
1432         opt_s.opf.violation;    %% 16 - OPF_VIOLATION
1433         CONSTR_TOL_X;           %% 17 - CONSTR_TOL_X
1434         CONSTR_TOL_F;           %% 18 - CONSTR_TOL_F
1435         CONSTR_MAX_IT;          %% 19 - CONSTR_MAX_IT
1436         old(20);                %% 20 - LPC_TOL_GRAD (removed)
1437         old(21);                %% 21 - LPC_TOL_X (removed)
1438         old(22);                %% 22 - LPC_MAX_IT (removed)
1439         old(23);                %% 23 - LPC_MAX_RESTART (removed)
1440         OPF_FLOW_LIM;           %% 24 - OPF_FLOW_LIM
1441         opt_s.opf.ignore_angle_lim; %% 25 - OPF_IGNORE_ANG_LIM
1442         OPF_ALG_DC;             %% 26 - OPF_ALG_DC
1443         0;                      %% 27 - RESERVED27
1444         0;                      %% 28 - RESERVED28
1445         0;                      %% 29 - RESERVED29
1446         0;                      %% 30 - RESERVED30
1447         
1448         %% output options
1449         opt_s.verbose;          %% 31 - VERBOSE
1450         opt_s.out.all;          %% 32 - OUT_ALL
1451         opt_s.out.sys_sum;      %% 33 - OUT_SYS_SUM
1452         opt_s.out.area_sum;     %% 34 - OUT_AREA_SUM
1453         opt_s.out.bus;          %% 35 - OUT_BUS
1454         opt_s.out.branch;       %% 36 - OUT_BRANCH
1455         opt_s.out.gen;          %% 37 - OUT_GEN
1456         opt_s.out.lim.all;      %% 38 - OUT_ALL_LIM
1457         opt_s.out.lim.v;        %% 39 - OUT_V_LIM
1458         opt_s.out.lim.line;     %% 40 - OUT_LINE_LIM
1459         opt_s.out.lim.pg;       %% 41 - OUT_PG_LIM
1460         opt_s.out.lim.qg;       %% 42 - OUT_QG_LIM
1461         0;                      %% 43 - RESERVED43 (was OUT_RAW)
1462         opt_s.out.force;        %% 44 - OUT_FORCE
1463         0;                      %% 45 - RESERVED45
1464         0;                      %% 46 - RESERVED46
1465         0;                      %% 47 - RESERVED47
1466         0;                      %% 48 - RESERVED48
1467         0;                      %% 49 - RESERVED49
1468         0;                      %% 50 - RESERVED50
1469         
1470         %% other options
1471         old(51);                %% 51 - SPARSE_QP (removed)
1472         opt_s.opf.return_raw_der;   %% 52 - RETURN_RAW_DER
1473         0;                      %% 53 - RESERVED53
1474         0;                      %% 54 - RESERVED54
1475         FMC_ALG;                %% 55 - FMC_ALG
1476         0;                      %% 56 - RESERVED56
1477         0;                      %% 57 - RESERVED57
1478         KNITRO_OPT;             %% 58 - KNITRO_OPT
1479         0;                      %% 59 - RESERVED59
1480         IPOPT_OPT;              %% 60 - IPOPT_OPT
1481         
1482         %% MINOPF options
1483         MINOPF_OPTS;            %% 61-73 - MNS_FEASTOL-MNS_MULT_PRICE
1484         0;                      %% 74 - RESERVED74
1485         0;                      %% 75 - RESERVED75
1486         0;                      %% 76 - RESERVED76
1487         0;                      %% 77 - RESERVED77
1488         0;                      %% 78 - RESERVED78
1489         0;                      %% 79 - RESERVED79
1490         FORCE_PC_EQ_P0;         %% 80 - FORCE_PC_EQ_P0, for c3sopf
1491         
1492         %% MIPS, PDIPM, SC-PDIPM, and TRALM options
1493         PDIPM_OPTS;             %% 81-86 - PDIPM_FEASTOL-SCPDIPM_RED_IT
1494         TRALM_OPTS;             %% 87-92 - TRALM_FEASTOL-TRALM_MINOR_IT
1495         SMOOTHING_RATIO;        %% 93 - SMOOTHING_RATIO
1496         0;                      %% 94 - RESERVED94
1497         
1498         %% CPLEX options
1499         CPLEX_OPTS;             %% 95-97 - CPLEX_LPMETHOD-CPLEX_OPT
1500         0;                      %% 98 - RESERVED98
1501         0;                      %% 99 - RESERVED99
1502         0;                      %% 100 - RESERVED100
1503         0;                      %% 101 - RESERVED101
1504         0;                      %% 102 - RESERVED102
1505         0;                      %% 103 - RESERVED103
1506         0;                      %% 104 - RESERVED104
1507         0;                      %% 105 - RESERVED105
1508         0;                      %% 106 - RESERVED106
1509         0;                      %% 107 - RESERVED107
1510         0;                      %% 108 - RESERVED108
1511         0;                      %% 109 - RESERVED109
1512         0;                      %% 110 - RESERVED110
1513 
1514         %% MOSEK options
1515         MOSEK_OPTS;             %% 111-116 - MOSEK_LP_ALG-MOSEK_OPT
1516         0;                      %% 117 - RESERVED117
1517         0;                      %% 118 - RESERVED118
1518         0;                      %% 119 - RESERVED119
1519         0;                      %% 120 - RESERVED120
1520 
1521         %% Gurobi options
1522         GUROBI_OPTS;            %% 121-124 - GRB_METHOD-GRB_OPT
1523     ];
1524 
1525 
1526 %%-------------------------------------------------------------------
1527 function optt = mpoption_default()
1528 if DEBUG, fprintf('mpoption_default()\n'); end
1529 persistent opt;             %% cache this for speed
1530 if ~isstruct(opt)
1531     opt = struct(...
1532         'v',                    mpoption_version, ...   %% version
1533         'model',                'AC', ...
1534         'pf',                   struct(...
1535             'alg',                  'NR', ...
1536             'current_balance',      0, ...
1537             'v_cartesian',          0, ...
1538             'tol',                  1e-8, ...
1539             'nr',                   struct(...
1540                 'max_it',               10, ...
1541                 'lin_solver',            '' ), ...
1542             'fd',                   struct(...
1543                 'max_it',               30  ), ...
1544             'gs',                   struct(...
1545                 'max_it',               1000    ), ...
1546             'radial',               struct(...
1547                 'max_it',               20   , ...
1548                 'vcorr',                 0  ), ...
1549             'enforce_q_lims',       0   ), ...
1550         'cpf',                  struct(...
1551             'parameterization',     3, ...
1552             'stop_at',              'NOSE', ...     %% 'NOSE', <lam val>, 'FULL'
1553             'enforce_p_lims',       0, ...
1554             'enforce_q_lims',       0, ...
1555             'enforce_v_lims',       0, ...
1556             'enforce_flow_lims',    0, ...
1557             'step',                 0.05, ...
1558             'step_min',             1e-4, ...
1559             'step_max',             0.2, ...
1560             'adapt_step',           0, ...
1561             'adapt_step_damping',   0.7, ...
1562             'adapt_step_tol',       1e-3, ...
1563             'target_lam_tol',       1e-5, ...
1564             'nose_tol',             1e-5, ...
1565             'p_lims_tol',           0.01, ...
1566             'q_lims_tol',           0.01, ...
1567             'v_lims_tol',           1e-4, ...
1568             'flow_lims_tol',        0.01, ...
1569             'plot',                 struct(...
1570                 'level',                0, ...
1571                 'bus',                  []  ), ...
1572             'user_callback',        ''  ), ...
1573         'opf',                  struct(...
1574             'ac',                   struct(...
1575                 'solver',               'DEFAULT'   ), ...
1576             'dc',                   struct(...
1577                 'solver',               'DEFAULT'   ), ...
1578             'current_balance',      0, ...
1579             'v_cartesian',          0, ...
1580             'violation',            5e-6, ...
1581             'use_vg',               0, ...
1582             'flow_lim',             'S', ...
1583             'ignore_angle_lim',     0, ...
1584             'softlims',             struct(...
1585                 'default',               1  ), ...
1586             'init_from_mpc',        -1, ...
1587             'start',                0, ...
1588             'return_raw_der',       0   ), ...
1589         'verbose',              1, ...
1590         'out',                  struct(...
1591             'all',                  -1, ...
1592             'sys_sum',              1, ...
1593             'area_sum',             0, ...
1594             'bus',                  1, ...
1595             'branch',               1, ...
1596             'gen',                  0, ...
1597             'lim',                  struct(...
1598                 'all',                  -1, ...
1599                 'v',                    1, ...
1600                 'line',                 1, ...
1601                 'pg',                   1, ...
1602                 'qg',                   1   ), ...
1603             'force',                0, ...
1604             'suppress_detail',      -1  ), ...
1605         'mips',                 struct(...  %% see mpoption_info_mips() for optional fields
1606             'step_control',         0, ...
1607             'linsolver',            '', ...
1608             'feastol',              0, ...
1609             'gradtol',              1e-6, ...
1610             'comptol',              1e-6, ...
1611             'costtol',              1e-6, ...
1612             'max_it',               150, ...
1613             'sc',                   struct(...
1614                 'red_it',               20  )), ...
1615         'exp',                  struct(... %% experimental options
1616             'sys_wide_zip_loads',   struct(...
1617                 'pw',                   [], ...
1618                 'qw',                   []  )) ...
1619     );
1620     opt_pkgs = mpoption_optional_pkgs();
1621     for k = 1:length(opt_pkgs)
1622         fname = ['mpoption_info_' opt_pkgs{k}];
1623         if exist(fname, 'file') == 2
1624             opt = nested_struct_copy(opt, feval(fname, 'D'));
1625         end
1626     end
1627 end
1628 optt = opt;
1629 
1630 %%-------------------------------------------------------------------
1631 function optt = mpoption_optional_fields()
1632 if DEBUG, fprintf('mpoption_optional_fields()\n'); end
1633 persistent opt;         %% cache this for speed
1634 if ~isstruct(opt)
1635     opt_pkgs = mpoption_optional_pkgs();
1636     opt = struct;
1637     for k = 1:length(opt_pkgs)
1638         fname = ['mpoption_info_' opt_pkgs{k}];
1639         if exist(fname, 'file') == 2
1640             opt = nested_struct_copy(opt, feval(fname, 'V'));
1641         end
1642     end
1643 end
1644 optt = opt;
1645 
1646 %% globals
1647 %%-------------------------------------------------------------------
1648 function v = mpoption_version
1649 v = 20;     %% version number of MATPOWER options struct
1650             %% (must be incremented every time structure is updated)
1651             %% v1   - first version based on struct (MATPOWER 5.0b1)
1652             %% v2   - added 'linprog' and 'quadprog' fields
1653             %% v3   - (forgot to increment v) added 'out.suppress_detail'
1654             %%        field
1655             %% v4   - (forgot to increment v) MIPS 1.1, added optional
1656             %%        fields to 'mips' options: xi, sigma, z0, alpha_min,
1657             %%        rho_min, rho_max, mu_threshold and max_stepsize
1658             %% v5   - (forgot to increment v) added 'opf.init_from_mpc'
1659             %%        field (MATPOWER 5.0)
1660             %% v6   - added 'clp' field
1661             %% v7   - added 'intlinprog' field
1662             %% v8   - MIPS 1.2, added 'linsolver' field to
1663             %%        'mips' options
1664             %% v9   - added 'exp' for experimental fields, specifically
1665             %%        'sys_wide_zip_loads.pw', 'sys_wide_zip_loads.qw'
1666             %% v10  - added 'most' field
1667             %% v11  - added 'cpf' options 'adapt_step_damping',
1668             %%        'enforce_p_lims', 'enforce_q_lims', 'target_lam_tol'
1669             %%        'nose_tol', 'p_lims_tol' and 'q_lims_tol',
1670             %%        removed option 'cpf.user_callback_args'
1671             %% v12  - added option 'opf.use_vg'
1672             %% v13  - added 'pf.radial.max_it', 'pf.radial.vcorr'
1673             %% v14  - added 'pf.nr.lin_solver'
1674             %% v15  - added 'cpf.enforce_v_lims', 'cpf.enforce_flow_lims',
1675             %%        'cpf.v_lims_tol', and 'cpf.flow_lims_tol'
1676             %% v16  - added 'opf.start' (deprecated 'opf.init_from_mpc')
1677             %% v17  - added 'knitro.maxit'
1678             %% v18  - added 'opf.current_balance' and 'opf.v_cartesian'
1679             %% v19  - added 'opf.softlims.default'
1680             %% v20  - added 'pf.current_balance' and 'pf.v_cartesian'
1681 
1682 %%-------------------------------------------------------------------
1683 function db_level = DEBUG
1684 db_level = 0;
1685 
1686 %%-------------------------------------------------------------------
1687 function pkgs = mpoption_optional_pkgs()
1688 pkgs = {...
1689     'clp', 'cplex', 'fmincon', 'gurobi', 'glpk', 'intlinprog', 'ipopt', ...
1690     'knitro', 'linprog', 'minopf', 'most', 'mosek', 'quadprog', 'sdp_pf', ...
1691     'sopf', 'tspopf', 'yalmip' ...
1692 };

Generated on Mon 24-Jun-2019 15:58:45 by m2html © 2005