Actual source code: gklanczos.c
slepc-3.15.2 2021-09-20
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: "lanczos"
13: Method: Explicitly restarted Lanczos
15: Algorithm:
17: Golub-Kahan-Lanczos bidiagonalization with explicit restart.
19: References:
21: [1] G.H. Golub and W. Kahan, "Calculating the singular values
22: and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
23: B 2:205-224, 1965.
25: [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
26: efficient parallel SVD solver based on restarted Lanczos
27: bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
28: 2008.
29: */
31: #include <slepc/private/svdimpl.h>
33: typedef struct {
34: PetscBool oneside;
35: } SVD_LANCZOS;
37: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
38: {
40: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
41: PetscInt N;
44: SVDCheckStandard(svd);
45: MatGetSize(svd->A,NULL,&N);
46: SVDSetDimensions_Default(svd);
47: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
48: if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
49: svd->leftbasis = PetscNot(lanczos->oneside);
50: SVDAllocateSolution(svd,1);
51: DSSetType(svd->ds,DSSVD);
52: DSSetCompact(svd->ds,PETSC_TRUE);
53: DSAllocate(svd->ds,svd->ncv);
54: return(0);
55: }
57: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt n)
58: {
60: PetscInt i;
61: Vec u,v;
64: BVGetColumn(svd->V,k,&v);
65: BVGetColumn(svd->U,k,&u);
66: MatMult(svd->A,v,u);
67: BVRestoreColumn(svd->V,k,&v);
68: BVRestoreColumn(svd->U,k,&u);
69: BVOrthogonalizeColumn(svd->U,k,NULL,alpha+k,NULL);
70: BVScaleColumn(U,k,1.0/alpha[k]);
72: for (i=k+1;i<n;i++) {
73: BVGetColumn(svd->V,i,&v);
74: BVGetColumn(svd->U,i-1,&u);
75: MatMult(svd->AT,u,v);
76: BVRestoreColumn(svd->V,i,&v);
77: BVRestoreColumn(svd->U,i-1,&u);
78: BVOrthogonalizeColumn(svd->V,i,NULL,beta+i-1,NULL);
79: BVScaleColumn(V,i,1.0/beta[i-1]);
81: BVGetColumn(svd->V,i,&v);
82: BVGetColumn(svd->U,i,&u);
83: MatMult(svd->A,v,u);
84: BVRestoreColumn(svd->V,i,&v);
85: BVRestoreColumn(svd->U,i,&u);
86: BVOrthogonalizeColumn(svd->U,i,NULL,alpha+i,NULL);
87: BVScaleColumn(U,i,1.0/alpha[i]);
88: }
90: BVGetColumn(svd->V,n,&v);
91: BVGetColumn(svd->U,n-1,&u);
92: MatMult(svd->AT,u,v);
93: BVRestoreColumn(svd->V,n,&v);
94: BVRestoreColumn(svd->U,n-1,&u);
95: BVOrthogonalizeColumn(svd->V,n,NULL,beta+n-1,NULL);
96: return(0);
97: }
99: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
100: {
102: PetscInt i,bvl,bvk;
103: PetscReal a,b;
104: Vec z,temp;
107: BVGetActiveColumns(V,&bvl,&bvk);
108: BVGetColumn(V,k,&z);
109: MatMult(svd->A,z,u);
110: BVRestoreColumn(V,k,&z);
112: for (i=k+1;i<n;i++) {
113: BVGetColumn(V,i,&z);
114: MatMult(svd->AT,u,z);
115: BVRestoreColumn(V,i,&z);
116: VecNormBegin(u,NORM_2,&a);
117: BVSetActiveColumns(V,0,i);
118: BVDotColumnBegin(V,i,work);
119: VecNormEnd(u,NORM_2,&a);
120: BVDotColumnEnd(V,i,work);
121: VecScale(u,1.0/a);
122: BVMultColumn(V,-1.0/a,1.0/a,i,work);
124: /* h = V^* z, z = z - V h */
125: BVDotColumn(V,i,work);
126: BVMultColumn(V,-1.0,1.0,i,work);
127: BVNormColumn(V,i,NORM_2,&b);
128: if (PetscAbsReal(b)<10*PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Recurrence generated a zero vector; use a two-sided variant");
129: BVScaleColumn(V,i,1.0/b);
131: BVGetColumn(V,i,&z);
132: MatMult(svd->A,z,u_1);
133: BVRestoreColumn(V,i,&z);
134: VecAXPY(u_1,-b,u);
135: alpha[i-1] = a;
136: beta[i-1] = b;
137: temp = u;
138: u = u_1;
139: u_1 = temp;
140: }
142: BVGetColumn(V,n,&z);
143: MatMult(svd->AT,u,z);
144: BVRestoreColumn(V,n,&z);
145: VecNormBegin(u,NORM_2,&a);
146: BVDotColumnBegin(V,n,work);
147: VecNormEnd(u,NORM_2,&a);
148: BVDotColumnEnd(V,n,work);
149: VecScale(u,1.0/a);
150: BVMultColumn(V,-1.0/a,1.0/a,n,work);
152: /* h = V^* z, z = z - V h */
153: BVDotColumn(V,n,work);
154: BVMultColumn(V,-1.0,1.0,n,work);
155: BVNormColumn(V,i,NORM_2,&b);
157: alpha[n-1] = a;
158: beta[n-1] = b;
159: BVSetActiveColumns(V,bvl,bvk);
160: return(0);
161: }
163: PetscErrorCode SVDSolve_Lanczos(SVD svd)
164: {
166: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
167: PetscReal *alpha,*beta,lastbeta,norm,resnorm;
168: PetscScalar *swork,*w,*Q,*PT;
169: PetscInt i,k,j,nv,ld;
170: Vec u=0,u_1=0;
171: Mat U,VT;
172: PetscBool conv;
175: /* allocate working space */
176: DSGetLeadingDimension(svd->ds,&ld);
177: PetscMalloc2(ld,&w,svd->ncv,&swork);
179: if (lanczos->oneside) {
180: MatCreateVecs(svd->A,NULL,&u);
181: MatCreateVecs(svd->A,NULL,&u_1);
182: }
184: /* normalize start vector */
185: if (!svd->nini) {
186: BVSetRandomColumn(svd->V,0);
187: BVNormColumn(svd->V,0,NORM_2,&norm);
188: BVScaleColumn(svd->V,0,1.0/norm);
189: }
191: while (svd->reason == SVD_CONVERGED_ITERATING) {
192: svd->its++;
194: /* inner loop */
195: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
196: BVSetActiveColumns(svd->V,svd->nconv,nv);
197: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
198: beta = alpha + ld;
199: if (lanczos->oneside) {
200: SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
201: } else {
202: BVSetActiveColumns(svd->U,svd->nconv,nv);
203: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,nv);
204: }
205: lastbeta = beta[nv-1];
206: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
208: /* compute SVD of bidiagonal matrix */
209: DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
210: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
211: DSSolve(svd->ds,w,NULL);
212: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
213: DSSynchronize(svd->ds,w,NULL);
215: /* compute error estimates */
216: k = 0;
217: conv = PETSC_TRUE;
218: DSGetArray(svd->ds,DS_MAT_U,&Q);
219: for (i=svd->nconv;i<nv;i++) {
220: svd->sigma[i] = PetscRealPart(w[i]);
221: resnorm = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
222: (*svd->converged)(svd,svd->sigma[i],resnorm,&svd->errest[i],svd->convergedctx);
223: if (conv) {
224: if (svd->errest[i] < svd->tol) k++;
225: else conv = PETSC_FALSE;
226: }
227: }
228: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
229: if (svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it) k = svd->nsv;
231: /* check convergence */
232: (*svd->stopping)(svd,svd->its,svd->max_it,svd->nconv+k,svd->nsv,&svd->reason,svd->stoppingctx);
234: /* compute restart vector */
235: DSGetArray(svd->ds,DS_MAT_VT,&PT);
236: if (svd->reason == SVD_CONVERGED_ITERATING) {
237: if (k<nv-svd->nconv) {
238: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
239: BVMultColumn(svd->V,1.0,0.0,nv,swork);
240: } else {
241: /* all approximations have converged, generate a new initial vector */
242: BVSetRandomColumn(svd->V,nv);
243: BVOrthogonalizeColumn(svd->V,nv,NULL,&norm,NULL);
244: BVScaleColumn(svd->V,nv,1.0/norm);
245: }
246: }
247: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
249: /* compute converged singular vectors */
250: DSGetMat(svd->ds,DS_MAT_VT,&VT);
251: BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k);
252: MatDestroy(&VT);
253: if (!lanczos->oneside) {
254: DSGetMat(svd->ds,DS_MAT_U,&U);
255: BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k);
256: MatDestroy(&U);
257: }
259: /* copy restart vector from the last column */
260: if (svd->reason == SVD_CONVERGED_ITERATING) {
261: BVCopyColumn(svd->V,nv,svd->nconv+k);
262: }
264: svd->nconv += k;
265: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
266: }
268: /* free working space */
269: VecDestroy(&u);
270: VecDestroy(&u_1);
271: PetscFree2(w,swork);
272: return(0);
273: }
275: PetscErrorCode SVDSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
276: {
278: PetscBool set,val;
279: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
282: PetscOptionsHead(PetscOptionsObject,"SVD Lanczos Options");
284: PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
285: if (set) { SVDLanczosSetOneSide(svd,val); }
287: PetscOptionsTail();
288: return(0);
289: }
291: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
292: {
293: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
296: if (lanczos->oneside != oneside) {
297: lanczos->oneside = oneside;
298: svd->state = SVD_STATE_INITIAL;
299: }
300: return(0);
301: }
303: /*@
304: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
305: to be used is one-sided or two-sided.
307: Logically Collective on svd
309: Input Parameters:
310: + svd - singular value solver
311: - oneside - boolean flag indicating if the method is one-sided or not
313: Options Database Key:
314: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
316: Note:
317: By default, a two-sided variant is selected, which is sometimes slightly
318: more robust. However, the one-sided variant is faster because it avoids
319: the orthogonalization associated to left singular vectors. It also saves
320: the memory required for storing such vectors.
322: Level: advanced
324: .seealso: SVDTRLanczosSetOneSide()
325: @*/
326: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
327: {
333: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
334: return(0);
335: }
337: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
338: {
339: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
342: *oneside = lanczos->oneside;
343: return(0);
344: }
346: /*@
347: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
348: to be used is one-sided or two-sided.
350: Not Collective
352: Input Parameters:
353: . svd - singular value solver
355: Output Parameters:
356: . oneside - boolean flag indicating if the method is one-sided or not
358: Level: advanced
360: .seealso: SVDLanczosSetOneSide()
361: @*/
362: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
363: {
369: PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
370: return(0);
371: }
373: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
374: {
378: PetscFree(svd->data);
379: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
380: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
381: return(0);
382: }
384: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
385: {
387: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
388: PetscBool isascii;
391: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
392: if (isascii) {
393: PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
394: }
395: return(0);
396: }
398: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
399: {
401: SVD_LANCZOS *ctx;
404: PetscNewLog(svd,&ctx);
405: svd->data = (void*)ctx;
407: svd->ops->setup = SVDSetUp_Lanczos;
408: svd->ops->solve = SVDSolve_Lanczos;
409: svd->ops->destroy = SVDDestroy_Lanczos;
410: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
411: svd->ops->view = SVDView_Lanczos;
412: svd->ops->computevectors = SVDComputeVectors_Left;
413: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
414: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
415: return(0);
416: }