예제 #1
0
파일: umfpack.c 프로젝트: MGKhKhD/cvxopt
static PyObject* symbolic(PyObject *self, PyObject *args)
{
    spmatrix *A;
    double info[UMFPACK_INFO];
    void *symbolic;

    if (!PyArg_ParseTuple(args, "O", &A)) return NULL;

    if (!SpMatrix_Check(A)) PY_ERR_TYPE("A must be a sparse matrix");
    if (SP_NCOLS(A) == 0 || SP_NROWS(A) == 0) {
        PyErr_SetString(PyExc_ValueError, "A must have at least one "
            "row and column");
        return NULL;
    }

    switch (SP_ID(A)){
        case DOUBLE:
            UMFD(symbolic)(SP_NROWS(A), SP_NCOLS(A), SP_COL(A),
                SP_ROW(A), SP_VAL(A), &symbolic, NULL, info);
            if (info[UMFPACK_STATUS] == UMFPACK_OK)
#if PY_MAJOR_VERSION >= 3
                return (PyObject *) PyCapsule_New( (void *) symbolic, 
                    "UMFPACK SYM D FACTOR", 
                    (PyCapsule_Destructor) &free_umfpack_d_symbolic);  
#else
                return (PyObject *) PyCObject_FromVoidPtrAndDesc(
                    (void *) symbolic, "UMFPACK SYM D FACTOR",
                    free_umfpack_d_symbolic);
#endif
            else
                UMFD(free_symbolic)(&symbolic);
            break;

        case COMPLEX:
            UMFZ(symbolic)(SP_NROWS(A), SP_NCOLS(A), SP_COL(A),
                SP_ROW(A), SP_VAL(A), NULL, &symbolic, NULL, info);
            if (info[UMFPACK_STATUS] == UMFPACK_OK)
#if PY_MAJOR_VERSION >= 3
                return (PyObject *) PyCapsule_New(
                    (void *) symbolic, "UMFPACK SYM Z FACTOR",
                    (PyCapsule_Destructor) &free_umfpack_z_symbolic);
#else
                return (PyObject *) PyCObject_FromVoidPtrAndDesc(
                    (void *) symbolic, "UMFPACK SYM Z FACTOR",
                    free_umfpack_z_symbolic);
#endif
            else
예제 #2
0
static cholmod_sparse * create_matrix(spmatrix *A)
{
    cholmod_sparse *B;

    if (!(B = CHOL(allocate_sparse)(SP_NROWS(A), SP_NCOLS(A), 0,
        1, 0, 0, (SP_ID(A) == DOUBLE ? CHOLMOD_REAL : CHOLMOD_COMPLEX),
        &Common))) return NULL;

    int i;
    for (i=0; i<SP_NCOLS(A); i++)
        ((int_t *)B->nz)[i] = SP_COL(A)[i+1] - SP_COL(A)[i];
    B->x = SP_VAL(A);
    B->i = SP_ROW(A);
    B->nzmax = SP_NNZ(A);
    memcpy(B->p, SP_COL(A), (SP_NCOLS(A)+1)*sizeof(int_t));
    return B;
}
예제 #3
0
파일: misc.c 프로젝트: cvxopt/smcp
static PyObject *phase1_sdp
(PyObject *self, PyObject *args, PyObject *kwrds) {

  matrix *u;
  spmatrix *Ai,*Ao;
  int_t k,i,j,n,m,nnz,nz,col;


  if(!PyArg_ParseTuple(args,"OO",&Ai,&u)) return NULL;
  n = (int_t) sqrt((double)SP_NROWS(Ai));
  m = SP_NCOLS(Ai) - 1;
  nnz = SP_NNZ(Ai) - SP_COL(Ai)[1] + 1 + m + n + 1;

  Ao = SpMatrix_New((n+2)*(n+2),m+2,nnz,DOUBLE);
  if (!Ao) return PyErr_NoMemory();

  // A_0
  SP_VALD(Ao)[0] = 1.0;
  SP_ROW(Ao)[0] = n*(n+2)+n;
  SP_COL(Ao)[0] = 0;
  SP_COL(Ao)[1] = 1;

  // A_i, i=1,..,m
  for (i=1;i<=m;i++){
    k = SP_COL(Ao)[i];
    nz = SP_COL(Ai)[i+1]-SP_COL(Ai)[i]; // nonzeros in Ai
    // copy Ai
    memcpy(SP_VALD(Ao)+k,SP_VALD(Ai)+SP_COL(Ai)[i],nz*sizeof(double));
    // insert -u[i]
    SP_VALD(Ao)[k+nz] = -MAT_BUFD(u)[i-1];
    // update row colptr
    SP_COL(Ao)[i+1] = SP_COL(Ao)[i]+nz+1;
    // generate row indices
    for (j=0;j<nz;j++) {
      col = SP_ROW(Ai)[SP_COL(Ai)[i]+j] / n;
      SP_ROW(Ao)[k+j] = SP_ROW(Ai)[SP_COL(Ai)[i]+j] + col*2;
    }
    SP_ROW(Ao)[k+nz] = n*(n+2)+n;
  }
  // last constraint
  k = SP_COL(Ao)[m+1];
  for (i=0;i<n;i++){
    SP_VALD(Ao)[k+i] = 1.0;
    SP_ROW(Ao)[k+i] = i*(n+2)+i;
  }
  SP_VALD(Ao)[k+n] = 1.0;
  SP_ROW(Ao)[k+n] = (n+2)*(n+2)-1;
  SP_COL(Ao)[m+2] = SP_COL(Ao)[m+1] + n + 1;

  return (PyObject*) Ao;
}
예제 #4
0
파일: misc.c 프로젝트: cvxopt/smcp
static PyObject *nzcolumns
(PyObject *self, PyObject *args)
{
  PyObject *A;
  matrix *Nz;
  int_t m,n,i,j,p,nnz,sum, *tmp;

  if (!PyArg_ParseTuple(args,"O",&A)) return NULL;

  n = (int_t) sqrt((double)SP_NROWS(A));
  m = SP_NCOLS(A)-1;

  Nz = Matrix_New(m,1,INT);
  if (!Nz) return PyErr_NoMemory();

  tmp = malloc(n*sizeof(int_t));
  //tmp = Matrix_New(n,1,INT);
  if (!tmp) return PyErr_NoMemory();

  // erase workspace
  for (i=0;i<n;i++) tmp[i] = 0;

  for (j=0;j<m;j++){
    p = SP_COL(A)[j+1];
    nnz = SP_COL(A)[j+2]-p;
    if (nnz) {
      // Find nonzero cols
      for (i=0;i<nnz;i++) {
	tmp[SP_ROW(A)[p+i] % n] += 1;
	tmp[SP_ROW(A)[p+i] / n] += 1;
      }
      // Count nonzero cols and reset workspace
      MAT_BUFI(Nz)[j] = 0;
      sum = 0;
#pragma omp parallel for shared(tmp,Nz,j,n) private(i) reduction(+:sum)
      for (i=0;i<n;i++) {
	if(tmp[i]) {
	  tmp[i] = 0;
	  sum++;
	}
      }
      MAT_BUFI(Nz)[j] = sum;
    }
  }

  free(tmp);

  return (PyObject*) Nz;
}
예제 #5
0
static PyObject* splinsolve(PyObject *self, PyObject *args,
    PyObject *kwrds)
{
    spmatrix *A, *B, *X;
    matrix *P=NULL;
    int n, nnz;
    cholmod_sparse *Ac=NULL, *Bc=NULL, *Xc=NULL;
    cholmod_factor *L=NULL;
#if PY_MAJOR_VERSION >= 3
    int uplo_='L';
#endif
    char uplo='L';
    char *kwlist[] = {"A", "B", "p", "uplo", NULL};

    if (!set_options()) return NULL;
#if PY_MAJOR_VERSION >= 3
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OO|OC", kwlist, &A,
        &B, &P, &uplo_)) return NULL;
    uplo = (char) uplo_;
#else
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OO|Oc", kwlist, &A,
        &B, &P, &uplo)) return NULL;
