Actual source code: trlanczos.c
slepc-3.5.2 2014-10-10
1: /*
3: SLEPc singular value solver: "trlanczos"
5: Method: Thick-restart Lanczos
7: Algorithm:
9: Golub-Kahan-Lanczos bidiagonalization with thick-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_TRLANCZOS;
50: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
51: {
53: PetscInt N;
56: SVDMatGetSize(svd,NULL,&N);
57: SVDSetDimensions_Default(svd);
58: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
59: if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
60: svd->leftbasis = PETSC_TRUE;
61: SVDAllocateSolution(svd,1);
62: DSSetType(svd->ds,DSSVD);
63: DSSetCompact(svd->ds,PETSC_TRUE);
64: DSAllocate(svd->ds,svd->ncv);
65: return(0);
66: }
70: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n)
71: {
73: PetscReal a,b;
74: PetscScalar gamma;
75: PetscInt i,k=nconv+l;
76: Vec ui,ui1,vi;
79: BVGetColumn(V,k,&vi);
80: BVGetColumn(U,k,&ui);
81: SVDMatMult(svd,PETSC_FALSE,vi,ui);
82: BVRestoreColumn(V,k,&vi);
83: BVRestoreColumn(U,k,&ui);
84: if (l>0) {
85: BVMultColumn(U,-1.0,1.0,k,&gamma);
86: beta[nconv] = PetscRealPart(gamma);
87: }
88: BVNormColumn(U,k,NORM_2,&a);
89: BVScaleColumn(U,k,1.0/a);
90: alpha[k] = a;
92: for (i=k+1;i<n;i++) {
93: BVGetColumn(V,i,&vi);
94: BVGetColumn(U,i-1,&ui1);
95: SVDMatMult(svd,PETSC_TRUE,ui1,vi);
96: BVRestoreColumn(V,i,&vi);
97: BVRestoreColumn(U,i-1,&ui1);
98: BVOrthogonalizeColumn(V,i,NULL,&b,NULL);
99: BVScaleColumn(V,i,1.0/b);
100: beta[i-1] = b;
102: BVGetColumn(V,i,&vi);
103: BVGetColumn(U,i,&ui);
104: SVDMatMult(svd,PETSC_FALSE,vi,ui);
105: BVRestoreColumn(V,i,&vi);
106: BVGetColumn(U,i-1,&ui1);
107: VecAXPY(ui,-b,ui1);
108: BVRestoreColumn(U,i-1,&ui1);
109: BVRestoreColumn(U,i,&ui);
110: BVNormColumn(U,i,NORM_2,&a);
111: BVScaleColumn(U,i,1.0/a);
112: alpha[i] = a;
113: }
115: BVGetColumn(V,n,&vi);
116: BVGetColumn(U,n-1,&ui1);
117: SVDMatMult(svd,PETSC_TRUE,ui1,vi);
118: BVRestoreColumn(V,n,&vi);
119: BVRestoreColumn(U,n-1,&ui1);
120: BVOrthogonalizeColumn(V,n,NULL,&b,NULL);
121: beta[n-1] = b;
122: return(0);
123: }
127: /*
128: Custom CGS orthogonalization, preprocess after first orthogonalization
129: */
130: static PetscErrorCode SVDOrthogonalizeCGS(BV V,PetscInt i,PetscScalar* h,PetscReal a,BVOrthogRefineType refine,PetscReal eta,PetscReal *norm)
131: {
133: PetscReal sum,onorm;
134: PetscScalar dot;
135: PetscInt j;
138: switch (refine) {
139: case BV_ORTHOG_REFINE_NEVER:
140: BVNormColumn(V,i,NORM_2,norm);
141: break;
142: case BV_ORTHOG_REFINE_ALWAYS:
143: BVSetActiveColumns(V,0,i);
144: BVDotColumn(V,i,h);
145: BVMultColumn(V,-1.0,1.0,i,h);
146: BVNormColumn(V,i,NORM_2,norm);
147: break;
148: case BV_ORTHOG_REFINE_IFNEEDED:
149: dot = h[i];
150: onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
151: sum = 0.0;
152: for (j=0;j<i;j++) {
153: sum += PetscRealPart(h[j] * PetscConj(h[j]));
154: }
155: *norm = PetscRealPart(dot)/(a*a) - sum;
156: if (*norm>0.0) *norm = PetscSqrtReal(*norm);
157: else {
158: BVNormColumn(V,i,NORM_2,norm);
159: }
160: if (*norm < eta*onorm) {
161: BVSetActiveColumns(V,0,i);
162: BVDotColumn(V,i,h);
163: BVMultColumn(V,-1.0,1.0,i,h);
164: BVNormColumn(V,i,NORM_2,norm);
165: }
166: break;
167: }
168: return(0);
169: }
173: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
174: {
175: PetscErrorCode ierr;
176: PetscReal a,b,eta;
177: PetscScalar gamma;
178: PetscInt i,j,k=nconv+l;
179: Vec ui,ui1,vi;
180: BVOrthogRefineType refine;
183: BVGetColumn(V,k,&vi);
184: BVGetColumn(U,k,&ui);
185: SVDMatMult(svd,PETSC_FALSE,vi,ui);
186: BVRestoreColumn(V,k,&vi);
187: BVRestoreColumn(U,k,&ui);
188: if (l>0) {
189: BVMultColumn(U,-1.0,1.0,k,&gamma);
190: beta[nconv] = PetscRealPart(gamma);
191: }
192: BVGetOrthogonalization(V,NULL,&refine,&eta);
194: for (i=k+1;i<n;i++) {
195: BVGetColumn(V,i,&vi);
196: BVGetColumn(U,i-1,&ui1);
197: SVDMatMult(svd,PETSC_TRUE,ui1,vi);
198: BVRestoreColumn(V,i,&vi);
199: BVRestoreColumn(U,i-1,&ui1);
200: BVNormColumn(U,i-1,NORM_2,&a);
201: if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
202: BVSetActiveColumns(V,0,i+1);
203: BVGetColumn(V,i,&vi);
204: BVDotVec(V,vi,work);
205: BVRestoreColumn(V,i,&vi);
206: BVSetActiveColumns(V,0,i);
207: } else {
208: BVSetActiveColumns(V,0,i);
209: BVDotColumn(V,i,work);
210: }
212: BVScaleColumn(U,i-1,1.0/a);
213: for (j=0;j<i;j++) work[j] = work[j] / a;
214: BVMultColumn(V,-1.0,1.0/a,i,work);
215: SVDOrthogonalizeCGS(V,i,work,a,refine,eta,&b);
216: BVScaleColumn(V,i,1.0/b);
218: BVGetColumn(V,i,&vi);
219: BVGetColumn(U,i,&ui);
220: BVGetColumn(U,i-1,&ui1);
221: SVDMatMult(svd,PETSC_FALSE,vi,ui);
222: VecAXPY(ui,-b,ui1);
223: BVRestoreColumn(V,i,&vi);
224: BVRestoreColumn(U,i,&ui);
225: BVRestoreColumn(U,i-1,&ui1);
227: alpha[i-1] = a;
228: beta[i-1] = b;
229: }
231: BVGetColumn(V,n,&vi);
232: BVGetColumn(U,n-1,&ui1);
233: SVDMatMult(svd,PETSC_TRUE,ui1,vi);
234: BVRestoreColumn(V,n,&vi);
235: BVRestoreColumn(U,n-1,&ui1);
237: BVNormColumn(svd->U,n-1,NORM_2,&a);
238: if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
239: BVSetActiveColumns(V,0,n+1);
240: BVGetColumn(V,n,&vi);
241: BVDotVec(V,vi,work);
242: BVRestoreColumn(V,n,&vi);
243: } else {
244: BVSetActiveColumns(V,0,n);
245: BVDotColumn(V,n,work);
246: }
248: BVScaleColumn(U,n-1,1.0/a);
249: for (j=0;j<n;j++) work[j] = work[j] / a;
250: BVMultColumn(V,-1.0,1.0/a,n,work);
251: SVDOrthogonalizeCGS(V,n,work,a,refine,eta,&b);
252: BVSetActiveColumns(V,nconv,n);
253: alpha[n-1] = a;
254: beta[n-1] = b;
255: return(0);
256: }
260: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
261: {
263: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
264: PetscReal *alpha,*beta,lastbeta,norm;
265: PetscScalar *Q,*swork=NULL,*w;
266: PetscInt i,k,l,nv,ld;
267: Mat U,VT;
268: PetscBool conv;
269: BVOrthogType orthog;
272: /* allocate working space */
273: DSGetLeadingDimension(svd->ds,&ld);
274: BVGetOrthogonalization(svd->V,&orthog,NULL,NULL);
275: PetscMalloc1(ld,&w);
276: if (lanczos->oneside && orthog == BV_ORTHOG_CGS) {
277: PetscMalloc1(svd->ncv+1,&swork);
278: }
280: /* normalize start vector */
281: if (!svd->nini) {
282: BVSetRandomColumn(svd->V,0,svd->rand);
283: BVNormColumn(svd->V,0,NORM_2,&norm);
284: BVScaleColumn(svd->V,0,1.0/norm);
285: }
287: l = 0;
288: while (svd->reason == SVD_CONVERGED_ITERATING) {
289: svd->its++;
291: /* inner loop */
292: nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
293: BVSetActiveColumns(svd->V,svd->nconv,nv);
294: BVSetActiveColumns(svd->U,svd->nconv,nv);
295: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
296: beta = alpha + ld;
297: if (lanczos->oneside) {
298: if (orthog == BV_ORTHOG_MGS) {
299: SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv);
300: } else {
301: SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
302: }
303: } else {
304: SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv+l,nv);
305: }
306: lastbeta = beta[nv-1];
307: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
308: BVScaleColumn(svd->V,nv,1.0/lastbeta);
310: /* compute SVD of general matrix */
311: DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);
312: if (l==0) {
313: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
314: } else {
315: DSSetState(svd->ds,DS_STATE_RAW);
316: }
317: DSSolve(svd->ds,w,NULL);
318: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
320: /* compute error estimates */
321: k = 0;
322: conv = PETSC_TRUE;
323: DSGetArray(svd->ds,DS_MAT_U,&Q);
324: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
325: beta = alpha + ld;
326: for (i=svd->nconv;i<nv;i++) {
327: svd->sigma[i] = PetscRealPart(w[i]);
328: beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta;
329: svd->errest[i] = PetscAbsScalar(beta[i]);
330: if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
331: if (conv) {
332: if (svd->errest[i] < svd->tol) k++;
333: else conv = PETSC_FALSE;
334: }
335: }
336: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
337: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
339: /* check convergence and update l */
340: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
341: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
342: if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
343: else l = PetscMax((nv-svd->nconv-k)/2,0);
345: /* compute converged singular vectors and restart vectors */
346: DSGetMat(svd->ds,DS_MAT_VT,&VT);
347: BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k+l);
348: MatDestroy(&VT);
349: DSGetMat(svd->ds,DS_MAT_U,&U);
350: BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k+l);
351: MatDestroy(&U);
353: /* copy the last vector to be the next initial vector */
354: if (svd->reason == SVD_CONVERGED_ITERATING) {
355: BVCopyColumn(svd->V,nv,svd->nconv+k+l);
356: }
358: svd->nconv += k;
359: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
360: }
362: /* orthonormalize U columns in one side method */
363: if (lanczos->oneside) {
364: for (i=0;i<svd->nconv;i++) {
365: BVOrthogonalizeColumn(svd->U,i,NULL,&norm,NULL);
366: BVScaleColumn(svd->U,i,1.0/norm);
367: }
368: }
370: /* free working space */
371: PetscFree(w);
372: if (swork) { PetscFree(swork); }
373: return(0);
374: }
378: PetscErrorCode SVDSetFromOptions_TRLanczos(SVD svd)
379: {
381: PetscBool set,val;
382: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
385: PetscOptionsHead("SVD TRLanczos Options");
386: PetscOptionsBool("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
387: if (set) {
388: SVDTRLanczosSetOneSide(svd,val);
389: }
390: PetscOptionsTail();
391: return(0);
392: }
396: static PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
397: {
398: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
401: lanczos->oneside = oneside;
402: return(0);
403: }
407: /*@
408: SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
409: to be used is one-sided or two-sided.
411: Logically Collective on SVD
413: Input Parameters:
414: + svd - singular value solver
415: - oneside - boolean flag indicating if the method is one-sided or not
417: Options Database Key:
418: . -svd_trlanczos_oneside <boolean> - Indicates the boolean flag
420: Note:
421: By default, a two-sided variant is selected, which is sometimes slightly
422: more robust. However, the one-sided variant is faster because it avoids
423: the orthogonalization associated to left singular vectors.
425: Level: advanced
427: .seealso: SVDLanczosSetOneSide()
428: @*/
429: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
430: {
436: PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
437: return(0);
438: }
442: /*@
443: SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
444: to be used is one-sided or two-sided.
446: Not Collective
448: Input Parameters:
449: . svd - singular value solver
451: Output Parameters:
452: . oneside - boolean flag indicating if the method is one-sided or not
454: Level: advanced
456: .seealso: SVDTRLanczosSetOneSide()
457: @*/
458: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
459: {
465: PetscTryMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
466: return(0);
467: }
471: static PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
472: {
473: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
476: *oneside = lanczos->oneside;
477: return(0);
478: }
482: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
483: {
487: PetscFree(svd->data);
488: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",NULL);
489: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",NULL);
490: return(0);
491: }
495: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
496: {
498: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
501: PetscViewerASCIIPrintf(viewer," TRLanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
502: return(0);
503: }
507: PETSC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD svd)
508: {
510: SVD_TRLANCZOS *ctx;
513: PetscNewLog(svd,&ctx);
514: svd->data = (void*)ctx;
516: svd->ops->setup = SVDSetUp_TRLanczos;
517: svd->ops->solve = SVDSolve_TRLanczos;
518: svd->ops->destroy = SVDDestroy_TRLanczos;
519: svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
520: svd->ops->view = SVDView_TRLanczos;
521: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",SVDTRLanczosSetOneSide_TRLanczos);
522: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",SVDTRLanczosGetOneSide_TRLanczos);
523: return(0);
524: }