/* ************************************************************ PROCEDURE vecsymPSD - Let y = (x+x')/2. INPUT x - length sum(K.s.^2) vector, to be symmetrized. rsdpN - number of real PSD blocks sdpN - total number of PSD blocks, length(K.s) sdpNL - K.s, length sdpN array listing the orders of the PSD blocks. OUTPUT y - length sum(K.s.^2) vector, y = vecsym(x,K). ************************************************************ */ void vecsymPSD(double *y, const double *x,const int rsdpN,const int sdpN, const double *sdpNL) { int k,nk,nksqr; /* ------------------------------------------------------------ Make real PSD blocks symmetric ------------------------------------------------------------ */ for(k = 0; k < rsdpN; k++){ nk = sdpNL[k]; nksqr = SQR(nk); symproj(y,x,nk); x += nksqr; y += nksqr; } /* ------------------------------------------------------------ Make complex PSD blocks Hermitian ------------------------------------------------------------ */ for(; k < sdpN; k++){ nk = sdpNL[k]; nksqr = SQR(nk); symproj(y,x,nk); x += nksqr; y += nksqr; skewproj(y,x,nk); x += nksqr; y += nksqr; } }
/* ************************************************************ PROCEDURE mexFunction - Entry for Matlab [lab,q] = eigK(x,K) Computes spectral coefficients of x w.r.t. K REMARK If this function is used internally by SeDuMi, then complex numbers are stored in a single real vector. To make it invokable from the Matlab command-line by the user, we also allow Matlab complex vector x. ************************************************************ */ void mexFunction(const int nlhs, mxArray *plhs[], const int nrhs, const mxArray *prhs[]) { mxArray *output_array[3], *Xk, *hXk; coneK cK; int k, nk, nksqr, lendiag,i,ii,nkp1, lenfull; double *lab,*q,*qpi,*labk,*xwork,*xpiwork; const double *x,*xpi; /* ------------------------------------------------------------ Check for proper number of arguments ------------------------------------------------------------ */ mxAssert(nrhs >= NPARIN, "eigK requires more input arguments"); mxAssert(nlhs <= NPAROUT, "eigK produces less output arguments"); /* ------------------------------------------------------------ Disassemble cone K structure ------------------------------------------------------------ */ conepars(K_IN, &cK); /* ------------------------------------------------------------ Compute statistics based on cone K structure ------------------------------------------------------------ */ lendiag = cK.lpN + 2 * (cK.lorN + cK.rconeN) + cK.rLen + cK.hLen; lenfull = cK.lpN + cK.qDim + cK.rDim + cK.hDim; if(cK.rconeN > 0) for(i = 0; i < cK.rconeN; i++) lenfull += cK.rconeNL[i]; /* ------------------------------------------------------------ Get input vector x ------------------------------------------------------------ */ mxAssert(mxGetM(X_IN) * mxGetN(X_IN) == lenfull, "Size mismatch x"); mxAssert(!mxIsSparse(X_IN), "x must be full (not sparse)."); x = mxGetPr(X_IN); if(mxIsComplex(X_IN)) xpi = mxGetPi(X_IN) + cK.lpN; /* ------------------------------------------------------------ Allocate output LAB(diag), eigvec Q(full for psd) ------------------------------------------------------------ */ LAB_OUT = mxCreateDoubleMatrix(lendiag, 1, mxREAL); lab = mxGetPr(LAB_OUT); if(nlhs > 1){ if(mxIsComplex(X_IN)){ Q_OUT = mxCreateDoubleMatrix(cK.rDim, 1, mxCOMPLEX); qpi = mxGetPi(Q_OUT); } else Q_OUT = mxCreateDoubleMatrix(cK.rDim + cK.hDim, 1, mxREAL); q = mxGetPr(Q_OUT); } /* ------------------------------------------------------------ Allocate working arrays: ------------------------------------------------------------ */ Xk = mxCreateDoubleMatrix(0,0,mxREAL); hXk = mxCreateDoubleMatrix(0,0,mxCOMPLEX); if(mxIsComplex(X_IN)){ xwork = (double *) mxCalloc(MAX(1,2 * SQR(cK.rMaxn)), sizeof(double)); xpiwork = xwork + SQR(cK.rMaxn); } else xwork =(double *) mxCalloc(MAX(1,SQR(cK.rMaxn)+2*SQR(cK.hMaxn)), sizeof(double)); /* ------------------------------------------------------------ The actual job is done here:. ------------------------------------------------------------ */ if(cK.lpN){ /* ------------------------------------------------------------ LP: lab = x ------------------------------------------------------------ */ memcpy(lab, x, cK.lpN * sizeof(double)); lab += cK.lpN; x += cK.lpN; } /* ------------------------------------------------------------ CONSIDER FIRST MATLAB-REAL-TYPE: ------------------------------------------------------------ */ if(!mxIsComplex(X_IN)){ /* Not Matlab-type complex */ /* ------------------------------------------------------------ LORENTZ: (I) lab = qeig(x) ------------------------------------------------------------ */ for(k = 0; k < cK.lorN; k++){ nk = cK.lorNL[k]; qeig(lab,x,nk); lab += 2; x += nk; } /* ------------------------------------------------------------ RCONE: LAB = eig(X) (Lorentz-Rcone's are not used internally) ------------------------------------------------------------ */ for(k = 0; k < cK.rconeN; k++){ nk = cK.rconeNL[k]; rconeeig(lab,x[0],x[1],realssqr(x+2,nk-2)); lab += 2; x += nk; } /* ------------------------------------------------------------ PSD: (I) LAB = eig(X) ------------------------------------------------------------ */ if(nlhs < 2){ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(1, output_array, 1, &Xk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); /* ------------------------------------------------------------ With mexCallMATLAB, we invoked the mexFunction "eig", which allocates a matrix struct *output_array[0], AND a block for the float data of that matrix. ==> mxDestroyArray() does not only free the float data, it also releases the matrix struct (and this is what we want). ------------------------------------------------------------ */ mxDestroyArray(output_array[0]); lab += nk; x += SQR(nk); } /* ------------------------------------------------------------ WARNING: Matlab's eig doesn't recognize Hermitian, hence VERY slow ------------------------------------------------------------ */ for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork + nksqr); mexCallMATLAB(1, output_array, 1, &hXk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); mxDestroyArray(output_array[0]); lab += nk; x += 2 * nksqr; } } else{ /* ------------------------------------------------------------ SDP: (II) (Q,LAB) = eig(X) ------------------------------------------------------------ */ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(2, output_array, 1, &Xk, "eig"); nksqr = SQR(nk); /* copy Q-matrix */ memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); lab += nk; x += nksqr; q += nksqr; } for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork+nksqr); mexCallMATLAB(2, output_array, 1, &hXk, "eig"); memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); q += nksqr; if(mxIsComplex(output_array[0])) /* if any imaginary part */ memcpy(q, mxGetPi(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); lab += nk; x += 2 * nksqr; q += nksqr; } } /* [lab,q] = eigK */ } /* !iscomplex */ else{ /* is MATLAB type complex */ /* ------------------------------------------------------------ LORENTZ: (I) lab = qeig(x) ------------------------------------------------------------ */ for(k = 0; k < cK.lorN; k++){ nk = cK.lorNL[k]; cxqeig(lab,x,xpi,nk); lab += 2; x += nk; xpi += nk; } /* ------------------------------------------------------------ RCONE: LAB = eig(X) (Lorentz-Rcone's are not used internally) ------------------------------------------------------------ */ for(k = 0; k < cK.rconeN; k++){ nk = cK.rconeNL[k]; rconeeig(lab,x[0],x[1], realssqr(x+2,nk-2) + realssqr(xpi+2,nk-2)); lab += 2; x += nk; xpi += nk; } /* ------------------------------------------------------------ PSD: (I) LAB = eig(X) ------------------------------------------------------------ */ for(k = 0; k < cK.sdpN; k++){ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xpiwork,xpi,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xpiwork); if(nlhs < 2){ mexCallMATLAB(1, output_array, 1, &hXk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); } else{ mexCallMATLAB(2, output_array, 1, &hXk, "eig"); memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); if(mxIsComplex(output_array[0])) /* if any imaginary part */ memcpy(qpi, mxGetPi(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[1]); q += nksqr; qpi += nksqr; } mxDestroyArray(output_array[0]); lab += nk; x += nksqr; xpi += nksqr; } } /* iscomplex */ /* ------------------------------------------------------------ Release PSD-working arrays. ------------------------------------------------------------ */ mxSetM(Xk,0); mxSetN(Xk,0); mxSetPr(Xk, (double *) NULL); mxDestroyArray(Xk); mxSetM(hXk,0); mxSetN(hXk,0); mxSetPr(hXk, (double *) NULL); mxSetPi(hXk, (double *) NULL); mxDestroyArray(hXk); mxFree(xwork); }
/* ************************************************************ PROCEDURE mexFunction - Entry for Matlab ************************************************************ */ void mexFunction(const int nlhs, mxArray *plhs[], const int nrhs, const mxArray *prhs[]) { mxArray *output_array[3], *Xk, *hXk; coneK cK; int k, nk, nkp1, nksqr, lendiag,lenud, lenfull, i,ii; double *lab,*q,*labk,*xwork; const double *x; /* ------------------------------------------------------------ Check for proper number of arguments ------------------------------------------------------------ */ mxAssert(nrhs >= NPARIN, "psdeig requires more input arguments"); mxAssert(nlhs <= NPAROUT, "psdeig produces less output arguments"); /* ------------------------------------------------------------ Disassemble cone K structure ------------------------------------------------------------ */ conepars(K_IN, &cK); /* ------------------------------------------------------------ Compute statistics based on cone K structure ------------------------------------------------------------ */ lendiag = cK.rLen + cK.hLen; lenud = cK.rDim + cK.hDim; lenfull = cK.lpN + cK.qDim + lenud; /* ------------------------------------------------------------ Get input vector x ------------------------------------------------------------ */ mxAssert(!mxIsSparse(X_IN), "x must be full (not sparse)."); x = mxGetPr(X_IN); if(mxGetM(X_IN) * mxGetN(X_IN) != lenud){ mxAssert(mxGetM(X_IN) * mxGetN(X_IN) == lenfull, "Size mismatch x"); x += cK.lpN + cK.qDim; /* point to PSD part */ } /* ------------------------------------------------------------ Allocate output LAB(diag), eigvec Q(full for psd) ------------------------------------------------------------ */ LAB_OUT = mxCreateDoubleMatrix(lendiag, 1, mxREAL); lab = mxGetPr(LAB_OUT); if(nlhs > 1){ Q_OUT = mxCreateDoubleMatrix(lenud, 1, mxREAL); q = mxGetPr(Q_OUT); } /* ------------------------------------------------------------ Allocate working arrays: ------------------------------------------------------------ */ Xk = mxCreateDoubleMatrix(0,0,mxREAL); hXk = mxCreateDoubleMatrix(0,0,mxCOMPLEX); xwork =(double *) mxCalloc(MAX(1,SQR(cK.rMaxn)+2*SQR(cK.hMaxn)), sizeof(double)); /* ------------------------------------------------------------ PSD: (I) LAB = eig(X) ------------------------------------------------------------ */ if(nlhs < 2){ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(1, output_array, 1, &Xk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); /* ------------------------------------------------------------ With mexCallMATLAB, we invoked the mexFunction "eig", which allocates a matrix struct *output_array[0], AND a block for the float data of that matrix. ==> mxDestroyArray() does not only free the float data, it also releases the matrix struct (and this is what we want). ------------------------------------------------------------ */ mxDestroyArray(output_array[0]); lab += nk; x += SQR(nk); } /* ------------------------------------------------------------ WARNING: Matlab's eig doesn't recognize Hermitian, hence VERY slow ------------------------------------------------------------ */ for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork + nksqr); mexCallMATLAB(1, output_array, 1, &hXk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); mxDestroyArray(output_array[0]); lab += nk; x += 2 * nksqr; } } /* nlhs < 2 */ else{ /* ------------------------------------------------------------ SDP: (II) (Q,LAB) = eig(X) ------------------------------------------------------------ */ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(2, output_array, 1, &Xk, "eig"); nksqr = SQR(nk); /* copy Q-matrix */ memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); lab += nk; x += nksqr; q += nksqr; } for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork+nksqr); #ifdef USE_SVD mexCallMATLAB(3, output_array, 1, &hXk, "svd"); #else mexCallMATLAB(2, output_array, 1, &hXk, "eig"); #endif memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); q += nksqr; if(mxIsComplex(output_array[0])) /* if any imaginary part */ memcpy(q, mxGetPi(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); #ifdef USE_SVD mxDestroyArray(output_array[2]); #endif lab += nk; x += 2 * nksqr; q += nksqr; } } /* [lab,q] = eigK */ /* ------------------------------------------------------------ Release PSD-working arrays. ------------------------------------------------------------ */ mxSetM(Xk,0); mxSetN(Xk,0); mxSetPr(Xk, (double *) NULL); mxDestroyArray(Xk); mxSetM(hXk,0); mxSetN(hXk,0); mxSetPr(hXk, (double *) NULL); mxSetPi(hXk, (double *) NULL); mxDestroyArray(hXk); mxFree(xwork); }
void mexFunction( const int nlhs, mxArray *plhs[], const int nrhs, const mxArray *prhs[] ) { mxArray *output_array[3], *Xk; coneK cK; mwSize nk, nkp1, nksqr, lendiag, lenud, lenfull, nmax; mwIndex k, i, ii; double *lab, *q, *labk; const double *x; /* Argument check */ mxAssert(nrhs >= NPARIN, "psdeig requires more input arguments"); mxAssert(nlhs <= NPAROUT, "psdeig produces less output arguments"); /* Disassemble cone structure and determine output sizes */ conepars(K_IN, &cK); lendiag = cK.rLen + cK.hLen; lenud = cK.rDim + cK.hDim; lenfull = cK.lpN + cK.qDim + lenud; /* Get input vector x and skip to the PSD terms */ mxAssert( !mxIsSparse(X_IN), "x must be full (not sparse)." ); x = mxGetPr(X_IN); if ( mxGetM(X_IN) * mxGetN(X_IN) != lenud ) { mxAssert( mxGetM(X_IN) * mxGetN(X_IN) == lenfull, "Size mismatch x" ); x += cK.lpN + cK.qDim; } /* Allocate the output arrays */ LAB_OUT = mxCreateDoubleMatrix( lendiag, (mwSize)1, mxREAL ); lab = mxGetPr( LAB_OUT ); if ( nlhs > 1 ) { Q_OUT = mxCreateDoubleMatrix( lenud, (mwSize)1, mxREAL ); q = mxGetPr( Q_OUT ); } /* * Real symmetric matrices */ if ( cK.rsdpN != 0 ) { nmax = 1; for ( k = 0 ; k != cK.rsdpN ; ++k ) { nk = (mwSize)cK.sdpNL[k]; if ( nmax < nk ) nmax = nk; } Xk = mxCreateDoubleMatrix( nmax, nmax, mxREAL ); for ( k = 0 ; k != cK.rsdpN ; ++k ) { nk = (mwSize)cK.sdpNL[k]; nksqr = SQR(nk); mxSetM( Xk, nk ); mxSetN( Xk, nk ); /* Symmetric projection onto the work array */ symproj( mxGetPr(Xk), x, nk ); if ( nlhs <= 1 ) { /* One argument only: lab */ mexCallMATLAB( 1, output_array, 1, &Xk, "eig" ); memcpy( lab, mxGetPr(output_array[0]), nk * sizeof(double) ); mxDestroyArray( output_array[0] ); } else { /* First argument: Q */ mexCallMATLAB( 2, output_array, 1, &Xk, "eig" ); memcpy( q, mxGetPr(output_array[0]), nksqr * sizeof(double) ); q += nksqr; /* Second argument: extract diag(Lab) */ nkp1 = nk + 1; labk = mxGetPr( output_array[1] ); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray( output_array[0] ); mxDestroyArray( output_array[1] ); } lab += nk; x += nksqr; } mxDestroyArray( Xk ); } /* Complex Hermitian matrices * Note that if a complex argument is offered, even the "real" matrices * must be passed through here, just to be safe. */ if ( cK.sdpN != cK.rsdpN ) { nmax = 1; for ( k = cK.rsdpN ; k != cK.sdpN ; ++k ) { nk = (mwSize)cK.sdpNL[k]; if ( nmax < nk ) nmax = nk; } Xk = mxCreateDoubleMatrix( nmax, nmax, mxCOMPLEX ); for ( k = cK.rsdpN ; k != cK.sdpN ; ++k ) { nk = (mwSize)cK.sdpNL[k]; nksqr = SQR(nk); mxSetM(Xk, nk); mxSetN(Xk, nk); /* Skew-symmetric projection onto the work matrix */ symproj( mxGetPr(Xk), x, nk ); skewproj( mxGetPi(Xk), x + nksqr, nk ); if ( nlhs <= 1 ) { /* One argument only: lab */ mexCallMATLAB( 1, output_array, 1, &Xk, "eig" ); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); mxDestroyArray(output_array[0]); } else { #ifdef USE_SVD mexCallMATLAB(3, output_array, 1, &Xk, "svd"); #else mexCallMATLAB(2, output_array, 1, &Xk, "eig"); #endif /* First argument: Q */ memcpy( q, mxGetPr(output_array[0]), nksqr * sizeof(double) ); q += nksqr; if( mxIsComplex(output_array[0]) ) /* if any imaginary part */ memcpy( q, mxGetPi(output_array[0]), nksqr * sizeof(double) ); q += nksqr; /* Second argument: extract diag(Lab) */ nkp1 = nk + 1; labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); #ifdef USE_SVD mxDestroyArray(output_array[2]); #endif } lab += nk; x += 2 * nksqr; } mxDestroyArray( Xk ); } }