Actual source code: gklanczos.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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> /*I "slepcsvd.h" I*/
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: SVDMatGetSize(svd,NULL,&N);
45: SVDSetDimensions_Default(svd);
46: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
47: if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
48: svd->leftbasis = PetscNot(lanczos->oneside);
49: SVDAllocateSolution(svd,1);
50: DSSetType(svd->ds,DSSVD);
51: DSSetCompact(svd->ds,PETSC_TRUE);
52: DSAllocate(svd->ds,svd->ncv);
53: return(0);
54: }
56: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt n)
57: {
59: PetscInt i;
60: Vec u,v;
63: BVGetColumn(svd->V,k,&v);
64: BVGetColumn(svd->U,k,&u);
65: SVDMatMult(svd,PETSC_FALSE,v,u);
66: BVRestoreColumn(svd->V,k,&v);
67: BVRestoreColumn(svd->U,k,&u);
68: BVOrthogonalizeColumn(svd->U,k,NULL,alpha+k,NULL);
69: BVScaleColumn(U,k,1.0/alpha[k]);
71: for (i=k+1;i<n;i++) {
72: BVGetColumn(svd->V,i,&v);
73: BVGetColumn(svd->U,i-1,&u);
74: SVDMatMult(svd,PETSC_TRUE,u,v);
75: BVRestoreColumn(svd->V,i,&v);
76: BVRestoreColumn(svd->U,i-1,&u);
77: BVOrthogonalizeColumn(svd->V,i,NULL,beta+i-1,NULL);
78: BVScaleColumn(V,i,1.0/beta[i-1]);
80: BVGetColumn(svd->V,i,&v);
81: BVGetColumn(svd->U,i,&u);
82: SVDMatMult(svd,PETSC_FALSE,v,u);
83: BVRestoreColumn(svd->V,i,&v);
84: BVRestoreColumn(svd->U,i,&u);
85: BVOrthogonalizeColumn(svd->U,i,NULL,alpha+i,NULL);
86: BVScaleColumn(U,i,1.0/alpha[i]);
87: }
89: BVGetColumn(svd->V,n,&v);
90: BVGetColumn(svd->U,n-1,&u);
91: SVDMatMult(svd,PETSC_TRUE,u,v);
92: BVRestoreColumn(svd->V,n,&v);
93: BVRestoreColumn(svd->U,n-1,&u);
94: BVOrthogonalizeColumn(svd->V,n,NULL,beta+n-1,NULL);
95: return(0);
96: }
98: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
99: {
101: PetscInt i,bvl,bvk;
102: PetscReal a,b;
103: Vec z,temp;
106: BVGetActiveColumns(V,&bvl,&bvk);
107: BVGetColumn(V,k,&z);
108: SVDMatMult(svd,PETSC_FALSE,z,u);
109: BVRestoreColumn(V,k,&z);
111: for (i=k+1;i<n;i++) {
112: BVGetColumn(V,i,&z);
113: SVDMatMult(svd,PETSC_TRUE,u,z);
114: BVRestoreColumn(V,i,&z);
115: VecNormBegin(u,NORM_2,&a);
116: BVSetActiveColumns(V,0,i);
117: BVDotColumnBegin(V,i,work);
118: VecNormEnd(u,NORM_2,&a);
119: BVDotColumnEnd(V,i,work);
120: VecScale(u,1.0/a);
121: BVMultColumn(V,-1.0/a,1.0/a,i,work);
123: /* h = V^* z, z = z - V h */
124: BVDotColumn(V,i,work);
125: BVMultColumn(V,-1.0,1.0,i,work);
126: BVNormColumn(V,i,NORM_2,&b);
127: if (PetscAbsReal(b)<10*PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Recurrence generated a zero vector; use a two-sided variant");
128: BVScaleColumn(V,i,1.0/b);
130: BVGetColumn(V,i,&z);
131: SVDMatMult(svd,PETSC_FALSE,z,u_1);
132: BVRestoreColumn(V,i,&z);
133: VecAXPY(u_1,-b,u);
134: alpha[i-1] = a;
135: beta[i-1] = b;
136: temp = u;
137: u = u_1;
138: u_1 = temp;
139: }
141: BVGetColumn(V,n,&z);
142: SVDMatMult(svd,PETSC_TRUE,u,z);
143: BVRestoreColumn(V,n,&z);
144: VecNormBegin(u,NORM_2,&a);
145: BVDotColumnBegin(V,n,work);
146: VecNormEnd(u,NORM_2,&a);
147: BVDotColumnEnd(V,n,work);
148: VecScale(u,1.0/a);
149: BVMultColumn(V,-1.0/a,1.0/a,n,work);
151: /* h = V^* z, z = z - V h */
152: BVDotColumn(V,n,work);
153: BVMultColumn(V,-1.0,1.0,n,work);
154: BVNormColumn(V,i,NORM_2,&b);
156: alpha[n-1] = a;
157: beta[n-1] = b;
158: BVSetActiveColumns(V,bvl,bvk);
159: return(0);
160: }
162: PetscErrorCode SVDSolve_Lanczos(SVD svd)
163: {
165: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
166: PetscReal *alpha,*beta,lastbeta,norm,resnorm;
167: PetscScalar *swork,*w,*Q,*PT;
168: PetscInt i,k,j,nv,ld;
169: Vec u=0,u_1=0;
170: Mat U,VT;
171: PetscBool conv;
174: /* allocate working space */
175: DSGetLeadingDimension(svd->ds,&ld);
176: PetscMalloc2(ld,&w,svd->ncv,&swork);
178: if (lanczos->oneside) {
179: SVDMatCreateVecs(svd,NULL,&u);
180: SVDMatCreateVecs(svd,NULL,&u_1);
181: }
183: /* normalize start vector */
184: if (!svd->nini) {
185: BVSetRandomColumn(svd->V,0);
186: BVNormColumn(svd->V,0,NORM_2,&norm);
187: BVScaleColumn(svd->V,0,1.0/norm);
188: }
190: while (svd->reason == SVD_CONVERGED_ITERATING) {
191: svd->its++;
193: /* inner loop */
194: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
195: BVSetActiveColumns(svd->V,svd->nconv,nv);
196: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
197: beta = alpha + ld;
198: if (lanczos->oneside) {
199: SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
200: } else {
201: BVSetActiveColumns(svd->U,svd->nconv,nv);
202: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,nv);
203: }
204: lastbeta = beta[nv-1];
205: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
207: /* compute SVD of bidiagonal matrix */
208: DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
209: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
210: DSSolve(svd->ds,w,NULL);
211: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
212: DSSynchronize(svd->ds,w,NULL);
214: /* compute error estimates */
215: k = 0;
216: conv = PETSC_TRUE;
217: DSGetArray(svd->ds,DS_MAT_U,&Q);
218: for (i=svd->nconv;i<nv;i++) {
219: svd->sigma[i] = PetscRealPart(w[i]);
220: resnorm = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
221: (*svd->converged)(svd,svd->sigma[i],resnorm,&svd->errest[i],svd->convergedctx);
222: if (conv) {
223: if (svd->errest[i] < svd->tol) k++;
224: else conv = PETSC_FALSE;
225: }
226: }
227: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
229: /* check convergence */
230: (*svd->stopping)(svd,svd->its,svd->max_it,svd->nconv+k,svd->nsv,&svd->reason,svd->stoppingctx);
232: /* compute restart vector */
233: DSGetArray(svd->ds,DS_MAT_VT,&PT);
234: if (svd->reason == SVD_CONVERGED_ITERATING) {
235: if (k<nv-svd->nconv) {
236: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
237: BVMultColumn(svd->V,1.0,0.0,nv,swork);
238: } else {
239: /* all approximations have converged, generate a new initial vector */
240: BVSetRandomColumn(svd->V,nv);
241: BVOrthogonalizeColumn(svd->V,nv,NULL,&norm,NULL);
242: BVScaleColumn(svd->V,nv,1.0/norm);
243: }
244: }
245: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
247: /* compute converged singular vectors */
248: DSGetMat(svd->ds,DS_MAT_VT,&VT);
249: BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k);
250: MatDestroy(&VT);
251: if (!lanczos->oneside) {
252: DSGetMat(svd->ds,DS_MAT_U,&U);
253: BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k);
254: MatDestroy(&U);
255: }
257: /* copy restart vector from the last column */
258: if (svd->reason == SVD_CONVERGED_ITERATING) {
259: BVCopyColumn(svd->V,nv,svd->nconv+k);
260: }
262: svd->nconv += k;
263: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
264: }
266: /* free working space */
267: VecDestroy(&u);
268: VecDestroy(&u_1);
269: PetscFree2(w,swork);
270: return(0);
271: }
273: PetscErrorCode SVDSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
274: {
276: PetscBool set,val;
277: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
280: PetscOptionsHead(PetscOptionsObject,"SVD Lanczos Options");
282: PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
283: if (set) { SVDLanczosSetOneSide(svd,val); }
285: PetscOptionsTail();
286: return(0);
287: }
289: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
290: {
291: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
294: if (lanczos->oneside != oneside) {
295: lanczos->oneside = oneside;
296: svd->state = SVD_STATE_INITIAL;
297: }
298: return(0);
299: }
301: /*@
302: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
303: to be used is one-sided or two-sided.
305: Logically Collective on SVD
307: Input Parameters:
308: + svd - singular value solver
309: - oneside - boolean flag indicating if the method is one-sided or not
311: Options Database Key:
312: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
314: Note:
315: By default, a two-sided variant is selected, which is sometimes slightly
316: more robust. However, the one-sided variant is faster because it avoids
317: the orthogonalization associated to left singular vectors. It also saves
318: the memory required for storing such vectors.
320: Level: advanced
322: .seealso: SVDTRLanczosSetOneSide()
323: @*/
324: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
325: {
331: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
332: return(0);
333: }
335: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
336: {
337: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
340: *oneside = lanczos->oneside;
341: return(0);
342: }
344: /*@
345: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
346: to be used is one-sided or two-sided.
348: Not Collective
350: Input Parameters:
351: . svd - singular value solver
353: Output Parameters:
354: . oneside - boolean flag indicating if the method is one-sided or not
356: Level: advanced
358: .seealso: SVDLanczosSetOneSide()
359: @*/
360: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
361: {
367: PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
368: return(0);
369: }
371: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
372: {
376: PetscFree(svd->data);
377: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
378: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
379: return(0);
380: }
382: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
383: {
385: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
386: PetscBool isascii;
389: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
390: if (isascii) {
391: PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
392: }
393: return(0);
394: }
396: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
397: {
399: SVD_LANCZOS *ctx;
402: PetscNewLog(svd,&ctx);
403: svd->data = (void*)ctx;
405: svd->ops->setup = SVDSetUp_Lanczos;
406: svd->ops->solve = SVDSolve_Lanczos;
407: svd->ops->destroy = SVDDestroy_Lanczos;
408: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
409: svd->ops->view = SVDView_Lanczos;
410: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
411: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
412: return(0);
413: }