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-2010, Universidad 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
30: #include private/ipimpl.h
31: #include slepcblaslapack.h
33: typedef struct {
34: PetscTruth oneside;
35: } SVD_LANCZOS;
39: PetscErrorCode SVDSetUp_LANCZOS(SVD svd)
40: {
41: PetscErrorCode ierr;
42: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
43: PetscInt i,N,nloc;
44: PetscScalar *pU;
47: SVDMatGetSize(svd,PETSC_NULL,&N);
48: if (svd->ncv) { /* ncv set */
49: if (svd->ncv<svd->nsv) SETERRQ(1,"The value of ncv must be at least nsv");
50: }
51: else if (svd->mpd) { /* mpd set */
52: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
53: }
54: else { /* neither set: defaults depend on nsv being small or large */
55: if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
56: else { svd->mpd = 500; svd->ncv = PetscMin(N,svd->nsv+svd->mpd); }
57: }
58: if (!svd->mpd) svd->mpd = svd->ncv;
59: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(1,"The value of ncv must not be larger than nev+mpd");
60: if (!svd->max_it)
61: svd->max_it = PetscMax(N/svd->ncv,100);
62: if (svd->U) {
63: VecGetArray(svd->U[0],&pU);
64: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
65: PetscFree(pU);
66: PetscFree(svd->U);
67: }
68: if (!lanczos->oneside) {
69: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
70: SVDMatGetLocalSize(svd,&nloc,PETSC_NULL);
71: PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pU);
72: for (i=0;i<svd->ncv;i++) {
73: VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pU+i*nloc,&svd->U[i]);
74: }
75: }
76: return(0);
77: }
81: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec *U,PetscInt k,PetscInt n,PetscScalar* work)
82: {
84: PetscInt i;
85:
87: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
88: IPOrthogonalize(svd->ip,0,PETSC_NULL,k,PETSC_NULL,U,U[k],work,alpha,PETSC_NULL);
89: VecScale(U[k],1.0/alpha[0]);
90: for (i=k+1;i<n;i++) {
91: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
92: IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,beta+i-k-1,PETSC_NULL);
93: VecScale(V[i],1.0/beta[i-k-1]);
95: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
96: IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,U,U[i],work,alpha+i-k,PETSC_NULL);
97: VecScale(U[i],1.0/alpha[i-k]);
98: }
99: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
100: IPOrthogonalize(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,beta+n-k-1,PETSC_NULL);
101: return(0);
102: }
106: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
107: {
109: PetscInt i;
110: PetscReal a,b;
111: Vec temp;
112:
114: SVDMatMult(svd,PETSC_FALSE,V[k],u);
115: for (i=k+1;i<n;i++) {
116: SVDMatMult(svd,PETSC_TRUE,u,V[i]);
117: IPNormBegin(svd->ip,u,&a);
118: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
119: IPNormEnd(svd->ip,u,&a);
120: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
121:
122: VecScale(u,1.0/a);
123: SlepcVecMAXPBY(V[i],1.0/a,-1.0/a,i,work,V);
125: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
126: VecScale(V[i],1.0/b);
127:
128: SVDMatMult(svd,PETSC_FALSE,V[i],u_1);
129: VecAXPY(u_1,-b,u);
131: alpha[i-k-1] = a;
132: beta[i-k-1] = b;
133: temp = u;
134: u = u_1;
135: u_1 = temp;
136: }
137: SVDMatMult(svd,PETSC_TRUE,u,v);
138: IPNormBegin(svd->ip,u,&a);
139: IPMInnerProductBegin(svd->ip,v,n,V,work);
140: IPNormEnd(svd->ip,u,&a);
141: IPMInnerProductEnd(svd->ip,v,n,V,work);
142:
143: VecScale(u,1.0/a);
144: SlepcVecMAXPBY(v,1.0/a,-1.0/a,n,work,V);
146: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
147:
148: alpha[n-k-1] = a;
149: beta[n-k-1] = b;
150: return(0);
151: }
155: PetscErrorCode SVDSolve_LANCZOS(SVD svd)
156: {
157: #if defined(SLEPC_MISSING_LAPACK_BDSDC)
159: SETERRQ(PETSC_ERR_SUP,"BDSDC - Lapack routine is unavailable.");
160: #else
162: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
163: PetscReal *alpha,*beta,norm,*work,*Q,*PT;
164: PetscScalar *swork;
165: PetscBLASInt n,info,*iwork;
166: PetscInt i,j,k,m,nv;
167: Vec v,u,u_1;
168: PetscTruth conv;
169:
171: /* allocate working space */
172: PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
173: PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
174: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&Q);
175: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&PT);
176: PetscMalloc(sizeof(PetscReal)*(3*svd->n+4)*svd->n,&work);
177: PetscMalloc(sizeof(PetscBLASInt)*8*svd->n,&iwork);
178: #if !defined(PETSC_USE_COMPLEX)
179: if (svd->which == SVD_SMALLEST) {
180: #endif
181: PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&swork);
182: #if !defined(PETSC_USE_COMPLEX)
183: } else {
184: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
185: }
186: #endif
188: VecDuplicate(svd->V[0],&v);
189: if (lanczos->oneside) {
190: SVDMatGetVecs(svd,PETSC_NULL,&u);
191: SVDMatGetVecs(svd,PETSC_NULL,&u_1);
192: }
193:
194: /* normalize start vector */
195: if (svd->nini==0) {
196: SlepcVecSetRandom(svd->V[0],svd->rand);
197: }
198: VecNormalize(svd->V[0],&norm);
199:
200: while (svd->reason == SVD_CONVERGED_ITERATING) {
201: svd->its++;
203: /* inner loop */
204: nv = PetscMin(svd->nconv+svd->mpd,svd->n);
205: if (lanczos->oneside) {
206: SVDOneSideLanczos(svd,alpha,beta,svd->V,v,u,u_1,svd->nconv,nv,swork);
207: } else {
208: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,nv,swork);
209: }
211: /* compute SVD of bidiagonal matrix */
212: n = nv - svd->nconv;
213: PetscMemzero(PT,sizeof(PetscReal)*n*n);
214: PetscMemzero(Q,sizeof(PetscReal)*n*n);
215: for (i=0;i<n;i++)
216: PT[i*n+i] = Q[i*n+i] = 1.0;
217: PetscLogEventBegin(SVD_Dense,0,0,0,0);
218: LAPACKbdsdc_("U","I",&n,alpha,beta,Q,&n,PT,&n,PETSC_NULL,PETSC_NULL,work,iwork,&info);
219: PetscLogEventEnd(SVD_Dense,0,0,0,0);
221: /* compute error estimates */
222: k = 0;
223: conv = PETSC_TRUE;
224: for (i=svd->nconv;i<nv;i++) {
225: if (svd->which == SVD_SMALLEST) j = n-i+svd->nconv-1;
226: else j = i-svd->nconv;
227: svd->sigma[i] = alpha[j];
228: svd->errest[i] = PetscAbsScalar(Q[j*n+n-1])*beta[n-1];
229: if (alpha[j] > svd->tol) svd->errest[i] /= alpha[j];
230: if (conv) {
231: if (svd->errest[i] < svd->tol) k++;
232: else conv = PETSC_FALSE;
233: }
234: }
235:
236: /* check convergence */
237: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
238: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
239:
240: /* compute restart vector */
241: if (svd->reason == SVD_CONVERGED_ITERATING) {
242: if (svd->which == SVD_SMALLEST) j = n-k-1;
243: else j = k;
244: for (m=0;m<n;m++) swork[m] = PT[m*n+j];
245: SlepcVecMAXPBY(v,0.0,1.0,n,swork,svd->V+svd->nconv);
246: }
247:
248: /* compute converged singular vectors */
249: #if !defined(PETSC_USE_COMPLEX)
250: if (svd->which == SVD_SMALLEST) {
251: #endif
252: for (i=0;i<k;i++) {
253: if (svd->which == SVD_SMALLEST) j = n-i-1;
254: else j = i;
255: for (m=0;m<n;m++) swork[i*n+m] = PT[m*n+j];
256: }
257: SlepcUpdateVectors(n,svd->V+svd->nconv,0,k,swork,n,PETSC_FALSE);
258: if (!lanczos->oneside) {
259: for (i=0;i<k;i++) {
260: if (svd->which == SVD_SMALLEST) j = n-i-1;
261: else j = i;
262: for (m=0;m<n;m++) swork[i*n+m] = Q[j*n+m];
263: }
264: SlepcUpdateVectors(n,svd->U+svd->nconv,0,k,swork,n,PETSC_FALSE);
265: }
266: #if !defined(PETSC_USE_COMPLEX)
267: } else {
268: SlepcUpdateVectors(n,svd->V+svd->nconv,0,k,PT,n,PETSC_TRUE);
269: if (!lanczos->oneside) {
270: SlepcUpdateVectors(n,svd->U+svd->nconv,0,k,Q,n,PETSC_FALSE);
271: }
272: }
273: #endif
274:
275: /* copy restart vector from temporary space */
276: if (svd->reason == SVD_CONVERGED_ITERATING) {
277: VecCopy(v,svd->V[svd->nconv+k]);
278: }
279:
280: svd->nconv += k;
281: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
282: }
283:
284: /* free working space */
285: VecDestroy(v);
286: if (lanczos->oneside) {
287: VecDestroy(u);
288: VecDestroy(u_1);
289: }
290: PetscFree(alpha);
291: PetscFree(beta);
292: PetscFree(Q);
293: PetscFree(PT);
294: PetscFree(work);
295: PetscFree(iwork);
296: PetscFree(swork);
297: return(0);
298: #endif
299: }
303: PetscErrorCode SVDSetFromOptions_LANCZOS(SVD svd)
304: {
306: PetscTruth set,val;
307: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
310: PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"LANCZOS Singular Value Solver Options","SVD");
311: PetscOptionsTruth("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
312: if (set) {
313: SVDLanczosSetOneSide(svd,val);
314: }
315: PetscOptionsEnd();
316: return(0);
317: }
322: PetscErrorCode SVDLanczosSetOneSide_LANCZOS(SVD svd,PetscTruth oneside)
323: {
324: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
327: if (lanczos->oneside != oneside) {
328: lanczos->oneside = oneside;
329: svd->setupcalled = 0;
330: }
331: return(0);
332: }
337: /*@
338: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
339: to be used is one-sided or two-sided.
341: Collective on SVD
343: Input Parameters:
344: + svd - singular value solver
345: - oneside - boolean flag indicating if the method is one-sided or not
347: Options Database Key:
348: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
350: Note:
351: By default, a two-sided variant is selected, which is sometimes slightly
352: more robust. However, the one-sided variant is faster because it avoids
353: the orthogonalization associated to left singular vectors. It also saves
354: the memory required for storing such vectors.
356: Level: advanced
358: .seealso: SVDTRLanczosSetOneSide()
359: @*/
360: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscTruth oneside)
361: {
362: PetscErrorCode ierr, (*f)(SVD,PetscTruth);
366: PetscObjectQueryFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",(void (**)())&f);
367: if (f) {
368: (*f)(svd,oneside);
369: }
370: return(0);
371: }
375: /*@
376: SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
377: to be used is one-sided or two-sided.
379: Collective on SVD
381: Input Parameters:
382: . svd - singular value solver
384: Output Parameters:
385: . oneside - boolean flag indicating if the method is one-sided or not
387: Level: advanced
389: .seealso: SVDLanczosSetOneSide()
390: @*/
391: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscTruth *oneside)
392: {
393: PetscErrorCode ierr, (*f)(SVD,PetscTruth*);
397: PetscObjectQueryFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",(void (**)())&f);
398: if (f) {
399: (*f)(svd,oneside);
400: }
401: return(0);
402: }
407: PetscErrorCode SVDLanczosGetOneSide_LANCZOS(SVD svd,PetscTruth *oneside)
408: {
409: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
413: *oneside = lanczos->oneside;
414: return(0);
415: }
420: PetscErrorCode SVDDestroy_LANCZOS(SVD svd)
421: {
426: SVDDestroy_Default(svd);
427: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","",PETSC_NULL);
428: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosGetOneSide_C","",PETSC_NULL);
429: return(0);
430: }
434: PetscErrorCode SVDView_LANCZOS(SVD svd,PetscViewer viewer)
435: {
437: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
440: PetscViewerASCIIPrintf(viewer,"Lanczos reorthogonalization: %s\n",lanczos->oneside ? "one-side" : "two-side");
441: return(0);
442: }
447: PetscErrorCode SVDCreate_LANCZOS(SVD svd)
448: {
450: SVD_LANCZOS *lanczos;
453: PetscNew(SVD_LANCZOS,&lanczos);
454: PetscLogObjectMemory(svd,sizeof(SVD_LANCZOS));
455: svd->data = (void *)lanczos;
456: svd->ops->setup = SVDSetUp_LANCZOS;
457: svd->ops->solve = SVDSolve_LANCZOS;
458: svd->ops->destroy = SVDDestroy_LANCZOS;
459: svd->ops->setfromoptions = SVDSetFromOptions_LANCZOS;
460: svd->ops->view = SVDView_LANCZOS;
461: lanczos->oneside = PETSC_FALSE;
462: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","SVDLanczosSetOneSide_LANCZOS",SVDLanczosSetOneSide_LANCZOS);
463: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosGetOneSide_C","SVDLanczosGetOneSide_LANCZOS",SVDLanczosGetOneSide_LANCZOS);
464: return(0);
465: }