/* ************************************************************ PROCEDURE cxqeig - computes the 2 spectral values w.r.t. Lorentz cone, complex version. INPUT: x,xpi - full n x 1, real and imaginary parts n - length of x OUTPUT: lab - 2*1, the two spectral values cxqeig(x). ************************************************************ */ void cxqeig(double *lab,const double *x,const double *xpi,const int n) { double x1, nx2; /* ------------------------------------------------------------ x1 = x(1), x2ssqr = norm( x(2:n) + i* xpi(2:n) ); labx = [x1 - nx2; x1 + nx2]/sqrt(2); ------------------------------------------------------------ */ x1 = x[0]; nx2 = sqrt(realssqr(x+1,n-1) + realssqr(xpi+1,n-1)); lab[0] = (x1 - nx2) / M_SQRT2; lab[1] = (x1 + nx2) / M_SQRT2; }
/* ************************************************************ PROCEDURE qdet - computes det(x) = lab1 * lab2 w.r.t. Lorentz cone INPUT: x - full n x 1 n - length of x RETURNS: determinant qdet(x). ************************************************************ */ double qdet(const double *x,const int n) { double y; /* ------------------------------------------------------------ qdet(x) = x'*J*x / 2 = (x_1^2 - |x_2|^2)/2 For stability, we evaluate it is (x_1+|x_2|)(x_1-|x_2|)/2 ------------------------------------------------------------ */ y = sqrt(realssqr(x+1,n-1)); return ((x[0]+y) * (x[0]-y)) / 2; }
/* ************************************************************ PROCEDURE: getada2 - Let ADA += ddota'*ddota. INPUT ada.{jc,ir} - sparsity structure of ada. ddota - sparse lorN x m matrix. perm, invperm - length(m) array, ordering in which ADA should be computed, and its inverse. We compute in order triu(ADA(perm,perm)), but store at original places. OPTIMAL PERM: sort(sum(spones(ddota))), i.e. start with sparsest. m - order of ADA, number of constraints. lorN - length(K.q), number of Lorentz blocks. UPDATED ada.pr - ada(i,j) += ddotai'*ddotaj. ONLY triu(ADA(perm,perm)) is updated. (So caller typically should symmetrize afterwards.) WORKING ARRAYS ddotaj - work vector, size lorN. ************************************************************ */ void getada2(jcir ada, jcir ddota, const mwIndex *perm, const mwIndex *invperm, const mwIndex m, const mwIndex lorN, double *ddotaj) { mwIndex i,j, knz,inz, permj; double adaij; /* ------------------------------------------------------------ Init ddotaj = all-0 (for Lorentz) ------------------------------------------------------------ */ fzeros(ddotaj, lorN); /* ============================================================ MAIN getada LOOP: loop over nodes perm(0:m-1) ============================================================ */ for(j = 0; j < m; j++){ permj = perm[j]; if(ddota.jc[permj] < ddota.jc[permj+1]){ /* Only work if nonempty */ /* ------------------------------------------------------------ Let ddotaj = ddota(:,j) in full ------------------------------------------------------------ */ for(i = ddota.jc[permj]; i < ddota.jc[permj+1]; i++) ddotaj[ddota.ir[i]] = ddota.pr[i]; /* ------------------------------------------------------------ For all i with invpermi < j: ada_ij += ddota_i'*ddotaj. ------------------------------------------------------------ */ for(inz = ada.jc[permj]; inz < ada.jc[permj+1]; inz++){ i = ada.ir[inz]; if(invperm[i] <= j){ adaij = ada.pr[inz]; if(invperm[i] < j) for(knz = ddota.jc[i]; knz < ddota.jc[i+1]; knz++) adaij += ddota.pr[knz] * ddotaj[ddota.ir[knz]]; else /* diag entry: += ||ddota(:,permj)||^2 */ adaij += realssqr(ddota.pr + ddota.jc[i], ddota.jc[i+1]-ddota.jc[i]); ada.pr[inz] = adaij; } } /* ------------------------------------------------------------ Re-initialize ddotaj = 0. ------------------------------------------------------------ */ for(i = ddota.jc[permj]; i < ddota.jc[permj+1]; i++) /* Lorentz */ ddotaj[ddota.ir[i]] = 0.0; } } /* j = 0:m-1 */ }
/* ************************************************************ PROCEDURE qrfac - QR factorization for nxn matrix. INPUT n - order of matrix to be factored UPDATED u - Full nxn. On input, u is matrix to be factored. On output, triu(u) = uppertriangular factor; tril(u,-1) = undefined. OUTPUT beta - length n vector. kth Householder reflection is Qk = I-qk*qk' / beta[k], where qk = q(k:n-1,k). q - n x (n-1) matrix; each column is a Householder reflection. ************************************************************ */ void qrfac(double *beta, double *q, double *u, const mwIndex n) { mwIndex i,k, kcol, nmink, icol; double dk, betak, qkui, qkk; for(k = 0, kcol = 0; k < n-1; k++, kcol += n+1){ /* ------------------------------------------------------------ kth Householder reflection: dk = sign(xkk) * ||xk(k:n)||, qk(k+1:n) = x(k+1:n); qkk = xkk+dk, betak = dk*qkk, ukk = -dk. ------------------------------------------------------------ */ qkk = u[kcol]; dk = SIGN(qkk) * sqrt(realssqr(u+kcol,n-k)); memcpy(q + kcol+1, u+kcol+1, (n-k-1) * sizeof(double)); qkk += dk; betak = dk * qkk; q[kcol] = qkk; if(betak == 0.0) /* If xk is all-0 then set beta = 1. */ betak = 1.0; beta[k] = betak; u[kcol] = -dk; /* ------------------------------------------------------------ Reflect columns k+1:n-1, i.e. xi -= (qk'*xi / betak) * qk, where xi = x(k:n-1, i). ------------------------------------------------------------ */ nmink = n-k; betak = -betak; for(i = k + 1, icol = kcol + n; i < n; i++, icol += n){ qkui = realdot(q+kcol, u+icol, nmink); addscalarmul(u+icol, qkui/betak, q+kcol, nmink); } } }
/* ************************************************************ 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 rotorder UPDATED u - full n x n matrix. On input, triu(u) is possibly unstable factor. On output, triu(u(:,perm)) is a stable factor. U_OUT = Q*U_IN, where Q is a sequence of givens rotations, given in g. OUTPUT perm - length n stable (column) pivot ordering. gjc - The givens rotations at step k are g[gjc[k]:gjc[k+1]-1]. The order in each column is bottom up. g - length gjc[n] <= n(n-1)/2 array of givens rotations. At worst we need n-1-k rotations in iter k=0:n-2. ************************************************************ */ void rotorder(int *perm, double *u, int *gjc, twodouble *g, double *d, const double maxusqr, const int n) { int i,j,k,inz, pivk, m; double *uj, *rowuk; double dk,y,nexty, h, uki,ukmax; twodouble gi; /* ------------------------------------------------------------ Initialize: Let perm = 1:n, inz = 0. (inz points into rotation list r) Let d(0) = 0, h = 1: this will let us compute all d's (since d(0)<h). ------------------------------------------------------------ */ for(j = 0; j < n; j++) perm[j] = j; inz = 0; d[0] = 0.0; h = 1.0; for(k = 0, rowuk = u; k < n-1; k++, rowuk++) { gjc[k] = inz; /* ------------------------------------------------------------ If current d's are not reliable then compute d(i) = sum(u(k:n-1,i).^2) from scratch. ------------------------------------------------------------ */ if(d[perm[k]] <= h) { for(j = k; j < n; j++) { i = perm[j]; d[i] = realssqr(rowuk + i*n,j+1-k); } h = d[perm[k]] * DRELTOL; } /* ------------------------------------------------------------ Let ukmax = max(U(k,perm(k+1:n)).^2) ------------------------------------------------------------ */ ukmax = 0.0; for(j = k + 1; j < n; j++) { uki = rowuk[perm[j] * n]; uki *= uki; ukmax = MAX(ukmax, uki); } /* ------------------------------------------------------------ If ukmax > maxusqr * d(k), then pivot k is unstable. If so, find best pivot: (pivk, dk) = max(perm(d(k:n))). ------------------------------------------------------------ */ if(ukmax > maxusqr * d[perm[k]]) { dk = 0.0; for(j = k+1; j < n; j++) if(d[perm[j]] > dk) { pivk = j; dk = d[perm[j]]; } /* ------------------------------------------------------------ Pivot on column pivk, and make U(:,perm) upper-triangular by pivk - k givens rotations on U(:,perm(k:n)). Givens at row i is {u(i,j), norm( u(i+1:pivk,j) )} for j=perm[pivk] and i = k:pivk-1. ------------------------------------------------------------ */ m = pivk - k; /* number of Givens rotations needed */ j = perm[pivk]; /* uj(1:m) should become 0 */ uj = rowuk + j * n; nexty = uj[m]; /* last nonzero in col uj */ y = SQR(nexty); for(i = m-1; i >= 0; i--) { gi.x = uj[i]; gi.y = nexty; y += SQR(gi.x); nexty = sqrt(y); gi.x /= nexty; /* Normalize to rotation [x,y; y,-x] */ gi.y /= nexty; g[i] = gi; } /* y == d[j] after loop */ uj[0] = nexty; /* New pivotal diagonal entry */ /* ------------------------------------------------------------ move pivot j=perm[pivk] to head of perm (shifting old k:pivk-1) ------------------------------------------------------------ */ memmove(perm+k+1, perm+k, m * sizeof(int)); /* move 1-> */ perm[k] = j; /* inserted at k */ /* ------------------------------------------------------------ Apply rotations to columns perm(k+1:n-1). Apply 1,2,...,m rotations on column k+1,..,k+m=pivk, and m rotations on cols pivk+1:n-1. ------------------------------------------------------------ */ for(i = 1; i <= m; i++) givensrotuj(rowuk + perm[k+i] * n, g,i); for(i += k; i < n; i++) givensrot(rowuk + perm[i] * n, g,m); inz += m; /* point to next avl. place */ g += m; /* ------------------------------------------------------------ Update d(perm(k+1:n)) -= u(k,perm(k+1:n)).^2. ------------------------------------------------------------ */ for(j = k + 1; j < n; j++) { i = perm[j]; d[i] -= SQR(rowuk[i * n]); } } } /* ------------------------------------------------------------ We have reordered n-1 columns of U using inz Givens-rotations. ------------------------------------------------------------ */ mxAssert(n > 0,""); gjc[n-1] = inz; }
/* ************************************************************ PROCEDURE prpiqrfac - QR factorization for nxn matrix. INPUT n - order of matrix to be factored UPDATED u - Full nxn. On input, u is matrix to be factored. On output, triu(u) = uppertriangular factor; tril(u,-1) = undefined. OUTPUT beta - length n vector. kth Householder reflection is Qk = I-qk*qk' / beta[k], where qk = q(k:n-1,k). q,qpi - n x n matrix; each of the first n-1 columns is a Householder reflection. The n-th column gives Qn = diag(q(:,n)), which is NOT Hermitian (viz. diag complex rotations). We have u_IN = Q_1*Q_2*.. *Q_n * triu(u_OUT), u_OUT = Q_n' * Q_{n-1}* .. * Q_2*Q_1*u_IN. ************************************************************ */ void prpiqrfac(double *beta, double *q, double *qpi, double *u, double *upi, const mwIndex n) { mwIndex i,j,k, kcol, nmink, icol; double betak, qkui,qkuiim, absxkk, normxk, xkk,xkkim; double *ui,*uipi, *qk, *qkpi; for(k = 0, kcol = 0; k < n-1; k++, kcol += n+1){ qk = q+kcol; qkpi = qpi + kcol; /* ------------------------------------------------------------ kth Householder reflection: Set absxkk = |xkk| and normxk = norm(xk(k:n)), then ukk = -sign(xkk) * normxk, qkk = xkk - ukk, qk(k+1:n) = xk(k+1:n). Remark: sign(xkk) := xkk/|xkk|, a complex number. ------------------------------------------------------------ */ xkk = u[kcol]; xkkim = upi[kcol]; absxkk = SQR(xkk) + SQR(xkkim); normxk = absxkk + realssqr(u+kcol+1,n-k-1) + realssqr(upi+kcol+1,n-k-1); memcpy(qk+1, u+kcol+1, (n-k-1) * sizeof(double)); /* real */ memcpy(qkpi+1, upi+kcol+1, (n-k-1) * sizeof(double)); /* imag */ absxkk = sqrt(absxkk); normxk = sqrt(normxk); if(absxkk > 0.0){ u[kcol] = -(xkk / absxkk) * normxk; /* ukk = -sign(xkk) * normxk */ upi[kcol] = -(xkkim / absxkk) * normxk; } else u[kcol] = -normxk; /* sign(0) := 1 */ qk[0] = xkk - u[kcol]; /* qkk = xkk - ukk */ qkpi[0] = xkkim - upi[kcol]; /* ------------------------------------------------------------ betak = normxk * (normxk + absxkk) EXCEPTION: if xk is all-0 then set beta = 1 (to avoid division by 0). ------------------------------------------------------------ */ if(normxk == 0.0) betak = 1.0; else betak = normxk * (normxk + absxkk); beta[k] = betak; /* ------------------------------------------------------------ Reflect columns k+1:n-1, i.e. xi -= (qk'*xi / betak) * qk, where xi = x(k:n-1, i). ------------------------------------------------------------ */ nmink = n-k; for(i = k + 1, icol = kcol + n; i < n; i++, icol += n){ ui = u + icol; uipi = upi+icol; qkui = realdot(qk, ui, nmink) + realdot(qkpi, uipi, nmink); qkuiim = realdot(qk, uipi, nmink) - realdot(qkpi, ui, nmink); qkui /= betak; qkuiim /= betak; /* for all j, we have x(j,i) -= (qkui + i * qkuiim) * (qk[j] + i*qkpi[j]) */ for(j = 0; j < nmink; j++){ ui[j] -= qkui * qk[j] - qkuiim * qkpi[j]; uipi[j] -= qkui * qkpi[j] + qkuiim * qk[j]; } } } /* k = 0:n-2 */ /* ------------------------------------------------------------ The Q*R decomposition is now done, but IM diag(u) may be nonzero. Therefore, we multiply each row with conj(sign(u_ii)) = conj(u_ii)/|u_ii|. Let q(1:n,n) = sign(diag(u)) ------------------------------------------------------------ */ mxAssert(n>=0,""); if(n > 0){ kcol = n * (n-1); /* sign column q(:,n) */ qk = q + kcol; qkpi = qpi + kcol; /* Let icol point to (i,i) entry */ for(i = 0, icol = 0; i < n; i++, icol += n+1){ xkk = u[icol]; /* get u(i,i) */ xkkim = upi[icol]; absxkk = sqrt(SQR(xkk) + SQR(xkkim)); qk[i] = xkk / absxkk; /* q(i) = sign(u(i,i)) */ qkpi[i] = xkkim / absxkk; u[icol] = absxkk; /* new u(i,i) = |uOLD(i,i)| */ upi[icol] = 0.0; } /* ------------------------------------------------------------ Let Unew = Q_n'*Uold, i.e. u(i,k) = conj(qn(i)) * uOLD(i,k), i < k ------------------------------------------------------------ */ for(k = 1, kcol = n; k < n; k++, kcol += n) isconjhadamul(u+kcol, upi+kcol, qk,qkpi,k); } }