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-2009, 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 slepcblaslapack.h
32: typedef struct {
33: PetscTruth oneside;
34: } SVD_LANCZOS;
38: PetscErrorCode SVDSetUp_LANCZOS(SVD svd)
39: {
40: PetscErrorCode ierr;
41: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
42: PetscInt i,N,nloc;
43: PetscScalar *pU;
46: SVDMatGetSize(svd,PETSC_NULL,&N);
47: if (svd->ncv) { /* ncv set */
48: if (svd->ncv<svd->nsv) SETERRQ(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(1,"The value of ncv must not be larger than nev+mpd");
59: if (!svd->max_it)
60: svd->max_it = PetscMax(N/svd->ncv,100);
61: if (svd->U) {
62: VecGetArray(svd->U[0],&pU);
63: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
64: PetscFree(pU);
65: PetscFree(svd->U);
66: }
67: if (!lanczos->oneside) {
68: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
69: SVDMatGetLocalSize(svd,&nloc,PETSC_NULL);
70: PetscMalloc(svd->ncv*nloc*sizeof(PetscScalar),&pU);
71: for (i=0;i<svd->ncv;i++) {
72: VecCreateMPIWithArray(((PetscObject)svd)->comm,nloc,PETSC_DECIDE,pU+i*nloc,&svd->U[i]);
73: }
74: }
75: return(0);
76: }
80: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec *U,PetscInt k,PetscInt n,PetscScalar* work,Vec wv,Vec wu)
81: {
83: PetscInt i;
84:
86: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
87: IPOrthogonalize(svd->ip,k,PETSC_NULL,U,U[k],work,alpha,PETSC_NULL,wu,PETSC_NULL);
88: VecScale(U[k],1.0/alpha[0]);
89: for (i=k+1;i<n;i++) {
90: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
91: IPOrthogonalize(svd->ip,i,PETSC_NULL,V,V[i],work,beta+i-k-1,PETSC_NULL,wv,PETSC_NULL);
92: VecScale(V[i],1.0/beta[i-k-1]);
94: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
95: IPOrthogonalize(svd->ip,i,PETSC_NULL,U,U[i],work,alpha+i-k,PETSC_NULL,wu,PETSC_NULL);
96: VecScale(U[i],1.0/alpha[i-k]);
97: }
98: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
99: IPOrthogonalize(svd->ip,n,PETSC_NULL,V,v,work,beta+n-k-1,PETSC_NULL,wv,PETSC_NULL);
100: return(0);
101: }
105: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work,Vec wv)
106: {
108: PetscInt i,j;
109: PetscReal a,b;
110: Vec temp;
111:
113: SVDMatMult(svd,PETSC_FALSE,V[k],u);
114: for (i=k+1;i<n;i++) {
115: SVDMatMult(svd,PETSC_TRUE,u,V[i]);
116: IPNormBegin(svd->ip,u,&a);
117: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
118: IPNormEnd(svd->ip,u,&a);
119: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
120:
121: VecScale(u,1.0/a);
122: VecScale(V[i],1.0/a);
123: for (j=0;j<i;j++) work[j] = - work[j] / a;
124: VecMAXPY(V[i],i,work,V);
126: IPOrthogonalizeCGS(svd->ip,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b,wv);
127: VecScale(V[i],1.0/b);
128:
129: SVDMatMult(svd,PETSC_FALSE,V[i],u_1);
130: VecAXPY(u_1,-b,u);
132: alpha[i-k-1] = a;
133: beta[i-k-1] = b;
134: temp = u;
135: u = u_1;
136: u_1 = temp;
137: }
138: SVDMatMult(svd,PETSC_TRUE,u,v);
139: IPNormBegin(svd->ip,u,&a);
140: IPMInnerProductBegin(svd->ip,v,n,V,work);
141: IPNormEnd(svd->ip,u,&a);
142: IPMInnerProductEnd(svd->ip,v,n,V,work);
143:
144: VecScale(u,1.0/a);
145: VecScale(v,1.0/a);
146: for (j=0;j<n;j++) work[j] = - work[j] / a;
147: VecMAXPY(v,n,work,V);
149: IPOrthogonalizeCGS(svd->ip,n,PETSC_NULL,V,v,work,PETSC_NULL,&b,wv);
150:
151: alpha[n-k-1] = a;
152: beta[n-k-1] = b;
153: return(0);
154: }
158: PetscErrorCode SVDSolve_LANCZOS(SVD svd)
159: {
160: #if defined(SLEPC_MISSING_LAPACK_BDSDC)
162: SETERRQ(PETSC_ERR_SUP,"BDSDC - Lapack routine is unavailable.");
163: #else
165: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
166: PetscReal *alpha,*beta,norm,*work,*Q,*PT;
167: PetscScalar *swork;
168: PetscBLASInt n,info,*iwork;
169: PetscInt i,j,k,m,nv;
170: Vec v,u,u_1,wv,wu;
171: PetscTruth conv;
172:
174: /* allocate working space */
175: PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
176: PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
177: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&Q);
178: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&PT);
179: PetscMalloc(sizeof(PetscReal)*(3*svd->n+4)*svd->n,&work);
180: PetscMalloc(sizeof(PetscBLASInt)*8*svd->n,&iwork);
181: #if !defined(PETSC_USE_COMPLEX)
182: if (svd->which == SVD_SMALLEST) {
183: #endif
184: PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&swork);
185: #if !defined(PETSC_USE_COMPLEX)
186: } else {
187: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
188: }
189: #endif
191: VecDuplicate(svd->V[0],&v);
192: VecDuplicate(svd->V[0],&wv);
193: if (lanczos->oneside) {
194: SVDMatGetVecs(svd,PETSC_NULL,&u);
195: SVDMatGetVecs(svd,PETSC_NULL,&u_1);
196: } else {
197: VecDuplicate(svd->U[0],&wu);
198: }
199:
200: /* normalize start vector */
201: VecCopy(svd->vec_initial,svd->V[0]);
202: VecNormalize(svd->V[0],&norm);
203:
204: while (svd->reason == SVD_CONVERGED_ITERATING) {
205: svd->its++;
207: /* inner loop */
208: nv = PetscMin(svd->nconv+svd->mpd,svd->n);
209: if (lanczos->oneside) {
210: SVDOneSideLanczos(svd,alpha,beta,svd->V,v,u,u_1,svd->nconv,nv,swork,wv);
211: } else {
212: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,nv,swork,wv,wu);
213: }
215: /* compute SVD of bidiagonal matrix */
216: n = nv - svd->nconv;
217: PetscMemzero(PT,sizeof(PetscReal)*n*n);
218: PetscMemzero(Q,sizeof(PetscReal)*n*n);
219: for (i=0;i<n;i++)
220: PT[i*n+i] = Q[i*n+i] = 1.0;
221: PetscLogEventBegin(SVD_Dense,0,0,0,0);
222: LAPACKbdsdc_("U","I",&n,alpha,beta,Q,&n,PT,&n,PETSC_NULL,PETSC_NULL,work,iwork,&info);
223: PetscLogEventEnd(SVD_Dense,0,0,0,0);
225: /* compute error estimates */
226: k = 0;
227: conv = PETSC_TRUE;
228: for (i=svd->nconv;i<nv;i++) {
229: if (svd->which == SVD_SMALLEST) j = n-i+svd->nconv-1;
230: else j = i-svd->nconv;
231: svd->sigma[i] = alpha[j];
232: svd->errest[i] = PetscAbsScalar(Q[j*n+n-1])*beta[n-1];
233: if (alpha[j] > svd->tol) svd->errest[i] /= alpha[j];
234: if (conv) {
235: if (svd->errest[i] < svd->tol) k++;
236: else conv = PETSC_FALSE;
237: }
238: }
239:
240: /* check convergence */
241: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
242: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
243:
244: /* compute restart vector */
245: if (svd->reason == SVD_CONVERGED_ITERATING) {
246: if (svd->which == SVD_SMALLEST) j = n-k-1;
247: else j = k;
248: VecSet(v,0.0);
249: for (m=0;m<n;m++) swork[m] = PT[m*n+j];
250: VecMAXPY(v,n,swork,svd->V+svd->nconv);
251: }
252:
253: /* compute converged singular vectors */
254: #if !defined(PETSC_USE_COMPLEX)
255: if (svd->which == SVD_SMALLEST) {
256: #endif
257: for (i=0;i<k;i++) {
258: if (svd->which == SVD_SMALLEST) j = n-i-1;
259: else j = i;
260: for (m=0;m<n;m++) swork[i*n+m] = PT[m*n+j];
261: }
262: SlepcUpdateVectors(n,svd->V+svd->nconv,0,k,swork,n,PETSC_FALSE);
263: if (!lanczos->oneside) {
264: for (i=0;i<k;i++) {
265: if (svd->which == SVD_SMALLEST) j = n-i-1;
266: else j = i;
267: for (m=0;m<n;m++) swork[i*n+m] = Q[j*n+m];
268: }
269: SlepcUpdateVectors(n,svd->U+svd->nconv,0,k,swork,n,PETSC_FALSE);
270: }
271: #if !defined(PETSC_USE_COMPLEX)
272: } else {
273: SlepcUpdateVectors(n,svd->V+svd->nconv,0,k,PT,n,PETSC_TRUE);
274: if (!lanczos->oneside) {
275: SlepcUpdateVectors(n,svd->U+svd->nconv,0,k,Q,n,PETSC_FALSE);
276: }
277: }
278: #endif
279:
280: /* copy restart vector from temporary space */
281: if (svd->reason == SVD_CONVERGED_ITERATING) {
282: VecCopy(v,svd->V[svd->nconv+k]);
283: }
284:
285: svd->nconv += k;
286: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
287: }
288:
289: /* free working space */
290: VecDestroy(v);
291: VecDestroy(wv);
292: if (lanczos->oneside) {
293: VecDestroy(u);
294: VecDestroy(u_1);
295: } else {
296: VecDestroy(wu);
297: }
298: PetscFree(alpha);
299: PetscFree(beta);
300: PetscFree(Q);
301: PetscFree(PT);
302: PetscFree(work);
303: PetscFree(iwork);
304: PetscFree(swork);
305: return(0);
306: #endif
307: }
311: PetscErrorCode SVDSetFromOptions_LANCZOS(SVD svd)
312: {
314: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
317: PetscOptionsBegin(((PetscObject)svd)->comm,((PetscObject)svd)->prefix,"LANCZOS Singular Value Solver Options","SVD");
318: PetscOptionsTruth("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",PETSC_FALSE,&lanczos->oneside,PETSC_NULL);
319: PetscOptionsEnd();
320: return(0);
321: }
326: PetscErrorCode SVDLanczosSetOneSide_LANCZOS(SVD svd,PetscTruth oneside)
327: {
328: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
331: if (lanczos->oneside != oneside) {
332: lanczos->oneside = oneside;
333: svd->setupcalled = 0;
334: }
335: return(0);
336: }
341: /*@
342: SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
343: to be used is one-sided or two-sided.
345: Collective on SVD
347: Input Parameters:
348: + svd - singular value solver
349: - oneside - boolean flag indicating if the method is one-sided or not
351: Options Database Key:
352: . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
354: Note:
355: By default, a two-sided variant is selected, which is sometimes slightly
356: more robust. However, the one-sided variant is faster because it avoids
357: the orthogonalization associated to left singular vectors. It also saves
358: the memory required for storing such vectors.
360: Level: advanced
362: .seealso: SVDTRLanczosSetOneSide()
363: @*/
364: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscTruth oneside)
365: {
366: PetscErrorCode ierr, (*f)(SVD,PetscTruth);
370: PetscObjectQueryFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",(void (**)())&f);
371: if (f) {
372: (*f)(svd,oneside);
373: }
374: return(0);
375: }
379: PetscErrorCode SVDView_LANCZOS(SVD svd,PetscViewer viewer)
380: {
382: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
385: PetscViewerASCIIPrintf(viewer,"Lanczos reorthogonalization: %s\n",lanczos->oneside ? "one-side" : "two-side");
386: return(0);
387: }
392: PetscErrorCode SVDCreate_LANCZOS(SVD svd)
393: {
395: SVD_LANCZOS *lanczos;
398: PetscNew(SVD_LANCZOS,&lanczos);
399: PetscLogObjectMemory(svd,sizeof(SVD_LANCZOS));
400: svd->data = (void *)lanczos;
401: svd->ops->setup = SVDSetUp_LANCZOS;
402: svd->ops->solve = SVDSolve_LANCZOS;
403: svd->ops->destroy = SVDDestroy_Default;
404: svd->ops->setfromoptions = SVDSetFromOptions_LANCZOS;
405: svd->ops->view = SVDView_LANCZOS;
406: lanczos->oneside = PETSC_FALSE;
407: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","SVDLanczosSetOneSide_LANCZOS",SVDLanczosSetOneSide_LANCZOS);
408: return(0);
409: }