ivan@137: /************************************************************************** ivan@137: * ivan@137: * File name: ompcore.c ivan@137: * ivan@137: * Ron Rubinstein ivan@137: * Computer Science Department ivan@137: * Technion, Haifa 32000 Israel ivan@137: * ronrubin@cs ivan@137: * ivan@137: * Last Updated: 25.8.2009 ivan@137: * ivan@137: *************************************************************************/ ivan@137: ivan@137: ivan@137: #include "ompcore.h" ivan@137: #include "omputils.h" ivan@137: #include "ompprof.h" ivan@137: #include "myblas.h" ivan@137: #include ivan@137: #include ivan@137: ivan@137: ivan@137: ivan@137: /****************************************************************************** ivan@137: * * ivan@137: * Batch-OMP Implementation * ivan@137: * * ivan@137: ******************************************************************************/ ivan@137: ivan@137: mxArray* ompcore(double D[], double x[], double DtX[], double XtX[], double G[], mwSize n, mwSize m, mwSize L, ivan@137: int T, double eps, int gamma_mode, int profile, double msg_delta, int erroromp) ivan@137: { ivan@137: ivan@137: profdata pd; ivan@137: mxArray *Gamma; ivan@137: mwIndex i, j, signum, pos, *ind, *gammaIr, *gammaJc, gamma_count; ivan@137: mwSize allocated_coefs, allocated_cols; ivan@137: int DtX_specified, XtX_specified, batchomp, standardomp, *selected_atoms; ivan@137: double *alpha, *r, *Lchol, *c, *Gsub, *Dsub, sum, *gammaPr, *tempvec1, *tempvec2; ivan@137: double eps2, resnorm, delta, deltaprev, secs_remain; ivan@137: int mins_remain, hrs_remain; ivan@137: clock_t lastprint_time, starttime; ivan@137: ivan@137: ivan@137: ivan@137: /*** status flags ***/ ivan@137: ivan@137: DtX_specified = (DtX!=0); /* indicates whether D'*x was provided */ ivan@137: XtX_specified = (XtX!=0); /* indicates whether sum(x.*x) was provided */ ivan@137: ivan@137: standardomp = (G==0); /* batch-omp or standard omp are selected depending on availability of G */ ivan@137: batchomp = !standardomp; ivan@137: ivan@137: ivan@137: ivan@137: /*** allocate output matrix ***/ ivan@137: ivan@137: ivan@137: if (gamma_mode == FULL_GAMMA) { ivan@137: ivan@137: /* allocate full matrix of size m X L */ ivan@137: ivan@137: Gamma = mxCreateDoubleMatrix(m, L, mxREAL); ivan@137: gammaPr = mxGetPr(Gamma); ivan@137: gammaIr = 0; ivan@137: gammaJc = 0; ivan@137: } ivan@137: else { ivan@137: ivan@137: /* allocate sparse matrix with room for allocated_coefs nonzeros */ ivan@137: ivan@137: /* for error-omp, begin with L*sqrt(n)/2 allocated nonzeros, otherwise allocate L*T nonzeros */ ivan@137: allocated_coefs = erroromp ? (mwSize)(ceil(L*sqrt((double)n)/2.0) + 1.01) : L*T; ivan@137: Gamma = mxCreateSparse(m, L, allocated_coefs, mxREAL); ivan@137: gammaPr = mxGetPr(Gamma); ivan@137: gammaIr = mxGetIr(Gamma); ivan@137: gammaJc = mxGetJc(Gamma); ivan@137: gamma_count = 0; ivan@137: gammaJc[0] = 0; ivan@137: } ivan@137: ivan@137: ivan@137: /*** helper arrays ***/ ivan@137: ivan@137: alpha = (double*)mxMalloc(m*sizeof(double)); /* contains D'*residual */ ivan@137: ind = (mwIndex*)mxMalloc(n*sizeof(mwIndex)); /* indices of selected atoms */ ivan@137: selected_atoms = (int*)mxMalloc(m*sizeof(int)); /* binary array with 1's for selected atoms */ ivan@137: c = (double*)mxMalloc(n*sizeof(double)); /* orthogonal projection result */ ivan@137: ivan@137: /* current number of columns in Dsub / Gsub / Lchol */ ivan@137: allocated_cols = erroromp ? (mwSize)(ceil(sqrt((double)n)/2.0) + 1.01) : T; ivan@137: ivan@137: /* Cholesky decomposition of D_I'*D_I */ ivan@137: Lchol = (double*)mxMalloc(n*allocated_cols*sizeof(double)); ivan@137: ivan@137: /* temporary vectors for various computations */ ivan@137: tempvec1 = (double*)mxMalloc(m*sizeof(double)); ivan@137: tempvec2 = (double*)mxMalloc(m*sizeof(double)); ivan@137: ivan@137: if (batchomp) { ivan@137: /* matrix containing G(:,ind) - the columns of G corresponding to the selected atoms, in order of selection */ ivan@137: Gsub = (double*)mxMalloc(m*allocated_cols*sizeof(double)); ivan@137: } ivan@137: else { ivan@137: /* matrix containing D(:,ind) - the selected atoms from D, in order of selection */ ivan@137: Dsub = (double*)mxMalloc(n*allocated_cols*sizeof(double)); ivan@137: ivan@137: /* stores the residual */ ivan@137: r = (double*)mxMalloc(n*sizeof(double)); ivan@137: } ivan@137: ivan@137: if (!DtX_specified) { ivan@137: /* contains D'*x for the current signal */ ivan@137: DtX = (double*)mxMalloc(m*sizeof(double)); ivan@137: } ivan@137: ivan@137: ivan@137: ivan@137: /*** initializations for error omp ***/ ivan@137: ivan@137: if (erroromp) { ivan@137: eps2 = eps*eps; /* compute eps^2 */ ivan@137: if (T<0 || T>n) { /* unspecified max atom num - set max atoms to n */ ivan@137: T = n; ivan@137: } ivan@137: } ivan@137: ivan@137: ivan@137: ivan@137: /*** initialize timers ***/ ivan@137: ivan@137: initprofdata(&pd); /* initialize profiling counters */ ivan@137: starttime = clock(); /* record starting time for eta computations */ ivan@137: lastprint_time = starttime; /* time of last status display */ ivan@137: ivan@137: ivan@137: ivan@137: /********************** perform omp for each signal **********************/ ivan@137: ivan@137: ivan@137: ivan@137: for (signum=0; signumeps2 && T>0) { ivan@137: ivan@137: /* compute DtX */ ivan@137: ivan@137: if (!DtX_specified) { ivan@137: matT_vec(1, D, x+n*signum, DtX, n, m); ivan@137: addproftime(&pd, DtX_TIME); ivan@137: } ivan@137: ivan@137: ivan@137: /* initialize alpha := DtX */ ivan@137: ivan@137: memcpy(alpha, DtX + m*signum*DtX_specified, m*sizeof(double)); ivan@137: ivan@137: ivan@137: /* mark all atoms as unselected */ ivan@137: ivan@137: for (i=0; ieps2 && i=allocated_cols) { ivan@137: ivan@137: allocated_cols = (mwSize)(ceil(allocated_cols*MAT_INC_FACTOR) + 1.01); ivan@137: ivan@137: Lchol = (double*)mxRealloc(Lchol,n*allocated_cols*sizeof(double)); ivan@137: ivan@137: batchomp ? (Gsub = (double*)mxRealloc(Gsub,m*allocated_cols*sizeof(double))) : ivan@137: (Dsub = (double*)mxRealloc(Dsub,n*allocated_cols*sizeof(double))) ; ivan@137: } ivan@137: ivan@137: ivan@137: /* append column to Gsub or Dsub */ ivan@137: ivan@137: if (batchomp) { ivan@137: memcpy(Gsub+i*m, G+pos*m, m*sizeof(double)); ivan@137: } ivan@137: else { ivan@137: memcpy(Dsub+i*n, D+pos*n, n*sizeof(double)); ivan@137: } ivan@137: ivan@137: ivan@137: /*** Cholesky update ***/ ivan@137: ivan@137: if (i==0) { ivan@137: *Lchol = 1; ivan@137: } ivan@137: else { ivan@137: ivan@137: /* incremental Cholesky decomposition: compute next row of Lchol */ ivan@137: ivan@137: if (standardomp) { ivan@137: matT_vec(1, Dsub, D+n*pos, tempvec1, n, i); /* compute tempvec1 := Dsub'*d where d is new atom */ ivan@137: addproftime(&pd, DtD_TIME); ivan@137: } ivan@137: else { ivan@137: vec_assign(tempvec1, Gsub+i*m, ind, i); /* extract tempvec1 := Gsub(ind,i) */ ivan@137: } ivan@137: backsubst('L', Lchol, tempvec1, tempvec2, n, i); /* compute tempvec2 = Lchol \ tempvec1 */ ivan@137: for (j=0; j selected atoms are dependent */ ivan@137: break; ivan@137: } ivan@137: Lchol[i*n+i] = sqrt(1-sum); ivan@137: } ivan@137: ivan@137: addproftime(&pd, LCHOL_TIME); ivan@137: ivan@137: i++; ivan@137: ivan@137: ivan@137: /* perform orthogonal projection and compute sparse coefficients */ ivan@137: ivan@137: vec_assign(tempvec1, DtX + m*signum*DtX_specified, ind, i); /* extract tempvec1 = DtX(ind) */ ivan@137: cholsolve('L', Lchol, tempvec1, c, n, i); /* solve LL'c = tempvec1 for c */ ivan@137: addproftime(&pd, COMPCOEF_TIME); ivan@137: ivan@137: ivan@137: /* update alpha = D'*residual */ ivan@137: ivan@137: if (standardomp) { ivan@137: mat_vec(-1, Dsub, c, r, n, i); /* compute r := -Dsub*c */ ivan@137: vec_sum(1, x+n*signum, r, n); /* compute r := x+r */ ivan@137: ivan@137: ivan@137: /*memcpy(r, x+n*signum, n*sizeof(double)); /* assign r := x */ ivan@137: /*mat_vec1(-1, Dsub, c, 1, r, n, i); /* compute r := r-Dsub*c */ ivan@137: ivan@137: addproftime(&pd, COMPRES_TIME); ivan@137: matT_vec(1, D, r, alpha, n, m); /* compute alpha := D'*r */ ivan@137: addproftime(&pd, DtR_TIME); ivan@137: ivan@137: /* update residual norm */ ivan@137: if (erroromp) { ivan@137: resnorm = dotprod(r, r, n); ivan@137: addproftime(&pd, UPDATE_RESNORM_TIME); ivan@137: } ivan@137: } ivan@137: else { ivan@137: mat_vec(1, Gsub, c, tempvec1, m, i); /* compute tempvec1 := Gsub*c */ ivan@137: memcpy(alpha, DtX + m*signum*DtX_specified, m*sizeof(double)); /* set alpha = D'*x */ ivan@137: vec_sum(-1, tempvec1, alpha, m); /* compute alpha := alpha - tempvec1 */ ivan@137: addproftime(&pd, UPDATE_DtR_TIME); ivan@137: ivan@137: /* update residual norm */ ivan@137: if (erroromp) { ivan@137: vec_assign(tempvec2, tempvec1, ind, i); /* assign tempvec2 := tempvec1(ind) */ ivan@137: delta = dotprod(c,tempvec2,i); /* compute c'*tempvec2 */ ivan@137: resnorm = resnorm - delta + deltaprev; /* residual norm update */ ivan@137: deltaprev = delta; ivan@137: addproftime(&pd, UPDATE_RESNORM_TIME); ivan@137: } ivan@137: } ivan@137: } ivan@137: ivan@137: ivan@137: /*** generate output vector gamma ***/ ivan@137: ivan@137: if (gamma_mode == FULL_GAMMA) { /* write the coefs in c to their correct positions in gamma */ ivan@137: for (j=0; j= allocated_coefs) { ivan@137: ivan@137: while(gamma_count+i >= allocated_coefs) { ivan@137: allocated_coefs = (mwSize)(ceil(GAMMA_INC_FACTOR*allocated_coefs) + 1.01); ivan@137: } ivan@137: ivan@137: mxSetNzmax(Gamma, allocated_coefs); ivan@137: mxSetPr(Gamma, mxRealloc(gammaPr, allocated_coefs*sizeof(double))); ivan@137: mxSetIr(Gamma, mxRealloc(gammaIr, allocated_coefs*sizeof(mwIndex))); ivan@137: ivan@137: gammaPr = mxGetPr(Gamma); ivan@137: gammaIr = mxGetIr(Gamma); ivan@137: } ivan@137: ivan@137: /* append coefs to gamma and update the indices */ ivan@137: for (j=0; j0 && (clock()-lastprint_time)/(double)CLOCKS_PER_SEC >= msg_delta) ivan@137: { ivan@137: lastprint_time = clock(); ivan@137: ivan@137: /* estimated remainig time */ ivan@137: secs2hms( ((L-signum-1)/(double)(signum+1)) * ((lastprint_time-starttime)/(double)CLOCKS_PER_SEC) , ivan@137: &hrs_remain, &mins_remain, &secs_remain); ivan@137: ivan@137: mexPrintf("omp: signal %d / %d, estimated remaining time: %02d:%02d:%05.2f\n", ivan@137: signum+1, L, hrs_remain, mins_remain, secs_remain); ivan@137: mexEvalString("drawnow;"); ivan@137: } ivan@137: ivan@137: } ivan@137: ivan@137: /* end omp */ ivan@137: ivan@137: ivan@137: ivan@137: /*** print final messages ***/ ivan@137: ivan@137: if (msg_delta>0) { ivan@137: mexPrintf("omp: signal %d / %d\n", signum, L); ivan@137: } ivan@137: ivan@137: if (profile) { ivan@137: printprofinfo(&pd, erroromp, batchomp, L); ivan@137: } ivan@137: ivan@137: ivan@137: ivan@137: /* free memory */ ivan@137: ivan@137: if (!DtX_specified) { ivan@137: mxFree(DtX); ivan@137: } ivan@137: if (standardomp) { ivan@137: mxFree(r); ivan@137: mxFree(Dsub); ivan@137: } ivan@137: else { ivan@137: mxFree(Gsub); ivan@137: } ivan@137: mxFree(tempvec2); ivan@137: mxFree(tempvec1); ivan@137: mxFree(Lchol); ivan@137: mxFree(c); ivan@137: mxFree(selected_atoms); ivan@137: mxFree(ind); ivan@137: mxFree(alpha); ivan@137: ivan@137: return Gamma; ivan@137: }