Interfaces

Matlab

There are two ways to call SNOPT and SQOPT from Matlab. The SNOPT routines utilize the snoptA interface.

snopt.m

% Calling sequence 1:
[ ] = snopt(x, xlow, xupp, xmul, xstate,...
            Flow, Fupp, Fmul, Fstate, userfun,...
            [options])

% Calling sequence 2:
[ ] = snopt(x, xlow, xupp, xmul, xstate,...
            Flow, Fupp, Fmul, Fstate, userfun,...
            ObjAdd, ObjRow, [options])

% Calling sequence 3:
[ ] = snopt(x, xlow, xupp, xmul, xstate,...
            Flow, Fupp, Fmul, Fstate, userfun,...
            ObjAdd, ObjRow,
            A, G, [options])

% Output from snopt:
[x,F,info,xmul,Fmul,xstate,Fstate,output] = snopt(...)

The vectors xlow and xupp and Flow and Fupp define the lower and upper bounds on the variables x and the constraints F. ObjRow indicates the row of F defining the objective (by default, ObjRow = 1). The vectors (x, xmul, xstate) define the initial point, multipliers, and state of the variables. Similarly, (F, Fmul, Fstate) for the constraints. ObjAdd is a constant term of the objective.

If derivative information is not provided, snJac is automatically called to determine the sparsity pattern of the Jacobian. Otherwise, input argument A and G are defined to be either a struct or a matrix (dense or sparse) that defines the constant and nonlinear elements in the Jacobian of F.

If A is a struct, the structure is provided in coordinate form with fields: A.row, A.col, and A.val. If i = A.row(k) and j = A.col(k), then A.val(k) is the constant (i,j)-th element in A.

Input argument G defines the nonlinear elements in the Jacobian. If G is a struct, the Jacobian structure is provided in coordinate form with fields: G.row and G.col. If G is a matrix, then nonzero elements in the matrix denote the nonzero elements of the nonlinear elements of the Jacobian.

userfun is the handle for a Matlab funcion that defines the elements of F and optionally, their nonlinear derivatives.

snsolve.m

snsolve.m matches the call sequence of Matlab’s fmincon function. The only exception is the options structure.

% Calling sequence for snsolve:
[ ] = snsolve(obj,x0,A,b)
[ ] = snsolve(obj,x0,A,b,options)

[ ] = snsolve(obj,x0,A,b,Aeq,beq)
[ ] = snsolve(obj,x0,A,b,Aeq,beq,options)

[ ] = snsolve(obj,x0,A,b,Aeq,beq,xlow,xupp)
[ ] = snsolve(obj,x0,A,b,Aeq,beq,xlow,xupp,options)

[ ] = snsolve(obj,x0,A,b,Aeq,beq,xlow,xupp,nonlcon)
[ ] = snsolve(obj,x0,A,b,Aeq,beq,xlow,xupp,nonlcon,options)

[ ] = snsolve(obj,x0,A,b,Aeq,beq,xlow,xupp,nonlcon,lambda,states)
[ ] = snsolve(obj,x0,A,b,Aeq,beq,xlow,xupp,nonlcon,lambda,states,options)

% Output from snsolve:
[x,fval,info,output,lambda,states] = snsolve(...)

sqopt.m

% Calling sequence 1:
[] = sqopt(H, c, x0, xl, xu, A, al, au, [options])

% Calling sequence 2:
[] = sqopt(H, c, x0, xl, xu, A, al, au, states, lambda, [options])

% Output from sqopt:
[x,obj,info,output,lambda,states] = sqopt(...)

The SQOPT Matlab interface solves the quadratic problem

\[ \begin{align} \min_x \quad & q(x) = c^T x + \frac{1}{2} x^T Hx \\ \text{subject to} \quad & x_l \le x \le x_u, \quad a_l \le Ax \le a_u \\ \end{align} \]

sqsolve.m

sqsolve.m matches the call sequence of Matlab’s quadprog function. The only exception is the options structure.

[] = sqsolve(H, f)
[] = sqsolve(H, f, A, b)
[] = sqsolve(H, f, A, b, Aeq, beq)
[] = sqsolve(H, f, A, b, Aeq, beq, lb, ub)
[] = sqsolve(H, f, A, b, Aeq, beq, lb, ub, x0)
[] = sqsolve(H, f, A, b, Aeq, beq, lb, ub, x0, options)
[] = sqsolve(H, f, A, b, Aeq, beq, lb, ub, x0, lambda, states, options)

% Output from sqsolve:
[x,fval,exitflag,output,lambda] = sqsolve(H, f, ...)

Options Structure

