Actual source code: gklanczos.c
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: Last update: Jun 2007
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: SLEPc - Scalable Library for Eigenvalue Problem Computations
26: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
28: This file is part of SLEPc.
30: SLEPc is free software: you can redistribute it and/or modify it under the
31: terms of version 3 of the GNU Lesser General Public License as published by
32: the Free Software Foundation.
34: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
35: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
36: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
37: more details.
39: You should have received a copy of the GNU Lesser General Public License
40: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42: */
44: #include <slepc-private/svdimpl.h> /*I "slepcsvd.h" I*/
45: #include <slepc-private/ipimpl.h>
46: #include <slepcblaslapack.h>
48: typedef struct {
49: PetscBool oneside;
50: } SVD_LANCZOS;
54: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
55: {
57: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
58: PetscInt N;
61: SVDMatGetSize(svd,NULL,&N);
62: if (svd->ncv) { /* ncv set */
63: if (svd->ncv<svd->nsv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must be at least nsv");
64: } else if (svd->mpd) { /* mpd set */
65: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
66: } else { /* neither set: defaults depend on nsv being small or large */
67: if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
68: else {
69: svd->mpd = 500;
70: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
71: }
72: }
73: if (!svd->mpd) svd->mpd = svd->ncv;
74: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
75: if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
76: if (!lanczos->oneside && svd->ncv != svd->n) {
77: VecDestroyVecs(svd->n,&svd->U);
78: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
79: PetscLogObjectParents(svd,svd->ncv,svd->U);
80: }
81: DSSetType(svd->ds,DSSVD);
82: DSSetCompact(svd->ds,PETSC_TRUE);
83: DSAllocate(svd->ds,svd->ncv);
84: return(0);
85: }
89: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec *U,PetscInt k,PetscInt n,PetscScalar* work)
90: {
92: PetscInt i;
95: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
96: IPOrthogonalize(svd->ip,0,NULL,k,NULL,U,U[k],work,alpha+k,NULL);
97: VecScale(U[k],1.0/alpha[k]);
98: for (i=k+1;i<n;i++) {
99: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
100: IPOrthogonalize(svd->ip,0,NULL,i,NULL,V,V[i],work,beta+i-1,NULL);
101: VecScale(V[i],1.0/beta[i-1]);
103: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
104: IPOrthogonalize(svd->ip,0,NULL,i,NULL,U,U[i],work,alpha+i,NULL);
105: VecScale(U[i],1.0/alpha[i]);
106: }
107: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
108: IPOrthogonalize(svd->ip,0,NULL,n,NULL,V,v,work,beta+n-1,NULL);
109: return(0);
110: }
114: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
115: {
117: PetscInt i;
118: PetscReal a,b;
119: Vec temp;
122: SVDMatMult(svd,PETSC_FALSE,V[k],u);
123: for (i=k+1;i<n;i++) {
124: SVDMatMult(svd,PETSC_TRUE,u,V[i]);
125: IPNormBegin(svd->ip,u,&a);
126: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
127: IPNormEnd(svd->ip,u,&a);
128: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
130: VecScale(u,1.0/a);
131: SlepcVecMAXPBY(V[i],1.0/a,-1.0/a,i,work,V);
133: IPOrthogonalizeCGS1(svd->ip,0,NULL,i,NULL,V,V[i],work,NULL,&b);
134: VecScale(V[i],1.0/b);
136: SVDMatMult(svd,PETSC_FALSE,V[i],u_1);
137: VecAXPY(u_1,-b,u);
139: alpha[i-1] = a;
140: beta[i-1] = b;
141: temp = u;
142: u = u_1;
143: u_1 = temp;
144: }
145: SVDMatMult(svd,PETSC_TRUE,u,v);
146: IPNormBegin(svd->ip,u,&a);
147: IPMInnerProductBegin(svd->ip,v,n,V,work);
148: IPNormEnd(svd->ip,u,&a);
149: IPMInnerProductEnd(svd->ip,v,n,V,work);
151: VecScale(u,1.0/a);
152: SlepcVecMAXPBY(v,1.0/a,-1.0/a,n,work,V);
154: IPOrthogonalizeCGS1(svd->ip,0,NULL,n,NULL,V,v,work,NULL,&b);
156: alpha[n-1] = a;
157: beta[n-1] = b;
158: return(0);
159: }
163: PetscErrorCode SVDSolve_Lanczos(SVD svd)
164: {
166: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
167: PetscReal *alpha,*beta,lastbeta;
168: PetscScalar *swork,*w,*Q,*PT;
169: PetscInt i,k,j,nv,ld,off;
170: Vec v,u=0,u_1=0;
171: PetscBool conv;
174: /* allocate working space */
175: DSGetLeadingDimension(svd->ds,&ld);
176: PetscMalloc(sizeof(PetscScalar)*ld,&w);
177: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
179: VecDuplicate(svd->V[0],&v);
180: if (lanczos->oneside) {
181: SVDMatGetVecs(svd,NULL,&u);
182: SVDMatGetVecs(svd,NULL,&u_1);
183: }
185: /* normalize start vector */
186: if (!svd->nini) {
187: SlepcVecSetRandom(svd->V[0],svd->rand);
188: }
189: VecNormalize(svd->V[0],NULL);
191: while (svd->reason == SVD_CONVERGED_ITERATING) {
192: svd->its++;
194: /* inner loop */
195: nv = PetscMin(svd->nconv+svd->mpd,svd->n);
196: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
197: beta = alpha + ld;
198: if (lanczos->oneside) {
199: SVDOneSideLanczos(svd,alpha,beta,svd->V,v,u,u_1,svd->nconv,nv,swork);
200: } else {
201: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,nv,swork);
202: }
203: lastbeta = beta[nv-1];
204: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
206: /* compute SVD of bidiagonal matrix */
207: DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
208: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
209: DSSolve(svd->ds,w,NULL);
210: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
212: /* compute error estimates */
213: k = 0;
214: conv = PETSC_TRUE;
215: DSGetArray(svd->ds,DS_MAT_U,&Q);
216: for (i=svd->nconv;i<nv;i++) {
217: svd->sigma[i] = PetscRealPart(w[i]);
218: svd->errest[i] = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
219: if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
220: if (conv) {
221: if (svd->errest[i] < svd->tol) k++;
222: else conv = PETSC_FALSE;
223: }
224: }
226: /* check convergence */
227: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
228: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
230: /* compute restart vector */
231: DSGetArray(svd->ds,DS_MAT_VT,&PT);
232: if (svd->reason == SVD_CONVERGED_ITERATING) {
233: for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
234: SlepcVecMAXPBY(v,0.0,1.0,nv-svd->nconv,swork,svd->V+svd->nconv);
235: }
237: /* compute converged singular vectors */
238: off = svd->nconv+svd->nconv*ld;
239: SlepcUpdateVectors(nv-svd->nconv,svd->V+svd->nconv,0,k,PT+off,ld,PETSC_TRUE);
240: if (!lanczos->oneside) {
241: SlepcUpdateVectors(nv-svd->nconv,svd->U+svd->nconv,0,k,Q+off,ld,PETSC_FALSE);
242: }
243: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
244: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
246: /* copy restart vector from temporary space */
247: if (svd->reason == SVD_CONVERGED_ITERATING) {
248: VecCopy(v,svd->V[svd->nconv+k]);
249: }
251: svd->nconv += k;
252: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
253: }
255: /* free working space */
256: VecDestroy(&v);
257: VecDestroy(&u);
258: VecDestroy(&u_1);
259: PetscFree(w);
260: PetscFree(swork);
261: return(0);
262: }
266: PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd)
267: {
269: PetscBool set,val;
270: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
273: PetscOptionsHead("SVD Lanczos Options");
274: PetscOptionsBool("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
275: if (set) {
276: SVDLanczosSetOneSide(svd,val);
277: }
278: PetscOptionsTail();
279: return(0);
280: }
284: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
285: {
286: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
289: if (lanczos->oneside != oneside) {
290: lanczos->oneside = oneside;
291: svd->setupcalled = 0;
292: }
293: return(0);
294: }
298: /*@
299: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
300: to be used is one-sided or two-sided.
302: Logically Collective on SVD
304: Input Parameters:
305: + svd - singular value solver
306: - oneside - boolean flag indicating if the method is one-sided or not
308: Options Database Key:
309: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
311: Note:
312: By default, a two-sided variant is selected, which is sometimes slightly
313: more robust. However, the one-sided variant is faster because it avoids
314: the orthogonalization associated to left singular vectors. It also saves
315: the memory required for storing such vectors.
317: Level: advanced
319: .seealso: SVDTRLanczosSetOneSide()
320: @*/
321: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
322: {
328: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
329: return(0);
330: }
334: /*@
335: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
336: to be used is one-sided or two-sided.
338: Not Collective
340: Input Parameters:
341: . svd - singular value solver
343: Output Parameters:
344: . oneside - boolean flag indicating if the method is one-sided or not
346: Level: advanced
348: .seealso: SVDLanczosSetOneSide()
349: @*/
350: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
351: {
357: PetscTryMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
358: return(0);
359: }
363: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
364: {
365: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
368: *oneside = lanczos->oneside;
369: return(0);
370: }
374: PetscErrorCode SVDReset_Lanczos(SVD svd)
375: {
379: VecDestroyVecs(svd->n,&svd->U);
380: return(0);
381: }
385: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
386: {
390: PetscFree(svd->data);
391: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
392: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
393: return(0);
394: }
398: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
399: {
401: SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
404: PetscViewerASCIIPrintf(viewer," Lanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
405: return(0);
406: }
410: PETSC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
411: {
415: PetscNewLog(svd,SVD_LANCZOS,&svd->data);
416: svd->ops->setup = SVDSetUp_Lanczos;
417: svd->ops->solve = SVDSolve_Lanczos;
418: svd->ops->destroy = SVDDestroy_Lanczos;
419: svd->ops->reset = SVDReset_Lanczos;
420: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
421: svd->ops->view = SVDView_Lanczos;
422: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
423: PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
424: return(0);
425: }