Actual source code: gklanczos.c
slepc-3.23.1 2025-05-01
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: "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: static PetscErrorCode SVDSetUp_Lanczos(SVD svd)
38: {
39: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
40: PetscInt N;
42: PetscFunctionBegin;
43: SVDCheckStandard(svd);
44: SVDCheckDefinite(svd);
45: PetscCall(MatGetSize(svd->A,NULL,&N));
46: PetscCall(SVDSetDimensions_Default(svd));
47: PetscCheck(svd->ncv<=svd->nsv+svd->mpd,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
48: if (svd->max_it==PETSC_DETERMINE) svd->max_it = PetscMax(N/svd->ncv,100)*((svd->stop==SVD_STOP_THRESHOLD)?10:1);
49: svd->leftbasis = PetscNot(lanczos->oneside);
50: PetscCall(SVDAllocateSolution(svd,1));
51: PetscCall(DSSetType(svd->ds,DSSVD));
52: PetscCall(DSSetCompact(svd->ds,PETSC_TRUE));
53: PetscCall(DSSetExtraRow(svd->ds,PETSC_TRUE));
54: PetscCall(DSAllocate(svd->ds,svd->ncv+1));
55: PetscFunctionReturn(PETSC_SUCCESS);
56: }
58: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
59: {
60: PetscInt i;
61: Vec u,v;
62: PetscBool lindep=PETSC_FALSE;
64: PetscFunctionBegin;
65: PetscCall(BVGetColumn(svd->V,k,&v));
66: PetscCall(BVGetColumn(svd->U,k,&u));
67: PetscCall(MatMult(svd->A,v,u));
68: PetscCall(BVRestoreColumn(svd->V,k,&v));
69: PetscCall(BVRestoreColumn(svd->U,k,&u));
70: PetscCall(BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep));
71: if (PetscUnlikely(lindep)) {
72: *n = k;
73: if (breakdown) *breakdown = lindep;
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: for (i=k+1;i<*n;i++) {
78: PetscCall(BVGetColumn(svd->V,i,&v));
79: PetscCall(BVGetColumn(svd->U,i-1,&u));
80: PetscCall(MatMult(svd->AT,u,v));
81: PetscCall(BVRestoreColumn(svd->V,i,&v));
82: PetscCall(BVRestoreColumn(svd->U,i-1,&u));
83: PetscCall(BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep));
84: if (PetscUnlikely(lindep)) {
85: *n = i;
86: break;
87: }
88: PetscCall(BVGetColumn(svd->V,i,&v));
89: PetscCall(BVGetColumn(svd->U,i,&u));
90: PetscCall(MatMult(svd->A,v,u));
91: PetscCall(BVRestoreColumn(svd->V,i,&v));
92: PetscCall(BVRestoreColumn(svd->U,i,&u));
93: PetscCall(BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep));
94: if (PetscUnlikely(lindep)) {
95: *n = i;
96: break;
97: }
98: }
100: if (!lindep) {
101: PetscCall(BVGetColumn(svd->V,*n,&v));
102: PetscCall(BVGetColumn(svd->U,*n-1,&u));
103: PetscCall(MatMult(svd->AT,u,v));
104: PetscCall(BVRestoreColumn(svd->V,*n,&v));
105: PetscCall(BVRestoreColumn(svd->U,*n-1,&u));
106: PetscCall(BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep));
107: }
108: if (breakdown) *breakdown = lindep;
109: PetscFunctionReturn(PETSC_SUCCESS);
110: }
112: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
113: {
114: PetscInt i,bvl,bvk;
115: PetscReal a,b;
116: Vec z,temp;
118: PetscFunctionBegin;
119: PetscCall(BVGetActiveColumns(V,&bvl,&bvk));
120: PetscCall(BVGetColumn(V,k,&z));
121: PetscCall(MatMult(svd->A,z,u));
122: PetscCall(BVRestoreColumn(V,k,&z));
124: for (i=k+1;i<n;i++) {
125: PetscCall(BVGetColumn(V,i,&z));
126: PetscCall(MatMult(svd->AT,u,z));
127: PetscCall(BVRestoreColumn(V,i,&z));
128: PetscCall(VecNormBegin(u,NORM_2,&a));
129: PetscCall(BVSetActiveColumns(V,0,i));
130: PetscCall(BVDotColumnBegin(V,i,work));
131: PetscCall(VecNormEnd(u,NORM_2,&a));
132: PetscCall(BVDotColumnEnd(V,i,work));
133: PetscCall(VecScale(u,1.0/a));
134: PetscCall(BVMultColumn(V,-1.0/a,1.0/a,i,work));
136: /* h = V^* z, z = z - V h */
137: PetscCall(BVDotColumn(V,i,work));
138: PetscCall(BVMultColumn(V,-1.0,1.0,i,work));
139: PetscCall(BVNormColumn(V,i,NORM_2,&b));
140: PetscCheck(PetscAbsReal(b)>10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_PLIB,"Recurrence generated a zero vector; use a two-sided variant");
141: PetscCall(BVScaleColumn(V,i,1.0/b));
143: PetscCall(BVGetColumn(V,i,&z));
144: PetscCall(MatMult(svd->A,z,u_1));
145: PetscCall(BVRestoreColumn(V,i,&z));
146: PetscCall(VecAXPY(u_1,-b,u));
147: alpha[i-1] = a;
148: beta[i-1] = b;
149: temp = u;
150: u = u_1;
151: u_1 = temp;
152: }
154: PetscCall(BVGetColumn(V,n,&z));
155: PetscCall(MatMult(svd->AT,u,z));
156: PetscCall(BVRestoreColumn(V,n,&z));
157: PetscCall(VecNormBegin(u,NORM_2,&a));
158: PetscCall(BVDotColumnBegin(V,n,work));
159: PetscCall(VecNormEnd(u,NORM_2,&a));
160: PetscCall(BVDotColumnEnd(V,n,work));
161: PetscCall(VecScale(u,1.0/a));
162: PetscCall(BVMultColumn(V,-1.0/a,1.0/a,n,work));
164: /* h = V^* z, z = z - V h */
165: PetscCall(BVDotColumn(V,n,work));
166: PetscCall(BVMultColumn(V,-1.0,1.0,n,work));
167: PetscCall(BVNormColumn(V,i,NORM_2,&b));
169: alpha[n-1] = a;
170: beta[n-1] = b;
171: PetscCall(BVSetActiveColumns(V,bvl,bvk));
172: PetscFunctionReturn(PETSC_SUCCESS);
173: }
175: /*
176: SVDKrylovConvergence - Implements the loop that checks for convergence
177: in Krylov methods.
179: Input Parameters:
180: svd - the solver; some error estimates are updated in svd->errest
181: getall - whether all residuals must be computed
182: kini - initial value of k (the loop variable)
183: nits - number of iterations of the loop
184: normr - norm of triangular factor of qr([A;B]), used only in GSVD
186: Output Parameter:
187: kout - the first index where the convergence test failed
188: */
189: PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal normr,PetscInt *kout)
190: {
191: PetscInt k,marker,ld;
192: PetscReal *alpha,*beta,*betah,resnorm;
193: PetscBool extra;
195: PetscFunctionBegin;
196: if (PetscUnlikely(svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it)) *kout = svd->nsv;
197: else {
198: PetscCall(DSGetLeadingDimension(svd->ds,&ld));
199: PetscCall(DSGetExtraRow(svd->ds,&extra));
200: PetscCheck(extra,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Only implemented for DS with extra row");
201: marker = -1;
202: if (svd->trackall) getall = PETSC_TRUE;
203: PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
204: beta = alpha + ld;
205: betah = alpha + 2*ld;
206: for (k=kini;k<kini+nits;k++) {
207: if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k])*normr;
208: else resnorm = PetscAbsReal(beta[k]);
209: PetscCall((*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx));
210: if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
211: if (marker!=-1 && !getall) break;
212: }
213: PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
214: if (marker!=-1) k = marker;
215: *kout = k;
216: }
217: PetscFunctionReturn(PETSC_SUCCESS);
218: }
220: static PetscErrorCode SVDSolve_Lanczos(SVD svd)
221: {
222: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
223: PetscReal *alpha,*beta;
224: PetscScalar *swork,*w,*P,*aux1,*aux2;
225: PetscInt i,k,j,nv,ld;
226: Vec u=NULL,u_1=NULL;
227: Mat U,V;
229: PetscFunctionBegin;
230: /* allocate working space */
231: PetscCall(DSGetLeadingDimension(svd->ds,&ld));
232: PetscCall(PetscMalloc2(ld,&w,svd->ncv,&swork));
234: if (lanczos->oneside) {
235: PetscCall(MatCreateVecs(svd->A,NULL,&u));
236: PetscCall(MatCreateVecs(svd->A,NULL,&u_1));
237: }
239: /* normalize start vector */
240: if (!svd->nini) {
241: PetscCall(BVSetRandomColumn(svd->V,0));
242: PetscCall(BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL));
243: }
245: while (svd->reason == SVD_CONVERGED_ITERATING) {
246: svd->its++;
248: /* inner loop */
249: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
250: PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
251: beta = alpha + ld;
252: if (lanczos->oneside) PetscCall(SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork));
253: else {
254: PetscCall(SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL));
255: PetscCall(BVSetActiveColumns(svd->U,svd->nconv,nv));
256: }
257: PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
258: PetscCall(BVSetActiveColumns(svd->V,svd->nconv,nv));
260: /* compute SVD of bidiagonal matrix */
261: PetscCall(DSSetDimensions(svd->ds,nv,svd->nconv,0));
262: PetscCall(DSSVDSetDimensions(svd->ds,nv));
263: PetscCall(DSSetState(svd->ds,DS_STATE_INTERMEDIATE));
264: PetscCall(DSSolve(svd->ds,w,NULL));
265: PetscCall(DSSort(svd->ds,w,NULL,NULL,NULL,NULL));
266: PetscCall(DSUpdateExtraRow(svd->ds));
267: PetscCall(DSSynchronize(svd->ds,w,NULL));
268: for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);
270: /* check convergence */
271: PetscCall(SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k));
272: SVDSetCtxThreshold(svd,svd->sigma,k);
273: PetscCall((*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx));
275: /* compute restart vector */
276: if (svd->reason == SVD_CONVERGED_ITERATING) {
277: if (k<nv) {
278: PetscCall(DSGetArray(svd->ds,DS_MAT_V,&P));
279: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
280: PetscCall(DSRestoreArray(svd->ds,DS_MAT_V,&P));
281: PetscCall(BVMultColumn(svd->V,1.0,0.0,nv,swork));
282: } else {
283: /* all approximations have converged, generate a new initial vector */
284: PetscCall(BVSetRandomColumn(svd->V,nv));
285: PetscCall(BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL));
286: }
287: }
289: /* compute converged singular vectors */
290: PetscCall(DSGetMat(svd->ds,DS_MAT_V,&V));
291: PetscCall(BVMultInPlace(svd->V,V,svd->nconv,k));
292: PetscCall(DSRestoreMat(svd->ds,DS_MAT_V,&V));
293: if (!lanczos->oneside) {
294: PetscCall(DSGetMat(svd->ds,DS_MAT_U,&U));
295: PetscCall(BVMultInPlace(svd->U,U,svd->nconv,k));
296: PetscCall(DSRestoreMat(svd->ds,DS_MAT_U,&U));
297: }
299: if (svd->reason == SVD_CONVERGED_ITERATING) {
300: PetscCall(BVCopyColumn(svd->V,nv,k)); /* copy restart vector from the last column */
301: if (svd->stop==SVD_STOP_THRESHOLD && nv-k<5) { /* reallocate */
302: svd->ncv = svd->mpd+k;
303: PetscCall(SVDReallocateSolution(svd,svd->ncv+1));
304: for (i=nv;i<svd->ncv;i++) svd->perm[i] = i;
305: PetscCall(DSReallocate(svd->ds,svd->ncv+1));
306: aux1 = w;
307: aux2 = swork;
308: PetscCall(PetscMalloc2(svd->ncv+1,&w,svd->ncv,&swork));
309: PetscCall(PetscArraycpy(w,aux1,ld));
310: PetscCall(PetscFree2(aux1,aux2));
311: PetscCall(DSGetLeadingDimension(svd->ds,&ld));
312: }
313: }
315: svd->nconv = k;
316: PetscCall(SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv));
317: }
319: /* free working space */
320: PetscCall(VecDestroy(&u));
321: PetscCall(VecDestroy(&u_1));
322: PetscCall(PetscFree2(w,swork));
323: PetscFunctionReturn(PETSC_SUCCESS);
324: }
326: static PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd,PetscOptionItems PetscOptionsObject)
327: {
328: PetscBool set,val;
329: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
331: PetscFunctionBegin;
332: PetscOptionsHeadBegin(PetscOptionsObject,"SVD Lanczos Options");
334: PetscCall(PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set));
335: if (set) PetscCall(SVDLanczosSetOneSide(svd,val));
337: PetscOptionsHeadEnd();
338: PetscFunctionReturn(PETSC_SUCCESS);
339: }
341: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
342: {
343: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
345: PetscFunctionBegin;
346: if (lanczos->oneside != oneside) {
347: lanczos->oneside = oneside;
348: svd->state = SVD_STATE_INITIAL;
349: }
350: PetscFunctionReturn(PETSC_SUCCESS);
351: }
353: /*@
354: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
355: to be used is one-sided or two-sided.
357: Logically Collective
359: Input Parameters:
360: + svd - singular value solver
361: - oneside - boolean flag indicating if the method is one-sided or not
363: Options Database Key:
364: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
366: Note:
367: By default, a two-sided variant is selected, which is sometimes slightly
368: more robust. However, the one-sided variant is faster because it avoids
369: the orthogonalization associated to left singular vectors. It also saves
370: the memory required for storing such vectors.
372: Level: advanced
374: .seealso: SVDTRLanczosSetOneSide()
375: @*/
376: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
377: {
378: PetscFunctionBegin;
381: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
382: PetscFunctionReturn(PETSC_SUCCESS);
383: }
385: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
386: {
387: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
389: PetscFunctionBegin;
390: *oneside = lanczos->oneside;
391: PetscFunctionReturn(PETSC_SUCCESS);
392: }
394: /*@
395: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
396: to be used is one-sided or two-sided.
398: Not Collective
400: Input Parameters:
401: . svd - singular value solver
403: Output Parameters:
404: . oneside - boolean flag indicating if the method is one-sided or not
406: Level: advanced
408: .seealso: SVDLanczosSetOneSide()
409: @*/
410: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
411: {
412: PetscFunctionBegin;
414: PetscAssertPointer(oneside,2);
415: PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
416: PetscFunctionReturn(PETSC_SUCCESS);
417: }
419: static PetscErrorCode SVDDestroy_Lanczos(SVD svd)
420: {
421: PetscFunctionBegin;
422: PetscCall(PetscFree(svd->data));
423: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL));
424: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL));
425: PetscFunctionReturn(PETSC_SUCCESS);
426: }
428: static PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
429: {
430: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
431: PetscBool isascii;
433: PetscFunctionBegin;
434: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
435: if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two"));
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
440: {
441: SVD_LANCZOS *ctx;
443: PetscFunctionBegin;
444: PetscCall(PetscNew(&ctx));
445: svd->data = (void*)ctx;
447: svd->ops->setup = SVDSetUp_Lanczos;
448: svd->ops->solve = SVDSolve_Lanczos;
449: svd->ops->destroy = SVDDestroy_Lanczos;
450: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
451: svd->ops->view = SVDView_Lanczos;
452: svd->ops->computevectors = SVDComputeVectors_Left;
453: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos));
454: PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos));
455: PetscFunctionReturn(PETSC_SUCCESS);
456: }