Actual source code: gklanczos.c
slepc-3.16.1 2021-11-17
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),PETSC_ERR_USER_INPUT,"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: DSSetExtraRow(svd->ds,PETSC_TRUE);
54: DSAllocate(svd->ds,svd->ncv+1);
55: return(0);
56: }
58: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
59: {
61: PetscInt i;
62: Vec u,v;
63: PetscBool lindep=PETSC_FALSE;
66: BVGetColumn(svd->V,k,&v);
67: BVGetColumn(svd->U,k,&u);
68: MatMult(svd->A,v,u);
69: BVRestoreColumn(svd->V,k,&v);
70: BVRestoreColumn(svd->U,k,&u);
71: BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep);
72: if (lindep) {
73: *n = k;
74: if (breakdown) *breakdown = lindep;
75: return(0);
76: }
78: for (i=k+1;i<*n;i++) {
79: BVGetColumn(svd->V,i,&v);
80: BVGetColumn(svd->U,i-1,&u);
81: MatMult(svd->AT,u,v);
82: BVRestoreColumn(svd->V,i,&v);
83: BVRestoreColumn(svd->U,i-1,&u);
84: BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep);
85: if (lindep) {
86: *n = i;
87: break;
88: }
89: BVGetColumn(svd->V,i,&v);
90: BVGetColumn(svd->U,i,&u);
91: MatMult(svd->A,v,u);
92: BVRestoreColumn(svd->V,i,&v);
93: BVRestoreColumn(svd->U,i,&u);
94: BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep);
95: if (lindep) {
96: *n = i;
97: break;
98: }
99: }
101: if (!lindep) {
102: BVGetColumn(svd->V,*n,&v);
103: BVGetColumn(svd->U,*n-1,&u);
104: MatMult(svd->AT,u,v);
105: BVRestoreColumn(svd->V,*n,&v);
106: BVRestoreColumn(svd->U,*n-1,&u);
107: BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep);
108: }
109: if (breakdown) *breakdown = lindep;
110: return(0);
111: }
113: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
114: {
116: PetscInt i,bvl,bvk;
117: PetscReal a,b;
118: Vec z,temp;
121: BVGetActiveColumns(V,&bvl,&bvk);
122: BVGetColumn(V,k,&z);
123: MatMult(svd->A,z,u);
124: BVRestoreColumn(V,k,&z);
126: for (i=k+1;i<n;i++) {
127: BVGetColumn(V,i,&z);
128: MatMult(svd->AT,u,z);
129: BVRestoreColumn(V,i,&z);
130: VecNormBegin(u,NORM_2,&a);
131: BVSetActiveColumns(V,0,i);
132: BVDotColumnBegin(V,i,work);
133: VecNormEnd(u,NORM_2,&a);
134: BVDotColumnEnd(V,i,work);
135: VecScale(u,1.0/a);
136: BVMultColumn(V,-1.0/a,1.0/a,i,work);
138: /* h = V^* z, z = z - V h */
139: BVDotColumn(V,i,work);
140: BVMultColumn(V,-1.0,1.0,i,work);
141: BVNormColumn(V,i,NORM_2,&b);
142: if (PetscAbsReal(b)<10*PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_PLIB,"Recurrence generated a zero vector; use a two-sided variant");
143: BVScaleColumn(V,i,1.0/b);
145: BVGetColumn(V,i,&z);
146: MatMult(svd->A,z,u_1);
147: BVRestoreColumn(V,i,&z);
148: VecAXPY(u_1,-b,u);
149: alpha[i-1] = a;
150: beta[i-1] = b;
151: temp = u;
152: u = u_1;
153: u_1 = temp;
154: }
156: BVGetColumn(V,n,&z);
157: MatMult(svd->AT,u,z);
158: BVRestoreColumn(V,n,&z);
159: VecNormBegin(u,NORM_2,&a);
160: BVDotColumnBegin(V,n,work);
161: VecNormEnd(u,NORM_2,&a);
162: BVDotColumnEnd(V,n,work);
163: VecScale(u,1.0/a);
164: BVMultColumn(V,-1.0/a,1.0/a,n,work);
166: /* h = V^* z, z = z - V h */
167: BVDotColumn(V,n,work);
168: BVMultColumn(V,-1.0,1.0,n,work);
169: BVNormColumn(V,i,NORM_2,&b);
171: alpha[n-1] = a;
172: beta[n-1] = b;
173: BVSetActiveColumns(V,bvl,bvk);
174: return(0);
175: }
177: /*
178: SVDKrylovConvergence - Implements the loop that checks for convergence
179: in Krylov methods.
181: Input Parameters:
182: svd - the solver; some error estimates are updated in svd->errest
183: getall - whether all residuals must be computed
184: kini - initial value of k (the loop variable)
185: nits - number of iterations of the loop
187: Output Parameter:
188: kout - the first index where the convergence test failed
189: */
190: PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscInt *kout)
191: {
193: PetscInt k,marker,ld;
194: PetscReal *alpha,*beta,*betah,resnorm;
195: PetscBool extra;
198: if (svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it) *kout = svd->nsv;
199: else {
200: DSGetLeadingDimension(svd->ds,&ld);
201: DSGetExtraRow(svd->ds,&extra);
202: if (!extra) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Only implemented for DS with extra row");
203: marker = -1;
204: if (svd->trackall) getall = PETSC_TRUE;
205: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
206: beta = alpha + ld;
207: betah = alpha + 2*ld;
208: for (k=kini;k<kini+nits;k++) {
209: if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k]);
210: else resnorm = PetscAbsReal(beta[k]);
211: (*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx);
212: if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
213: if (marker!=-1 && !getall) break;
214: }
215: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
216: if (marker!=-1) k = marker;
217: *kout = k;
218: }
219: return(0);
220: }
222: PetscErrorCode SVDSolve_Lanczos(SVD svd)
223: {
225: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
226: PetscReal *alpha,*beta;
227: PetscScalar *swork,*w,*P;
228: PetscInt i,k,j,nv,ld;
229: Vec u=0,u_1=0;
230: Mat U,V;
233: /* allocate working space */
234: DSGetLeadingDimension(svd->ds,&ld);
235: PetscMalloc2(ld,&w,svd->ncv,&swork);
237: if (lanczos->oneside) {
238: MatCreateVecs(svd->A,NULL,&u);
239: MatCreateVecs(svd->A,NULL,&u_1);
240: }
242: /* normalize start vector */
243: if (!svd->nini) {
244: BVSetRandomColumn(svd->V,0);
245: BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL);
246: }
248: while (svd->reason == SVD_CONVERGED_ITERATING) {
249: svd->its++;
251: /* inner loop */
252: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
253: BVSetActiveColumns(svd->V,svd->nconv,nv);
254: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
255: beta = alpha + ld;
256: if (lanczos->oneside) {
257: SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
258: } else {
259: BVSetActiveColumns(svd->U,svd->nconv,nv);
260: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL);
261: }
262: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
264: /* compute SVD of bidiagonal matrix */
265: DSSetDimensions(svd->ds,nv,svd->nconv,0);
266: DSSVDSetDimensions(svd->ds,nv);
267: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
268: DSSolve(svd->ds,w,NULL);
269: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
270: DSUpdateExtraRow(svd->ds);
271: DSSynchronize(svd->ds,w,NULL);
272: for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);
274: /* check convergence */
275: SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,&k);
276: (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);
278: /* compute restart vector */
279: if (svd->reason == SVD_CONVERGED_ITERATING) {
280: if (k<nv) {
281: DSGetArray(svd->ds,DS_MAT_V,&P);
282: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
283: DSRestoreArray(svd->ds,DS_MAT_V,&P);
284: BVMultColumn(svd->V,1.0,0.0,nv,swork);
285: } else {
286: /* all approximations have converged, generate a new initial vector */
287: BVSetRandomColumn(svd->V,nv);
288: BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL);
289: }
290: }
292: /* compute converged singular vectors */
293: DSGetMat(svd->ds,DS_MAT_V,&V);
294: BVMultInPlace(svd->V,V,svd->nconv,k);
295: MatDestroy(&V);
296: if (!lanczos->oneside) {
297: DSGetMat(svd->ds,DS_MAT_U,&U);
298: BVMultInPlace(svd->U,U,svd->nconv,k);
299: MatDestroy(&U);
300: }
302: /* copy restart vector from the last column */
303: if (svd->reason == SVD_CONVERGED_ITERATING) {
304: BVCopyColumn(svd->V,nv,k);
305: }
307: svd->nconv = k;
308: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
309: }
311: /* free working space */
312: VecDestroy(&u);
313: VecDestroy(&u_1);
314: PetscFree2(w,swork);
315: return(0);
316: }
318: PetscErrorCode SVDSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
319: {
321: PetscBool set,val;
322: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
325: PetscOptionsHead(PetscOptionsObject,"SVD Lanczos Options");
327: PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
328: if (set) { SVDLanczosSetOneSide(svd,val); }
330: PetscOptionsTail();
331: return(0);
332: }
334: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
335: {
336: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
339: if (lanczos->oneside != oneside) {
340: lanczos->oneside = oneside;
341: svd->state = SVD_STATE_INITIAL;
342: }
343: return(0);
344: }
346: /*@
347: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
348: to be used is one-sided or two-sided.
350: Logically Collective on svd
352: Input Parameters:
353: + svd - singular value solver
354: - oneside - boolean flag indicating if the method is one-sided or not
356: Options Database Key:
357: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
359: Note:
360: By default, a two-sided variant is selected, which is sometimes slightly
361: more robust. However, the one-sided variant is faster because it avoids
362: the orthogonalization associated to left singular vectors. It also saves
363: the memory required for storing such vectors.
365: Level: advanced
367: .seealso: SVDTRLanczosSetOneSide()
368: @*/
369: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
370: {
376: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
377: return(0);
378: }
380: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
381: {
382: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
385: *oneside = lanczos->oneside;
386: return(0);
387: }
389: /*@
390: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
391: to be used is one-sided or two-sided.
393: Not Collective
395: Input Parameters:
396: . svd - singular value solver
398: Output Parameters:
399: . oneside - boolean flag indicating if the method is one-sided or not
401: Level: advanced
403: .seealso: SVDLanczosSetOneSide()
404: @*/
405: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
406: {
412: PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
413: return(0);
414: }
416: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
417: {
421: PetscFree(svd->data);
422: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
423: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
424: return(0);
425: }
427: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
428: {
430: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
431: PetscBool isascii;
434: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
435: if (isascii) {
436: PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
437: }
438: return(0);
439: }
441: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
442: {
444: SVD_LANCZOS *ctx;
447: PetscNewLog(svd,&ctx);
448: svd->data = (void*)ctx;
450: svd->ops->setup = SVDSetUp_Lanczos;
451: svd->ops->solve = SVDSolve_Lanczos;
452: svd->ops->destroy = SVDDestroy_Lanczos;
453: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
454: svd->ops->view = SVDView_Lanczos;
455: svd->ops->computevectors = SVDComputeVectors_Left;
456: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
457: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
458: return(0);
459: }