Actual source code: rsvd.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc singular value solver: "randomized"

 13:    Method: RSVD

 15:    Algorithm:

 17:        Randomized singular value decomposition.

 19:    References:

 21:        [1] N. Halko, P.-G. Martinsson, and J. A. Tropp, "Finding
 22:            structure with randomness: Probabilistic algorithms for
 23:            constructing approximate matrix decompositions", SIAM Rev.,
 24:            53(2):217-288, 2011.
 25: */

 27: #include <slepc/private/svdimpl.h>

 29: PetscErrorCode SVDSetUp_Randomized(SVD svd)
 30: {
 31:   PetscInt       N;

 33:   SVDCheckStandard(svd);
 34:   SVDCheckDefinite(svd);
 36:   MatGetSize(svd->A,NULL,&N);
 37:   SVDSetDimensions_Default(svd);
 39:   if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
 40:   svd->leftbasis = PETSC_TRUE;
 41:   svd->mpd = svd->ncv;
 42:   SVDAllocateSolution(svd,0);
 43:   DSSetType(svd->ds,DSSVD);
 44:   DSAllocate(svd->ds,svd->ncv);
 45:   SVDSetWorkVecs(svd,1,1);
 46:   return 0;
 47: }

 49: static PetscErrorCode SVDRandomizedResidualNorm(SVD svd,PetscInt i,PetscScalar sigma,PetscReal *res)
 50: {
 51:   PetscReal      norm1,norm2;
 52:   Vec            u,v,wu,wv;

 54:   *res = 1.0;
 55:   if (svd->conv!=SVD_CONV_MAXIT) {
 56:     wu = svd->swapped? svd->workr[0]: svd->workl[0];
 57:     wv = svd->swapped? svd->workl[0]: svd->workr[0];
 58:     BVGetColumn(svd->V,i,&v);
 59:     BVGetColumn(svd->U,i,&u);
 60:     /* norm1 = ||A*v-sigma*u||_2 */
 61:     MatMult(svd->A,v,wu);
 62:     VecAXPY(wu,-sigma,u);
 63:     VecNorm(wu,NORM_2,&norm1);
 64:     /* norm2 = ||A^T*u-sigma*v||_2 */
 65:     MatMult(svd->AT,u,wv);
 66:     VecAXPY(wv,-sigma,v);
 67:     VecNorm(wv,NORM_2,&norm2);
 68:     BVRestoreColumn(svd->V,i,&v);
 69:     BVRestoreColumn(svd->U,i,&u);
 70:     *res = SlepcAbs(norm1,norm2);
 71:   }
 72:   return 0;
 73: }

 75: /* If A is a virtual Hermitian transpose, then BVMatMult will fail if PRODUCT_AhB is not implemented */
 76: static PetscErrorCode BlockMatMult(BV V,Mat A,BV Y,Mat AT)
 77: {
 78:   if (!PetscDefined(USE_COMPLEX)) BVMatMult(V,A,Y);
 79:   else {
 80:     PetscBool flg=PETSC_FALSE;
 81:     PetscObjectTypeCompare((PetscObject)A,MATHERMITIANTRANSPOSEVIRTUAL,&flg);
 82:     if (flg) BVMatMultHermitianTranspose(V,AT,Y);
 83:     else BVMatMult(V,A,Y);
 84:   }
 85:   return 0;
 86: }

 88: PetscErrorCode SVDSolve_Randomized(SVD svd)
 89: {
 90:   PetscScalar    *w;
 91:   PetscReal      res=1.0;
 92:   PetscInt       i,k=0;
 93:   Mat            A,U,V;

 95:   /* Form random matrix, G. Complete the initial basis with random vectors */
 96:   BVSetActiveColumns(svd->V,svd->nini,svd->ncv);
 97:   BVSetRandomNormal(svd->V);
 98:   PetscCalloc1(svd->ncv,&w);

100:   /* Subspace Iteration */
101:   do {
102:     svd->its++;
103:     BVSetActiveColumns(svd->V,svd->nconv,svd->ncv);
104:     BVSetActiveColumns(svd->U,svd->nconv,svd->ncv);
105:     /* Form AG */
106:     BlockMatMult(svd->V,svd->A,svd->U,svd->AT);
107:     /* Orthogonalization Q=qr(AG)*/
108:     BVOrthogonalize(svd->U,NULL);
109:     /* Form B^*= AQ */
110:     BlockMatMult(svd->U,svd->AT,svd->V,svd->A);

112:     DSSetDimensions(svd->ds,svd->ncv,svd->nconv,svd->ncv);
113:     DSSVDSetDimensions(svd->ds,svd->ncv);
114:     DSGetMat(svd->ds,DS_MAT_A,&A);
115:     MatZeroEntries(A);
116:     BVOrthogonalize(svd->V,A);
117:     DSRestoreMat(svd->ds,DS_MAT_A,&A);
118:     DSSetState(svd->ds,DS_STATE_RAW);
119:     DSSolve(svd->ds,w,NULL);
120:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
121:     DSSynchronize(svd->ds,w,NULL);
122:     DSGetMat(svd->ds,DS_MAT_U,&U);
123:     DSGetMat(svd->ds,DS_MAT_V,&V);
124:     BVMultInPlace(svd->U,V,svd->nconv,svd->ncv);
125:     BVMultInPlace(svd->V,U,svd->nconv,svd->ncv);
126:     DSRestoreMat(svd->ds,DS_MAT_U,&U);
127:     DSRestoreMat(svd->ds,DS_MAT_V,&V);
128:     /* Check convergence */
129:     k = 0;
130:     for (i=svd->nconv;i<svd->ncv;i++) {
131:       SVDRandomizedResidualNorm(svd,i,w[i],&res);
132:       svd->sigma[i] = PetscRealPart(w[i]);
133:       (*svd->converged)(svd,svd->sigma[i],res,&svd->errest[i],svd->convergedctx);
134:       if (svd->errest[i] < svd->tol) k++;
135:       else break;
136:     }
137:     if (svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it) {
138:       k = svd->nsv;
139:       for (i=0;i<svd->ncv;i++) svd->sigma[i] = PetscRealPart(w[i]);
140:     }
141:     (*svd->stopping)(svd,svd->its,svd->max_it,svd->nconv+k,svd->nsv,&svd->reason,svd->stoppingctx);
142:     svd->nconv += k;
143:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);
144:   } while (svd->reason == SVD_CONVERGED_ITERATING);
145:   PetscFree(w);
146:   return 0;
147: }

149: SLEPC_EXTERN PetscErrorCode SVDCreate_Randomized(SVD svd)
150: {
151:   svd->ops->setup          = SVDSetUp_Randomized;
152:   svd->ops->solve          = SVDSolve_Randomized;
153:   return 0;
154: }