Actual source code: rsvd.c
slepc-3.17.2 2022-08-09
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;
34: MatGetSize(svd->A,NULL,&N);
35: SVDSetDimensions_Default(svd);
37: if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
38: svd->leftbasis = PETSC_TRUE;
39: svd->mpd = svd->ncv;
40: SVDAllocateSolution(svd,0);
41: DSSetType(svd->ds,DSSVD);
42: DSAllocate(svd->ds,svd->ncv);
43: SVDSetWorkVecs(svd,1,1);
44: PetscFunctionReturn(0);
45: }
47: static PetscErrorCode SVDRandomizedResidualNorm(SVD svd,PetscInt i,PetscScalar sigma,PetscReal *res)
48: {
49: PetscReal norm1,norm2;
50: Vec u,v,wu,wv;
52: wu = svd->swapped? svd->workr[0]: svd->workl[0];
53: wv = svd->swapped? svd->workl[0]: svd->workr[0];
54: if (svd->conv!=SVD_CONV_MAXIT) {
55: BVGetColumn(svd->V,i,&v);
56: BVGetColumn(svd->U,i,&u);
57: /* norm1 = ||A*v-sigma*u||_2 */
58: MatMult(svd->A,v,wu);
59: VecAXPY(wu,-sigma,u);
60: VecNorm(wu,NORM_2,&norm1);
61: /* norm2 = ||A^T*u-sigma*v||_2 */
62: MatMult(svd->AT,u,wv);
63: VecAXPY(wv,-sigma,v);
64: VecNorm(wv,NORM_2,&norm2);
65: BVRestoreColumn(svd->V,i,&v);
66: BVRestoreColumn(svd->U,i,&u);
67: *res = PetscSqrtReal(norm1*norm1+norm2*norm2);
68: } else {
69: *res = 1.0;
70: }
71: PetscFunctionReturn(0);
72: }
74: PetscErrorCode SVDSolve_Randomized(SVD svd)
75: {
76: PetscScalar *w;
77: PetscReal res=1.0;
78: PetscInt i,k=0;
79: Mat A,U,V;
81: /* Form random matrix, G. Complete the initial basis with random vectors */
82: BVSetActiveColumns(svd->V,svd->nini,svd->ncv);
83: BVSetRandomNormal(svd->V);
84: PetscCalloc1(svd->ncv,&w);
86: /* Subspace Iteration */
87: do {
88: svd->its++;
89: BVSetActiveColumns(svd->V,svd->nconv,svd->ncv);
90: BVSetActiveColumns(svd->U,svd->nconv,svd->ncv);
91: /* Form AG */
92: BVMatMult(svd->V,svd->A,svd->U);
93: /* Orthogonalization Q=qr(AG)*/
94: BVOrthogonalize(svd->U,NULL);
95: /* Form B^*= AQ */
96: BVMatMult(svd->U,svd->AT,svd->V);
98: DSSetDimensions(svd->ds,svd->ncv,svd->nconv,svd->ncv);
99: DSSVDSetDimensions(svd->ds,svd->ncv);
100: DSGetMat(svd->ds,DS_MAT_A,&A);
101: MatZeroEntries(A);
102: BVOrthogonalize(svd->V,A);
103: DSRestoreMat(svd->ds,DS_MAT_A,&A);
104: DSSetState(svd->ds,DS_STATE_RAW);
105: DSSolve(svd->ds,w,NULL);
106: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
107: DSSynchronize(svd->ds,w,NULL);
108: DSGetMat(svd->ds,DS_MAT_U,&U);
109: DSGetMat(svd->ds,DS_MAT_V,&V);
110: BVMultInPlace(svd->U,V,svd->nconv,svd->ncv);
111: BVMultInPlace(svd->V,U,svd->nconv,svd->ncv);
112: MatDestroy(&U);
113: MatDestroy(&V);
114: /* Check convergence */
115: k = 0;
116: for (i=svd->nconv;i<svd->ncv;i++) {
117: SVDRandomizedResidualNorm(svd,i,w[i],&res);
118: svd->sigma[i] = PetscRealPart(w[i]);
119: (*svd->converged)(svd,svd->sigma[i],res,&svd->errest[i],svd->convergedctx);
120: if (svd->errest[i] < svd->tol) k++;
121: else break;
122: }
123: if (svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it) {
124: k = svd->nsv;
125: for (i=0;i<svd->ncv;i++) svd->sigma[i] = PetscRealPart(w[i]);
126: }
127: (*svd->stopping)(svd,svd->its,svd->max_it,svd->nconv+k,svd->nsv,&svd->reason,svd->stoppingctx);
128: svd->nconv += k;
129: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->ncv);
130: } while (svd->reason == SVD_CONVERGED_ITERATING);
131: PetscFree(w);
132: PetscFunctionReturn(0);
133: }
135: SLEPC_EXTERN PetscErrorCode SVDCreate_Randomized(SVD svd)
136: {
137: svd->ops->setup = SVDSetUp_Randomized;
138: svd->ops->solve = SVDSolve_Randomized;
139: PetscFunctionReturn(0);
140: }