Actual source code: gklanczos.c
1: /*
3: SLEPc singular value solver: "lanczos"
5: Method: Golub-Kahan-Lanczos bidiagonalization
7: Last update: Jun 2007
9: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
10: SLEPc - Scalable Library for Eigenvalue Problem Computations
11: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
13: This file is part of SLEPc.
14:
15: SLEPc is free software: you can redistribute it and/or modify it under the
16: terms of version 3 of the GNU Lesser General Public License as published by
17: the Free Software Foundation.
19: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
20: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
21: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
22: more details.
24: You should have received a copy of the GNU Lesser General Public License
25: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
26: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27: */
29: #include <private/svdimpl.h> /*I "slepcsvd.h" I*/
30: #include <private/ipimpl.h> /*I "slepcip.h" I*/
31: #include <slepcblaslapack.h>
33: typedef struct {
34: PetscBool oneside;
35: } SVD_LANCZOS;
39: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
40: {
42: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
43: PetscInt N;
46: SVDMatGetSize(svd,PETSC_NULL,&N);
47: if (svd->ncv) { /* ncv set */
48: if (svd->ncv<svd->nsv) SETERRQ(((PetscObject)svd)->comm,1,"The value of ncv must be at least nsv");
49: }
50: else if (svd->mpd) { /* mpd set */
51: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
52: }
53: else { /* neither set: defaults depend on nsv being small or large */
54: if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
55: else { svd->mpd = 500; svd->ncv = PetscMin(N,svd->nsv+svd->mpd); }
56: }
57: if (!svd->mpd) svd->mpd = svd->ncv;
58: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(((PetscObject)svd)->comm,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: if (!lanczos->oneside && svd->ncv != svd->n) {
61: if (svd->n) { VecDestroyVecs(svd->n,&svd->U); }
62: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
63: }
64: return(0);
65: }
69: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec *U,PetscInt k,PetscInt n,PetscScalar* work)
70: {
72: PetscInt i;
73:
75: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
76: IPOrthogonalize(svd->ip,0,PETSC_NULL,k,PETSC_NULL,U,U[k],work,alpha,PETSC_NULL);
77: VecScale(U[k],1.0/alpha[0]);
78: for (i=k+1;i<n;i++) {
79: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
80: IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,beta+i-k-1,PETSC_NULL);
81: VecScale(V[i],1.0/beta[i-k-1]);
83: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
84: IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,U,U[i],work,alpha+i-k,PETSC_NULL);
85: VecScale(U[i],1.0/alpha[i-k]);
86: }
87: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
88: IPOrthogonalize(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,beta+n-k-1,PETSC_NULL);
89: return(0);
90: }
94: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
95: {
97: PetscInt i;
98: PetscReal a,b;
99: Vec temp;
100:
102: SVDMatMult(svd,PETSC_FALSE,V[k],u);
103: for (i=k+1;i<n;i++) {
104: SVDMatMult(svd,PETSC_TRUE,u,V[i]);
105: IPNormBegin(svd->ip,u,&a);
106: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
107: IPNormEnd(svd->ip,u,&a);
108: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
109:
110: VecScale(u,1.0/a);
111: SlepcVecMAXPBY(V[i],1.0/a,-1.0/a,i,work,V);
113: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
114: VecScale(V[i],1.0/b);
115:
116: SVDMatMult(svd,PETSC_FALSE,V[i],u_1);
117: VecAXPY(u_1,-b,u);
119: alpha[i-k-1] = a;
120: beta[i-k-1] = b;
121: temp = u;
122: u = u_1;
123: u_1 = temp;
124: }
125: SVDMatMult(svd,PETSC_TRUE,u,v);
126: IPNormBegin(svd->ip,u,&a);
127: IPMInnerProductBegin(svd->ip,v,n,V,work);
128: IPNormEnd(svd->ip,u,&a);
129: IPMInnerProductEnd(svd->ip,v,n,V,work);
130:
131: VecScale(u,1.0/a);
132: SlepcVecMAXPBY(v,1.0/a,-1.0/a,n,work,V);
134: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
135:
136: alpha[n-k-1] = a;
137: beta[n-k-1] = b;
138: return(0);
139: }
143: PetscErrorCode SVDSolve_Lanczos(SVD svd)
144: {
145: #if defined(SLEPC_MISSING_LAPACK_BDSDC)
147: SETERRQ(((PetscObject)svd)->comm,PETSC_ERR_SUP,"BDSDC - Lapack routine is unavailable.");
148: #else
150: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
151: PetscReal *alpha,*beta,norm,*work,*Q,*PT;
152: PetscScalar *swork;
153: PetscBLASInt n,info,*iwork;
154: PetscInt i,j,k,m,nv;
155: Vec v,u=0,u_1=0;
156: PetscBool conv;
157:
159: /* allocate working space */
160: PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
161: PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
162: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&Q);
163: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&PT);
164: PetscMalloc(sizeof(PetscReal)*(3*svd->n+4)*svd->n,&work);
165: PetscMalloc(sizeof(PetscBLASInt)*8*svd->n,&iwork);
166: #if !defined(PETSC_USE_COMPLEX)
167: if (svd->which == SVD_SMALLEST) {
168: #endif
169: PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&swork);
170: #if !defined(PETSC_USE_COMPLEX)
171: } else {
172: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
173: }
174: #endif
176: VecDuplicate(svd->V[0],&v);
177: if (lanczos->oneside) {
178: SVDMatGetVecs(svd,PETSC_NULL,&u);
179: SVDMatGetVecs(svd,PETSC_NULL,&u_1);
180: }
181:
182: /* normalize start vector */
183: if (svd->nini==0) {
184: SlepcVecSetRandom(svd->V[0],svd->rand);
185: }
186: VecNormalize(svd->V[0],&norm);
187:
188: while (svd->reason == SVD_CONVERGED_ITERATING) {
189: svd->its++;
191: /* inner loop */
192: nv = PetscMin(svd->nconv+svd->mpd,svd->n);
193: if (lanczos->oneside) {
194: SVDOneSideLanczos(svd,alpha,beta,svd->V,v,u,u_1,svd->nconv,nv,swork);
195: } else {
196: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,nv,swork);
197: }
199: /* compute SVD of bidiagonal matrix */
200: n = nv - svd->nconv;
201: PetscMemzero(PT,sizeof(PetscReal)*n*n);
202: PetscMemzero(Q,sizeof(PetscReal)*n*n);
203: for (i=0;i<n;i++)
204: PT[i*n+i] = Q[i*n+i] = 1.0;
205: PetscLogEventBegin(SVD_Dense,0,0,0,0);
206: LAPACKbdsdc_("U","I",&n,alpha,beta,Q,&n,PT,&n,PETSC_NULL,PETSC_NULL,work,iwork,&info);
207: PetscLogEventEnd(SVD_Dense,0,0,0,0);
209: /* compute error estimates */
210: k = 0;
211: conv = PETSC_TRUE;
212: for (i=svd->nconv;i<nv;i++) {
213: if (svd->which == SVD_SMALLEST) j = n-i+svd->nconv-1;
214: else j = i-svd->nconv;
215: svd->sigma[i] = alpha[j];
216: svd->errest[i] = PetscAbsScalar(Q[j*n+n-1])*beta[n-1];
217: if (alpha[j] > svd->tol) svd->errest[i] /= alpha[j];
218: if (conv) {
219: if (svd->errest[i] < svd->tol) k++;
220: else conv = PETSC_FALSE;
221: }
222: }
223:
224: /* check convergence */
225: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
226: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
227:
228: /* compute restart vector */
229: if (svd->reason == SVD_CONVERGED_ITERATING) {
230: if (svd->which == SVD_SMALLEST) j = n-k-1;
231: else j = k;
232: for (m=0;m<n;m++) swork[m] = PT[m*n+j];
233: SlepcVecMAXPBY(v,0.0,1.0,n,swork,svd->V+svd->nconv);
234: }
235:
236: /* compute converged singular vectors */
237: #if !defined(PETSC_USE_COMPLEX)
238: if (svd->which == SVD_SMALLEST) {
239: #endif
240: for (i=0;i<k;i++) {
241: if (svd->which == SVD_SMALLEST) j = n-i-1;
242: else j = i;
243: for (m=0;m<n;m++) swork[i*n+m] = PT[m*n+j];
244: }
245: SlepcUpdateVectors(n,svd->V+svd->nconv,0,k,swork,n,PETSC_FALSE);
246: if (!lanczos->oneside) {
247: for (i=0;i<k;i++) {
248: if (svd->which == SVD_SMALLEST) j = n-i-1;
249: else j = i;
250: for (m=0;m<n;m++) swork[i*n+m] = Q[j*n+m];
251: }
252: SlepcUpdateVectors(n,svd->U+svd->nconv,0,k,swork,n,PETSC_FALSE);
253: }
254: #if !defined(PETSC_USE_COMPLEX)
255: } else {
256: SlepcUpdateVectors(n,svd->V+svd->nconv,0,k,PT,n,PETSC_TRUE);
257: if (!lanczos->oneside) {
258: SlepcUpdateVectors(n,svd->U+svd->nconv,0,k,Q,n,PETSC_FALSE);
259: }
260: }
261: #endif
262:
263: /* copy restart vector from temporary space */
264: if (svd->reason == SVD_CONVERGED_ITERATING) {
265: VecCopy(v,svd->V[svd->nconv+k]);
266: }
267:
268: svd->nconv += k;
269: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
270: }
271:
272: /* free working space */
273: VecDestroy(&v);
274: VecDestroy(&u);
275: VecDestroy(&u_1);
276: PetscFree(alpha);
277: PetscFree(beta);
278: PetscFree(Q);
279: PetscFree(PT);
280: PetscFree(work);
281: PetscFree(iwork);
282: PetscFree(swork);
283: return(0);
284: #endif
285: }
289: PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd)
290: {
292: PetscBool set,val;
293: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
296: PetscOptionsHead("SVD Lanczos Options");
297: PetscOptionsBool("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
298: if (set) {
299: SVDLanczosSetOneSide(svd,val);
300: }
301: PetscOptionsTail();
302: return(0);
303: }
308: PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
309: {
310: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
313: if (lanczos->oneside != oneside) {
314: lanczos->oneside = oneside;
315: svd->setupcalled = 0;
316: }
317: return(0);
318: }
323: /*@
324: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
325: to be used is one-sided or two-sided.
327: Logically Collective on SVD
329: Input Parameters:
330: + svd - singular value solver
331: - oneside - boolean flag indicating if the method is one-sided or not
333: Options Database Key:
334: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
336: Note:
337: By default, a two-sided variant is selected, which is sometimes slightly
338: more robust. However, the one-sided variant is faster because it avoids
339: the orthogonalization associated to left singular vectors. It also saves
340: the memory required for storing such vectors.
342: Level: advanced
344: .seealso: SVDTRLanczosSetOneSide()
345: @*/
346: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
347: {
353: PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
354: return(0);
355: }
359: /*@
360: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
361: to be used is one-sided or two-sided.
363: Not Collective
365: Input Parameters:
366: . svd - singular value solver
368: Output Parameters:
369: . oneside - boolean flag indicating if the method is one-sided or not
371: Level: advanced
373: .seealso: SVDLanczosSetOneSide()
374: @*/
375: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
376: {
382: PetscTryMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
383: return(0);
384: }
389: PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
390: {
391: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
394: *oneside = lanczos->oneside;
395: return(0);
396: }
401: PetscErrorCode SVDReset_Lanczos(SVD svd)
402: {
406: VecDestroyVecs(svd->n,&svd->U);
407: return(0);
408: }
412: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
413: {
417: PetscFree(svd->data);
418: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","",PETSC_NULL);
419: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosGetOneSide_C","",PETSC_NULL);
420: return(0);
421: }
425: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
426: {
428: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
431: PetscViewerASCIIPrintf(viewer," Lanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
432: return(0);
433: }
438: PetscErrorCode SVDCreate_Lanczos(SVD svd)
439: {
443: PetscNewLog(svd,SVD_LANCZOS,&svd->data);
444: svd->ops->setup = SVDSetUp_Lanczos;
445: svd->ops->solve = SVDSolve_Lanczos;
446: svd->ops->destroy = SVDDestroy_Lanczos;
447: svd->ops->reset = SVDReset_Lanczos;
448: svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
449: svd->ops->view = SVDView_Lanczos;
450: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","SVDLanczosSetOneSide_Lanczos",SVDLanczosSetOneSide_Lanczos);
451: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosGetOneSide_C","SVDLanczosGetOneSide_Lanczos",SVDLanczosGetOneSide_Lanczos);
452: if (!svd->ip) { SVDGetIP(svd,&svd->ip); }
453: return(0);
454: }