Actual source code: gklanczos.c
slepc-3.6.1 2015-09-03
1: /*
3: SLEPc singular value solver: "lanczos"
5: Method: Explicitly restarted Lanczos
7: Algorithm:
9: Golub-Kahan-Lanczos bidiagonalization with explicit restart.
11: References:
13: [1] G.H. Golub and W. Kahan, "Calculating the singular values
14: and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
15: B 2:205-224, 1965.
17: [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
18: efficient parallel SVD solver based on restarted Lanczos
19: bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
20: 2008.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: SLEPc - Scalable Library for Eigenvalue Problem Computations
24: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
26: This file is part of SLEPc.
28: SLEPc is free software: you can redistribute it and/or modify it under the
29: terms of version 3 of the GNU Lesser General Public License as published by
30: the Free Software Foundation.
32: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
33: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
34: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
35: more details.
37: You should have received a copy of the GNU Lesser General Public License
38: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
39: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: */
42: #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
44: typedef struct {
45: PetscBool oneside;
46: } SVD_LANCZOS;
50: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
51: {
53: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
54: PetscInt N;
57: SVDMatGetSize(svd,NULL,&N);
58: SVDSetDimensions_Default(svd);
59: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
60: if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
61: svd->leftbasis = (lanczos->oneside)? PETSC_FALSE: PETSC_TRUE;
62: SVDAllocateSolution(svd,1);
63: DSSetType(svd->ds,DSSVD);
64: DSSetCompact(svd->ds,PETSC_TRUE);
65: DSAllocate(svd->ds,svd->ncv);
66: return(0);
67: }
71: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt n)
72: {
74: PetscInt i;
75: Vec u,v;
78: BVGetColumn(svd->V,k,&v);
79: BVGetColumn(svd->U,k,&u);
80: SVDMatMult(svd,PETSC_FALSE,v,u);
81: BVRestoreColumn(svd->V,k,&v);
82: BVRestoreColumn(svd->U,k,&u);
83: BVOrthogonalizeColumn(svd->U,k,NULL,alpha+k,NULL);
84: BVScaleColumn(U,k,1.0/alpha[k]);
86: for (i=k+1;i<n;i++) {
87: BVGetColumn(svd->V,i,&v);
88: BVGetColumn(svd->U,i-1,&u);
89: SVDMatMult(svd,PETSC_TRUE,u,v);
90: BVRestoreColumn(svd->V,i,&v);
91: BVRestoreColumn(svd->U,i-1,&u);
92: BVOrthogonalizeColumn(svd->V,i,NULL,beta+i-1,NULL);
93: BVScaleColumn(V,i,1.0/beta[i-1]);
95: BVGetColumn(svd->V,i,&v);
96: BVGetColumn(svd->U,i,&u);
97: SVDMatMult(svd,PETSC_FALSE,v,u);
98: BVRestoreColumn(svd->V,i,&v);
99: BVRestoreColumn(svd->U,i,&u);
100: BVOrthogonalizeColumn(svd->U,i,NULL,alpha+i,NULL);
101: BVScaleColumn(U,i,1.0/alpha[i]);
102: }
104: BVGetColumn(svd->V,n,&v);
105: BVGetColumn(svd->U,n-1,&u);
106: SVDMatMult(svd,PETSC_TRUE,u,v);
107: BVRestoreColumn(svd->V,n,&v);
108: BVRestoreColumn(svd->U,n-1,&u);
109: BVOrthogonalizeColumn(svd->V,n,NULL,beta+n-1,NULL);
110: return(0);
111: }
115: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
116: {
118: PetscInt i,bvl,bvk;
119: PetscReal a,b;
120: Vec z,temp;
123: BVGetActiveColumns(V,&bvl,&bvk);
124: BVGetColumn(V,k,&z);
125: SVDMatMult(svd,PETSC_FALSE,z,u);
126: BVRestoreColumn(V,k,&z);
128: for (i=k+1;i<n;i++) {
129: BVGetColumn(V,i,&z);
130: SVDMatMult(svd,PETSC_TRUE,u,z);
131: BVRestoreColumn(V,i,&z);
132: VecNormBegin(u,NORM_2,&a);
133: BVSetActiveColumns(V,0,i);
134: BVDotColumnBegin(V,i,work);
135: VecNormEnd(u,NORM_2,&a);
136: BVDotColumnEnd(V,i,work);
137: VecScale(u,1.0/a);
138: BVMultColumn(V,-1.0/a,1.0/a,i,work);
140: /* h = V^* z, z = z - V h */
141: BVDotColumn(V,i,work);
142: BVMultColumn(V,-1.0,1.0,i,work);
143: BVNormColumn(V,i,NORM_2,&b);
144: BVScaleColumn(V,i,1.0/b);
146: BVGetColumn(V,i,&z);
147: SVDMatMult(svd,PETSC_FALSE,z,u_1);
148: BVRestoreColumn(V,i,&z);
149: VecAXPY(u_1,-b,u);
150: alpha[i-1] = a;
151: beta[i-1] = b;
152: temp = u;
153: u = u_1;
154: u_1 = temp;
155: }
157: BVGetColumn(V,n,&z);
158: SVDMatMult(svd,PETSC_TRUE,u,z);
159: BVRestoreColumn(V,n,&z);
160: VecNormBegin(u,NORM_2,&a);
161: BVDotColumnBegin(V,n,work);
162: VecNormEnd(u,NORM_2,&a);
163: BVDotColumnEnd(V,n,work);
164: VecScale(u,1.0/a);
165: BVMultColumn(V,-1.0/a,1.0/a,n,work);
167: /* h = V^* z, z = z - V h */
168: BVDotColumn(V,n,work);
169: BVMultColumn(V,-1.0,1.0,n,work);
170: BVNormColumn(V,i,NORM_2,&b);
172: alpha[n-1] = a;
173: beta[n-1] = b;
174: BVSetActiveColumns(V,bvl,bvk);
175: return(0);
176: }
180: PetscErrorCode SVDSolve_Lanczos(SVD svd)
181: {
183: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
184: PetscReal *alpha,*beta,lastbeta,norm;
185: PetscScalar *swork,*w,*Q,*PT;
186: PetscInt i,k,j,nv,ld;
187: Vec u=0,u_1=0;
188: Mat U,VT;
189: PetscBool conv;
192: /* allocate working space */
193: DSGetLeadingDimension(svd->ds,&ld);
194: PetscMalloc2(ld,&w,svd->ncv,&swork);
196: if (lanczos->oneside) {
197: SVDMatCreateVecs(svd,NULL,&u);
198: SVDMatCreateVecs(svd,NULL,&u_1);
199: }
201: /* normalize start vector */
202: if (!svd->nini) {
203: BVSetRandomColumn(svd->V,0,svd->rand);
204: BVNormColumn(svd->V,0,NORM_2,&norm);
205: BVScaleColumn(svd->V,0,1.0/norm);
206: }
208: while (svd->reason == SVD_CONVERGED_ITERATING) {
209: svd->its++;
211: /* inner loop */
212: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
213: BVSetActiveColumns(svd->V,svd->nconv,nv);
214: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
215: beta = alpha + ld;
216: if (lanczos->oneside) {
217: SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
218: } else {
219: BVSetActiveColumns(svd->U,svd->nconv,nv);
220: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,nv);
221: }
222: lastbeta = beta[nv-1];
223: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
225: /* compute SVD of bidiagonal matrix */
226: DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
227: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
228: DSSolve(svd->ds,w,NULL);
229: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
231: /* compute error estimates */
232: k = 0;
233: conv = PETSC_TRUE;
234: DSGetArray(svd->ds,DS_MAT_U,&Q);
235: for (i=svd->nconv;i<nv;i++) {
236: svd->sigma[i] = PetscRealPart(w[i]);
237: svd->errest[i] = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
238: if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
239: if (conv) {
240: if (svd->errest[i] < svd->tol) k++;
241: else conv = PETSC_FALSE;
242: }
243: }
244: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
246: /* check convergence */
247: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
248: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
250: /* compute restart vector */
251: DSGetArray(svd->ds,DS_MAT_VT,&PT);
252: if (svd->reason == SVD_CONVERGED_ITERATING) {
253: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
254: BVMultColumn(svd->V,1.0,0.0,nv,swork);
255: }
256: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
258: /* compute converged singular vectors */
259: DSGetMat(svd->ds,DS_MAT_VT,&VT);
260: BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k);
261: MatDestroy(&VT);
262: if (!lanczos->oneside) {
263: DSGetMat(svd->ds,DS_MAT_U,&U);
264: BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k);
265: MatDestroy(&U);
266: }
268: /* copy restart vector from the last column */
269: if (svd->reason == SVD_CONVERGED_ITERATING) {
270: BVCopyColumn(svd->V,nv,svd->nconv+k);
271: }
273: svd->nconv += k;
274: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
275: }
277: /* free working space */
278: VecDestroy(&u);
279: VecDestroy(&u_1);
280: PetscFree2(w,swork);
281: return(0);
282: }
286: PetscErrorCode SVDSetFromOptions_Lanczos(PetscOptions *PetscOptionsObject,SVD svd)
287: {
289: PetscBool set,val;
290: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
293: PetscOptionsHead(PetscOptionsObject,"SVD Lanczos Options");
294: PetscOptionsBool("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
295: if (set) {
296: SVDLanczosSetOneSide(svd,val);
297: }
298: PetscOptionsTail();
299: return(0);
300: }
304: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
305: {
306: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
309: if (lanczos->oneside != oneside) {
310: lanczos->oneside = oneside;
311: svd->state = SVD_STATE_INITIAL;
312: }
313: return(0);
314: }
318: /*@
319: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
320: to be used is one-sided or two-sided.
322: Logically Collective on SVD
324: Input Parameters:
325: + svd - singular value solver
326: - oneside - boolean flag indicating if the method is one-sided or not
328: Options Database Key:
329: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
331: Note:
332: By default, a two-sided variant is selected, which is sometimes slightly
333: more robust. However, the one-sided variant is faster because it avoids
334: the orthogonalization associated to left singular vectors. It also saves
335: the memory required for storing such vectors.
337: Level: advanced
339: .seealso: SVDTRLanczosSetOneSide()
340: @*/
341: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
342: {
348: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
349: return(0);
350: }
354: /*@
355: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
356: to be used is one-sided or two-sided.
358: Not Collective
360: Input Parameters:
361: . svd - singular value solver
363: Output Parameters:
364: . oneside - boolean flag indicating if the method is one-sided or not
366: Level: advanced
368: .seealso: SVDLanczosSetOneSide()
369: @*/
370: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
371: {
377: PetscTryMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
378: return(0);
379: }
383: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
384: {
385: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
388: *oneside = lanczos->oneside;
389: return(0);
390: }
394: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
395: {
399: PetscFree(svd->data);
400: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
401: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
402: return(0);
403: }
407: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
408: {
410: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
413: PetscViewerASCIIPrintf(viewer," Lanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
414: return(0);
415: }
419: PETSC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
420: {
422: SVD_LANCZOS *ctx;
425: PetscNewLog(svd,&ctx);
426: svd->data = (void*)ctx;
428: svd->ops->setup = SVDSetUp_Lanczos;
429: svd->ops->solve = SVDSolve_Lanczos;
430: svd->ops->destroy = SVDDestroy_Lanczos;
431: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
432: svd->ops->view = SVDView_Lanczos;
433: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
434: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
435: return(0);
436: }