static PyObject *SCMcolumn (PyObject *self, PyObject *args) { PyObject *H,*Av,*Ip,*Jp,*V; int_t m,n,K; int_t i,j,k,l,p,q,r,c; if(!PyArg_ParseTuple(args,"OOOOOn",&H,&Av,&V,&Ip,&Jp,&j)) return NULL; m = MAT_NCOLS(H); n = MAT_NROWS(V); K = MAT_NCOLS(V)/2; //#pragma omp parallel for shared(m,n,K,Av,Ip,Jp,H,V,j) private(i,l,k,r,c,q,p) for (i=j;i<m;i++) { p = SP_COL(Av)[i]; MAT_BUFD(H)[j*m+i] = 0; for (l=0;l<SP_COL(Av)[i+1]-p;l++) { q = SP_ROW(Av)[p+l]; r = MAT_BUFI(Ip)[q]; c = MAT_BUFI(Jp)[q]; for (k=0;k<K;k++) { MAT_BUFD(H)[j*m+i] += SP_VALD(Av)[p+l]* MAT_BUFD(V)[k*n+r]*MAT_BUFD(V)[(K+k)*n+c]; if (r != c) MAT_BUFD(H)[j*m+i] += SP_VALD(Av)[p+l]* MAT_BUFD(V)[k*n+c]*MAT_BUFD(V)[(K+k)*n+r]; } } } Py_RETURN_NONE; }
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; }
static PyObject *Av_to_spmatrix (PyObject *self, PyObject *args, PyObject *kwrds) { PyObject *scale = Py_False; spmatrix *Av,*Ip,*Jp; int_t i,j,n,nnz,c,ci,p,q; char *kwlist[] = {"Av","Ip","Jp","j","n","scale",NULL}; if (!PyArg_ParseTupleAndKeywords(args, kwrds, "OOOnn|O", kwlist, &Av,&Ip,&Jp,&j,&n, &scale)) return NULL; p = SP_COL(Av)[j]; nnz = SP_COL(Av)[j+1]-p; spmatrix *Aj = SpMatrix_New(n,n,nnz,DOUBLE); if (!Aj) return PyErr_NoMemory(); // Generate col-ptr and row index SP_COL(Aj)[0] = 0; if (scale==Py_False) { for (ci=0,i=0;i<nnz;i++) { q = SP_ROW(Av)[p+i]; SP_ROW(Aj)[i] = MAT_BUFI(Ip)[q]; c = MAT_BUFI(Jp)[q]; SP_VALD(Aj)[i] = SP_VALD(Av)[p+i]; while (ci < c) SP_COL(Aj)[++ci] = i; } while (ci < n) SP_COL(Aj)[++ci] = nnz; } else { for (ci=0,i=0;i<nnz;i++) { q = SP_ROW(Av)[p+i]; SP_ROW(Aj)[i] = MAT_BUFI(Ip)[q]; c = MAT_BUFI(Jp)[q]; SP_VALD(Aj)[i] = SP_VALD(Av)[p+i]; if (c == SP_ROW(Aj)[i]) SP_VALD(Aj)[i] *= 0.5; // scale diag. element while (ci < c) SP_COL(Aj)[++ci] = i; } while (ci < n) SP_COL(Aj)[++ci] = nnz; } return (PyObject *) Aj; }
static PyObject *sub2ind (PyObject *self, PyObject *args) { matrix *Im,*Jm; PyObject *siz; int_t i; int_t m,n; if (!PyArg_ParseTuple(args, "OOO", &siz, &Im, &Jm)) return NULL; if (!PyArg_ParseTuple(siz, "nn", &m, &n)) return NULL; matrix *Ind = Matrix_New(MAT_NROWS(Im),1,INT); if (!Ind) return PyErr_NoMemory(); for (i=0;i< MAT_NROWS(Im) ;i++) { // Add data check: MAT_BUFI(Ind)[i] = MAT_BUFI(Im)[i] + m*MAT_BUFI(Jm)[i]; } return Py_BuildValue("N", Ind); }
static PyObject *ind2sub (PyObject *self, PyObject *args) { matrix *Im; int_t i; int_t n; if (!PyArg_ParseTuple(args, "nO", &n, &Im)) return NULL; matrix *Il = Matrix_New(MAT_NROWS(Im),1,INT); if (!Il) return PyErr_NoMemory(); matrix *Jl = Matrix_New(MAT_NROWS(Im),1,INT); if (!Il) return PyErr_NoMemory(); for (i=0;i< MAT_NROWS(Im);i++) { MAT_BUFI(Il)[i] = MAT_BUFI(Im)[i] % n; MAT_BUFI(Jl)[i] = MAT_BUFI(Im)[i] / n; } return Py_BuildValue("NN", Il, Jl); }
static PyObject *matperm (PyObject *self, PyObject *args) { PyObject *nzc; matrix *pm; int_t Ns,Nd,m,i,Nmax; if (!PyArg_ParseTuple(args,"On",&nzc,&Nmax)) return NULL; m = MAT_NROWS(nzc); pm = Matrix_New(m,1,INT); if (!pm) return PyErr_NoMemory(); // Check Nmax if (Nmax<0) Nmax = 0; Ns = 0; Nd = 0; for (i=0;i<m;i++){ if(MAT_BUFI(nzc)[i] > Nmax) MAT_BUFI(pm)[Nd++] = i; else MAT_BUFI(pm)[m-1-Ns++] = i; } return Py_BuildValue("Nn",pm,Ns); }
static PyObject *scal_diag (PyObject *self, PyObject *args) { PyObject *Vp,*Id; int_t i,n; double t = 0.5; if(!PyArg_ParseTuple(args,"OO|d",&Vp,&Id,&t)) return NULL; n = MAT_NROWS(Id); for (i=0;i<n;i++){ SP_VALD(Vp)[MAT_BUFI(Id)[i]] *= t; } Py_RETURN_NONE; }
static PyObject *SCMcolumn2 (PyObject *self, PyObject *args) { PyObject *H,*Av,*Ip,*Jp,*V,*Kl; int_t m,n; int_t i,j,k,p,q,r,c,r1,c1,pj,pi; double alpha,beta; if(!PyArg_ParseTuple(args,"OOOOOOn",&H,&Av,&V,&Ip,&Jp,&Kl,&j)) return NULL; m = MAT_NCOLS(H); n = MAT_NROWS(V); for (i=j;i<m;i++) MAT_BUFD(H)[j*m+i] = 0; //#pragma omp parallel for shared(m,n,K,Av,Ip,Jp,H,V,j) private(i,l,k,r,c,q,p) pj = SP_COL(Av)[j]; for (p=0;p<SP_COL(Av)[j+1]-pj;p++) { alpha = SP_VALD(Av)[pj+p]; k = SP_ROW(Av)[pj+p]; r = MAT_BUFI(Ip)[k]; c = MAT_BUFI(Jp)[k]; if (r!=c) alpha*=2; // look up columns in V r = MAT_BUFI(Kl)[r]; c = MAT_BUFI(Kl)[c]; for(i=j;i<m;i++) { pi = SP_COL(Av)[i]; for (q=0;q<SP_COL(Av)[i+1]-pi;q++) { beta = SP_VALD(Av)[pi+q]; k = SP_ROW(Av)[pi+q]; r1 = MAT_BUFI(Ip)[k]; c1 = MAT_BUFI(Jp)[k]; MAT_BUFD(H)[j*m+i] += alpha*beta*MAT_BUFD(V)[n*r+r1]*MAT_BUFD(V)[n*c+c1]; if (r1!=c1) MAT_BUFD(H)[j*m+i] += alpha*beta*MAT_BUFD(V)[n*r+c1]*MAT_BUFD(V)[n*c+r1]; } } } Py_RETURN_NONE; }
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; }
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(""); }
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 }
static PyObject* sdpa_readhead (PyObject *self, PyObject *args) { int i,j,t; int_t m=0,n=0,nblocks=0; matrix *bstruct = NULL; PyObject *f; char buf[2048]; // buffer char *info; if (!PyArg_ParseTuple(args,"O",&f)) return NULL; #if PY_MAJOR_VERSION >= 3 if (PyUnicode_Check(f)) { const char* fname = PyUnicode_AsUTF8AndSize(f,NULL); #else if (PyString_Check(f)) { const char* fname = PyString_AsString(f); #endif FILE *fp = fopen(fname,"r"); if (!fp) { return NULL; } /* Skip comments and read m */ while (1) { info = fgets(buf,1024,fp); if (buf[0] != '*' && buf[0] != '"') { sscanf(buf,"%d",&i); break; } } m = (int_t) i; /* read nblocks */ j = fscanf(fp,"%d",&i); nblocks = (int_t) i; /* read blockstruct and compute block offsets*/ bstruct = Matrix_New(nblocks,1,INT); if (!bstruct) return PyErr_NoMemory(); n = 0; for (i=0; i<nblocks; i++) { j = fscanf(fp,"%*[^0-9+-]%d",&t); MAT_BUFI(bstruct)[i] = (int_t) t; n += (int_t) labs(MAT_BUFI(bstruct)[i]); } fclose(fp); } return Py_BuildValue("iiN",n,m,bstruct); } static char doc_sdpa_read[] = "Reads sparse SDPA data file (dat-s).\n" "\n" "A,b,bstruct = sdpa_read(f[,neg=False])\n" "\n" "PURPOSE\n" "Reads problem data from sparse SDPA data file for\n" "the semidefinite programs:\n" "\n" " (P) minimize <A0,X>\n" " subject to <Ai,X> = bi, i = 1,...,m\n" " X >= 0\n" "\n" " (D) maximize b'*y\n" " subject to sum_i Ai*yi + S = A0\n" " S >= 0\n" "\n" "Here '>=' means that X and S must be positive semidefinite.\n" "The matrices A0,A1,...Am are symmetric and of order n.\n" "If the optional argument 'neg' is True, the negative of the\n" "problem data is returned.\n" "\n" "ARGUMENTS\n" "f Python file object\n" "\n" "neg Python boolean (optional)\n" "\n" "RETURNS\n" "A CVXOPT sparse matrix of doubles with columns Ai[:]\n" " (Only lower trianglular elements of Ai are stored.)\n" "\n" "b CVXOPT matrix\n" "\n" "bstruct CVXOPT integer matrix\n"; static PyObject* sdpa_read (PyObject *self, PyObject *args, PyObject *kwrds) { int i,j,mno,bno,ii,jj,t; int_t k,m,n,nblocks,nlines; double v; long fpos; PyObject *f; PyObject *neg = Py_False; char *info; const char* fname; int_t* boff; // block offset char buf[2048]; // buffer char *kwlist[] = {"f","neg",NULL}; if (!PyArg_ParseTupleAndKeywords(args,kwrds,"O|O",kwlist,&f,&neg)) return NULL; #if PY_MAJOR_VERSION >= 3 if (PyUnicode_Check(f)) fname = PyUnicode_AsUTF8AndSize(f,NULL); #elif PY_MAJOR_VERSION == 2 if (PyString_Check(f)) fname = PyString_AsString(f); #endif FILE *fp = fopen(fname,"r"); if (!fp) { return NULL; } /* Skip comments and read m */ while (1) { info = fgets(buf,1024,fp); if (buf[0] != '*' && buf[0] != '"') { sscanf(buf,"%d",&i); break; } } m = (int_t) i; /* read nblocks */ j = fscanf(fp,"%d",&i); nblocks = (int_t) i; /* read blockstruct and compute block offsets*/ matrix *bstruct = Matrix_New(nblocks,1,INT); if (!bstruct) return PyErr_NoMemory(); boff = malloc(sizeof(int_t)*(nblocks+1)); if(!boff) return PyErr_NoMemory(); boff[0] = 0; n = 0; for (i=0; i<nblocks; i++) { j = fscanf(fp,"%*[^0-9+-]%d",&t); MAT_BUFI(bstruct)[i] = (int_t) t; n += (int_t) labs(MAT_BUFI(bstruct)[i]); boff[i+1] = n; } /* read vector b */ matrix *b = Matrix_New(m,1,DOUBLE); if (!b) return PyErr_NoMemory(); for (i=0;i<m;i++) { j = fscanf(fp,"%*[^0-9+-]%lf",&MAT_BUFD(b)[i]); if (neg == Py_True) MAT_BUFD(b)[i] *= -1; } /* count remaining lines */ fpos = ftell(fp); for (nlines = 0; fgets(buf, 1023, fp) != NULL; nlines++); //nlines--; fseek(fp,fpos,SEEK_SET); /* Create data matrix A */ spmatrix *A = SpMatrix_New(n*n,m+1,nlines,DOUBLE); if (!A) return PyErr_NoMemory(); // read data matrices fseek(fp,fpos,SEEK_SET); for (i=0,j=-1,k=0;k<nlines;k++){ if (fscanf(fp,"%*[^0-9+-]%d",&mno) <=0 ) break; if (fscanf(fp,"%*[^0-9+-]%d",&bno) <=0 ) break; if (fscanf(fp,"%*[^0-9+-]%d",&ii) <=0 ) break; if (fscanf(fp,"%*[^0-9+-]%d",&jj) <=0 ) break; if (fscanf(fp,"%*[^0-9+-]%lf",&v) <=0 ) break; // check that value is nonzero if (v != 0) { // add block offset ii += boff[bno-1]; jj += boff[bno-1]; // insert index and value SP_ROW(A)[i] = (int_t) ((ii-1)*n + (jj-1)); if (neg == Py_True) SP_VALD(A)[i] = -v; else SP_VALD(A)[i] = v; // update col. ptr. while (mno > j) SP_COL(A)[++j] = i; i++; } } // update last element(s) of col. ptr. while (m+1 > j) SP_COL(A)[++j] = i; fclose(fp); // free temp. memory free(boff); return Py_BuildValue("NNN",A,b,bstruct); }
static PyObject* sdpa_write (PyObject *self, PyObject *args, PyObject *kwrds) { int i,Il,Jl,Bl,Ml; int_t n; spmatrix *A; matrix *b,*bstruct; PyObject *f; PyObject *neg = Py_False; char *kwlist[] = {"f","A","b","bstruct","neg",NULL}; const char* fname; double v; if (!PyArg_ParseTupleAndKeywords(args,kwrds, "OOOO|O", kwlist, &f, &A, &b, &bstruct,&neg)) return NULL; #if PY_MAJOR_VERSION >= 3 if (PyUnicode_Check(f)) fname = PyUnicode_AsUTF8AndSize(f,NULL); #elif PY_MAJOR_VERSION == 2 if (PyString_Check(f)) fname = PyString_AsString(f); #endif FILE *fp = fopen(fname,"r"); if (!fp) { Py_DECREF(f); return NULL; } fprintf(fp,"* sparse SDPA data file (created by SMCP)\n"); fprintf(fp,"%i = m\n",(int) MAT_NROWS(b)); fprintf(fp,"%i = nBlocks\n", (int) MAT_NROWS(bstruct)); // compute n and write blockstruct n = 0; for (i=0;i<MAT_NROWS(bstruct);i++) { fprintf(fp,"%i ", (int) MAT_BUFI(bstruct)[i]); n += (int_t) labs(MAT_BUFI(bstruct)[i]); } fprintf(fp,"\n"); // write vector b if (neg == Py_True) { for (i=0;i<MAT_NROWS(b);i++) fprintf(fp,"%.12g ",-MAT_BUFD(b)[i]); } else { for (i=0;i<MAT_NROWS(b);i++) fprintf(fp,"%.12g ",MAT_BUFD(b)[i]); } fprintf(fp,"\n"); // Write data matrices A0,A1,A2,...,Am for (Ml=0;Ml<=MAT_NROWS(b);Ml++) { for (i=0;i<SP_COL(A)[Ml+1]-SP_COL(A)[Ml];i++){ Jl = 1 + SP_ROW(A)[SP_COL(A)[Ml]+i] / n; Il = 1 + SP_ROW(A)[SP_COL(A)[Ml]+i] % n; // Skip if element is in strict upper triangle if (Jl > Il) PyErr_Warn(PyExc_Warning,"Ignored strictly upper triangular element."); Bl = 1; while ((Il > labs(MAT_BUFI(bstruct)[Bl-1])) && (Jl > labs(MAT_BUFI(bstruct)[Bl-1]))) { Il -= (int_t) labs(MAT_BUFI(bstruct)[Bl-1]); Jl -= (int_t) labs(MAT_BUFI(bstruct)[Bl-1]); Bl += 1; } /* Error check */ if ((Il > labs(MAT_BUFI(bstruct)[Bl-1])) || (Jl > labs(MAT_BUFI(bstruct)[Bl-1]))) printf("Error: Matrix contains elements outside blocks!\n"); // print upper triangle entries: // <matno> <blkno> <i> <j> <entry> v = SP_VALD(A)[SP_COL(A)[Ml]+i]; if ( v != 0.0) { if (neg == Py_True) fprintf(fp,"%i %i %i %i %.12g\n", (int) Ml,(int) Bl,(int) Jl,(int) Il, -v); else fprintf(fp,"%i %i %i %i %.12g\n", (int) Ml,(int) Bl,(int) Jl,(int) Il, v); } } } fclose(fp); Py_DECREF(f); Py_RETURN_NONE; }
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),¶m_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, ¶m_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; }