Actual source code: gklanczos.c
slepc-3.5.2 2014-10-10
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-2014, 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: VecNorm(u,NORM_2,&a);
133: BVSetActiveColumns(V,0,i);
134: BVDotColumn(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: BVScaleColumn(V,i,1.0/b);
144: BVGetColumn(V,i,&z);
145: SVDMatMult(svd,PETSC_FALSE,z,u_1);
146: BVRestoreColumn(V,i,&z);
147: VecAXPY(u_1,-b,u);
148: alpha[i-1] = a;
149: beta[i-1] = b;
150: temp = u;
151: u = u_1;
152: u_1 = temp;
153: }
155: BVGetColumn(V,n,&z);
156: SVDMatMult(svd,PETSC_TRUE,u,z);
157: BVRestoreColumn(V,n,&z);
158: VecNorm(u,NORM_2,&a);
159: BVDotColumn(V,n,work);
160: VecScale(u,1.0/a);
161: BVMultColumn(V,-1.0/a,1.0/a,n,work);
163: /* h = V^* z, z = z - V h */
164: BVDotColumn(V,n,work);
165: BVMultColumn(V,-1.0,1.0,n,work);
166: BVNormColumn(V,i,NORM_2,&b);
168: alpha[n-1] = a;
169: beta[n-1] = b;
170: BVSetActiveColumns(V,bvl,bvk);
171: return(0);
172: }
176: PetscErrorCode SVDSolve_Lanczos(SVD svd)
177: {
179: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
180: PetscReal *alpha,*beta,lastbeta,norm;
181: PetscScalar *swork,*w,*Q,*PT;
182: PetscInt i,k,j,nv,ld;
183: Vec u=0,u_1=0;
184: Mat U,VT;
185: PetscBool conv;
188: /* allocate working space */
189: DSGetLeadingDimension(svd->ds,&ld);
190: PetscMalloc2(ld,&w,svd->ncv,&swork);
192: if (lanczos->oneside) {
193: SVDMatGetVecs(svd,NULL,&u);
194: SVDMatGetVecs(svd,NULL,&u_1);
195: }
197: /* normalize start vector */
198: if (!svd->nini) {
199: BVSetRandomColumn(svd->V,0,svd->rand);
200: BVNormColumn(svd->V,0,NORM_2,&norm);
201: BVScaleColumn(svd->V,0,1.0/norm);
202: }
204: while (svd->reason == SVD_CONVERGED_ITERATING) {
205: svd->its++;
207: /* inner loop */
208: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
209: BVSetActiveColumns(svd->V,svd->nconv,nv);
210: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
211: beta = alpha + ld;
212: if (lanczos->oneside) {
213: SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
214: } else {
215: BVSetActiveColumns(svd->U,svd->nconv,nv);
216: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,nv);
217: }
218: lastbeta = beta[nv-1];
219: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
221: /* compute SVD of bidiagonal matrix */
222: DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
223: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
224: DSSolve(svd->ds,w,NULL);
225: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
227: /* compute error estimates */
228: k = 0;
229: conv = PETSC_TRUE;
230: DSGetArray(svd->ds,DS_MAT_U,&Q);
231: for (i=svd->nconv;i<nv;i++) {
232: svd->sigma[i] = PetscRealPart(w[i]);
233: svd->errest[i] = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
234: if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
235: if (conv) {
236: if (svd->errest[i] < svd->tol) k++;
237: else conv = PETSC_FALSE;
238: }
239: }
240: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
242: /* check convergence */
243: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
244: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
246: /* compute restart vector */
247: DSGetArray(svd->ds,DS_MAT_VT,&PT);
248: if (svd->reason == SVD_CONVERGED_ITERATING) {
249: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
250: BVMultColumn(svd->V,1.0,0.0,nv,swork);
251: }
252: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
254: /* compute converged singular vectors */
255: DSGetMat(svd->ds,DS_MAT_VT,&VT);
256: BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k);
257: MatDestroy(&VT);
258: if (!lanczos->oneside) {
259: DSGetMat(svd->ds,DS_MAT_U,&U);
260: BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k);
261: MatDestroy(&U);
262: }
264: /* copy restart vector from the last column */
265: if (svd->reason == SVD_CONVERGED_ITERATING) {
266: BVCopyColumn(svd->V,nv,svd->nconv+k);
267: }
269: svd->nconv += k;
270: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
271: }
273: /* free working space */
274: VecDestroy(&u);
275: VecDestroy(&u_1);
276: PetscFree2(w,swork);
277: return(0);
278: }
282: PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd)
283: {
285: PetscBool set,val;
286: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
289: PetscOptionsHead("SVD Lanczos Options");
290: PetscOptionsBool("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
291: if (set) {
292: SVDLanczosSetOneSide(svd,val);
293: }
294: PetscOptionsTail();
295: return(0);
296: }
300: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
301: {
302: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
305: if (lanczos->oneside != oneside) {
306: lanczos->oneside = oneside;
307: svd->setupcalled = 0;
308: }
309: return(0);
310: }
314: /*@
315: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
316: to be used is one-sided or two-sided.
318: Logically Collective on SVD
320: Input Parameters:
321: + svd - singular value solver
322: - oneside - boolean flag indicating if the method is one-sided or not
324: Options Database Key:
325: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
327: Note:
328: By default, a two-sided variant is selected, which is sometimes slightly
329: more robust. However, the one-sided variant is faster because it avoids
330: the orthogonalization associated to left singular vectors. It also saves
331: the memory required for storing such vectors.
333: Level: advanced
335: .seealso: SVDTRLanczosSetOneSide()
336: @*/
337: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
338: {
344: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
345: return(0);
346: }
350: /*@
351: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
352: to be used is one-sided or two-sided.
354: Not Collective
356: Input Parameters:
357: . svd - singular value solver
359: Output Parameters:
360: . oneside - boolean flag indicating if the method is one-sided or not
362: Level: advanced
364: .seealso: SVDLanczosSetOneSide()
365: @*/
366: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
367: {
373: PetscTryMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
374: return(0);
375: }
379: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
380: {
381: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
384: *oneside = lanczos->oneside;
385: return(0);
386: }
390: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
391: {
395: PetscFree(svd->data);
396: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
397: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
398: return(0);
399: }
403: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
404: {
406: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
409: PetscViewerASCIIPrintf(viewer," Lanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
410: return(0);
411: }
415: PETSC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
416: {
418: SVD_LANCZOS *ctx;
421: PetscNewLog(svd,&ctx);
422: svd->data = (void*)ctx;
424: svd->ops->setup = SVDSetUp_Lanczos;
425: svd->ops->solve = SVDSolve_Lanczos;
426: svd->ops->destroy = SVDDestroy_Lanczos;
427: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
428: svd->ops->view = SVDView_Lanczos;
429: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
430: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
431: return(0);
432: }