#endif

    if (!SpMatrix_Check(A) || SP_NROWS(A) != SP_NCOLS(A))
        PY_ERR_TYPE("A is not a square sparse matrix");
    n = SP_NROWS(A);
    nnz = SP_NNZ(A);

    if (!SpMatrix_Check(B) || SP_ID(A) != SP_ID(B))
        PY_ERR_TYPE("B must be a sparse matrix of the same type as A");
    if (SP_NROWS(B) != n)
        PY_ERR(PyExc_ValueError, "incompatible dimensions for B");

    if (P) {
        if (!Matrix_Check(P) || MAT_ID(P) != INT) err_int_mtrx("p");
        if (MAT_LGT(P) != n) err_buf_len("p");
        if (!CHOL(check_perm)(P->buffer, n, n, &Common))
            PY_ERR(PyExc_ValueError, "not a valid permutation");
    }

    if (uplo != 'U' && uplo != 'L') err_char("uplo", "'L', 'U'");
    if (!(Ac = pack(A, uplo))) return PyErr_NoMemory();

    L = CHOL(analyze_p) (Ac, P ? MAT_BUFI(P): NULL, NULL, 0, &Common);
    if (Common.status != CHOLMOD_OK){
        CHOL(free_factor)(&L, &Common);
        CHOL(free_sparse)(&Ac, &Common);
        if (Common.status == CHOLMOD_OUT_OF_MEMORY)
            return PyErr_NoMemory();
        else {
            PyErr_SetString(PyExc_ValueError, "symbolic factorization "
                "failed");
            return NULL;
        }
    }

    CHOL(factorize) (Ac, L, &Common);
    CHOL(free_sparse)(&Ac, &Common);
    if (Common.status > 0) switch (Common.status) {
        case CHOLMOD_NOT_POSDEF:
            PyErr_SetObject(PyExc_ArithmeticError, Py_BuildValue("i",
                L->minor));
            CHOL(free_factor)(&L, &Common);
            return NULL;
            break;

        case CHOLMOD_DSMALL:
            /* This never happens unless we change the default value
             * of Common.dbound (0.0).  */
            if (L->is_ll)
                PyErr_Warn(PyExc_RuntimeWarning, "tiny diagonal "
                    "elements in L");
            else
                PyErr_Warn(PyExc_RuntimeWarning, "tiny diagonal "
                    "elements in D");
            break;

        default:
            PyErr_Warn(PyExc_UserWarning, "");
    }

    if (L->minor<n) {
        CHOL(free_factor)(&L, &Common);
        PY_ERR(PyExc_ArithmeticError, "singular matrix");
    }
    if (!(Bc = create_matrix(B))) {
      CHOL(free_factor)(&L, &Common);
      return PyErr_NoMemory();
    }

    Xc = CHOL(spsolve)(0, L, Bc, &Common);
    free_matrix(Bc);
    CHOL(free_factor)(&L, &Common);
    if (Common.status != CHOLMOD_OK){
        CHOL(free_sparse)(&Xc, &Common);
        if (Common.status == CHOLMOD_OUT_OF_MEMORY)
            return PyErr_NoMemory();
        else
            PY_ERR(PyExc_ValueError, "solve step failed");
    }

    if (!(X = SpMatrix_New(Xc->nrow, Xc->ncol,
        ((int_t*)Xc->p)[Xc->ncol], SP_ID(A)))) {
        CHOL(free_sparse)(&Xc, &Common);
        return PyErr_NoMemory();
    }
    memcpy(SP_COL(X), (int_t *) Xc->p, (Xc->ncol+1)*sizeof(int_t));
    memcpy(SP_ROW(X), (int_t *) Xc->i,
        ((int_t *) Xc->p)[Xc->ncol]*sizeof(int_t));
    memcpy(SP_VAL(X), (double *) Xc->x,
        ((int_t *) Xc->p)[Xc->ncol]*E_SIZE[SP_ID(X)]);
    CHOL(free_sparse)(&Xc, &Common);
    return (PyObject *) X;
}
예제 #6
0
static PyObject* linsolve(PyObject *self, PyObject *args,
    PyObject *kwrds)
{
    spmatrix *A;
    matrix *B, *P=NULL;
    int i, n, nnz, oB=0, ldB=0, nrhs=-1;
    cholmod_sparse *Ac=NULL;
    cholmod_factor *L=NULL;
    cholmod_dense *x=NULL, *b=NULL;
    void *b_old;
#if PY_MAJOR_VERSION >= 3
    int uplo_ = 'L';
#endif
    char uplo='L';
    char *kwlist[] = {"A", "B", "p", "uplo", "nrhs", "ldB", "offsetB",
        NULL};

    if (!set_options()) return NULL;
#if PY_MAJOR_VERSION >= 3
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OO|OCiii", kwlist,
        &A,  &B, &P, &uplo_, &nrhs, &ldB, &oB)) return NULL;
    uplo = (char) uplo_;
#else
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OO|Ociii", kwlist,
        &A,  &B, &P, &uplo, &nrhs, &ldB, &oB)) return NULL;
