Actual source code: gklanczos.c
slepc-3.17.0 2022-03-31
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: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
38: {
39: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
40: PetscInt N;
42: SVDCheckStandard(svd);
43: MatGetSize(svd->A,NULL,&N);
44: SVDSetDimensions_Default(svd);
46: if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
47: svd->leftbasis = PetscNot(lanczos->oneside);
48: SVDAllocateSolution(svd,1);
49: DSSetType(svd->ds,DSSVD);
50: DSSetCompact(svd->ds,PETSC_TRUE);
51: DSSetExtraRow(svd->ds,PETSC_TRUE);
52: DSAllocate(svd->ds,svd->ncv+1);
53: PetscFunctionReturn(0);
54: }
56: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
57: {
58: PetscInt i;
59: Vec u,v;
60: PetscBool lindep=PETSC_FALSE;
62: BVGetColumn(svd->V,k,&v);
63: BVGetColumn(svd->U,k,&u);
64: MatMult(svd->A,v,u);
65: BVRestoreColumn(svd->V,k,&v);
66: BVRestoreColumn(svd->U,k,&u);
67: BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep);
68: if (PetscUnlikely(lindep)) {
69: *n = k;
70: if (breakdown) *breakdown = lindep;
71: PetscFunctionReturn(0);
72: }
74: for (i=k+1;i<*n;i++) {
75: BVGetColumn(svd->V,i,&v);
76: BVGetColumn(svd->U,i-1,&u);
77: MatMult(svd->AT,u,v);
78: BVRestoreColumn(svd->V,i,&v);
79: BVRestoreColumn(svd->U,i-1,&u);
80: BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep);
81: if (PetscUnlikely(lindep)) {
82: *n = i;
83: break;
84: }
85: BVGetColumn(svd->V,i,&v);
86: BVGetColumn(svd->U,i,&u);
87: MatMult(svd->A,v,u);
88: BVRestoreColumn(svd->V,i,&v);
89: BVRestoreColumn(svd->U,i,&u);
90: BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep);
91: if (PetscUnlikely(lindep)) {
92: *n = i;
93: break;
94: }
95: }
97: if (!lindep) {
98: BVGetColumn(svd->V,*n,&v);
99: BVGetColumn(svd->U,*n-1,&u);
100: MatMult(svd->AT,u,v);
101: BVRestoreColumn(svd->V,*n,&v);
102: BVRestoreColumn(svd->U,*n-1,&u);
103: BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep);
104: }
105: if (PetscUnlikely(breakdown)) *breakdown = lindep;
106: PetscFunctionReturn(0);
107: }
109: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
110: {
111: PetscInt i,bvl,bvk;
112: PetscReal a,b;
113: Vec z,temp;
115: BVGetActiveColumns(V,&bvl,&bvk);
116: BVGetColumn(V,k,&z);
117: MatMult(svd->A,z,u);
118: BVRestoreColumn(V,k,&z);
120: for (i=k+1;i<n;i++) {
121: BVGetColumn(V,i,&z);
122: MatMult(svd->AT,u,z);
123: BVRestoreColumn(V,i,&z);
124: VecNormBegin(u,NORM_2,&a);
125: BVSetActiveColumns(V,0,i);
126: BVDotColumnBegin(V,i,work);
127: VecNormEnd(u,NORM_2,&a);
128: BVDotColumnEnd(V,i,work);
129: VecScale(u,1.0/a);
130: BVMultColumn(V,-1.0/a,1.0/a,i,work);
132: /* h = V^* z, z = z - V h */
133: BVDotColumn(V,i,work);
134: BVMultColumn(V,-1.0,1.0,i,work);
135: BVNormColumn(V,i,NORM_2,&b);
137: BVScaleColumn(V,i,1.0/b);
139: BVGetColumn(V,i,&z);
140: MatMult(svd->A,z,u_1);
141: BVRestoreColumn(V,i,&z);
142: VecAXPY(u_1,-b,u);
143: alpha[i-1] = a;
144: beta[i-1] = b;
145: temp = u;
146: u = u_1;
147: u_1 = temp;
148: }
150: BVGetColumn(V,n,&z);
151: MatMult(svd->AT,u,z);
152: BVRestoreColumn(V,n,&z);
153: VecNormBegin(u,NORM_2,&a);
154: BVDotColumnBegin(V,n,work);
155: VecNormEnd(u,NORM_2,&a);
156: BVDotColumnEnd(V,n,work);
157: VecScale(u,1.0/a);
158: BVMultColumn(V,-1.0/a,1.0/a,n,work);
160: /* h = V^* z, z = z - V h */
161: BVDotColumn(V,n,work);
162: BVMultColumn(V,-1.0,1.0,n,work);
163: BVNormColumn(V,i,NORM_2,&b);
165: alpha[n-1] = a;
166: beta[n-1] = b;
167: BVSetActiveColumns(V,bvl,bvk);
168: PetscFunctionReturn(0);
169: }
171: /*
172: SVDKrylovConvergence - Implements the loop that checks for convergence
173: in Krylov methods.
175: Input Parameters:
176: svd - the solver; some error estimates are updated in svd->errest
177: getall - whether all residuals must be computed
178: kini - initial value of k (the loop variable)
179: nits - number of iterations of the loop
180: normr - norm of triangular factor of qr([A;B]), used only in GSVD
182: Output Parameter:
183: kout - the first index where the convergence test failed
184: */
185: PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal normr,PetscInt *kout)
186: {
187: PetscInt k,marker,ld;
188: PetscReal *alpha,*beta,*betah,resnorm;
189: PetscBool extra;
191: if (PetscUnlikely(svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it)) *kout = svd->nsv;
192: else {
193: DSGetLeadingDimension(svd->ds,&ld);
194: DSGetExtraRow(svd->ds,&extra);
196: marker = -1;
197: if (svd->trackall) getall = PETSC_TRUE;
198: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
199: beta = alpha + ld;
200: betah = alpha + 2*ld;
201: for (k=kini;k<kini+nits;k++) {
202: if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k])*normr;
203: else resnorm = PetscAbsReal(beta[k]);
204: (*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx);
205: if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
206: if (marker!=-1 && !getall) break;
207: }
208: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
209: if (marker!=-1) k = marker;
210: *kout = k;
211: }
212: PetscFunctionReturn(0);
213: }
215: PetscErrorCode SVDSolve_Lanczos(SVD svd)
216: {
217: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
218: PetscReal *alpha,*beta;
219: PetscScalar *swork,*w,*P;
220: PetscInt i,k,j,nv,ld;
221: Vec u=0,u_1=0;
222: Mat U,V;
224: /* allocate working space */
225: DSGetLeadingDimension(svd->ds,&ld);
226: PetscMalloc2(ld,&w,svd->ncv,&swork);
228: if (lanczos->oneside) {
229: MatCreateVecs(svd->A,NULL,&u);
230: MatCreateVecs(svd->A,NULL,&u_1);
231: }
233: /* normalize start vector */
234: if (!svd->nini) {
235: BVSetRandomColumn(svd->V,0);
236: BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL);
237: }
239: while (svd->reason == SVD_CONVERGED_ITERATING) {
240: svd->its++;
242: /* inner loop */
243: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
244: BVSetActiveColumns(svd->V,svd->nconv,nv);
245: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
246: beta = alpha + ld;
247: if (lanczos->oneside) SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
248: else {
249: BVSetActiveColumns(svd->U,svd->nconv,nv);
250: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL);
251: }
252: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
254: /* compute SVD of bidiagonal matrix */
255: DSSetDimensions(svd->ds,nv,svd->nconv,0);
256: DSSVDSetDimensions(svd->ds,nv);
257: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
258: DSSolve(svd->ds,w,NULL);
259: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
260: DSUpdateExtraRow(svd->ds);
261: DSSynchronize(svd->ds,w,NULL);
262: for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);
264: /* check convergence */
265: SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k);
266: (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);
268: /* compute restart vector */
269: if (svd->reason == SVD_CONVERGED_ITERATING) {
270: if (k<nv) {
271: DSGetArray(svd->ds,DS_MAT_V,&P);
272: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
273: DSRestoreArray(svd->ds,DS_MAT_V,&P);
274: BVMultColumn(svd->V,1.0,0.0,nv,swork);
275: } else {
276: /* all approximations have converged, generate a new initial vector */
277: BVSetRandomColumn(svd->V,nv);
278: BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL);
279: }
280: }
282: /* compute converged singular vectors */
283: DSGetMat(svd->ds,DS_MAT_V,&V);
284: BVMultInPlace(svd->V,V,svd->nconv,k);
285: MatDestroy(&V);
286: if (!lanczos->oneside) {
287: DSGetMat(svd->ds,DS_MAT_U,&U);
288: BVMultInPlace(svd->U,U,svd->nconv,k);
289: MatDestroy(&U);
290: }
292: /* copy restart vector from the last column */
293: if (svd->reason == SVD_CONVERGED_ITERATING) BVCopyColumn(svd->V,nv,k);
295: svd->nconv = k;
296: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
297: }
299: /* free working space */
300: VecDestroy(&u);
301: VecDestroy(&u_1);
302: PetscFree2(w,swork);
303: PetscFunctionReturn(0);
304: }
306: PetscErrorCode SVDSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
307: {
308: PetscBool set,val;
309: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
311: PetscOptionsHead(PetscOptionsObject,"SVD Lanczos Options");
313: PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
314: if (set) SVDLanczosSetOneSide(svd,val);
316: PetscOptionsTail();
317: PetscFunctionReturn(0);
318: }
320: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
321: {
322: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
324: if (lanczos->oneside != oneside) {
325: lanczos->oneside = oneside;
326: svd->state = SVD_STATE_INITIAL;
327: }
328: PetscFunctionReturn(0);
329: }
331: /*@
332: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
333: to be used is one-sided or two-sided.
335: Logically Collective on svd
337: Input Parameters:
338: + svd - singular value solver
339: - oneside - boolean flag indicating if the method is one-sided or not
341: Options Database Key:
342: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
344: Note:
345: By default, a two-sided variant is selected, which is sometimes slightly
346: more robust. However, the one-sided variant is faster because it avoids
347: the orthogonalization associated to left singular vectors. It also saves
348: the memory required for storing such vectors.
350: Level: advanced
352: .seealso: SVDTRLanczosSetOneSide()
353: @*/
354: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
355: {
358: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
359: PetscFunctionReturn(0);
360: }
362: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
363: {
364: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
366: *oneside = lanczos->oneside;
367: PetscFunctionReturn(0);
368: }
370: /*@
371: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
372: to be used is one-sided or two-sided.
374: Not Collective
376: Input Parameters:
377: . svd - singular value solver
379: Output Parameters:
380: . oneside - boolean flag indicating if the method is one-sided or not
382: Level: advanced
384: .seealso: SVDLanczosSetOneSide()
385: @*/
386: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
387: {
390: PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
391: PetscFunctionReturn(0);
392: }
394: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
395: {
396: PetscFree(svd->data);
397: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
398: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
399: PetscFunctionReturn(0);
400: }
402: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
403: {
404: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
405: PetscBool isascii;
407: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
408: if (isascii) PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
409: PetscFunctionReturn(0);
410: }
412: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
413: {
414: SVD_LANCZOS *ctx;
416: PetscNewLog(svd,&ctx);
417: svd->data = (void*)ctx;
419: svd->ops->setup = SVDSetUp_Lanczos;
420: svd->ops->solve = SVDSolve_Lanczos;
421: svd->ops->destroy = SVDDestroy_Lanczos;
422: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
423: svd->ops->view = SVDView_Lanczos;
424: svd->ops->computevectors = SVDComputeVectors_Left;
425: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
426: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
427: PetscFunctionReturn(0);
428: }