The options structure can be used to set options for SNOPT or SQOPT by creating an entry with a field name equal to the options keyword with spaces replaced by underscores ‘_’. For example, options.iterations_limit = 250. (See SNOPT Options for more on SNOPT/SQOPT options.

Additional keywords include:

Field name

Description

options.name

is the problem name

options.start

‘Cold’, ‘Warm’

options.screen

is a string set to ‘on’ or ‘off’. Summary to the screen is controlled by this option. (default ‘on’)

options.printfile

is a string denoting the print file. By default, no print file is created. Not setting this option or setting it to ‘’ turns off print output.

options.specsfile

is a string denoting the options filename.

options.stop

is the snSTOP function called at every major iteration. (See snSTOP.)

options.iwork

is an integer defining the integer SNOPT workspace length

options.rwork

is an integer defining the real SNOPT workspace length.

C

The SNOPT C interface consists of several C subroutines that

typedef struct {
  char   *name;

  int     memCalled;
  int     initCalled;
  int     userWork;

  isnSTOP snSTOP;
  isnLog  snLog;
  isnLog2 snLog2;
  isqLog  sqLog;

  int     lenrw, leniw;
  int    *iw;
  double *rw;

  int     lenru, leniu;
  int    *iu;
  double *ru;

} snProblem;

Initializing/Finalizing

void snInitX( snProblem* prob, char* name,
              char* prtfile, int iprint, char* sumfile, int isumm );
void snInit( snProblem* prob, char* name, char* prtfile, int summOn );

void snInitXW( snProblem* prob, char* name,
               char* prtfile, int iprint, char* sumfile, int isumm,
              int *iw, int leniw, double *rw, int lenrw );
void snInitW( snProblem* prob, char* name, char* prtfile, int summOn,
              int *iw, int leniw, double *rw, int lenrw );

void deleteSNOPT( snProblem* prob );

Setting and Getting Options

void setPrintfile   ( snProblem* prob, char* prtname );
int  setSpecsfile   ( snProblem* prob, char* spcname );

void setPrintfileX  ( snProblem* prob, char* prtname, int iprint );
int  setSpecsfileX  ( snProblem* prob, char* spcname, int ispecs );

int setParameter    ( snProblem* prob, char stropt[] );
int setIntParameter ( snProblem* prob, char stropt[], int opt );
int getIntParameter ( snProblem* prob, char stropt[], int *opt );
int setRealParameter( snProblem* prob, char stropt[], double opt );
int getRealParameter( snProblem* prob, char stropt[], double *opt );

Advanced Features in C

void setUserI       ( snProblem* prob, int *iu, int leniu );
void setUserR       ( snProblem* prob, double *ru, int lenru );
void setUserspace   ( snProblem* prob, int *iu, int leniu, double *ru, int lenru );

void setLog         ( snProblem* prob, isnLog snLog, isnLog2 snLog2, isqLog sqLog );
void setSTOP        ( snProblem* prob, isnSTOP snSTOP );

void setWorkspace   ( snProblem* prob, int m, int n, int ne,
                      int negCon, int nnCon, int nnObj, int nnJac);
void setWorkspaceA  ( snProblem* prob, int nF, int n, int neA, int neG);

Solving with SNOPT in C

int snJac( snProblem* prob,
           int nF, int n, snFunA usrfun,
           double *x, double *xlow, double *xupp,
           int *neA, int *iAfun, int *jAvar, double *A,
           int *neG, int *iGfun, int *jGvar );

 int snoptA( snProblem* prob, int start,
             int nF, int n, double ObjAdd, int ObjRow,
             snFunA usrfun,
             int neA, int *iAfun, int *jAvar, double *A,
             int neG, int *iGfun, int *jGvar,
             double *xlow, double *xupp, double *Flow, double *Fupp,
             double *x, int *xstate, double *xmul,
             double *F, int *Fstate, double *Fmul,
             int* nS, int* nInf, double* sInf );

 int snoptB( snProblem* prob, int start, int m, int n, int ne,
             int nnCon, int nnObj, int nnJac, int iObj, double ObjAdd,
             snConB funcon, snObjB funobj,
             double *valJ, int *indJ, int *locJ,
             double *bl, double *bu, int *hs, double *x,
             double *pi, double *rc, double* objective,
             int* nS, int* nInf, double* sInf );

 int snoptC( snProblem* prob, int start, int m, int n, int ne,
             int nnCon, int nnObj, int nnJac, int iObj, double ObjAdd,
             snFunC usrfun,
             double *valJ, int *indJ, int *locJ,
             double *bl, double *bu, int *hs, double *x,
             double *pi, double *rc, double* objective,
             int* nS, int* nInf, double* sInf );

C++

The SNOPT C++ interface is defined by the snoptProblem class.

snoptProblem Class

Constructors for SNOPTA

SNOPTA()
SNOPTA(const std::string &name);
SNOPTA(Int aleniw, Int alenrw);
SNOPTA(Int *aiw, Int aleniw, Real *arw, Int alenrw);
SNOPTA(const std::string &name, Int aleniw, Int alenrw);
SNOPTA(const std::string &name, Int *iw, Int aleniw, Real *rw, Int alenrw);
~SNOPTA();

Constructors for SNOPTB and SNOPTC

SNOPTBC();
SNOPTBC(const std::string &name);
SNOPTBC(Int aleniw, Int alenrw);
SNOPTBC(Int *aiw, Int aleniw, Real *arw, Int alenrw);
SNOPTBC(const std::string &name, Int aleniw, Int alenrw);
SNOPTBC(const std::string &name, Int *iw, Int aleniw, Real *rw, Int alenrw);
~SNOPTBC();

Initializing/Finalizing

void initialize(const std::string &prtfile, Int summOn);
void initialize(const std::string &prtfile, Int iprInt,
                const std::string &sumfile, Int isumm);

Setting/Getting Options

void setProbName(const std::string &Prob);
void setPrintFile(const std::string &prtname);
void setPrintFile(const std::string &prtname, Int iprInt);

Int setSpecsFile(const std::string &specname);
Int setSpecsFile(const std::string &specname, Int ispecs);

Int setParameter(const std::string &stroptin);
Int setParameter(const std::string &stropt, Int opt);
Int setParameter(const std::string &stropt, Real opt);

Int getParameter(const std::string &stropt, Int &opt);
Int getParameter(const std::string &stropt, Real &opt);

Advanced Features in C++

void setUserI(Int  *iu, Int leniu);
void setUserR(Real *ru, Int lenru);
void setUserspace(Int *iu, Int leniu, Real *ru, Int lenru);

void setWorkspace(Int neF, Int n, Int neA, Int neG);

SNOPTA in C++

Int  computeJac(Int neF, Int n, snFunA usrfunA,
                Real *x, Real *xlow, Real*xupp,
                Int *&iAfun, Int *&jAvar, Real *&A, Int &neA,
                Int *&iGfun, Int *&jGvar, Int &neG);

Int  solve(Int starttype, Int nF, Int n, Real ObjAdd,
           Int ObjRow, snFunA usrfunA,
           Real *xlow, Real *xupp, Real *Flow, Real *Fupp,
           Real *x, Int *xstate, Real *xmul,
           Real *F, Int *Fstate, Real *Fmul,
           Int &nS, Int &nInf, Real &sInf);
Int  solve(Int starttype, Int nF, Int n, Real ObjAdd,
           Int ObjRow, snFunA usrfunA,
           Real *xlow, Real *xupp, Real *Flow, Real *Fupp,
           Real *x, Int *xstate, Real *xmul,
           Real *F, Int *Fstate, Real *Fmul,
           Int &nS, Int &nInf, Real &sInf,
           void* usrdat);
Int  solve(Int starttype, Int nF, Int n, Real ObjAdd,
           Int ObjRow, snFunA usrfunA,
           Int *iAfun, Int *jAvar, Real *A, Int neA,
           Int *iGfun, Int *jGvar, Int neG,
           Real *xlow, Real *xupp, Real *Flow, Real *Fupp,
           Real *x, Int *xstate, Real *xmul,
           Real *F, Int *Fstate, Real *Fmul,
           Int &nS, Int &nInf, Real &sInf);
Int  solve(Int starttype, Int nF, Int n, Real ObjAdd,
           Int ObjRow, snFunA usrfunA,
           Int *iAfun, Int *jAvar, Real *A, Int neA,
           Int *iGfun, Int *jGvar, Int neG,
           Real *xlow, Real *xupp, Real *Flow, Real *Fupp,
           Real *x, Int *xstate, Real *xmul,
           Real *F, Int *Fstate, Real *Fmul,
           Int &nS, Int &nInf, Real &sInf,
           void* usrdat);

SNOPTB/SNOPTC in C++

void setWorkspace(Int m, Int n, Int ne, Int negCon, Int nnCon, Int nnObj, Int nnJac);

// SNOPTB
Int solve(Int starttype, Int m, Int n, Int ne, Int negCon,
          Int nnCon, Int nnObj, Int nnJac,
          Int iObj, Real ObjAdd,
          snConB funcon, snObjB funobj,
          Real *Jval, Int *indJ, Int *locJ,
          Real *bl, Real *bu, Int *hs,
          Real *x, Real *pi, Real *rc,
          Int &nS, Int &nInf, Real &sInf, Real &objective);
Int solve(Int starttype, Int m, Int n, Int ne, Int negCon,
          Int nnCon, Int nnObj, Int nnJac,
          Int iObj, Real ObjAdd,
          snConB funcon, snObjB funobj,
          Real *Jval, Int *indJ, Int *locJ,
          Real *bl, Real *bu, Int *hs,
          Real *x, Real *pi, Real *rc,
          Int &nS, Int &nInf, Real &sInf, Real &objective,
          void* usrdat);

// SNOPTC
Int solve(Int starttype, Int m, Int n, Int ne, Int negCon,
          Int nnCon, Int nnObj, Int nnJac,
          Int iObj, Real ObjAdd, snFunC usrfunC,
          Real *Jval, Int *indJ, Int *locJ,
          Real *bl, Real *bu, Int *hs,
          Real *x, Real *pi, Real *rc,
          Int &nS, Int &nInf, Real &sInf, Real &objective);

Int solve(Int starttype, Int m, Int n, Int ne, Int negCon,
          Int nnCon, Int nnObj, Int nnJac,
          Int iObj, Real ObjAdd, snFunC usrfunC,
          Real *Jval, Int *indJ, Int *locJ,
          Real *bl, Real *bu, Int *hs,
          Real *x, Real *pi, Real *rc,
          Int &nS, Int &nInf, Real &sInf, Real &objective,
          void* usrdat);