#endif

    if (!SpMatrix_Check(A) || SP_NROWS(A) != SP_NCOLS(A))
        PY_ERR_TYPE("A is not a sparse matrix");
    n = SP_NROWS(A);
    nnz = SP_NNZ(A);

    if (!Matrix_Check(B) || MAT_ID(B) != SP_ID(A))
        PY_ERR_TYPE("B must be a dense matrix of the same numerical "
            "type as A");
    if (nrhs < 0) nrhs = MAT_NCOLS(B);
    if (n == 0 || nrhs == 0) return Py_BuildValue("");
    if (ldB == 0) ldB = MAX(1,MAT_NROWS(B));
    if (ldB < MAX(1,n)) err_ld("ldB");
    if (oB < 0) err_nn_int("offsetB");
    if (oB + (nrhs-1)*ldB + n > MAT_LGT(B)) err_buf_len("B");

    if (P) {
        if (!Matrix_Check(P) || MAT_ID(P) != INT) err_int_mtrx("p");
        if (MAT_LGT(P) != n) err_buf_len("p");
        if (!CHOL(check_perm)(P->buffer, n, n, &Common))
            PY_ERR(PyExc_ValueError, "not a valid permutation");
    }
    if (uplo != 'U' && uplo != 'L') err_char("uplo", "'L', 'U'");

    if (!(Ac = pack(A, uplo))) return PyErr_NoMemory();
    L = CHOL(analyze_p)(Ac, P ? MAT_BUFI(P): NULL, NULL, 0, &Common);
    if (Common.status != CHOLMOD_OK){
        free_matrix(Ac);
        CHOL(free_sparse)(&Ac, &Common);
        CHOL(free_factor)(&L, &Common);
        if (Common.status == CHOLMOD_OUT_OF_MEMORY)
            return PyErr_NoMemory();
        else {
            PyErr_SetString(PyExc_ValueError, "symbolic factorization "
                "failed");
            return NULL;
        }
    }

    CHOL(factorize) (Ac, L, &Common);
    CHOL(free_sparse)(&Ac, &Common);
    if (Common.status < 0) {
        CHOL(free_factor)(&L, &Common);
        switch (Common.status) {
            case CHOLMOD_OUT_OF_MEMORY:
                return PyErr_NoMemory();

            default:
                PyErr_SetString(PyExc_ValueError, "factorization "
                    "failed");
                return NULL;
        }
    }
    if (Common.status > 0) switch (Common.status) {
        case CHOLMOD_NOT_POSDEF:
            PyErr_SetObject(PyExc_ArithmeticError,
                Py_BuildValue("i", L->minor));
            CHOL(free_factor)(&L, &Common);
            return NULL;
            break;

        case CHOLMOD_DSMALL:
            /* This never happens unless we change the default value
             * of Common.dbound (0.0).  */
            if (L->is_ll)
                PyErr_Warn(PyExc_RuntimeWarning, "tiny diagonal "
                    "elements in L");
            else
                PyErr_Warn(PyExc_RuntimeWarning, "tiny diagonal "
                    "elements in D");
            break;

        default:
            PyErr_Warn(PyExc_UserWarning, "");
    }

    if (L->minor<n) {
        CHOL(free_factor)(&L, &Common);
        PY_ERR(PyExc_ArithmeticError, "singular matrix");
    }
    b = CHOL(allocate_dense)(n, 1, n, (MAT_ID(B) == DOUBLE ?
        CHOLMOD_REAL : CHOLMOD_COMPLEX) , &Common);
    if (Common.status == CHOLMOD_OUT_OF_MEMORY) {
        CHOL(free_factor)(&L, &Common);
        CHOL(free_dense)(&b, &Common);
        return PyErr_NoMemory();
    }
    b_old = b->x;
    for (i=0; i<nrhs; i++) {
        b->x = MAT_BUF(B) + (i*ldB + oB)*E_SIZE[MAT_ID(B)];
        x = CHOL(solve) (CHOLMOD_A, L, b, &Common);
        if (Common.status != CHOLMOD_OK){
            PyErr_SetString(PyExc_ValueError, "solve step failed");
            CHOL(free_factor)(&L, &Common);
            b->x = b_old;
            CHOL(free_dense)(&b, &Common);
            CHOL(free_dense)(&x, &Common);
            return NULL;
        }
        memcpy(b->x, x->x, SP_NROWS(A)*E_SIZE[MAT_ID(B)]);
        CHOL(free_dense)(&x, &Common);
    }
    b->x = b_old;
    CHOL(free_dense)(&b, &Common);
    CHOL(free_factor)(&L, &Common);
    return Py_BuildValue("");
}
예제 #7
0
static PyObject* numeric(PyObject *self, PyObject *args)
{
    spmatrix *A;
    PyObject *F;
    cholmod_factor *Lc;
    cholmod_sparse *Ac = NULL;
    char uplo;
#if PY_MAJOR_VERSION >= 3
    const char *descr; 
#else
    char *descr; 
#endif

    if (!set_options()) return NULL;

    if (!PyArg_ParseTuple(args, "OO", &A, &F)) return NULL;

    if (!SpMatrix_Check(A) || SP_NROWS(A) != SP_NCOLS(A))
        PY_ERR_TYPE("A is not a sparse matrix");

#if PY_MAJOR_VERSION >= 3
    if (!PyCapsule_CheckExact(F) || !(descr = PyCapsule_GetName(F)))
        err_CO("F");
#else
    if (!PyCObject_Check(F)) err_CO("F");
    descr = PyCObject_GetDesc(F);
    if (!descr) PY_ERR_TYPE("F is not a CHOLMOD factor");
#endif
    if (SP_ID(A) == DOUBLE){
        if (!strcmp(descr, "CHOLMOD FACTOR D L"))
	    uplo = 'L';
	else if (!strcmp(descr, "CHOLMOD FACTOR D U"))
	    uplo = 'U';
        else
	    PY_ERR_TYPE("F is not the CHOLMOD factor of a 'd' matrix");
    } else {
        if (!strcmp(descr, "CHOLMOD FACTOR Z L"))
	    uplo = 'L';
	else if (!strcmp(descr, "CHOLMOD FACTOR Z U"))
	    uplo = 'U';
        else
	    PY_ERR_TYPE("F is not the CHOLMOD factor of a 'z' matrix");
    }
#if PY_MAJOR_VERSION >= 3
    Lc = (cholmod_factor *) PyCapsule_GetPointer(F, descr);
#else
    Lc = (cholmod_factor *) PyCObject_AsVoidPtr(F);
#endif
    if (!(Ac = pack(A, uplo))) return PyErr_NoMemory();
    CHOL(factorize) (Ac, Lc, &Common);
    CHOL(free_sparse)(&Ac, &Common);

    if (Common.status < 0) switch (Common.status) {
        case CHOLMOD_OUT_OF_MEMORY:
            return PyErr_NoMemory();

        default:
            PyErr_SetString(PyExc_ValueError, "factorization failed");
            return NULL;
    }

    if (Common.status > 0) switch (Common.status) {
        case CHOLMOD_NOT_POSDEF:
            PyErr_SetObject(PyExc_ArithmeticError, Py_BuildValue("i",
                Lc->minor));
            return NULL;
            break;

        case CHOLMOD_DSMALL:
            /* This never happens unless we change the default value
             * of Common.dbound (0.0).  */
            if (Lc->is_ll)
                PyErr_Warn(PyExc_RuntimeWarning, "tiny diagonal "\
                    "elements in L");
            else
                PyErr_Warn(PyExc_RuntimeWarning, "tiny diagonal "\
                    "elements in D");
            break;

        default:
            PyErr_Warn(PyExc_UserWarning, "");
    }

    return Py_BuildValue("");
}
예제 #8
0
static PyObject* symbolic(PyObject *self, PyObject *args,
    PyObject *kwrds)
{
    spmatrix *A;
    cholmod_sparse *Ac = NULL;
    cholmod_factor *L;
    matrix *P=NULL;
#if PY_MAJOR_VERSION >= 3
    int uplo_='L';
#endif
    char uplo='L';
    int n;
    char *kwlist[] = {"A", "p", "uplo", NULL};

    if (!set_options()) return NULL;

#if PY_MAJOR_VERSION >= 3
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "O|OC", kwlist, &A,
        &P, &uplo_)) return NULL;
    uplo = (char) uplo_;
#else
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "O|Oc", kwlist, &A,
        &P, &uplo)) return NULL;
#endif
    if (!SpMatrix_Check(A) || SP_NROWS(A) != SP_NCOLS(A))
        PY_ERR_TYPE("A is not a square sparse matrix");
    n = SP_NROWS(A);

    if (P) {
        if (!Matrix_Check(P) || MAT_ID(P) != INT) err_int_mtrx("p");
        if (MAT_LGT(P) != n) err_buf_len("p");
        if (!CHOL(check_perm)(P->buffer, n, n, &Common))
            PY_ERR(PyExc_ValueError, "p is not a valid permutation");
    }
    if (uplo != 'U' && uplo != 'L') err_char("uplo", "'L', 'U'");
    if (!(Ac = pack(A, uplo))) return PyErr_NoMemory();
    L = CHOL(analyze_p)(Ac, P ? MAT_BUFI(P): NULL, NULL, 0, &Common);
    CHOL(free_sparse)(&Ac, &Common);

    if (Common.status != CHOLMOD_OK){
        if (Common.status == CHOLMOD_OUT_OF_MEMORY)
            return PyErr_NoMemory();
        else{
            PyErr_SetString(PyExc_ValueError, "symbolic factorization "
                "failed");
            return NULL;
        }
    }
#if PY_MAJOR_VERSION >= 3
    return (PyObject *) PyCapsule_New((void *) L, SP_ID(A)==DOUBLE ?  
        (uplo == 'L' ?  "CHOLMOD FACTOR D L" : "CHOLMOD FACTOR D U") :
        (uplo == 'L' ?  "CHOLMOD FACTOR Z L" : "CHOLMOD FACTOR Z U"),
        (PyCapsule_Destructor) &cvxopt_free_cholmod_factor); 
#else
    return (PyObject *) PyCObject_FromVoidPtrAndDesc((void *) L,
        SP_ID(A)==DOUBLE ?  
        (uplo == 'L' ?  "CHOLMOD FACTOR D L" : "CHOLMOD FACTOR D U") :
        (uplo == 'L' ?  "CHOLMOD FACTOR Z L" : "CHOLMOD FACTOR Z U"),
	cvxopt_free_cholmod_factor);
