Actual source code: rsvd.c

slepc-3.16.0 2021-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: {
 32:   PetscInt       N;

 35:   if (svd->which!=SVD_LARGEST) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"This solver supports only largest singular values");
 36:   MatGetSize(svd->A,NULL,&N);
 37:   SVDSetDimensions_Default(svd);
 38:   if (svd->ncv<svd->nsv) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The value of ncv must not be smaller than nsv");
 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: {
 52:   PetscReal      norm1,norm2;
 53:   Vec            u,v,wu,wv;

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

 78: PetscErrorCode SVDSolve_Randomized(SVD svd)
 79: {
 81:   PetscScalar    *w;
 82:   PetscReal      res=1.0;
 83:   PetscInt       i,k=0;
 84:   Mat            A,U,V;

 87:   /* Form random matrix, G. Complete the initial basis with random vectors */
 88:   BVSetActiveColumns(svd->V,svd->nini,svd->ncv);
 89:   BVSetRandomNormal(svd->V);
 90:   PetscCalloc1(svd->ncv,&w);

 92:   /* Subspace Iteration */
 93:   do {
 94:     svd->its++;
 95:     BVSetActiveColumns(svd->V,svd->nconv,svd->ncv);
 96:     BVSetActiveColumns(svd->U,svd->nconv,svd->ncv);
 97:     /* Form AG */
 98:     BVMatMult(svd->V,svd->A,svd->U);
 99:     /* Orthogonalization Q=qr(AG)*/
100:     BVOrthogonalize(svd->U,NULL);
101:     /* Form B^*= AQ */
102:     BVMatMult(svd->U,svd->AT,svd->V);

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

141: SLEPC_EXTERN PetscErrorCode SVDCreate_Randomized(SVD svd)
142: {
144:   svd->ops->setup          = SVDSetUp_Randomized;
145:   svd->ops->solve          = SVDSolve_Randomized;
146:   return(0);
147: }