5.12. pyopus.problems.cuteritf — CUTEr problem binary interace source code

C source code of the CUTEr interface

This module is independent of PyOPUS. It can be taken as is and used as a module in some other package.

pyopus.problems.cuteritf.itf_c_source = '\n/* CUTEr2 interface to Python and NumPy */\n/* (c)2011 Arpad Buermen */\n/* Licensed under GPL V3 */\n\n/* Note that in Windows we do not use Debug compile because we don\'t have the debug version\n of Python libraries and interpreter. We use Release version instead where optimizations\n are disabled. Such a Release version can be debugged. \n */\n \n/* Unused CUTEr tools - sparse finite element matrices and banded matrices\n cdimse\n\t ceh\n\t csgreh\n\t ubandh\n\t udimse\n\t ueh\n ugreh\n\t \n CUTEr tools that are not used because they duplicate functionality or are obsolete\n cnames\t\t... used pbname, varnames, connames\n\t udimen\t\t... used cdimen \n\t ufn\t\t... used uofg\n\t ugr\t\t... used uofg\n\t unames\t\t... used pbname, varnames\n cscfg\t\t... obsolete\n\t cscifg\t\t... obsolete\n*/\n\n#define NPY_NO_DEPRECATED_API NPY_1_7_API_VERSION\n \n#include "Python.h"\n#include "cuter.h"\n#include "numpy/arrayobject.h"\n#include <math.h>\n#include <stdio.h>\n\n/* Debug switch - uncomment to enable debug messages */\n/* #undefine PYDEBUG */\n\n/* Debug file */\n#define df stdout\n\n#ifndef WIN32\n#define __declspec(a) \n#endif\n\n/* Safeguard against C++ symbol mangling */\n#ifdef __cplusplus\nextern "C" {\n#endif\n\n/* NumPy type for CUTEr integer \n INTEGER type in FORTRAN is by default 4 bytes wide (tested with gfortran \n on AMD64 Linux). The C prototypes of FORTRAN-based CUTEr toos use \n integer C type which is typedefed to long int in cuter.h. \n On IA32 Linux long int is 4 bytes wide. Unfortunately on AMD64 it is 8 \n bytes wide. This causes problems if we are passing arguments to\n FORTRAN functions. A serious problem appears when an input or output \n argument is an array because the elements of the array are twice the size\n the FORTRAN expects them to be. \n Until this bug is fixed in CUTEr the solution is to use npy_integer for\n variables and npy_integer_type_num for creating NumPy arrays. To get rid \n of warnings integer arguments are explicitly converted to (integer *) when \n they are passed to FORTRAN functions. This is harmless since FORTRAN \n arguments are passed by reference (as pointers) and the type conversion does \n not affect the pointer\'s value. \n The same goes for logical typedef of CUTEr. This interface uses npy_logical \n which is actually int, but type casts npy_logical* pointers to logical*.\n*/\nstatic int npy_integer_type_num=NPY_INT;\ntypedef int npy_integer; \ntypedef int npy_logical; \n\n\n/* Prototypes */\nstatic PyObject *cuter__dims(PyObject *self, PyObject *args);\nstatic PyObject *cuter__setup(PyObject *self, PyObject *args);\nstatic PyObject *cuter__varnames(PyObject *self, PyObject *args);\nstatic PyObject *cuter__connames(PyObject *self, PyObject *args);\nstatic PyObject *cuter_objcons(PyObject *self, PyObject *args);\nstatic PyObject *cuter_obj(PyObject *self, PyObject *args);\nstatic PyObject *cuter_cons(PyObject *self, PyObject *args);\nstatic PyObject *cuter_lagjac(PyObject *self, PyObject *args);\nstatic PyObject *cuter_jprod(PyObject *self, PyObject *args);\nstatic PyObject *cuter_hess(PyObject *self, PyObject *args);\nstatic PyObject *cuter_ihess(PyObject *self, PyObject *args);\nstatic PyObject *cuter_hprod(PyObject *self, PyObject *args);\nstatic PyObject *cuter_gradhess(PyObject *self, PyObject *args);\nstatic PyObject *cuter__scons(PyObject *self, PyObject *args); \nstatic PyObject *cuter__slagjac(PyObject *self, PyObject *args); \nstatic PyObject *cuter__sphess(PyObject *self, PyObject *args); \nstatic PyObject *cuter__isphess(PyObject *self, PyObject *args); \nstatic PyObject *cuter__gradsphess(PyObject *self, PyObject *args); \nstatic PyObject *cuter_report(PyObject *self, PyObject *args); \n\n/* Persistent data */\n#define STR_LEN 10\nstatic npy_integer CUTEr_nvar = 0;\t\t/* number of variables */\nstatic npy_integer CUTEr_ncon = 0;\t\t/* number of constraints */\nstatic npy_integer CUTEr_nnzj = 0;\t\t/* nnz in Jacobian */\nstatic npy_integer CUTEr_nnzh = 0;\t\t/* nnz in upper triangular Hessian */\nstatic char CUTEr_probName[STR_LEN+1];\t/* problem name */\nstatic char setupCalled = 0;\t\t\t/* Flag to indicate if setup was called */\nstatic char dataFileOpen = 0;\t\t\t/* Flag to indicate if OUTSDIf is open */\n\nstatic npy_integer funit = 42;\t\t\t/* FORTRAN unit number for OUTSDIF.d */\nstatic npy_integer iout = 6;\t\t\t/* FORTRAN unit number for error output */\nstatic char fName[] = "OUTSDIF.d"; \t/* Data file name */\n/* Logical constants for FORTRAN calls */\nstatic logical somethingFalse = FALSE_, somethingTrue = TRUE_;\n\n\n/* Helper functions */\n\n/* Open data file, return 0 on error. */\nint open_datafile(void) {\n\tnpy_integer ioErr;\t\t\t\t\t/* Exit flag from OPEN and CLOSE */\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: opening data file\\n");\n#endif\n\tioErr = 0;\n\tif (! dataFileOpen) \n\t\tFORTRAN_OPEN((integer *)&funit, fName, (integer *)&ioErr);\n\tif (ioErr) {\n\t\tPyErr_SetString(PyExc_Exception, "Failed to open data file");\n\t\treturn 0;\n\t}\n\tdataFileOpen = 1;\n\treturn 1;\n}\n\n/* Close data file, return 0 on error. */\nint close_datafile(void) {\n\tnpy_integer ioErr;\t\t\t\t\t/* Exit flag from OPEN and CLOSE */\n\tioErr = 0;\n\tFORTRAN_CLOSE((integer *)&funit, (integer *)&ioErr);\n\tif (ioErr) {\n\t\tPyErr_SetString(PyExc_Exception, "Error closing data file");\n\t\treturn 0;\n\t}\n\tdataFileOpen = 0;\n\treturn 1; \n}\n\n/* Check if the problem is set up, return 0 if it is not. */\nint check_setup(void) {\n\tif (!setupCalled) {\n\t\tPyErr_SetString(PyExc_Exception, "Problem is not set up");\n\t\treturn 0;\n\t}\n\treturn 1;\n}\n\n/* Trim trailing spaces from a string starting at index n. */\nvoid trim_string(char *s, int n) {\n\tint i;\n\t\n\tfor(i=n;i>=0;i--) {\n\t\tif (s[i]!=\' \')\n\t\t\tbreak;\n\t}\n\ts[i+1]=0;\n}\n\n/* Decrese reference counf for newly created dictionary members */\nPyObject *decRefDict(PyObject *dict) {\n\tPyObject *key, *value;\n\tPy_ssize_t pos; \n\tpos=0;\n\twhile (PyDict_Next(dict, &pos, &key, &value)) {\n\t\tPy_XDECREF(value);\n\t}\n\treturn dict;\n}\n\n/* Decrease reference count for newly created tuple members */\nPyObject *decRefTuple(PyObject *tuple) {\n\tPy_ssize_t pos; \n\tfor(pos=0;pos<PyTuple_Size(tuple);pos++) {\n\t\tPy_XDECREF(PyTuple_GetItem(tuple, pos));\n\t}\n\treturn tuple;\n}\n\n/* Extract sparse gradient and Jacobian in form of NumPy arrays */\nvoid extract_sparse_gradient_jacobian(npy_integer nnzjplusno, npy_integer *sji, npy_integer *sjfi, npy_double *sjv, \n\t\tPyArrayObject **Mgi, PyArrayObject **Mgv, PyArrayObject **MJi, PyArrayObject **MJfi, PyArrayObject **MJv) {\n\tnpy_double *gv, *Jv;\n\tnpy_integer *gi, *Ji, *Jfi, nnzg, i, jg, jj; \n\tnpy_intp dims[1]; \n\t\n\t/* Get number of nonzeros in gradient vector */\n\tnnzg=0;\n\tfor(i=0;i<nnzjplusno;i++) {\n\t\tif (sjfi[i]==0)\n\t\t\tnnzg++;\n\t}\n\t\n\t/* Alocate and fill objective/Lagrangian gradient data and Jacobian data, \n\t convert indices from FORTRAN to C. */\n\tdims[0]=nnzg; \n\t*Mgi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t*Mgv=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\tgi=(npy_integer *)PyArray_DATA(*Mgi);\n\tgv=(npy_double *)PyArray_DATA(*Mgv);\n\tdims[0]=nnzjplusno-nnzg; \n\t*MJi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t*MJfi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t*MJv=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\tJi=(npy_integer *)PyArray_DATA(*MJi);\n\tJfi=(npy_integer *)PyArray_DATA(*MJfi);\n\tJv=(npy_double *)PyArray_DATA(*MJv);\n\tjg=0;\n\tjj=0;\n\tfor(i=0;i<nnzjplusno;i++) {\n\t\tif (sjfi[i]==0) {\n\t\t\tgi[jg]=sji[i]-1;\n\t\t\tgv[jg]=sjv[i];\n\t\t\tjg++;\n\t\t} else {\n\t\t\tJi[jj]=sji[i]-1;\n\t\t\tJfi[jj]=sjfi[i]-1;\n\t\t\tJv[jj]=sjv[i];\n\t\t\tjj++;\n\t\t}\n\t}\n}\n\n/* Extract sparse Hessian in form of NumPy arrays \n from sparse ijv format of Hessian\'s upper triangle + diagonal. \n Add elements to lower triangle. */\nvoid extract_sparse_hessian(npy_integer nnzho, npy_integer *si, npy_integer *sj, npy_double *sv, \n\t\t\t\t\tPyArrayObject **MHi, PyArrayObject **MHj, PyArrayObject **MHv) {\n\tnpy_integer *Hi, *Hj, nnzdiag, i, j;\n\tnpy_double *Hv; \n\tnpy_intp dims[1];\n\t\n\t\n\t/* Get number of nonzeros on the diagonal */\n\tnnzdiag=0;\n\tfor(i=0;i<nnzho;i++) {\n\t\tif (si[i]==sj[i])\n\t\t\tnnzdiag++;\n\t}\n\t\n\t/* Alocate and fill objective/Lagrangian gradient data and Jacobian data, \n\t convert indices from FORTRAN to C, fill lower triangle. */\n\tdims[0]=2*nnzho-nnzdiag; /* Do not duplicate diagonal elements */\n\t*MHi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t*MHj=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t*MHv=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\tHi=(npy_integer *)PyArray_DATA(*MHi);\n\tHj=(npy_integer *)PyArray_DATA(*MHj);\n\tHv=(npy_double *)PyArray_DATA(*MHv);\n\tj=0;\n\tfor(i=0;i<nnzho;i++) {\n\t\tHi[j]=si[i]-1;\n\t\tHj[j]=sj[i]-1;\n\t\tHv[j]=sv[i];\n\t\tj++;\n\t\tif (si[i]!=sj[i]) {\n\t\t\tHi[j]=sj[i]-1;\n\t\t\tHj[j]=si[i]-1;\n\t\t\tHv[j]=sv[i];\n\t\t\tj++;\n\t\t}\n\t}\n}\n\n\n/* Functions */\n\nstatic char cuter__dims_doc[]=\n"Returns the dimension of the problem and the number of constraints.\\n"\n"\\n"\n"(n, m)=_dims()\\n"\n"\\n"\n"Output\\n"\n"n -- number of variables\\n"\n"m -- number of constraints\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"__init__.py script when the test function interface is loaded.\\n"\n"If you decide to call it anyway, the working directory at the time of call\\n"\n"must be the one where the file OUTSIF.d can be found.\\n"\n"\\n"\n"CUTEr tools used: CDIMEN\\n";\n\nstatic PyObject *cuter__dims(PyObject *self, PyObject *args) {\n\tif (PyObject_Length(args)!=0)\n\t\tPyErr_SetString(PyExc_Exception, "_dims() takes no arguments");\n\n\tif (!open_datafile())\n\t\treturn NULL;\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Calling CDIMEN\\n");\n#endif\n\tCDIMEN((integer *)&funit, (integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon);\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: n = %-d, m = %-d\\n", CUTEr_nvar, CUTEr_ncon);\n#endif\n\n\treturn decRefTuple(PyTuple_Pack(2, PyLong_FromLong((long)CUTEr_nvar), PyLong_FromLong((long)CUTEr_ncon))); \n}\n\n\nstatic char cuter__setup_doc[]=\n"Sets up the problem.\\n"\n"\\n"\n"data=_setup(efirst, lfirst, nvfirst)\\n"\n"\\n"\n"Input\\n"\n"efirst -- if True, equation constraints are ordered before inequations.\\n"\n" Defaults to False.\\n"\n"lfirst -- if True, linear constraints are ordered before nonlinear ones.\\n"\n" Defaults to False.\\n"\n"nvfirst -- if True, nonlinear variables are ordered before linear ones.\\n"\n" Defaults to False.\\n"\n"\\n"\n"Setting both efirst and lfirst to True results in the following ordering:\\n"\n"linear equations, followed by linear inequations, nonlinear equations,\\n"\n"and finally nonlinear inequations.\\n"\n"\\n"\n"Output\\n"\n"data -- dictionary with the summary of test function\'s properties\\n"\n"\\n"\n"The problem data dictionary has the following members:\\n"\n"name -- problem name\\n"\n"n -- number of variables\\n"\n"m -- number of constraints (excluding bounds)\\n"\n"x -- initial point (1D array of length n)\\n"\n"bl -- vector of lower bounds on variables (1D array of length n)\\n"\n"bu -- vector of upper bounds on variables (1D array of length n)\\n"\n"nnzh -- number of nonzero elements in the diagonal and upper triangle of\\n"\n" sparse Hessian\\n"\n"vartype -- 1D integer array of length n storing variable type\\n"\n" 0=real, 1=boolean (0 or 1), 2=integer\\n"\n"\\n"\n"For constrained problems the following additional members are available\\n"\n"nnzj -- number of nonzero elements in sparse Jacobian of constraints\\n"\n"v -- initial value of Lagrange multipliers (1D array of length m)\\n"\n"cl -- lower bounds on constraint functions (1D array of length m)\\n"\n"cu -- upper bounds on constraint functions (1D array of length m)\\n"\n"equatn -- 1D boolean array of length m indicating whether a constraint\\n"\n" is an equation constraint\\n"\n"linear -- 1D boolean array of length m indicating whether a constraint\\n"\n" is a linear constraint\\n"\n"\\n"\n"-1e+20 and 1e+20 in bl, bu, cl, and cu stand for -infinity and +infinity.\\n"\n"\\n"\n"This function must be called before any other CUTEr function is called.\\n"\n"The only exception is the _dims() function.\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"__init__.py script when the test function interface is loaded.\\n"\n"If you decide to call it anyway, the working directory at the time of call\\n"\n"must be the one where the file OUTSIF.d can be found.\\n"\n"\\n"\n"CUTEr tools used: CDIMEN, CSETUP, USETUP, CVARTY, UVARTY, \\n"\n" CDIMSH, UDIMSH, CDIMSJ, PBNAME\\n";\n\nstatic PyObject *cuter__setup(PyObject *self, PyObject *args) {\n\tnpy_logical efirst = FALSE_, lfirst = FALSE_, nvfrst = FALSE_;\n\tint eFirst, lFirst, nvFirst;\n\tPyObject *dict;\n\tPyArrayObject *Mx, *Mbl, *Mbu, *Mv=NULL, *Mcl=NULL, *Mcu=NULL, *Meq=NULL, *Mlinear=NULL;\n\tPyArrayObject *Mvt; \n\tdoublereal *x, *bl, *bu, *v=NULL, *cl=NULL, *cu=NULL;\n\tnpy_integer *vartypes; \n\tnpy_logical *equatn=NULL, *linear=NULL;\n\tnpy_intp dims[1];\n\tint i;\n\t\n\tif (PyObject_Length(args)!=0 && PyObject_Length(args)!=3) {\n\t\tPyErr_SetString(PyExc_Exception, "_setup() takes 0 or 3 arguments");\n\t\treturn NULL; \n\t}\n\t\n\tif (PyObject_Length(args)==3) {\n\t\tif (!PyArg_ParseTuple(args, "iii", &eFirst, &lFirst, &nvFirst)) {\n\t\t\treturn NULL; \n\t\t}\n\n\t\tefirst = eFirst ? TRUE_ : FALSE_;\n\t\tlfirst = lFirst ? TRUE_ : FALSE_;\n\t\tnvfrst = nvFirst ? TRUE_ : FALSE_;\n\t}\n\t\n\tif (!open_datafile())\n\t\treturn NULL; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Calling CDIMEN\\n");\n#endif\n\tCDIMEN((integer *)&funit, (integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon);\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: n = %-d, m = %-d\\n", CUTEr_nvar, CUTEr_ncon);\n#endif\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Allocating space\\n");\n#endif\n\t/* Numpy arrays */\n\tdims[0]=CUTEr_nvar;\n\tMx=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tMbl=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tMbu=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tMvt=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_INT); \n\tif (CUTEr_ncon>0) {\n\t\tdims[0]=CUTEr_ncon;\n\t\tMv=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\tMcl=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\tMcu=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\tMeq=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_BOOL);\n\t\tMlinear=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_BOOL);\n\t}\n\t\n\t/* Get internal data buffers */\n\t/* Assume that npy_double is equivalent to double and npy_int is equivalent to integer */\n\tx = (npy_double *)PyArray_DATA(Mx); \n\tbl = (npy_double *)PyArray_DATA(Mbl); \n\tbu = (npy_double *)PyArray_DATA(Mbu); \n\tif (CUTEr_ncon>0) {\n\t\tv = (npy_double *)PyArray_DATA(Mv); \n\t\tcl = (npy_double *)PyArray_DATA(Mcl); \n\t\tcu = (npy_double *)PyArray_DATA(Mcu); \n\t\n\t\t/* Create temporary CUTEr logical arrays */\n\t\tequatn = (npy_logical *)malloc(CUTEr_ncon*sizeof(npy_logical));\n\t\tlinear = (npy_logical *)malloc(CUTEr_ncon*sizeof(npy_logical));\n\t}\n\tvartypes=(npy_integer *)malloc(CUTEr_nvar*sizeof(npy_integer)); \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Calling [UC]SETUP\\n");\n#endif\n\tif (CUTEr_ncon > 0)\n\t\tCSETUP((integer *)&funit, (integer *)&iout, (integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, bl, bu,\n\t\t\t\t(integer *)&CUTEr_nvar, (logical *)equatn, (logical *)linear, v, cl, cu, (integer *)&CUTEr_ncon,\n\t\t\t\t(logical *)&efirst, (logical *)&lfirst, (logical *)&nvfrst);\n\telse\n\t\tUSETUP((integer *)&funit, (integer *)&iout, (integer *)&CUTEr_nvar, x, bl, bu, (integer *)&CUTEr_nvar);\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: n = %-d, m = %-d\\n", CUTEr_nvar, CUTEr_ncon);\n#endif\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Calling [UC]VARTY\\n");\n#endif\n\tif (CUTEr_ncon > 0)\n\t\tCVARTY((integer *)&CUTEr_nvar, (integer *)vartypes);\n\telse\n\t\tUVARTY((integer *)&CUTEr_nvar, (integer *)vartypes);\n\t\t\n\t/* Copy logical values to NumPy bool arrays and free temporary storage */\n\tif (CUTEr_ncon > 0) {\n\t\tfor(i=0; i<CUTEr_ncon; i++) {\n\t\t\t*((npy_bool*)(PyArray_GETPTR1(Meq, i)))=(equatn[i]==TRUE_ ? NPY_TRUE : NPY_FALSE);\n\t\t\t*((npy_bool*)(PyArray_GETPTR1(Mlinear, i)))=(linear[i]==TRUE_ ? NPY_TRUE : NPY_FALSE);\n\t\t}\n\t\tfree(equatn);\n\t\tfree(linear);\n\t} \n\t\n\t/* Copy variable types to NumPy integer arrays and free temporary storage */\n\tfor(i=0; i<CUTEr_nvar; i++) {\n\t\t*((npy_int*)PyArray_GETPTR1(Mvt, i))=vartypes[i];\n\t}\n\tfree(vartypes);\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Calling [CU]DIMSH\\n");\n#endif\n\tif (CUTEr_ncon>0)\n\t\tCDIMSH((integer *)&CUTEr_nnzh);\n\telse \n\t\tUDIMSH((integer *)&CUTEr_nnzh);\n\t\t\n\tif (CUTEr_ncon > 0) {\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: Calling CDIMSJ\\n");\n#endif\n\t\tCDIMSJ((integer *)&CUTEr_nnzj);\n\t\tCUTEr_nnzj -= CUTEr_nvar;\n\t} \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: nnzh = %-d, nnzj = %-d\\n", CUTEr_nnzh, CUTEr_nnzj);\n\tfprintf(df, "PyCUTEr: Finding out problem name\\n");\n#endif\n\tPBNAME((integer *)&CUTEr_nvar, CUTEr_probName);\n\ttrim_string(CUTEr_probName, STR_LEN-1);\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: %-s\\n", CUTEr_probName);\n\tfprintf(df, "PyCUTEr: Closing data file\\n");\n#endif\n\tclose_datafile();\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: Building structure\\n");\n#endif\n\tdict=PyDict_New();\n\tPyDict_SetItemString(dict, "n", PyLong_FromLong((long)CUTEr_nvar)); \n\tPyDict_SetItemString(dict, "m", PyLong_FromLong((long)CUTEr_ncon)); \n\tPyDict_SetItemString(dict, "nnzh", PyLong_FromLong((long)CUTEr_nnzh)); \n\tPyDict_SetItemString(dict, "x", (PyObject *)Mx); \n\tPyDict_SetItemString(dict, "bl", (PyObject *)Mbl); \n\tPyDict_SetItemString(dict, "bu", (PyObject *)Mbu); \n\tPyDict_SetItemString(dict, "name", PyUnicode_FromString(CUTEr_probName)); \n\tPyDict_SetItemString(dict, "vartype", (PyObject *)Mvt); \n\tif (CUTEr_ncon > 0) {\n\t\tPyDict_SetItemString(dict, "nnzj", PyLong_FromLong((long)CUTEr_nnzj)); \n\t\tPyDict_SetItemString(dict, "v", (PyObject*)Mv); \n\t\tPyDict_SetItemString(dict, "cl", (PyObject*)Mcl); \n\t\tPyDict_SetItemString(dict, "cu", (PyObject*)Mcu); \n\t\tPyDict_SetItemString(dict, "equatn", (PyObject*)Meq); \n\t\tPyDict_SetItemString(dict, "linear", (PyObject*)Mlinear); \n\t}\n\t\n\tsetupCalled = 1;\n\t\n\treturn decRefDict(dict); \n}\n\n\nstatic char cuter__varnames_doc[]=\n"Returns the names of variables in the problem.\\n"\n"\\n"\n"namelist=_varnames()\\n"\n"\\n"\n"Output\\n"\n"namelist -- list of length n holding strings holding names of variables\\n"\n"\\n"\n"The list reflects the ordering imposed by the nvfirst argument to _setup().\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"__init__.py script when the test function interface is loaded.\\n"\n"\\n"\n"CUTEr tools used: VARNAMES\\n";\n\nstatic PyObject *cuter__varnames(PyObject *self, PyObject *args) {\n\tchar *Fvnames, Fvname[STR_LEN+1], *ptr;\n\tPyObject *list; \n\tint i, j; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\t\n\tif (PyObject_Length(args)!=0) {\n\t\tPyErr_SetString(PyExc_Exception, "_varnames() takes 0 arguments");\n\t\treturn NULL; \n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: allocating space\\n");\n#endif\n\tFvnames=(char *)malloc(CUTEr_nvar*STR_LEN*sizeof(char));\n\tlist=PyList_New(0);\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling VARNAMES\\n");\n#endif\n\tVARNAMES((integer *)&CUTEr_nvar, Fvnames);\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: building results\\n");\n#endif\n\tfor(i=0;i<CUTEr_nvar;i++) {\n\t\tptr=Fvnames+i*STR_LEN; \n\t\tfor(j=0;j<STR_LEN;j++) {\n\t\t\tFvname[j]=*ptr;\n\t\t\tptr++;\n\t\t}\n\t\ttrim_string(Fvname, STR_LEN-1);\n\t\tPyList_Append(list, PyUnicode_FromString(Fvname)); \n\t}\n\t\n\tfree(Fvnames); \n\t\n\treturn list;\n}\n\n\nstatic char cuter__connames_doc[]=\n"Returns the names of constraints in the problem.\\n"\n"\\n"\n"namelist=_connames()\\n"\n"\\n"\n"Output\\n"\n"namelist -- list of length m holding strings holding names of constraints\\n"\n"\\n"\n"The list is ordered in the way specified by efirst and lfirst arguments to\\n"\n"_setup().\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"__init__.py script when the test function interface is loaded.\\n"\n"\\n"\n"CUTEr tools used: CONNAMES\\n";\n\nstatic PyObject *cuter__connames(PyObject *self, PyObject *args) {\n\tchar *Fcnames, Fcname[STR_LEN+1], *ptr;\n\tPyObject *list; \n\tint i, j; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\t\n\tif (PyObject_Length(args)!=0) {\n\t\tPyErr_SetString(PyExc_Exception, "_connames() takes 0 arguments");\n\t\treturn NULL; \n\t}\n\t\n\tlist=PyList_New(0);\n\t\n\tif (CUTEr_ncon>0) {\n\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: allocating space\\n");\n#endif\n\t\tFcnames=(char *)malloc(CUTEr_ncon*STR_LEN*sizeof(char));\n\t\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: calling CONNAMES\\n");\n#endif\n\t\tCONNAMES((integer *)&CUTEr_ncon, Fcnames);\n\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: building results\\n");\n#endif\n\t\tfor(i=0;i<CUTEr_ncon;i++) {\n\t\t\tptr=Fcnames+i*STR_LEN; \n\t\t\tfor(j=0;j<STR_LEN;j++) {\n\t\t\t\tFcname[j]=*ptr;\n\t\t\t\tptr++;\n\t\t\t}\n\t\t\ttrim_string(Fcname, STR_LEN-1);\n\t\t\tPyList_Append(list, PyUnicode_FromString(Fcname)); \n\t\t}\n\t\n\t\tfree(Fcnames); \n\t}\n\t\n\treturn list;\n}\n\n\nstatic char cuter_objcons_doc[]=\n"Returns the value of objective and constraints at x.\\n"\n"\\n"\n"(f, c)=objcons(x)\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"\\n"\n"Output\\n"\n"f -- 1D array of length 1 holding the value of the function at x\\n"\n"c -- 1D array of length m holding the values of constraints at x\\n"\n"\\n"\n"CUTEr tools used: CFN\\n";\n\nstatic PyObject *cuter_objcons(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *Mf, *Mc; \n\tdoublereal *x, *f, *c; \n\tnpy_intp dims[1]; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\tif (!PyArg_ParseTuple(args, "O", &arg1)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct length and shape */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1)&& PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tdims[0]=1;\n\tMf=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tf=(npy_double *)PyArray_DATA(Mf);\n\tdims[0]=CUTEr_ncon;\n\tMc=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tc=(npy_double *)PyArray_DATA(Mc);\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CFN\\n");\n#endif\n\tCFN((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, f, (integer *)&CUTEr_ncon, c);\n\t\n\treturn decRefTuple(PyTuple_Pack(2, Mf, Mc)); \n}\n\n\nstatic char cuter_obj_doc[]=\n"Returns the value of objective and its gradient at x.\\n"\n"\\n"\n"f=obj(x)\\n"\n"(f, g)=obj(x, gradFlag)\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"gradFlag -- if given the function returns f and g; can be anything\\n"\n"\\n"\n"Output\\n"\n"f -- 1D array of length 1 holding the value of the function at x\\n"\n"g -- 1D array of length n holding the value of the gradient of f at x\\n"\n"\\n"\n"CUTEr tools used: UOFG, COFG\\n";\n\nstatic PyObject *cuter_obj(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1; \n\tPyObject *arg2;\n\tPyArrayObject *Mf, *Mg=NULL; \n\tdoublereal *x, *f, *g=NULL; \n\tnpy_intp dims[1];\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\tif (!PyArg_ParseTuple(args, "O|O", &arg1, &arg2)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct length and shape */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tdims[0]=1;\n\tMf=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tf=(npy_double *)PyArray_DATA(Mf);\n\tif (PyObject_Length(args)>1) {\n\t\tdims[0]=CUTEr_nvar;\n\t\tMg=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\tg=(npy_double *)PyArray_DATA(Mg);\n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling [UC]OFG\\n");\n#endif\n\tif (CUTEr_ncon == 0) {\n\t\tif (PyObject_Length(args)==1) {\n\t\t\tUOFG((integer *)&CUTEr_nvar, x, f, NULL, &somethingFalse);\n\t\t\treturn (PyObject *)Mf; \n\t\t} else {\n\t\t\tUOFG((integer *)&CUTEr_nvar, x, f, g, &somethingTrue);\n\t\t\treturn decRefTuple(PyTuple_Pack(2, Mf, Mg)); \n\t\t}\n\t} else {\n\t\tif (PyObject_Length(args)==1) {\n\t\t\tCOFG((integer *)&CUTEr_nvar, x, f, NULL, &somethingFalse);\n\t\t\treturn (PyObject *)Mf; \n\t\t} else {\n\t\t\tCOFG((integer *)&CUTEr_nvar, x, f, g, &somethingTrue);\n\t\t\treturn decRefTuple(PyTuple_Pack(2, Mf, Mg)); \n\t\t}\n\t}\n}\n\n\nstatic char cuter_cons_doc[]=\n"Returns the value of constraints and the Jacobian of constraints at x.\\n"\n"\\n"\n"c=cons(x) -- constraints\\n"\n"ci=cons(x, False, i) -- i-th constraint\\n"\n"(c, J)=cons(x, True) -- Jacobian of constraints\\n"\n"(ci, Ji)=cons(x, True, i) -- i-th constraint and its gradient\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"i -- integer index of constraint (between 0 and m-1)\\n"\n"\\n"\n"Output\\n"\n"c -- 1D array of length m holding the values of constraints at x\\n"\n"ci -- 1D array of length 1 holding the value of i-th constraint at x\\n"\n"J -- 2D array with m rows of n columns holding Jacobian of constraints at x\\n"\n"Ji -- 1D array of length n holding the gradient of i-th constraintn"\n" (also equal to the i-th row of Jacobian)\\n"\n"\\n"\n"CUTEr tools used: CCFG, CCIFG\\n";\n\nstatic PyObject *cuter_cons(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *Mc, *MJ; \n\tPyObject *arg2;\n\tdoublereal *x, *c, *J; \n\tint derivs, index, wantSingle;\n\tnpy_integer icon; \n\tnpy_integer zero = 0;\n\tnpy_intp dims[2];\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\tif (!PyArg_ParseTuple(args, "O|Oi", &arg1, &arg2, &index)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct length and shape */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\t/* Do we want derivatives */\n\tif (arg2!=NULL && arg2==Py_True)\n\t\tderivs=1;\n\telse\n\t\tderivs=0;\n\t\n\t/* Do we want a particular derivative */\n\tif (PyObject_Length(args)==3) {\n\t\t/* Check index */\n\t\tif (index<0 || index>=CUTEr_ncon) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 3 must be an integer between 0 and ncon-1");\n\t\t\treturn NULL; \n\t\t}\n\t\ticon=index+1;\n\t\twantSingle=1;\n\t} else {\n\t\twantSingle=0;\n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tif (!wantSingle) {\n\t\tdims[0]=CUTEr_ncon;\n\t\tMc=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\tc=(npy_double *)PyArray_DATA(Mc);\n\t\tif (derivs) {\n\t\t\tdims[0]=CUTEr_ncon;\n\t\t\tdims[1]=CUTEr_nvar;\n\t\t\t/* Create a FORTRAN style array (first index stride is 1) */\n\t\t\tMJ=(PyArrayObject *)PyArray_New(&PyArray_Type, 2, dims, NPY_DOUBLE, NULL, NULL, 0, NPY_ARRAY_F_CONTIGUOUS, NULL);\n\t\t\tJ=(npy_double *)PyArray_DATA(MJ);\n\t\t}\n\t} else {\n\t\tdims[0]=1;\n\t\tMc=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\tc=(npy_double *)PyArray_DATA(Mc);\n\t\tif (derivs) {\n\t\t\tdims[0]=CUTEr_nvar;\n\t\t\tMJ=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\t\t\tJ=(npy_double *)PyArray_DATA(MJ);\n\t\t}\n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CCFG/CCIFG\\n");\n#endif\n\tif (!wantSingle) {\n\t\tif (!derivs) {\n\t\t\tCCFG((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, (integer *)&CUTEr_ncon, c,\n\t\t\t\t\t&somethingFalse, (integer *)&zero, (integer *)&zero, NULL, &somethingFalse);\n\t\t\treturn (PyObject *)Mc;\n\t\t} else {\n\t\t\tCCFG((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, (integer *)&CUTEr_ncon, c,\n\t\t\t\t\t&somethingFalse, (integer *)&CUTEr_ncon, (integer *)&CUTEr_nvar, J,\n\t\t\t\t\t&somethingTrue);\n\t\t\treturn decRefTuple(PyTuple_Pack(2, Mc, MJ)); \n\t\t}\n\t} else {\n\t\tif (!derivs) {\n\t\t\tCCIFG((integer *)&CUTEr_nvar, (integer *)&icon, x, c, NULL, &somethingFalse);\n\t\t\treturn (PyObject *)Mc; \n\t\t} else {\n\t\t\tCCIFG((integer *)&CUTEr_nvar, (integer *)&icon, x, c, J, &somethingTrue);\n\t\t\treturn decRefTuple(PyTuple_Pack(2, Mc, MJ)); \n\t\t}\n\t}\n}\n\n\nstatic char cuter_lagjac_doc[]=\n"Returns the gradient of the objective or Lagrangian, and the Jacobian of\\n"\n"constraints at x. The gradient is the gradient with respect to the problem\'s\\n"\n"variables (has n components).\\n"\n"\\n"\n"(g, J)=lagjac(x) -- objective gradient and the Jacobian of constraints\\n"\n"(g, J)=lagjac(x, v) -- Lagrangian gradient and the Jacobian of constraints\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"v -- 1D array of length m with the Lagrange multipliers\\n"\n"\\n"\n"Output\\n"\n"g -- 1D array of length n holding the gradient at x\\n"\n"J -- 2D array with m rows of n columns holding Jacobian of constraints at x\\n"\n"\\n"\n"CUTEr tools used: CGR\\n";\n\nstatic PyObject *cuter_lagjac(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *Mg, *MJ; \n\tdoublereal *x, *v=NULL, *g, *J; \n\tint lagrangian; \n\tnpy_intp dims[2];\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\tif (!PyArg_ParseTuple(args, "O|O", &arg1, &arg2)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct length and shape */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check if v is double and of correct length and shape. */\n\tif (arg2!=NULL) {\n\t\tif (!(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length ncon");\n\t\t\treturn NULL; \n\t\t}\n\t\tlagrangian=1;\n\t} else {\n\t\tlagrangian=0;\n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tif (lagrangian) \n\t\tv=(npy_double *)PyArray_DATA(arg2);\n\tdims[0]=CUTEr_nvar;\n\tMg=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tg=(npy_double *)PyArray_DATA(Mg);\n\tdims[0]=CUTEr_ncon;\n\tdims[1]=CUTEr_nvar;\n\t/* Create a FORTRAN style array (first index stride is 1) */\n\tMJ=(PyArrayObject *)PyArray_New(&PyArray_Type, 2, dims, NPY_DOUBLE, NULL, NULL, 0, NPY_ARRAY_F_CONTIGUOUS, NULL);\n\tJ=(npy_double *)PyArray_DATA(MJ);\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CGR\\n");\n#endif\n\tif (!lagrangian) {\n\t\tCGR((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, &somethingFalse, (integer *)&CUTEr_ncon,\n\t\t\tNULL, g, &somethingFalse, (integer *)&CUTEr_ncon, (integer *)&CUTEr_nvar, J);\n\t} else {\n\t\tCGR((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, &somethingTrue, (integer *)&CUTEr_ncon,\n\t\t\tv, g, &somethingFalse, (integer *)&CUTEr_ncon, (integer *)&CUTEr_nvar, J);\n\t}\n\t\n\treturn decRefTuple(PyTuple_Pack(2, Mg, MJ)); \n}\n\n\nstatic char cuter_jprod_doc[]=\n"Returns the product of constraints Jacobian at x with vector p\\n"\n"\\n"\n"r=jprod(transpose, p, x) -- computes Jacobian at x before product calculation\\n"\n"r=jprod(transpose, p) -- uses last computed Jacobian\\n"\n"\\n"\n"Input\\n"\n"transpose -- boolean flag indicating that the Jacobian should be trasposed\\n"\n" before the product is calculated\\n"\n"p -- the vector that will be multiplied with the Jacobian\\n"\n" 1D array of length n (m) if transpose if False (True)\\n"\n"x -- 1D array of length n holding the values of variables used in the\\n"\n" evaluation of the constraints Jacobian\\n"\n"\\n"\n"Output\\n"\n"r -- 1D array of length m if transpose=False (or n if transpose=True)\\n"\n" with the result\\n"\n"\\n"\n"CUTEr tools used: CJPROD\\n";\n\nstatic PyObject *cuter_jprod(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg2, *arg3, *Mr; \n\tPyObject *arg1;\n\tdoublereal *p, *x, *r; \n\tint transpose; \n\tnpy_intp dims[1];\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ3=NULL;\n\tif (!PyArg_ParseTuple(args, "OO|O", &arg1, &arg2, &arg3)) \n\t\treturn NULL; \n\t\n\t/* Check if arg1 is True */\n\tif (arg1==Py_True) \n\t\ttranspose=1;\n\telse\n\t\ttranspose=0; \n\t\n\t/* Check if p is double and of correct dimension */\n\tif (!(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check length of p when J is not transposed */\n\tif (!transpose && !(PyArray_DIM(arg2, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be of length nvar (J is not transposed)");\n\t\treturn NULL;\n\t}\n\t\n\t/* Check length of p when J is transposed */\n\tif (transpose && !(PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be of length ncon (J is transposed)");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check if x is double and of correct length and shape. */\n\tif (arg3!=NULL) {\n\t\tif (!(arg3!=NULL && PyArray_Check(arg3) && PyArray_ISFLOAT(arg3) && PyArray_TYPE(arg3)==NPY_DOUBLE && PyArray_NDIM(arg3)==1 && PyArray_DIM(arg3, 0)==CUTEr_nvar)) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 3 must be a 1D double array of length nvar");\n\t\t\treturn NULL; \n\t\t}\n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tp=(npy_double *)PyArray_DATA(arg2);\n\tif (arg3!=NULL) \n\t\tx=(npy_double *)PyArray_DATA(arg3);\n\telse\n\t\tx=NULL;\n\tif (!transpose) {\n\t\tdims[0]=CUTEr_ncon;\n\t} else {\n\t\tdims[0]=CUTEr_nvar;\n\t}\n\tMr=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE);\n\tr=(npy_double *)PyArray_DATA(Mr);\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CJPROD\\n");\n#endif\n\tif (!transpose) {\n\t\tif (arg3==NULL) {\n\t\t\tCJPROD((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingTrue,\n\t\t\t\t\t&somethingFalse, NULL, p, (integer *)&CUTEr_nvar, r, (integer *)&CUTEr_ncon);\n\t\t} else {\n\t\t\tCJPROD((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingFalse,\n\t\t\t\t\t&somethingFalse, x, p, (integer *)&CUTEr_nvar, r, (integer *)&CUTEr_ncon);\n\t\t}\n\t} else {\n\t\tif (arg3==NULL) {\n\t\t\tCJPROD((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingTrue,\n\t\t\t\t\t&somethingTrue, NULL, p, (integer *)&CUTEr_ncon, r, (integer *)&CUTEr_nvar);\n\t\t} else {\n\t\t\tCJPROD((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingFalse,\n\t\t\t\t\t&somethingTrue, x, p, (integer *)&CUTEr_ncon, r, (integer *)&CUTEr_nvar);\n\t\t}\n\t}\n\t\n\treturn (PyObject *)Mr;\n}\n\n\nstatic char cuter_hess_doc[]=\n"Returns the Hessian of the objective (for unconstrained problems) or the\\n"\n"Hessian of the Lagrangian (for constrained problems) at x.\\n" \n"\\n"\n"H=hess(x) -- Hessian of objective at x for unconstrained problems\\n"\n"H=hess(x, v) -- Hessian of Lagrangian at (x, v) for constrained problems\\n"\n"\\n"\n"The first form can only be used for unconstrained problems. The second one\\n"\n"can only be used for constrained problems. For obtaining the Hessian of the\\n"\n"objective in case of a constrained problem use ihess().\\n"\n"\\n"\n"The Hessian is meant with respect to problem variables (has dimension n).\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n holding the values of variables\\n"\n"v -- 1D array of length m holding the values of Lagrange multipliers\\n"\n"\\n"\n"Output\\n"\n"H -- 2D array with n rows of n columns holding the Hessian at x (or (x, v))\\n"\n"\\n"\n"CUTEr tools used: CDH, UDH\\n";\n\nstatic PyObject *cuter_hess(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *MH; \n\tdoublereal *x, *v=NULL, *H; \n\tnpy_intp dims[2];\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\tif (!PyArg_ParseTuple(args, "O|O", &arg1, &arg2)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\tif (CUTEr_ncon>0) {\n\t\t/* Check if v is double and of correct dimension */\n\t\tif (arg2!=NULL) {\n\t\t\tif (!(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length ncon");\n\t\t\t\treturn NULL; \n\t\t\t}\n\t\t} else {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be specified for constrained problems. Use ihess().");\n\t\t\treturn NULL; \n\t\t}\n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tif (CUTEr_ncon>0) \n\t\tv=(npy_double *)PyArray_DATA(arg2);\n\tdims[0]=CUTEr_nvar;\n\tdims[1]=CUTEr_nvar;\n\t/* Create a FORTRAN style array (first index stride is 1) */\n\tMH=(PyArrayObject *)PyArray_New(&PyArray_Type, 2, dims, NPY_DOUBLE, NULL, NULL, 0, NPY_ARRAY_F_CONTIGUOUS, NULL);\n\tH=(npy_double *)PyArray_DATA(MH);\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling [CU]DH\\n");\n#endif\n\tif (CUTEr_ncon>0) {\n\t\tCDH((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, (integer *)&CUTEr_ncon, v, (integer *)&CUTEr_nvar, H);\n\t} else {\n\t\tUDH((integer *)&CUTEr_nvar, x, (integer *)&CUTEr_nvar, H);\n\t}\n\t\n\treturn (PyObject *)MH; \n}\n\n\nstatic char cuter_ihess_doc[]=\n"Returns the Hessian of the objective or the Hessian of i-th constraint at x.\\n"\n"\\n"\n"H=ihess(x) -- Hessian of the objective\\n"\n"H=ihess(x, i) -- Hessian of i-th constraint\\n"\n"\\n"\n"The first form can only be used for unconstrained problems. The second one\\n"\n"can only be used for constrained problems. For obtaining the Hessian of the\\n"\n"objective in case of a constrained problem use ihess().\\n"\n"\\n"\n"The Hessian is meant with respect to problem variables (has dimension n).\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n holding the values of variables\\n"\n"i -- integer holding the index of the constraint (between 0 and m-1)\\n"\n"\\n"\n"Output\\n"\n"H -- 2D array with n rows of n columns holding the Hessian at x\\n"\n"\\n"\n"CUTEr tools used: CIDH, UDH\\n";\n\nstatic PyObject *cuter_ihess(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *MH; \n\tdoublereal *x, *H; \n\tnpy_intp dims[2];\n\tint i; \n\tnpy_integer icon; \n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\tif (!PyArg_ParseTuple(args, "O|i", &arg1, &i)) \n\t\treturn NULL; \n\t\n\tif (PyObject_Length(args)>1) {\n\t\ticon=i+1;\n\t\tif (i<0 || i>=CUTEr_ncon) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be between 0 and ncon-1");\n\t\t\treturn NULL; \n\t\t}\n\t} else {\n\t\ticon=0; \n\t}\n\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tdims[0]=CUTEr_nvar;\n\tdims[1]=CUTEr_nvar;\n\t/* Create a FORTRAN style array (first index stride is 1) */\n\tMH=(PyArrayObject *)PyArray_New(&PyArray_Type, 2, dims, NPY_DOUBLE, NULL, NULL, 0, NPY_ARRAY_F_CONTIGUOUS, NULL);\n\tH=(npy_double *)PyArray_DATA(MH);\t\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CIDH/UDH\\n");\n#endif\n\tif (CUTEr_ncon>0) {\n\t\tCIDH((integer *)&CUTEr_nvar, x, (integer *)&icon, (integer *)&CUTEr_nvar, H);\n\t} else {\n\t\tUDH((integer *)&CUTEr_nvar, x, (integer *)&CUTEr_nvar, H);\n\t}\n\t\n\treturn (PyObject *)MH;\n}\n\n\nstatic char cuter_hprod_doc[]=\n"Returns the product of Hessian at x and vector p.\\n"\n"The Hessian is either the Hessian of objective or the Hessian of Lagrangian.\\n" \n"\\n"\n"r=hprod(p, x, v) -- use Hessian of Lagrangian at x (constrained problem)\\n"\n"r=hprod(p, x) -- use Hessian of objective at x (unconstrained problem)\\n"\n"r=hprod(p) -- use last computed Hessian\\n"\n"\\n"\n"The first form can only be used for constrained problems. The second one\\n"\n"can only be used for unconstrained problems.\\n"\n"\\n"\n"The Hessian is meant with respect to problem variables (has dimension n).\\n"\n"\\n"\n"Input\\n"\n"p -- 1D array of length n holding the components of the vector\\n"\n"x -- 1D array of length n holding the values of variables\\n"\n"v -- 1D array of length m holding the values of Lagrange multipliers\\n"\n"\\n"\n"Output\\n"\n"r -- 1D array of length n holding the result\\n"\n"\\n"\n"CUTEr tools used: CPROD, UPROD\\n";\n\nstatic PyObject *cuter_hprod(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *arg3, *Mr; \n\tdoublereal *p, *x=NULL, *v=NULL, *r; \n\tnpy_intp dims[1];\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=arg3=NULL; \n\tif (!PyArg_ParseTuple(args, "O|OO", &arg1, &arg2, &arg3)) \n\t\treturn NULL; \n\t\n\tif (CUTEr_ncon>0) {\n\t\tif (PyObject_Length(args)==2) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Need 1 or 3 arguments for constrained problems");\n\t\t\treturn NULL; \n\t\t}\n\t} else {\n\t\tif (PyObject_Length(args)==3) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Need 1 or 2 arguments for unconstrained problems");\n\t\t\treturn NULL; \n\t\t}\n\t}\n\t\n\t/* Check if p is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check if x is double and of correct dimension */\n\tif (arg2!=NULL && !(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check if v is double and of correct dimension */\n\tif (arg3!=NULL && !(PyArray_Check(arg3) && PyArray_ISFLOAT(arg3) && PyArray_TYPE(arg3)==NPY_DOUBLE && PyArray_NDIM(arg3)==1 && PyArray_DIM(arg3, 0)==CUTEr_ncon)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 3 must be a 1D double array of length ncon");\n\t\treturn NULL; \n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tp=(npy_double *)PyArray_DATA(arg1);\n\tif (arg2!=NULL)\n\t\tx=(npy_double *)PyArray_DATA(arg2);\n\tif (arg3!=NULL)\n\t\tv=(npy_double *)PyArray_DATA(arg3);\n\tdims[0]=CUTEr_nvar;\n\tMr=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\tr=(npy_double *)PyArray_DATA(Mr);\t\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling [CU]PROD\\n");\n#endif\n\tif (CUTEr_ncon>0) {\n\t\tif (arg2==NULL)\n\t\t\tCPROD((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingTrue, NULL, (integer *)&CUTEr_ncon, NULL, p, r);\n\t\telse \n\t\t\tCPROD((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingFalse, x, (integer *)&CUTEr_ncon, v, p, r);\n\t} else {\n\t\tif (arg2==NULL)\n\t\t\tUPROD((integer *)&CUTEr_nvar, &somethingTrue, NULL, p, r);\n\t\telse \n\t\t\tUPROD((integer *)&CUTEr_nvar, &somethingFalse, x, p, r);\n\t}\n\t\n\treturn (PyObject *)Mr;\n}\n\n\nstatic char cuter_gradhess_doc[]=\n"Returns the Hessian of the Lagrangian, the Jacobian of constraints, and the\\n"\n"gradient of the objective or the gradient of the Lagrangian at.\\n" \n"\\n"\n"(g, H)=gradhess(x) -- for unconstrained problems\\n"\n"(g, J, H)=gradhess(x, v, gradl) -- for constrained problems\\n"\n"\\n"\n"The first form can only be used for unconstrained problems. The second one\\n"\n"can only be used for constrained problems.\\n"\n"\\n"\n"The Hessian is meant with respect to problem variables (has dimension n).\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n holding the values of variables\\n"\n"v -- 1D array of length m holding the values of Lagrange multipliers\\n"\n"gradl -- boolean flag. If False the gradient of the objective is returned, \\n"\n" if True the gradient of the Lagrangian is returned.\\n"\n" Default is False.\\n"\n"\\n"\n"Output\\n"\n"g -- 1D array of length n holding\\n"\n" the gradient of objective at x (for unconstrained problems) or\\n"\n" the gradient of Lagrangian at (x, v) (for constrained problems)\\n"\n"J -- 2D array with m rows and n columns holding the Jacobian of constraints\\n"\n"H -- 2D array with n rows and n columns holding\\n"\n" the Hessian of the objective (for unconstrained problems) or\\n"\n" the Hessian of the Lagrangian (for constrained problems)\\n"\n"\\n"\n"CUTEr tools used: CGRDH, UGRDH\\n";\n\nstatic PyObject *cuter_gradhess(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *Mg, *MH, *MJ; \n\tPyObject *arg3; \n\tdoublereal *x, *v=NULL, *g, *H, *J; \n\tnpy_logical grlagf; \n\tnpy_intp dims[2];\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\targ3=NULL; \n\tif (!PyArg_ParseTuple(args, "O|OO", &arg1, &arg2, &arg3)) \n\t\treturn NULL; \n\t\n\tif (CUTEr_ncon>0) {\n\t\tif (PyObject_Length(args)<2) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Need at least 2 arguments for constrained problems");\n\t\t\treturn NULL; \n\t\t}\n\t} else {\n\t\tif (PyObject_Length(args)!=1) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Need 1 argument for unconstrained problems");\n\t\t\treturn NULL; \n\t\t}\n\t}\n\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check if v is double and of correct dimension */\n\tif (arg2!=NULL && !(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length ncon");\n\t\treturn NULL; \n\t}\n\t\n\t/* Are we computing the gradient of the Lagrangian */\n\tif (arg3!=NULL && arg3==Py_True) {\n\t\tgrlagf=TRUE_;\n\t} else {\n\t\tgrlagf=FALSE_; \n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tif (arg2!=NULL)\n\t\tv=(npy_double *)PyArray_DATA(arg2);\n\tdims[0]=CUTEr_nvar;\n\tMg=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\tg=(npy_double *)PyArray_DATA(Mg);\t\n\tdims[0]=CUTEr_nvar;\n\tdims[1]=CUTEr_nvar;\n\t/* Create a FORTRAN style array (first index stride is 1) */\n\tMH=(PyArrayObject *)PyArray_New(&PyArray_Type, 2, dims, NPY_DOUBLE, NULL, NULL, 0, NPY_ARRAY_F_CONTIGUOUS, NULL);\n\tH=(npy_double *)PyArray_DATA(MH);\n\tdims[0]=CUTEr_ncon;\n\tdims[1]=CUTEr_nvar;\n\t/* Create a FORTRAN style array (first index stride is 1) */\n\tMJ=(PyArrayObject *)PyArray_New(&PyArray_Type, 2, dims, NPY_DOUBLE, NULL, NULL, 0, NPY_ARRAY_F_CONTIGUOUS, NULL);\n\tJ=(npy_double *)PyArray_DATA(MJ);\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling [CU]GRDH\\n");\n#endif\n\tif (CUTEr_ncon>0) {\n\t\tCGRDH((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, (logical *)&grlagf, (integer *)&CUTEr_ncon, v,\n\t\t\t\tg, &somethingFalse, (integer *)&CUTEr_ncon, (integer *)&CUTEr_nvar, J, (integer *)&CUTEr_nvar, H);\n\t\treturn decRefTuple(PyTuple_Pack(3, Mg, MJ, MH)); \n\t} else {\n\t\tUGRDH((integer *)&CUTEr_nvar, x, g, (integer *)&CUTEr_nvar, H);\n\t\treturn decRefTuple(PyTuple_Pack(2, Mg, MH)); \n\t}\n}\n\n\nstatic char cuter__scons_doc[]=\n"Returns the value of constraints and the sparse Jacobian of constraints at x.\\n"\n"\\n"\n"(c, Jvi, Jfi, Jv)=_scons(x) -- Jacobian of constraints\\n"\n"(ci, gi, gv)=_scons(x, i) -- i-th constraint and its gradient\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"i -- integer index of constraint (between 0 and m-1)\\n"\n"\\n"\n"Output\\n"\n"c -- 1D array of length m holding the values of constraints at x\\n"\n"Jvi -- 1D array of integers holding the column indices (0 .. n-1)\\n"\n" of nozero elements in sparse Jacobian of constraints\\n"\n"Jfi -- 1D array of integers holding the row indices (0 .. m-1)\\n"\n" of nozero elements in sparse Jacobian of constraints\\n"\n"Jv -- 1D array holding the values of nonzero elements in the sparse Jacobian\\n"\n" of constraints at x. Has the same length as Jvi and Jfi.\\n"\n"ci -- 1D array of length 1 with the value of i-th constraint at x\\n"\n"gi -- 1D array of integers holding the indices (0 .. n-1) of nonzero\\n"\n" elements in the sparse gradient vector of i-th constraint\\n"\n"gv -- 1D array holding the values of nonzero elements in the sparse gradient\\n"\n" vector representing the gradient of i-th constraint at x.\\n"\n" Has the same length as gi. gi and gv corespond to the i-th row of\\n"\n" constraints Jacobian at x.\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"wrapper function scons().\\n"\n"\\n"\n"CUTEr tools used: CCFSG, CCIFSG\\n";\n\nstatic PyObject *cuter__scons(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *Mc, *MJi, *MJfi, *MJv, *Mgi, *Mgv; \n\tdoublereal *c, *Jv, *gv, *x, *sv; \n\tnpy_integer *Ji, *Jfi, *gi, *si; \n\tnpy_integer index, nnzsgc; \n\tint i; \n\tnpy_intp dims[1]; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\tif (!PyArg_ParseTuple(args, "O|i", &arg1, &i)) \n\t\treturn NULL; \n\t\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\tif (PyObject_Length(args)==2) {\n\t\tif (i<0 || i>=CUTEr_ncon) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be an integer between 0 and ncon-1");\n\t\t\treturn NULL; \n\t\t}\n\t\tindex=i+1;\n\t}\n\t\n\tif (PyObject_Length(args)==1) {\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\t\tx=(npy_double *)PyArray_DATA(arg1);\n\t\tdims[0]=CUTEr_nnzj;\n\t\tMJi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t\tMJfi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t\tMJv=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\t\tJi=(npy_integer *)PyArray_DATA(MJi);\t\n\t\tJfi=(npy_integer *)PyArray_DATA(MJfi);\t\n\t\tJv=(npy_double *)PyArray_DATA(MJv);\t\n\t\tdims[0]=CUTEr_ncon; \n\t\tMc=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\t\tc=(npy_double *)PyArray_DATA(Mc);\t\n\t\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: calling CCFSG\\n");\n#endif\n\t\tCCFSG((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, (integer *)&CUTEr_ncon, c, (integer *)&CUTEr_nnzj,\n (integer *)&CUTEr_nnzj, Jv, (integer *)Ji, (integer *)Jfi, &somethingTrue);\n\t\t\n\t\t/* Convert FORTRAN indices to C indices */\n\t\tfor(i=0;i<CUTEr_nnzj;i++) {\n\t\t\tJi[i]--;\n\t\t\tJfi[i]--;\n\t\t}\n\t\t\n\t\treturn decRefTuple(PyTuple_Pack(4, Mc, MJi, MJfi, MJv));\n\t} else {\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\t\tx=(npy_double *)PyArray_DATA(arg1);\n\t\tsi=(npy_integer *)malloc(CUTEr_nvar*sizeof(npy_integer));\n\t\tsv=(npy_double *)malloc(CUTEr_nvar*sizeof(npy_double));\n\t\tdims[0]=1;\n\t\tMc=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\t\tc=(npy_double *)PyArray_DATA(Mc);\t\n\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: calling CCIFSG\\n");\n#endif\n\t\tCCIFSG((integer *)&CUTEr_nvar, (integer *)&index, x, c, (integer *)&nnzsgc, (integer *)&CUTEr_nvar, sv, (integer *)si, &somethingTrue);\n\t\t\n\t\t/* Allocate and copy results, convert indices from FORTRAN to C, free storage */\n\t\tdims[0]=nnzsgc; \n\t\tMgi=(PyArrayObject *)PyArray_SimpleNew(1, dims, npy_integer_type_num); \n\t\tMgv=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\t\tgi=(npy_integer *)PyArray_DATA(Mgi);\n\t\tgv=(npy_double *)PyArray_DATA(Mgv);\t\n\t\tfor (i=0;i<nnzsgc;i++) {\n\t\t\tgi[i]=si[i]-1;\n\t\t\tgv[i]=sv[i];\n\t\t}\n\t\tfree(si);\n\t\tfree(sv);\n\t\t\n\t\treturn decRefTuple(PyTuple_Pack(3, Mc, Mgi, Mgv)); \n\t}\n}\n\n\nstatic char cuter__slagjac_doc[]=\n"Returns the sparse gradient of objective at x or Lagrangian at (x, v), \\n"\n"and the sparse Jacobian of constraints at x.\\n"\n"\\n"\n"(gi, gv, Jvi, Jfi, Jv)=_slagjac(x) -- objective gradient and Jacobian\\n"\n"(gi, gv, Jvi, Jfi, Jv)=_slagjac(x, v) -- Lagrangian gradient and Jacobian\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"v -- 1D array of length m with the values of Lagrange multipliers\\n"\n"\\n"\n"Output\\n"\n"gi -- 1D array of integers holding the indices (0 .. n-1) of nonzero\\n"\n" elements in the sparse gradient vector\\n"\n"gv -- 1D array holding the values of nonzero elements in the sparse gradient\\n"\n" vector. Has the same length as gi.\\n"\n"Jvi -- 1D array of integers holding the column indices (0 .. n-1)\\n"\n" of nozero elements in sparse Jacobian of constraints\\n"\n"Jfi -- 1D array of integers holding the row indices (0 .. m-1)\\n"\n" of nozero elements in sparse Jacobian of constraints\\n"\n"Jv -- 1D array holding the values of nonzero elements in the sparse Jacobian\\n"\n" of constraints at x. Has the same length as Jvi and Jfi.\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"wrapper function slagjac().\\n"\n"\\n"\n"CUTEr tools used: CSGR\\n";\n\nstatic PyObject *cuter__slagjac(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *Mgi, *Mgv, *MJi, *MJfi, *MJv; \n\tdoublereal *x, *v=NULL, *sv; \n\tnpy_integer *si, *sfi;\n\tnpy_integer nnzjplusn, nnzjplusno; \n\tint lagrangian; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\tif (!PyArg_ParseTuple(args, "O|O", &arg1, &arg2)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct length and shape */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\t/* Check if v is double and of correct length and shape. */\n\tif (arg2!=NULL) {\n\t\tif (!(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length ncon");\n\t\t\treturn NULL; \n\t\t}\n\t\tlagrangian=1;\n\t} else {\n\t\tlagrangian=0;\n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tif (lagrangian) \n\t\tv=(npy_double *)PyArray_DATA(arg2);\n\tnnzjplusn=CUTEr_nnzj+CUTEr_nvar; \n\tsi=(npy_integer *)malloc(nnzjplusn*sizeof(npy_integer));\n\tsfi=(npy_integer *)malloc(nnzjplusn*sizeof(npy_integer));\n\tsv=(npy_double *)malloc(nnzjplusn*sizeof(npy_double)); \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CSGR\\n");\n#endif\n\t/* Must use different variable for output NNZJ and input LCJAC */\n\tif (!lagrangian) {\n\t\tCSGR((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingFalse, (integer *)&CUTEr_ncon,\n\t\t\t\tNULL, x, (integer *)&nnzjplusno, (integer *)&nnzjplusn, sv, (integer *)si, (integer *)sfi);\n\t} else {\n\t\tCSGR((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, &somethingTrue, (integer *)&CUTEr_ncon,\n\t\t\t\tv, x, (integer *)&nnzjplusno, (integer *)&nnzjplusn, sv, (integer *)si, (integer *)sfi);\n\t}\n\t\n\textract_sparse_gradient_jacobian(nnzjplusno, si, sfi, sv, (PyArrayObject **)&Mgi, (PyArrayObject **)&Mgv, (PyArrayObject **)&MJi, (PyArrayObject **)&MJfi, (PyArrayObject **)&MJv); \n\t\n\t/* Free temporary storage */\n\tfree(si);\n\tfree(sfi);\n\tfree(sv);\n\t\n\treturn decRefTuple(PyTuple_Pack(5, Mgi, Mgv, MJi, MJfi, MJv)); \n}\n\n\nstatic char cuter__sphess_doc[]=\n"Returns the sparse Hessian of the objective at x (unconstrained problems) or\\n"\n"the sparse Hessian of the Lagrangian (constrained problems) at (x, v).\\n"\n"\\n"\n"(Hi, Hj, Hv)=_sphess(x) -- Hessian of objective (unconstrained problems)\\n"\n"(Hi, Hj, Hv)=_sphess(x, v) -- Hessian of Lagrangian (constrained problems)\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"v -- 1D array of length m with the values of Lagrange multipliers\\n"\n"\\n"\n"Output\\n"\n"Hi -- 1D array of integers holding the row indices (0 .. n-1)\\n"\n" of nozero elements in sparse Hessian\\n"\n"Hj -- 1D array of integers holding the column indices (0 .. n-1)\\n"\n" of nozero elements in sparse Hessian\\n"\n"Hv -- 1D array holding the values of nonzero elements in the sparse Hessian\\n"\n" Has the same length as Hi and Hj.\\n"\n"\\n"\n"Hi, Hj, and Hv represent the full Hessian and not only the diagonal and the\\n"\n"upper triangle. To obtain the Hessian of the objective of constrained\\n"\n"problems use _isphess().\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"wrapper function sphess().\\n"\n"\\n"\n"CUTEr tools used: CSH, USH\\n";\n\nstatic PyObject *cuter__sphess(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *MHi, *MHj, *MHv; \n\tdoublereal *x, *v=NULL, *sv; \n\tnpy_integer *si, *sj, nnzho; \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\tif (!PyArg_ParseTuple(args, "O|O", &arg1, &arg2)) \n\t\treturn NULL; \n\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\tif (CUTEr_ncon>0) {\n\t\t/* Check if v is double and of correct dimension */\n\t\tif (arg2!=NULL) {\n\t\t\tif (!(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length ncon");\n\t\t\t\treturn NULL; \n\t\t\t}\n\t\t} else {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be specified for constrained problems. Use _isphess().");\n\t\t\treturn NULL; \n\t\t}\n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tif (CUTEr_ncon>0) \n\t\tv=(npy_double *)PyArray_DATA(arg2);\n\tsi=(npy_integer *)malloc(CUTEr_nnzh*sizeof(npy_integer));\n\tsj=(npy_integer *)malloc(CUTEr_nnzh*sizeof(npy_integer));\n\tsv=(npy_double *)malloc(CUTEr_nnzh*sizeof(npy_double)); \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling [CU]SH\\n");\n#endif\n\tif (CUTEr_ncon>0) {\n\t\tCSH((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, (integer *)&CUTEr_ncon, v, (integer *)&nnzho, (integer *)&CUTEr_nnzh, \n\t\t\tsv, (integer *)si, (integer *)sj);\n\t} else {\n\t\tUSH((integer *)&CUTEr_nvar, x, (integer *)&nnzho, (integer *)&CUTEr_nnzh, \n\t\t\tsv, (integer *)si, (integer *)sj);\n\t}\n\t\n\textract_sparse_hessian(nnzho, si, sj, sv, (PyArrayObject **)&MHi, (PyArrayObject **)&MHj, (PyArrayObject **)&MHv); \n\t\n\t/* Free temporary storage */\n\tfree(si);\n\tfree(sj);\n\tfree(sv);\n\t\n\treturn decRefTuple(PyTuple_Pack(3, MHi, MHj, MHv)); \n}\n\n\nstatic char cuter__isphess_doc[]=\n"Returns the sparse Hessian of the objective or the sparse Hessian of i-th\\n"\n"constraint at x.\\n"\n"\\n"\n"(Hi, Hj, Hv)=_isphess(x) -- Hessian of objective\\n"\n"(Hi, Hj, Hv)=_isphess(x, i) -- Hessian of i-th constraint\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"i -- integer holding the index of constraint (between 0 and m-1)\\n"\n"\\n"\n"Output\\n"\n"Hi -- 1D array of integers holding the row indices (0 .. n-1)\\n"\n" of nozero elements in sparse Hessian\\n"\n"Hj -- 1D array of integers holding the column indices (0 .. n-1)\\n"\n" of nozero elements in sparse Hessian\\n"\n"Hv -- 1D array holding the values of nonzero elements in the sparse Hessian\\n"\n" Has the same length as Hi and Hj.\\n"\n"\\n"\n"Hi, Hj, and Hv represent the full Hessian and not only the diagonal and the\\n"\n"upper triangle.\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"wrapper function isphess().\\n"\n"\\n"\n"CUTEr tools used: CISH, USH\\n";\n\nstatic PyObject *cuter__isphess(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *MHi, *MHj, *MHv; \n\tdoublereal *x, *sv; \n\tnpy_integer *si, *sj, nnzho, i; \n\tnpy_integer icon; \n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\tif (!PyArg_ParseTuple(args, "O|i", &arg1, &i)) \n\t\treturn NULL; \n\t\n\tif (PyObject_Length(args)>1) {\n\t\ticon=i+1;\n\t\tif (i<0 || i>=CUTEr_ncon) {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be between 0 and ncon-1");\n\t\t\treturn NULL; \n\t\t}\n\t} else {\n\t\ticon=0; \n\t}\n\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tsi=(npy_integer *)malloc(CUTEr_nnzh*sizeof(npy_integer));\n\tsj=(npy_integer *)malloc(CUTEr_nnzh*sizeof(npy_integer));\n\tsv=(npy_double *)malloc(CUTEr_nnzh*sizeof(npy_double)); \n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: calling CISH/USH\\n");\n#endif\n\tif (CUTEr_ncon>0) {\n\t\tCISH((integer *)&CUTEr_nvar, x, (integer *)&icon, (integer *)&nnzho, (integer *)&CUTEr_nnzh, sv, (integer *)si, (integer *)sj);\n\t} else {\n\t\tUSH((integer *)&CUTEr_nvar, x, (integer *)&nnzho, (integer *)&CUTEr_nnzh, sv, (integer *)si, (integer *)sj);\n\t}\n\t\n\textract_sparse_hessian(nnzho, si, sj, sv, (PyArrayObject **)&MHi, (PyArrayObject **)&MHj, (PyArrayObject **)&MHv); \n\t\n\t/* Free temporary storage */\n\tfree(si);\n\tfree(sj);\n\tfree(sv);\n\t\n\treturn decRefTuple(PyTuple_Pack(3, MHi, MHj, MHv)); \n}\n\n\nstatic char cuter__gradsphess_doc[]=\n"Returns the sparse Hessian of the Lagrangian, the sparse Jacobian of\\n"\n"constraints, and the gradient of the objective or Lagrangian.\\n" \n"\\n"\n"(g, Hi, Hj, Hv)=_gradsphess(x) -- unconstrained problems\\n"\n"(gi, gv, Jvi, Jfi, Jv, Hi, Hj, Hv)=_gradsphess(x, v, gradl)\\n"\n" -- constrained problems\\n"\n"\\n"\n"Input\\n"\n"x -- 1D array of length n with the values of variables\\n"\n"v -- 1D array of length m holding the values of Lagrange multipliers\\n"\n"gradl -- boolean flag. If False the gradient of the objective is returned, \\n"\n" if True the gradient of the Lagrangian is returned. Default is False.\\n"\n"\\n"\n"Output\\n"\n"g -- 1D array of length n with the gradient of objective or Lagrangian\\n"\n"Hi -- 1D array of integers holding the row indices (0 .. n-1)\\n"\n" of nozero elements in sparse Hessian\\n"\n"Hj -- 1D array of integers holding the column indices (0 .. n-1)\\n"\n" of nozero elements in sparse Hessian\\n"\n"Hv -- 1D array holding the values of nonzero elements in the sparse Hessian\\n"\n" Has the same length as Hi and Hj.\\n"\n"gi -- 1D array of integers holding the indices (0 .. n-1) of nonzero\\n"\n" elements in the sparse gradient vector\\n"\n"gv -- 1D array holding the values of nonzero elements in the sparse gradient\\n"\n" vector. Has the same length as gi.\\n"\n"Jvi -- 1D array of integers holding the column indices (0 .. n-1)\\n"\n" of nozero elements in sparse Jacobian of constraints\\n"\n"Jfi -- 1D array of integers holding the row indices (0 .. m-1)\\n"\n" of nozero elements in sparse Jacobian of constraints\\n"\n"Jv -- 1D array holding the values of nonzero elements in the sparse Jacobian\\n"\n" of constraints at x. Has the same length as Jvi and Jfi.\\n"\n"\\n"\n"For constrained problems the gradient is returned in sparse format.\\n"\n"\\n"\n"Hi, Hj, and Hv represent the full Hessian and not only the diagonal and the\\n"\n"upper triangle.\\n"\n"\\n"\n"This function is not supposed to be called by the user. It is called by the\\n"\n"wrapper function gradsphess().\\n"\n"\\n"\n"CUTEr tools used: CSGRSH, UGRSH\\n";\n\nstatic PyObject *cuter__gradsphess(PyObject *self, PyObject *args) {\n\tPyArrayObject *arg1, *arg2, *Mg=NULL, *Mgi, *Mgv, *MJi, *MJfi, *MJv, *MHi, *MHj, *MHv; \n\tPyObject *arg3;\n\tdoublereal *x, *v, *g, *sv, *sjv; \n\tnpy_integer lagrangian;\n\tnpy_integer *si, *sj, *sji, *sjfi, nnzho, nnzjplusn, nnzjplusno; \n\tnpy_intp dims[1]; \n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\n\targ2=NULL;\n\targ3=NULL;\n\tif (!PyArg_ParseTuple(args, "O|OO", &arg1, &arg2, &arg3)) \n\t\treturn NULL; \t\n\t\n\t/* Check bool argument */\n\tif (arg3!=NULL && arg3==Py_True)\n\t\tlagrangian=1;\n\telse\n\t\tlagrangian=0;\n\t\n\t/* Check if x is double and of correct dimension */\n\tif (!(PyArray_Check(arg1) && PyArray_ISFLOAT(arg1) && PyArray_TYPE(arg1)==NPY_DOUBLE && PyArray_NDIM(arg1)==1 && PyArray_DIM(arg1, 0)==CUTEr_nvar)) {\n\t\tPyErr_SetString(PyExc_Exception, "Argument 1 must be a 1D double array of length nvar");\n\t\treturn NULL; \n\t}\n\t\n\tif (CUTEr_ncon>0) {\n\t\t/* Check if v is double and of correct dimension */\n\t\tif (arg2!=NULL) {\n\t\t\t/* Check if v is double and of correct dimension */\n\t\t\tif (!(PyArray_Check(arg2) && PyArray_ISFLOAT(arg2) && PyArray_TYPE(arg2)==NPY_DOUBLE && PyArray_NDIM(arg2)==1 && PyArray_DIM(arg2, 0)==CUTEr_ncon)) {\n\t\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be a 1D double array of length ncon");\n\t\t\t\treturn NULL; \n\t\t\t}\n\t\t} else {\n\t\t\tPyErr_SetString(PyExc_Exception, "Argument 2 must be specified for constrained problems.");\n\t\t\treturn NULL; \n\t\t}\n\t}\n\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: preparing for evaluation\\n");\n#endif\n\tx=(npy_double *)PyArray_DATA(arg1);\n\tsi=(npy_integer *)malloc(CUTEr_nnzh*sizeof(npy_integer));\n\tsj=(npy_integer *)malloc(CUTEr_nnzh*sizeof(npy_integer));\n\tsv=(npy_double *)malloc(CUTEr_nnzh*sizeof(npy_double)); \n\t\t\n\tif (CUTEr_ncon>0) {\n\t\tv=(npy_double *)PyArray_DATA(arg2);\n\t\tnnzjplusn=CUTEr_nnzj+CUTEr_nvar;\n\t\tsji=(npy_integer *)malloc(nnzjplusn*sizeof(npy_integer));\n\t\tsjfi=(npy_integer *)malloc(nnzjplusn*sizeof(npy_integer));\n\t\tsjv=(npy_double *)malloc(nnzjplusn*sizeof(npy_double)); \n\t\t\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: calling CSGRSH\\n");\n#endif\n\t\tif (lagrangian) {\n\t\t\tCSGRSH((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, &somethingTrue, (integer *)&CUTEr_ncon, v,\n\t\t\t\t\t(integer *)&nnzjplusno, (integer *)&nnzjplusn, sjv, (integer *)sji, (integer *)sjfi, \n\t\t\t\t\t(integer *)&nnzho, (integer *)&CUTEr_nnzh, sv, (integer *)si, (integer *)sj);\n\t\t} else { \n\t\t\tCSGRSH((integer *)&CUTEr_nvar, (integer *)&CUTEr_ncon, x, &somethingFalse, (integer *)&CUTEr_ncon, v,\n\t\t\t\t\t(integer *)&nnzjplusno, (integer *)&nnzjplusn, sjv, (integer *)sji, (integer *)sjfi, \n\t\t\t\t\t(integer *)&nnzho, (integer *)&CUTEr_nnzh, sv, (integer *)si, (integer *)sj);\n\t\t}\n\t\t\n\t\textract_sparse_gradient_jacobian(nnzjplusno, sji, sjfi, sjv, (PyArrayObject **)&Mgi, (PyArrayObject **)&Mgv, (PyArrayObject **)&MJi, (PyArrayObject **)&MJfi, (PyArrayObject **)&MJv); \n\t\t\n\t\t/* Free temporary storage - Jacobian */\n\t\tfree(sji);\n\t\tfree(sjfi);\n\t\tfree(sjv);\n\t} else {\n\t\tdims[0]=CUTEr_nvar; \n\t\tMg=(PyArrayObject *)PyArray_SimpleNew(1, dims, NPY_DOUBLE); \n\t\tg=(npy_double *)PyArray_DATA(Mg);\n\t\t\n#ifdef PYDEBUG\n\t\tfprintf(df, "PyCUTEr: calling UGRSH\\n");\n#endif\n\t\tUGRSH((integer *)&CUTEr_nvar, x, g, (integer *)&nnzho, (integer *)&CUTEr_nnzh, sv, (integer *)si, (integer *)sj);\n\t}\n\t\n\textract_sparse_hessian(nnzho, si, sj, sv, (PyArrayObject **)&MHi, (PyArrayObject **)&MHj, (PyArrayObject **)&MHv); \n\t\n\t/* Free temporary storage - Hessian */\n\tfree(si);\n\tfree(sj);\n\tfree(sv);\n\t\t\n\tif (CUTEr_ncon>0) {\n\t\treturn decRefTuple(PyTuple_Pack(8, Mgi, Mgv, MJi, MJfi, MJv, MHi, MHj, MHv)); \n\t} else {\n\t\treturn decRefTuple(PyTuple_Pack(4, Mg, MHi, MHj, MHv)); \n\t}\n}\n\n\nstatic char cuter_report_doc[]=\n"Reports usage statistics.\\n"\n"\\n"\n"stat=report()\\n"\n"\\n"\n"Output\\n"\n"stat -- dictionary with the usage statistics\\n"\n"\\n"\n"The usage statistics dictionary has the following members:\\n"\n"f -- number of objective evaluations\\n"\n"g -- number of objective gradient evaluations\\n"\n"H -- number of objective Hessian evaluations\\n"\n"Hprod -- number of Hessian multiplications with a vector\\n"\n"tsetup -- CPU time used in setup\\n"\n"trun -- CPU time used in run\\n"\n"\\n"\n"For constrained problems the following additional members are available\\n"\n"c -- number of constraint evaluations\\n"\n"cg -- number of constraint gradient evaluations\\n"\n"cH -- number of constraint Hessian evaluations\\n"\n"\\n"\n"CUTEr tools used: CREPRT, UREPRT\\n";\n\nstatic PyObject *cuter_report(PyObject *self, PyObject *args) {\n\treal calls[7], time[2];\n\tPyObject *dict;\n\t\n#ifdef PYDEBUG\n\tfprintf(df, "PyCUTEr: checking arguments\\n");\n#endif\n\tif (!check_setup())\n\t\treturn NULL; \n\t\t\n\tif (PyObject_Length(args)!=0) {\n\t\tPyErr_SetString(PyExc_Exception, "report() takes no arguments");\n\t\treturn NULL; \n\t}\n\t\n\tif (CUTEr_ncon>0) \n\t\tCREPRT(calls, time); \n\telse\n\t\tUREPRT(calls, time); \n\t\t\n\tdict=PyDict_New();\n\tPyDict_SetItemString(dict, "f", PyLong_FromLong((long)(calls[0]))); \n\tPyDict_SetItemString(dict, "g", PyLong_FromLong((long)(calls[1]))); \n\tPyDict_SetItemString(dict, "H", PyLong_FromLong((long)(calls[2]))); \n\tPyDict_SetItemString(dict, "Hprod", PyLong_FromLong((long)(calls[3]))); \n\tif (CUTEr_ncon>0) {\n\t\tPyDict_SetItemString(dict, "c", PyLong_FromLong((long)(calls[4]))); \n\t\tPyDict_SetItemString(dict, "cg", PyLong_FromLong((long)(calls[5]))); \n\t\tPyDict_SetItemString(dict, "cH", PyLong_FromLong((long)(calls[6]))); \n\t}\n\tPyDict_SetItemString(dict, "tsetup", PyFloat_FromDouble((long)(time[0]))); \n\tPyDict_SetItemString(dict, "trun", PyFloat_FromDouble((long)(time[1]))); \n\t\n\treturn decRefDict(dict); \n}\n\n/* Methods table */\nstatic PyMethodDef module_methods[] = {\n\t{"_dims", cuter__dims, METH_VARARGS, cuter__dims_doc},\n\t{"_setup", cuter__setup, METH_VARARGS, cuter__setup_doc},\n\t{"_varnames", cuter__varnames, METH_VARARGS, cuter__varnames_doc},\n\t{"_connames", cuter__connames, METH_VARARGS, cuter__connames_doc},\n\t{"objcons", cuter_objcons, METH_VARARGS, cuter_objcons_doc},\n\t{"obj", cuter_obj, METH_VARARGS, cuter_obj_doc},\n\t{"cons", cuter_cons, METH_VARARGS, cuter_cons_doc},\n\t{"lagjac", cuter_lagjac, METH_VARARGS, cuter_lagjac_doc},\n\t{"jprod", cuter_jprod, METH_VARARGS, cuter_jprod_doc},\n\t{"hess", cuter_hess, METH_VARARGS, cuter_hess_doc},\n\t{"ihess", cuter_ihess, METH_VARARGS, cuter_ihess_doc},\n\t{"hprod", cuter_hprod, METH_VARARGS, cuter_hprod_doc},\n\t{"gradhess", cuter_gradhess, METH_VARARGS, cuter_gradhess_doc},\n\t{"_scons", cuter__scons, METH_VARARGS, cuter__scons_doc},\n\t{"_slagjac", cuter__slagjac, METH_VARARGS, cuter__slagjac_doc},\n\t{"_sphess", cuter__sphess, METH_VARARGS, cuter__sphess_doc},\n\t{"_isphess", cuter__isphess, METH_VARARGS, cuter__isphess_doc},\n\t{"_gradsphess", cuter__gradsphess, METH_VARARGS, cuter__gradsphess_doc},\n\t{"report", cuter_report, METH_VARARGS, cuter_report_doc}, \n\t{NULL, NULL, 0, NULL} /* Marks the end of this structure */\n};\n\nstatic struct PyModuleDef module = {\n PyModuleDef_HEAD_INIT,\n "_pycuteritf", /* name of module */\n "Interface to a CUTEr function.\\n", /* module documentation, may be NULL */\n /* TODO: in future set this to 0 so that this module will work with sub-interpreters */\n -1, /* size of per-interpreter state of the module, or -1 if the module keeps state in global variables. */\n module_methods,\n NULL, /* Slots for multi phase initialization */\n NULL, /* Traversal function for GC */\n NULL, /* Clear function for clearing the module */\n NULL, /* Function for deallocating the module */\n};\n\n/* Module initialization \n Module name must be _pycuteritf in compile and link */\n\n/* Prototype */\nPyMODINIT_FUNC PyInit__pycuteritf(void);\n\n/* Definition */\nPyMODINIT_FUNC PyInit__pycuteritf() {\n\tPyObject *m=PyModule_Create(&module);\n\t\n\t/* For using NumPy */\n\timport_array(); \n\t\n\treturn m;\n}\n\n\n#ifdef __cplusplus\n}\n#endif\n'

CUTEr problem binary iterface source code.