#endif
}
예제 #9
0
파일: misc.c 프로젝트: cvxopt/smcp
static PyObject *robustLS_to_sdp
(PyObject *self, PyObject *args, PyObject *kwrds) {

  PyObject *Alist,*bt, *Ai;
  spmatrix *A,*b;
  int_t m,n,mp,np,pt,i,j,k,N,nnz=0,ri=0;
  char *kwlist[] = {"Alist","bt",NULL};

  if(!PyArg_ParseTupleAndKeywords(args,kwrds,"OO",kwlist,&Alist,&bt)) return NULL;

  if(!PyList_Check(Alist)) {
    PyErr_SetString(PyExc_TypeError,"Alist must be a list of matrices");
    return NULL;
  }

  // get pt = p + 1
  pt = PyList_Size(Alist);

  // get size of bt
  if(Matrix_Check(bt)){
    m = MAT_NROWS(bt);
    np = MAT_NCOLS(bt);
  }
  else if (SpMatrix_Check(bt)){
    m = SP_NROWS(bt);
    np = SP_NCOLS(bt);
  }
  else {
    PyErr_SetString(PyExc_TypeError,"b must be a vector");
    return NULL;
  }
  if (np!=1) {
    PyErr_SetString(PyExc_TypeError,"b must be a vector");
    return NULL;
  }

  // get n and check A0
  if (!(Ai = PyList_GetItem(Alist,0))) return NULL;
  if (Matrix_Check(Ai)) {
    n = MAT_NCOLS(Ai);
    nnz += m*n;
  }
  else if (SpMatrix_Check(Ai)) {
    n = SP_NCOLS(Ai);
    nnz += SP_NNZ(Ai);
  }
  else {
    PyErr_SetString(PyExc_TypeError,"only spmatrix and matrix types allowed");
    return NULL;
  }

  // check remaining matrices in Alist
  for (i=1;i<pt;i++) {
    if (!(Ai = PyList_GetItem(Alist,i))) return NULL;
    if (Matrix_Check(Ai)) {
      mp = MAT_NROWS(Ai);
      np = MAT_NCOLS(Ai);
      nnz += m*n;
    }
    else if (SpMatrix_Check(Ai)) {
      mp = SP_NROWS(Ai);
      np = SP_NCOLS(Ai);
      nnz += SP_NNZ(Ai);
    }
    else {
      PyErr_SetString(PyExc_TypeError,"only spmatrix and matrix types allowed");
      return NULL;
    }
    if (!(mp==m && np==n)){
      PyErr_SetString(PyExc_TypeError,"matrices in Alist must have same size");
      return NULL;
    }
  }
  nnz += 2*m + pt;

  // generate b
  b = SpMatrix_New(n+2,1,2,DOUBLE);
  if (!b) return PyErr_NoMemory();
  SP_COL(b)[0] = 0;
  SP_VALD(b)[0] = -1;
  SP_ROW(b)[0] = 0;
  SP_VALD(b)[1] = -1;
  SP_ROW(b)[1] = 1;
  SP_COL(b)[1] = 2;

  // generate A
  N = m+pt;
  A = SpMatrix_New(N*N,n+3,nnz,DOUBLE);
  if (!A) return PyErr_NoMemory();

  // build A0
  SP_COL(A)[0] = ri;
  for(i=0;i<m;i++){
    if(SpMatrix_Check(bt)){
      SP_VALD(A)[ri] = -SP_VALD(bt)[i];
      SP_ROW(A)[ri++] = pt+i;
    }
    else{
      SP_VALD(A)[ri] = -MAT_BUFD(bt)[i];
      SP_ROW(A)[ri++] = pt+i;
    }
  }
  for(i=0;i<m;i++) {
    SP_VALD(A)[ri] = 1;
    SP_ROW(A)[ri++] = (N+1)*pt + i*N+i;
  }

  // build A1
  SP_COL(A)[1] = ri;
  for(i=0;i<pt-1;i++){
    SP_VALD(A)[ri] = -1;
    SP_ROW(A)[ri++] = N+1 + i*N+i;
  }

  // build A2
  SP_COL(A)[2] = ri;
  SP_VALD(A)[ri] = -1;
  SP_ROW(A)[ri++] = 0;
  SP_COL(A)[3] = ri;

  // build A3,...
  for(j=0;j<n;j++){
    // generate col. i
    for(i=0;i<pt;i++){
      Ai = PyList_GetItem(Alist,i);
      if(SpMatrix_Check(Ai)) {
	nnz = SP_COL(Ai)[j+1]-SP_COL(Ai)[j];
	for(k=0;k<nnz;k++) {
	  SP_VALD(A)[ri] = -SP_VALD(Ai)[SP_COL(Ai)[j]+k];
	  SP_ROW(A)[ri++] = pt+i*N + SP_ROW(Ai)[SP_COL(Ai)[j]+k];
	}
      }
      else {
	for (k=0;k<m;k++) {
	  SP_VALD(A)[ri] = -MAT_BUFD(Ai)[j*m+k];
	  SP_ROW(A)[ri++] = pt+i*N + k;
	}
      }
    }
    SP_COL(A)[j+4] = ri;
  }

  return Py_BuildValue("NN",A,b);
}
예제 #10
0
파일: umfpack.c 프로젝트: MGKhKhD/cvxopt
static PyObject* linsolve(PyObject *self, PyObject *args,
    PyObject *kwrds)
{
    spmatrix *A;
    matrix *B;
#if PY_MAJOR_VERSION >= 3
    int trans_ = 'N';
#endif
    char trans='N';
    double info[UMFPACK_INFO];
    int oB=0, n, nrhs=-1, ldB=0, k;
    void *symbolic, *numeric, *x;
    char *kwlist[] = {"A", "B", "trans", "nrhs", "ldB", "offsetB",
        NULL};

#if PY_MAJOR_VERSION >= 3
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OO|Ciii", kwlist,
        &A, &B, &trans_, &nrhs, &ldB, &oB)) return NULL;
    trans = (char) trans_;
#else
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OO|ciii", kwlist,
        &A, &B, &trans, &nrhs, &ldB, &oB)) return NULL;
#endif

    if (!SpMatrix_Check(A) || SP_NROWS(A) != SP_NCOLS(A))
        PY_ERR_TYPE("A must be a square sparse matrix");
    n = SP_NROWS(A);
    if (!Matrix_Check(B) || MAT_ID(B) != SP_ID(A))
        PY_ERR_TYPE("B must a dense matrix of the same numeric type "
            "as A");

    if (nrhs < 0) nrhs = B->ncols;
    if (n == 0 || nrhs == 0) return Py_BuildValue("i", 0);
    if (ldB == 0) ldB = MAX(1,B->nrows);
    if (ldB < MAX(1,n)) err_ld("ldB");
    if (oB < 0) err_nn_int("offsetB");
    if (oB + (nrhs-1)*ldB + n > MAT_LGT(B)) err_buf_len("B");

    if (trans != 'N' && trans != 'T' && trans != 'C')
        err_char("trans", "'N', 'T', 'C'");

    if (SP_ID(A) == DOUBLE)
        UMFD(symbolic)(n, n, SP_COL(A), SP_ROW(A), SP_VAL(A), &symbolic,
            NULL, info);
    else
        UMFZ(symbolic)(n, n, SP_COL(A), SP_ROW(A), SP_VAL(A), NULL,
            &symbolic, NULL, info);

    if (info[UMFPACK_STATUS] != UMFPACK_OK){
        if (SP_ID(A) == DOUBLE)
            UMFD(free_symbolic)(&symbolic);
        else
            UMFZ(free_symbolic)(&symbolic);
        if (info[UMFPACK_STATUS] == UMFPACK_ERROR_out_of_memory)
            return PyErr_NoMemory();
        else {
            snprintf(umfpack_error,20,"%s %i","UMFPACK ERROR",
                (int) info[UMFPACK_STATUS]);
            PyErr_SetString(PyExc_ValueError, umfpack_error);
            return NULL;
        }
    }

    if (SP_ID(A) == DOUBLE) {
        UMFD(numeric)(SP_COL(A), SP_ROW(A), SP_VAL(A), symbolic,
            &numeric, NULL, info);
        UMFD(free_symbolic)(&symbolic);
    } else {
        UMFZ(numeric)(SP_COL(A), SP_ROW(A), SP_VAL(A), NULL, symbolic,
            &numeric, NULL, info);
        UMFZ(free_symbolic)(&symbolic);
    }
    if (info[UMFPACK_STATUS] != UMFPACK_OK){
        if (SP_ID(A) == DOUBLE)
            UMFD(free_numeric)(&numeric);
        else
            UMFZ(free_numeric)(&numeric);
        if (info[UMFPACK_STATUS] == UMFPACK_ERROR_out_of_memory)
            return PyErr_NoMemory();
        else {
            if (info[UMFPACK_STATUS] == UMFPACK_WARNING_singular_matrix)
                PyErr_SetString(PyExc_ArithmeticError, "singular "
                    "matrix");
            else {
                snprintf(umfpack_error,20,"%s %i","UMFPACK ERROR",
                    (int) info[UMFPACK_STATUS]);
                PyErr_SetString(PyExc_ValueError, umfpack_error);
            }
            return NULL;
        }
    }

    if (!(x = malloc(n*E_SIZE[SP_ID(A)]))) {
        if (SP_ID(A) == DOUBLE)
            UMFD(free_numeric)(&numeric);
        else
            UMFZ(free_numeric)(&numeric);
        return PyErr_NoMemory();
    }
    for (k=0; k<nrhs; k++){
        if (SP_ID(A) == DOUBLE)
            UMFD(solve)(trans == 'N' ? UMFPACK_A: UMFPACK_Aat,
                SP_COL(A), SP_ROW(A), SP_VAL(A), x,
                MAT_BUFD(B) + k*ldB + oB, numeric, NULL, info);
        else
            UMFZ(solve)(trans == 'N' ? UMFPACK_A: trans == 'C' ?
                UMFPACK_At : UMFPACK_Aat, SP_COL(A), SP_ROW(A),
                SP_VAL(A), NULL, x, NULL,
                (double *)(MAT_BUFZ(B) + k*ldB + oB), NULL, numeric,
                NULL, info);
        if (info[UMFPACK_STATUS] == UMFPACK_OK)
            memcpy(B->buffer + (k*ldB + oB)*E_SIZE[SP_ID(A)], x,
                n*E_SIZE[SP_ID(A)]);
        else
	    break;
    }
    free(x);
    if (SP_ID(A) == DOUBLE)
        UMFD(free_numeric)(&numeric);
    else
        UMFZ(free_numeric)(&numeric);

    if (info[UMFPACK_STATUS] != UMFPACK_OK){
        if (info[UMFPACK_STATUS] == UMFPACK_ERROR_out_of_memory)
            return PyErr_NoMemory();
        else {
            if (info[UMFPACK_STATUS] == UMFPACK_WARNING_singular_matrix)
                PyErr_SetString(PyExc_ArithmeticError, "singular "
                    "matrix");
            else {
                snprintf(umfpack_error,20,"%s %i","UMFPACK ERROR",
                    (int) info[UMFPACK_STATUS]);
                PyErr_SetString(PyExc_ValueError, umfpack_error);
            }
        return NULL;
        }
    }
    return Py_BuildValue("");
}
예제 #11
0
파일: glpk.c 프로젝트: ugonj/cvxopt
static PyObject *integer(PyObject *self, PyObject *args,
    PyObject *kwrds)
{
    matrix *c, *h, *b=NULL, *x=NULL;
    PyObject *G, *A=NULL, *IntSet=NULL, *BinSet = NULL;
    PyObject *t=NULL;
    pyiocp *iocpParm = NULL;;
    glp_iocp *options = NULL;
    glp_prob *lp;
    int m, n, p, i, j, k, nnz, nnzmax, *rn=NULL, *cn=NULL;
    double *a=NULL, val;
    char *kwlist[] = {"c", "G", "h", "A", "b", "I", "B","iocp", NULL};

    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OOO|OOOOO!", kwlist, &c,
	    &G, &h, &A, &b, &IntSet, &BinSet,iocp_t,&iocpParm)) return NULL;

    if(!iocpParm) 
    {
      iocpParm = (pyiocp*)malloc(sizeof(*iocpParm));
      glp_init_iocp(&(iocpParm->obj));
    }
    if(iocpParm) 
    {
      Py_INCREF(iocpParm);
      options = &iocpParm->obj;
      options->presolve = 1;
    }

    if ((Matrix_Check(G) && MAT_ID(G) != DOUBLE) ||
        (SpMatrix_Check(G) && SP_ID(G) != DOUBLE) ||
        (!Matrix_Check(G) && !SpMatrix_Check(G))){
        PyErr_SetString(PyExc_TypeError, "G must be a 'd' matrix");
        return NULL;
    }
    if ((m = Matrix_Check(G) ? MAT_NROWS(G) : SP_NROWS(G)) <= 0)
        err_p_int("m");
    if ((n = Matrix_Check(G) ? MAT_NCOLS(G) : SP_NCOLS(G)) <= 0)
        err_p_int("n");

    if (!Matrix_Check(h) || h->id != DOUBLE) err_dbl_mtrx("h");
    if (h->nrows != m || h->ncols != 1){
        PyErr_SetString(PyExc_ValueError, "incompatible dimensions");
        return NULL;
    }

    if (A){
        if ((Matrix_Check(A) && MAT_ID(A) != DOUBLE) ||
            (SpMatrix_Check(A) && SP_ID(A) != DOUBLE) ||
            (!Matrix_Check(A) && !SpMatrix_Check(A))){
                PyErr_SetString(PyExc_ValueError, "A must be a dense "
                    "'d' matrix or a general sparse matrix");
                return NULL;
	}
        if ((p = Matrix_Check(A) ? MAT_NROWS(A) : SP_NROWS(A)) < 0)
            err_p_int("p");
        if ((Matrix_Check(A) ? MAT_NCOLS(A) : SP_NCOLS(A)) != n){
            PyErr_SetString(PyExc_ValueError, "incompatible "
                "dimensions");
            return NULL;
	}
    }
    else p = 0;

    if (b && (!Matrix_Check(b) || b->id != DOUBLE)) err_dbl_mtrx("b");
    if ((b && (b->nrows != p || b->ncols != 1)) || (!b && p !=0 )){
        PyErr_SetString(PyExc_ValueError, "incompatible dimensions");
        return NULL;
    }

    if ((IntSet) && (!PyAnySet_Check(IntSet)))
      PY_ERR_TYPE("invalid integer index set");

    if ((BinSet) && (!PyAnySet_Check(BinSet)))
      PY_ERR_TYPE("invalid binary index set");

    lp = glp_create_prob();
    glp_add_rows(lp, m+p);
    glp_add_cols(lp, n);

    for (i=0; i<n; i++){
        glp_set_obj_coef(lp, i+1, MAT_BUFD(c)[i]);
        glp_set_col_bnds(lp, i+1, GLP_FR, 0.0, 0.0);
    }
    for (i=0; i<m; i++)
        glp_set_row_bnds(lp, i+1, GLP_UP, 0.0, MAT_BUFD(h)[i]);
    for (i=0; i<p; i++)
        glp_set_row_bnds(lp, i+m+1, GLP_FX, MAT_BUFD(b)[i],
            MAT_BUFD(b)[i]);

    nnzmax = (SpMatrix_Check(G) ? SP_NNZ(G) : m*n ) +
        ((A && SpMatrix_Check(A)) ? SP_NNZ(A) : p*n);
    a = (double *) calloc(nnzmax+1, sizeof(double));
    rn = (int *) calloc(nnzmax+1, sizeof(int));
    cn = (int *) calloc(nnzmax+1, sizeof(int));
    if (!a || !rn || !cn){
        free(a);  free(rn);  free(cn);  glp_delete_prob(lp);
        return PyErr_NoMemory();
    }

    nnz = 0;
    if (SpMatrix_Check(G)) {
        for (j=0; j<n; j++) for (k=SP_COL(G)[j]; k<SP_COL(G)[j+1]; k++)
            if ((val = SP_VALD(G)[k]) != 0.0){
                a[1+nnz] = val;
                rn[1+nnz] = SP_ROW(G)[k]+1;
                cn[1+nnz] = j+1;
                nnz++;
            }
    }
    else for (j=0; j<n; j++) for (i=0; i<m; i++)
        if ((val = MAT_BUFD(G)[i+j*m]) != 0.0){
            a[1+nnz] = val;
            rn[1+nnz] = i+1;
            cn[1+nnz] = j+1;
            nnz++;
        }

    if (A && SpMatrix_Check(A)){
        for (j=0; j<n; j++) for (k=SP_COL(A)[j]; k<SP_COL(A)[j+1]; k++)
            if ((val = SP_VALD(A)[k]) != 0.0){
                a[1+nnz] = val;
                rn[1+nnz] = m+SP_ROW(A)[k]+1;
                cn[1+nnz] = j+1;
                nnz++;
            }
    }
    else for (j=0; j<n; j++) for (i=0; i<p; i++)
        if ((val = MAT_BUFD(A)[i+j*p]) != 0.0){
            a[1+nnz] = val;
            rn[1+nnz] = m+i+1;
            cn[1+nnz] = j+1;
            nnz++;
        }

    glp_load_matrix(lp, nnz, rn, cn, a);
    free(rn);  free(cn);  free(a);

    if (!(t = PyTuple_New(2))) {
        glp_delete_prob(lp);
        return PyErr_NoMemory();
    }

    if (IntSet) {
      PyObject *iter = PySequence_Fast(IntSet, "Critical error: not sequence");

      for (i=0; i<PySet_GET_SIZE(IntSet); i++) {

	PyObject *tmp = PySequence_Fast_GET_ITEM(iter, i);
#if PY_MAJOR_VERSION >= 3
	if (!PyLong_Check(tmp)) {
#else
	if (!PyInt_Check(tmp)) {
#endif
	  glp_delete_prob(lp);
	  Py_DECREF(iter);
	  PY_ERR_TYPE("non-integer element in I");
	}
#if PY_MAJOR_VERSION >= 3
	int k = PyLong_AS_LONG(tmp);
#else
	int k = PyInt_AS_LONG(tmp);
#endif
	if ((k < 0) || (k >= n)) {
	  glp_delete_prob(lp);
	  Py_DECREF(iter);
	  PY_ERR(PyExc_IndexError, "index element out of range in I");
	}
	glp_set_col_kind(lp, k+1, GLP_IV);
      }

      Py_DECREF(iter);
    }

    if (BinSet) {
      PyObject *iter = PySequence_Fast(BinSet, "Critical error: not sequence");

      for (i=0; i<PySet_GET_SIZE(BinSet); i++) {

	PyObject *tmp = PySequence_Fast_GET_ITEM(iter, i);
#if PY_MAJOR_VERSION >= 3
	if (!PyLong_Check(tmp)) {
#else
	if (!PyInt_Check(tmp)) {
#endif
	  glp_delete_prob(lp);
	  Py_DECREF(iter);
	  PY_ERR_TYPE("non-binary element in I");
	}
#if PY_MAJOR_VERSION >= 3
	int k = PyLong_AS_LONG(tmp);
#else
	int k = PyInt_AS_LONG(tmp);
#endif
	if ((k < 0) || (k >= n)) {
	  glp_delete_prob(lp);
	  Py_DECREF(iter);
	  PY_ERR(PyExc_IndexError, "index element out of range in B");
	}
	glp_set_col_kind(lp, k+1, GLP_BV);
      }

      Py_DECREF(iter);

    }


      switch (glp_intopt(lp,options)){

          case 0:

              x = (matrix *) Matrix_New(n,1,DOUBLE);
              if (!x) {
                  Py_XDECREF(iocpParm);
                  Py_XDECREF(t);
                  glp_delete_prob(lp);
                  return PyErr_NoMemory();
              }
              set_output_string(t,"optimal");
              set_output_string(t,"optimal");

              for (i=0; i<n; i++)
                  MAT_BUFD(x)[i] = glp_mip_col_val(lp, i+1);
              PyTuple_SET_ITEM(t, 1, (PyObject *) x);

              Py_XDECREF(iocpParm);
              glp_delete_prob(lp);
              return (PyObject *) t;

          case GLP_ETMLIM:

              x = (matrix *) Matrix_New(n,1,DOUBLE);
              if (!x) {
                  Py_XDECREF(t);
                  Py_XDECREF(iocpParm);
                  glp_delete_prob(lp);
                  return PyErr_NoMemory();
              }
              set_output_string(t,"time limit exceeded");

              for (i=0; i<n; i++)
                  MAT_BUFD(x)[i] = glp_mip_col_val(lp, i+1);
              PyTuple_SET_ITEM(t, 1, (PyObject *) x);

              Py_XDECREF(iocpParm);
              glp_delete_prob(lp);
              return (PyObject *) t;


          case GLP_EBOUND:
              set_output_string(t,"incorrect bounds");
              break;
          case GLP_EFAIL:
              set_output_string(t,"invalid MIP formulation");
              break;

          case GLP_ENOPFS:
              set_output_string(t,"primal infeasible");
              break;

          case GLP_ENODFS:
              set_output_string(t,"dual infeasible");
              break;

          case GLP_EMIPGAP:
              set_output_string(t,"Relative mip gap tolerance reached");
              break;

              /*case LPX_E_ITLIM:

                set_output_string(t,"maxiters exceeded");
                break;*/

              /*case LPX_E_SING:

                set_output_string(t,"singular or ill-conditioned basis");
                break;*/


          default:

              set_output_string(t,"unknown");
      }

      Py_XDECREF(iocpParm);
    glp_delete_prob(lp);

    PyTuple_SET_ITEM(t, 1, Py_BuildValue(""));
    return (PyObject *) t;
}


static PyMethodDef glpk_functions[] = {
    {"lp", (PyCFunction) simplex, METH_VARARGS|METH_KEYWORDS,
        doc_simplex},
    {"ilp", (PyCFunction) integer, METH_VARARGS|METH_KEYWORDS,
        doc_integer},
    {NULL}  /* Sentinel */
};

#if PY_MAJOR_VERSION >= 3

static PyModuleDef glpk_module_def = {
    PyModuleDef_HEAD_INIT,
    "glpk",
    glpk__doc__,
    -1,
    glpk_functions,
    NULL, NULL, NULL, NULL
};

void addglpkConstants (void)
{
  PyModule_AddIntMacro(glpk_module, GLP_ON);
  PyModule_AddIntMacro(glpk_module,GLP_OFF);

  /* reason codes: */
  PyModule_AddIntMacro(glpk_module,GLP_IROWGEN);
  PyModule_AddIntMacro(glpk_module,GLP_IBINGO);
  PyModule_AddIntMacro(glpk_module,GLP_IHEUR);
  PyModule_AddIntMacro(glpk_module,GLP_ICUTGEN);
  PyModule_AddIntMacro(glpk_module,GLP_IBRANCH);
  PyModule_AddIntMacro(glpk_module,GLP_ISELECT);
  PyModule_AddIntMacro(glpk_module,GLP_IPREPRO);

  /* branch selection indicator: */
  PyModule_AddIntMacro(glpk_module,GLP_NO_BRNCH);
  PyModule_AddIntMacro(glpk_module,GLP_DN_BRNCH);
  PyModule_AddIntMacro(glpk_module,GLP_UP_BRNCH);

  /* return codes: */
  PyModule_AddIntMacro(glpk_module,GLP_EBADB);
  PyModule_AddIntMacro(glpk_module,GLP_ESING);
  PyModule_AddIntMacro(glpk_module,GLP_ECOND);
  PyModule_AddIntMacro(glpk_module,GLP_EBOUND);
  PyModule_AddIntMacro(glpk_module,GLP_EFAIL);
  PyModule_AddIntMacro(glpk_module,GLP_EOBJLL);
  PyModule_AddIntMacro(glpk_module,GLP_EOBJUL);
  PyModule_AddIntMacro(glpk_module,GLP_EITLIM);
  PyModule_AddIntMacro(glpk_module,GLP_ETMLIM);
  PyModule_AddIntMacro(glpk_module,GLP_ENOPFS);
  PyModule_AddIntMacro(glpk_module,GLP_ENODFS);
  PyModule_AddIntMacro(glpk_module,GLP_EROOT);
  PyModule_AddIntMacro(glpk_module,GLP_ESTOP);
  PyModule_AddIntMacro(glpk_module,GLP_EMIPGAP);
  PyModule_AddIntMacro(glpk_module,GLP_ENOFEAS);
  PyModule_AddIntMacro(glpk_module,GLP_ENOCVG);
  PyModule_AddIntMacro(glpk_module,GLP_EINSTAB);
  PyModule_AddIntMacro(glpk_module,GLP_EDATA);
  PyModule_AddIntMacro(glpk_module,GLP_ERANGE);

  /* condition indicator: */
  PyModule_AddIntMacro(glpk_module,GLP_KKT_PE);
  PyModule_AddIntMacro(glpk_module,GLP_KKT_PB);
  PyModule_AddIntMacro(glpk_module,GLP_KKT_DE);
  PyModule_AddIntMacro(glpk_module,GLP_KKT_DB);
  PyModule_AddIntMacro(glpk_module,GLP_KKT_CS);

  /* MPS file format: */
  PyModule_AddIntMacro(glpk_module,GLP_MPS_DECK);
  PyModule_AddIntMacro(glpk_module,GLP_MPS_FILE);

  /* simplex method control parameters */
  /* message level: */
  PyModule_AddIntMacro(glpk_module,GLP_MSG_OFF);
  PyModule_AddIntMacro(glpk_module,GLP_MSG_ERR);
  PyModule_AddIntMacro(glpk_module,GLP_MSG_ON);
  PyModule_AddIntMacro(glpk_module,GLP_MSG_ALL);
  PyModule_AddIntMacro(glpk_module,GLP_MSG_DBG);
  /* simplex method option: */
  PyModule_AddIntMacro(glpk_module,GLP_PRIMAL);
  PyModule_AddIntMacro(glpk_module,GLP_DUALP);
  PyModule_AddIntMacro(glpk_module,GLP_DUAL);
  /* pricing technique: */
  PyModule_AddIntMacro(glpk_module,GLP_PT_STD);
  PyModule_AddIntMacro(glpk_module,GLP_PT_PSE);
  /* ratio test technique: */
  PyModule_AddIntMacro(glpk_module,GLP_RT_STD);
  PyModule_AddIntMacro(glpk_module,GLP_RT_HAR);

  /* interior-point solver control parameters */
  /* ordering algorithm: */
  PyModule_AddIntMacro(glpk_module,GLP_ORD_NONE);
  PyModule_AddIntMacro(glpk_module,GLP_ORD_QMD);
  PyModule_AddIntMacro(glpk_module,GLP_ORD_AMD);
  PyModule_AddIntMacro(glpk_module,GLP_ORD_SYMAMD);
}

PyMODINIT_FUNC PyInit_glpk(void)
{
  if (!(glpk_module = PyModule_Create(&glpk_module_def))) return NULL;
  if (PyType_Ready(&iocp_t) < 0 || (PyType_Ready(&smcp_t) < 0)) return NULL;
  /*  Adding macros */
  addglpkConstants();
  /* Adding  option lists as objects */
  Py_INCREF(&smcp_t);
  PyModule_AddObject(glpk_module,"smcp",(PyObject*)&smcp_t);
  Py_INCREF(&iocp_t);
  PyModule_AddObject(glpk_module,"iocp",(PyObject*)&iocp_t);
  if (import_cvxopt() < 0) return NULL;
  return glpk_module;
}

#else

PyMODINIT_FUNC initglpk(void)
{
    glpk_module = Py_InitModule3("cvxopt.glpk", glpk_functions, 
            glpk__doc__);
    if (PyType_Ready(&iocp_t) < 0 || (PyType_Ready(&smcp_t) < 0)) return NULL;
    addglpkConstants();
    Py_INCREF(&smcp_t);
    PyModule_AddObject(glpk_module,"smcp",(PyObject*)&smcp_t);
    Py_INCREF(&iocp_t);
    PyModule_AddObject(glpk_module,"iocp",(PyObject*)&iocp_t);
    if (import_cvxopt() < 0) return;
}
예제 #12
0
파일: glpk.c 프로젝트: ugonj/cvxopt
static PyObject *simplex(PyObject *self, PyObject *args,
    PyObject *kwrds)
{
    matrix *c, *h, *b=NULL, *x=NULL, *z=NULL, *y=NULL;
    PyObject *G, *A=NULL, *t=NULL;
    glp_prob *lp;
    glp_smcp *options = NULL;
    pysmcp *smcpParm = NULL;
    int m, n, p, i, j, k, nnz, nnzmax, *rn=NULL, *cn=NULL;
    double *a=NULL, val;
    char *kwlist[] = {"c", "G", "h", "A", "b","options", NULL};

    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OOO|OOO!", kwlist, &c,
        &G, &h, &A, &b,&smcp_t,&smcpParm)) return NULL;

    if ((Matrix_Check(G) && MAT_ID(G) != DOUBLE) ||
        (SpMatrix_Check(G) && SP_ID(G) != DOUBLE) ||
        (!Matrix_Check(G) && !SpMatrix_Check(G))){
        PyErr_SetString(PyExc_TypeError, "G must be a 'd' matrix");
        return NULL;
    }
    if ((m = Matrix_Check(G) ? MAT_NROWS(G) : SP_NROWS(G)) <= 0)
        err_p_int("m");
    if ((n = Matrix_Check(G) ? MAT_NCOLS(G) : SP_NCOLS(G)) <= 0)
        err_p_int("n");

    if (!Matrix_Check(h) || h->id != DOUBLE) err_dbl_mtrx("h");
    if (h->nrows != m || h->ncols != 1){
        PyErr_SetString(PyExc_ValueError, "incompatible dimensions");
        return NULL;
    }

    if (A){
        if ((Matrix_Check(A) && MAT_ID(A) != DOUBLE) ||
            (SpMatrix_Check(A) && SP_ID(A) != DOUBLE) ||
            (!Matrix_Check(A) && !SpMatrix_Check(A))){
                PyErr_SetString(PyExc_ValueError, "A must be a dense "
                    "'d' matrix or a general sparse matrix");
                return NULL;
	}
        if ((p = Matrix_Check(A) ? MAT_NROWS(A) : SP_NROWS(A)) < 0)
            err_p_int("p");
        if ((Matrix_Check(A) ? MAT_NCOLS(A) : SP_NCOLS(A)) != n){
            PyErr_SetString(PyExc_ValueError, "incompatible "
                "dimensions");
            return NULL;
	}
    }
    else p = 0;

    if (b && (!Matrix_Check(b) || b->id != DOUBLE)) err_dbl_mtrx("b");
    if ((b && (b->nrows != p || b->ncols != 1)) || (!b && p !=0 )){
        PyErr_SetString(PyExc_ValueError, "incompatible dimensions");
        return NULL;
    }
    if(!smcpParm) 
    {
      smcpParm = (pysmcp*)malloc(sizeof(*smcpParm));
      glp_init_smcp(&(smcpParm->obj));
    }
    if(smcpParm) 
    {
      Py_INCREF(smcpParm);
      options = &smcpParm->obj;
      options->presolve = 1;
    }

    lp = glp_create_prob();
    glp_add_rows(lp, m+p);
    glp_add_cols(lp, n);

    for (i=0; i<n; i++){
        glp_set_obj_coef(lp, i+1, MAT_BUFD(c)[i]);
        glp_set_col_bnds(lp, i+1, GLP_FR, 0.0, 0.0);
    }
    for (i=0; i<m; i++)
        glp_set_row_bnds(lp, i+1, GLP_UP, 0.0, MAT_BUFD(h)[i]);
    for (i=0; i<p; i++)
        glp_set_row_bnds(lp, i+m+1, GLP_FX, MAT_BUFD(b)[i],
            MAT_BUFD(b)[i]);

    nnzmax = (SpMatrix_Check(G) ? SP_NNZ(G) : m*n ) +
        ((A && SpMatrix_Check(A)) ? SP_NNZ(A) : p*n);
    a = (double *) calloc(nnzmax+1, sizeof(double));
    rn = (int *) calloc(nnzmax+1, sizeof(int));
    cn = (int *) calloc(nnzmax+1, sizeof(int));
    if (!a || !rn || !cn){
        free(a);  free(rn);  free(cn);  glp_delete_prob(lp);
        return PyErr_NoMemory();
    }

    nnz = 0;
    if (SpMatrix_Check(G)) {
        for (j=0; j<n; j++) for (k=SP_COL(G)[j]; k<SP_COL(G)[j+1]; k++)
            if ((val = SP_VALD(G)[k]) != 0.0){
                a[1+nnz] = val;
                rn[1+nnz] = SP_ROW(G)[k]+1;
                cn[1+nnz] = j+1;
                nnz++;
            }
    }
    else for (j=0; j<n; j++) for (i=0; i<m; i++)
        if ((val = MAT_BUFD(G)[i+j*m]) != 0.0){
            a[1+nnz] = val;
            rn[1+nnz] = i+1;
            cn[1+nnz] = j+1;
            nnz++;
        }

    if (A && SpMatrix_Check(A)){
        for (j=0; j<n; j++) for (k=SP_COL(A)[j]; k<SP_COL(A)[j+1]; k++)
            if ((val = SP_VALD(A)[k]) != 0.0){
                a[1+nnz] = val;
                rn[1+nnz] = m+SP_ROW(A)[k]+1;
                cn[1+nnz] = j+1;
                nnz++;
            }
    }
    else for (j=0; j<n; j++) for (i=0; i<p; i++)
        if ((val = MAT_BUFD(A)[i+j*p]) != 0.0){
            a[1+nnz] = val;
            rn[1+nnz] = m+i+1;
            cn[1+nnz] = j+1;
            nnz++;
        }

    glp_load_matrix(lp, nnz, rn, cn, a);
    free(rn);  free(cn);  free(a);

    if (!(t = PyTuple_New(A ? 4 : 3))){
        glp_delete_prob(lp);
        return PyErr_NoMemory();
    }


    switch (glp_simplex(lp,options)){

        case 0:

            x = (matrix *) Matrix_New(n,1,DOUBLE);
            z = (matrix *) Matrix_New(m,1,DOUBLE);
            if (A) y = (matrix *) Matrix_New(p,1,DOUBLE);
            if (!x || !z || (A && !y)){
                Py_XDECREF(x);
                Py_XDECREF(z);
                Py_XDECREF(y);
                Py_XDECREF(t);
                Py_XDECREF(smcpParm);
                glp_delete_prob(lp);
                return PyErr_NoMemory();
            }

            set_output_string(t,"optimal");

            for (i=0; i<n; i++)
                MAT_BUFD(x)[i] = glp_get_col_prim(lp, i+1);
            PyTuple_SET_ITEM(t, 1, (PyObject *) x);

            for (i=0; i<m; i++)
                MAT_BUFD(z)[i] = -glp_get_row_dual(lp, i+1);
            PyTuple_SET_ITEM(t, 2, (PyObject *) z);

            if (A){
                for (i=0; i<p; i++)
                    MAT_BUFD(y)[i] = -glp_get_row_dual(lp, m+i+1);
                PyTuple_SET_ITEM(t, 3, (PyObject *) y);
            }

            Py_XDECREF(smcpParm);
            glp_delete_prob(lp);
            return (PyObject *) t;
        case GLP_EBADB:
            set_output_string(t,"incorrect initial basis");
            break;
        case GLP_ESING:
            set_output_string(t,"singular initial basis matrix");
            break;
        case GLP_ECOND:
            set_output_string(t,"ill-conditioned initial basis matrix");
            break;
        case GLP_EBOUND:
            set_output_string(t,"incorrect bounds");
            break;
        case GLP_EFAIL:
            set_output_string(t,"solver failure");
            break;
        case GLP_EOBJLL:
            set_output_string(t,"objective function reached lower limit");
            break;
        case GLP_EOBJUL:
            set_output_string(t,"objective function reached upper limit");
            break;
        case GLP_EITLIM:
            set_output_string(t,"iteration limit exceeded");
            break;
        case GLP_ETMLIM:
            set_output_string(t,"time limit exceeded");
            break;
        case GLP_ENOPFS:
            set_output_string(t,"primal infeasible");
            break;
        case GLP_ENODFS:
            set_output_string(t,"dual infeasible");
            break;
        default:
            set_output_string(t,"unknown");
            break;
    }

    Py_XDECREF(smcpParm);
    glp_delete_prob(lp);

    PyTuple_SET_ITEM(t, 1, Py_BuildValue(""));
    PyTuple_SET_ITEM(t, 2, Py_BuildValue(""));
    if (A) PyTuple_SET_ITEM(t, 3, Py_BuildValue(""));

    return (PyObject *) t;
}
예제 #13
0
파일: amd.c 프로젝트: ChiahungTai/cvxopt
static int set_defaults(double *control)
{
    int_t pos=0;
    int param_id;
    PyObject *param, *key, *value;
#if PY_MAJOR_VERSION < 3
    char *keystr; 
#endif
    char err_str[100];

    amd_defaults(control);

    if (!(param = PyObject_GetAttrString(amd_module, "options")) ||
        !PyDict_Check(param)){
        PyErr_SetString(PyExc_AttributeError, "missing amd.options"
            "dictionary");
        return 0;
    }
    while (PyDict_Next(param, &pos, &key, &value))
#if PY_MAJOR_VERSION >= 3
        if ((PyUnicode_Check(key)) && 
            get_param_idx(_PyUnicode_AsString(key),&param_id)) {
            if (!PyLong_Check(value) && !PyFloat_Check(value)){
                sprintf(err_str, "invalid value for AMD parameter: %-.20s",
                    _PyUnicode_AsString(key));
#else
        if ((keystr = PyString_AsString(key)) && get_param_idx(keystr,
            &param_id)) {
            if (!PyInt_Check(value) && !PyFloat_Check(value)){
                sprintf(err_str, "invalid value for AMD parameter: "
                    "%-.20s", keystr);
#endif
                PyErr_SetString(PyExc_ValueError, err_str);
                Py_DECREF(param);
                return 0;
            }
            control[param_id] = PyFloat_AsDouble(value);
        }
    Py_DECREF(param);
    return 1;
}


static char doc_order[] =
    "Computes the approximate minimum degree ordering of a square "
    "matrix.\n\n"
    "p = order(A, uplo='L')\n\n"
    "PURPOSE\n"
    "Computes a permutation p that reduces fill-in in the Cholesky\n"
    "factorization of A[p,p].\n\n"
    "ARGUMENTS\n"
    "A         square sparse matrix\n\n"
    "uplo      'L' or 'U'.  If uplo is 'L', the lower triangular part\n"
    "          of A is used and the upper triangular is ignored.  If\n"
    "          uplo is 'U', the upper triangular part is used and the\n"
    "          lower triangular part is ignored.\n\n"
    "p         'i' matrix of length equal to the order of A";


static PyObject* order_c(PyObject *self, PyObject *args, PyObject *kwrds)
{
    spmatrix *A;
    matrix *perm;
#if PY_MAJOR_VERSION >= 3
    int uplo_ = 'L';
#endif
    char uplo = 'L';
    int j, k, n, nnz, alloc=0, info;
    int_t *rowind=NULL, *colptr=NULL;
    double control[AMD_CONTROL];
    char *kwlist[] = {"A", "uplo", NULL};

#if PY_MAJOR_VERSION >= 3
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "O|C", kwlist, &A,
        &uplo_)) return NULL;
    uplo = (char) uplo_;
#else
    if (!PyArg_ParseTupleAndKeywords(args, kwrds, "O|c", kwlist, &A,
        &uplo)) return NULL;
#endif
    if (!set_defaults(control)) return NULL;
    if (!SpMatrix_Check(A) || SP_NROWS(A) != SP_NCOLS(A)){
        PyErr_SetString(PyExc_TypeError, "A must be a square sparse "
            "matrix");
        return NULL;
    }
    if (uplo != 'U' && uplo != 'L') err_char("uplo", "'L', 'U'");
    if (!(perm = (matrix *) Matrix_New((int)SP_NROWS(A),1,INT)))
        return PyErr_NoMemory();
    n = SP_NROWS(A);
    for (nnz=0, j=0; j<n; j++) {
        if (uplo == 'L'){
            for (k=SP_COL(A)[j]; k<SP_COL(A)[j+1] && SP_ROW(A)[k]<j; k++);
            nnz += SP_COL(A)[j+1] - k;
        }
        else {
            for (k=SP_COL(A)[j]; k<SP_COL(A)[j+1] && SP_ROW(A)[k] <= j;
                k++);
            nnz += k - SP_COL(A)[j];
        }
    }
    if (nnz == SP_NNZ(A)){
        colptr = (int_t *) SP_COL(A);
        rowind = (int_t *) SP_ROW(A);
    }
    else {
        alloc = 1;
        colptr = (int_t *) calloc(n+1, sizeof(int_t));
        rowind = (int_t *) calloc(nnz, sizeof(int_t));
        if (!colptr || !rowind) {
            Py_XDECREF(perm);  free(colptr);  free(rowind);
            return PyErr_NoMemory();
        }
        colptr[0] = 0;
        for (j=0; j<n; j++) {
            if (uplo == 'L'){
                for (k=SP_COL(A)[j]; k<SP_COL(A)[j+1] && SP_ROW(A)[k] < j; 
                    k++);
                nnz = SP_COL(A)[j+1] - k;
                colptr[j+1] = colptr[j] + nnz;
                memcpy(rowind + colptr[j], (int_t *) SP_ROW(A) + k,
                    nnz*sizeof(int_t));
            }
            else {
                for (k=SP_COL(A)[j]; k<SP_COL(A)[j+1] && SP_ROW(A)[k] <= j;
                    k++);
                nnz = k - SP_COL(A)[j];
                colptr[j+1] = colptr[j] + nnz;
                memcpy(rowind + colptr[j], (int_t *) (SP_ROW(A) +
                    SP_COL(A)[j]), nnz*sizeof(int_t));
            }
        }
    }
    info = amd_order(n, colptr, rowind, MAT_BUFI(perm), control, NULL);
    if (alloc){
        free(colptr);
        free(rowind);
    }
    switch (info) {
        case AMD_OUT_OF_MEMORY:
            Py_XDECREF(perm);
            return PyErr_NoMemory();

        case AMD_INVALID:
            Py_XDECREF(perm);
            return Py_BuildValue("");

        case AMD_OK:
            return (PyObject *) perm;
    }
    return Py_BuildValue("");
}

static PyMethodDef amd_functions[] = {
    {"order", (PyCFunction) order_c, METH_VARARGS|METH_KEYWORDS, doc_order},
    {NULL}  /* Sentinel */
};

#if PY_MAJOR_VERSION >= 3

static PyModuleDef amd_module_def = {
    PyModuleDef_HEAD_INIT,
    "amd",
    amd__doc__,
    -1,
    amd_functions,
    NULL, NULL, NULL, NULL
};

PyMODINIT_FUNC PyInit_amd(void)
{
    if (!(amd_module = PyModule_Create(&amd_module_def))) return NULL;
    PyModule_AddObject(amd_module, "options", PyDict_New());
    if (import_cvxopt() < 0) return NULL;
    return amd_module;
}

#else
PyMODINIT_FUNC initamd(void)
{
    amd_module = Py_InitModule3("cvxopt.amd", amd_functions, amd__doc__);
    PyModule_AddObject(amd_module, "options", PyDict_New());
    if (import_cvxopt() < 0) return;
}