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-2007, Universidad Politecnica de Valencia, Spain
13: This file is part of SLEPc. See the README file for conditions of use
14: and additional information.
15: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
16: */
18: #include src/svd/svdimpl.h
19: #include slepcblaslapack.h
21: typedef struct {
22: PetscTruth oneside;
23: } SVD_LANCZOS;
27: PetscErrorCode SVDSetUp_LANCZOS(SVD svd)
28: {
29: PetscErrorCode ierr;
30: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
31: PetscInt N;
32: int i;
35: SVDMatGetSize(svd,PETSC_NULL,&N);
36: if (svd->ncv == PETSC_DECIDE)
37: svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
38: if (svd->max_it == PETSC_DECIDE)
39: svd->max_it = PetscMax(N/svd->ncv,100);
40: if (svd->U) {
41: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
42: PetscFree(svd->U);
43: }
44: if (!lanczos->oneside) {
45: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
46: for (i=0;i<svd->ncv;i++) { SVDMatGetVecs(svd,PETSC_NULL,svd->U+i); }
47: }
48: return(0);
49: }
53: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec *U,int k,int n,PetscScalar* work,Vec wv,Vec wu)
54: {
56: int i;
57:
59: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
60: IPOrthogonalize(svd->ip,k,PETSC_NULL,U,U[k],work,alpha,PETSC_NULL,wu);
61: VecScale(U[k],1.0/alpha[0]);
62: for (i=k+1;i<n;i++) {
63: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
64: IPOrthogonalize(svd->ip,i,PETSC_NULL,V,V[i],work,beta+i-k-1,PETSC_NULL,wv);
65: VecScale(V[i],1.0/beta[i-k-1]);
67: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
68: IPOrthogonalize(svd->ip,i,PETSC_NULL,U,U[i],work,alpha+i-k,PETSC_NULL,wu);
69: VecScale(U[i],1.0/alpha[i-k]);
70: }
71: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
72: IPOrthogonalize(svd->ip,n,PETSC_NULL,V,v,work,beta+n-k-1,PETSC_NULL,wv);
73: return(0);
74: }
78: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec u,Vec u_1,int k,int n,PetscScalar* work,Vec wv)
79: {
81: int i,j;
82: PetscReal a,b;
83: Vec temp;
84:
86: SVDMatMult(svd,PETSC_FALSE,V[k],u);
87: for (i=k+1;i<n;i++) {
88: SVDMatMult(svd,PETSC_TRUE,u,V[i]);
89: IPNormBegin(svd->ip,u,&a);
90: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
91: IPNormEnd(svd->ip,u,&a);
92: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
93:
94: VecScale(u,1.0/a);
95: VecScale(V[i],1.0/a);
96: for (j=0;j<i;j++) work[j] = - work[j] / a;
97: VecMAXPY(V[i],i,work,V);
99: IPOrthogonalizeCGS(svd->ip,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b,wv);
100: VecScale(V[i],1.0/b);
101:
102: SVDMatMult(svd,PETSC_FALSE,V[i],u_1);
103: VecAXPY(u_1,-b,u);
105: alpha[i-k-1] = a;
106: beta[i-k-1] = b;
107: temp = u;
108: u = u_1;
109: u_1 = temp;
110: }
111: SVDMatMult(svd,PETSC_TRUE,u,v);
112: IPNormBegin(svd->ip,u,&a);
113: IPMInnerProductBegin(svd->ip,v,n,V,work);
114: IPNormEnd(svd->ip,u,&a);
115: IPMInnerProductEnd(svd->ip,v,n,V,work);
116:
117: VecScale(u,1.0/a);
118: VecScale(v,1.0/a);
119: for (j=0;j<n;j++) work[j] = - work[j] / a;
120: VecMAXPY(v,n,work,V);
122: IPOrthogonalizeCGS(svd->ip,n,PETSC_NULL,V,v,work,PETSC_NULL,&b,wv);
123:
124: alpha[n-k-1] = a;
125: beta[n-k-1] = b;
126: return(0);
127: }
131: PetscErrorCode SVDSolve_LANCZOS(SVD svd)
132: {
133: #if defined(SLEPC_MISSING_LAPACK_BDSDC)
135: SETERRQ(PETSC_ERR_SUP,"BDSDC - Lapack routine is unavailable.");
136: #else
138: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
139: PetscReal *alpha,*beta,norm,*work,*Q,*PT;
140: PetscScalar *swork;
141: PetscInt *perm;
142: int i,j,k,m,n,info,nwork=0,*iwork;
143: Vec v,u,u_1,wv,wu,*workV,*workU,*permV,*permU;
144: PetscTruth conv;
145:
147: /* allocate working space */
148: PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
149: PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
150: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&Q);
151: PetscMalloc(sizeof(PetscReal)*svd->n*svd->n,&PT);
152: PetscMalloc(sizeof(PetscReal)*(3*svd->n+4)*svd->n,&work);
153: PetscMalloc(sizeof(int)*8*svd->n,&iwork);
154: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
155: VecDuplicate(svd->V[0],&v);
156: VecDuplicate(svd->V[0],&wv);
157: PetscMalloc(sizeof(Vec)*svd->n,&workV);
158: if (lanczos->oneside) {
159: SVDMatGetVecs(svd,PETSC_NULL,&u);
160: SVDMatGetVecs(svd,PETSC_NULL,&u_1);
161: } else {
162: VecDuplicate(svd->U[0],&wu);
163: PetscMalloc(sizeof(Vec)*svd->n,&workU);
164: }
165:
166: /* normalize start vector */
167: VecCopy(svd->vec_initial,svd->V[0]);
168: VecNormalize(svd->V[0],&norm);
169:
170: while (svd->reason == SVD_CONVERGED_ITERATING) {
171: svd->its++;
173: /* inner loop */
174: if (lanczos->oneside) {
175: SVDOneSideLanczos(svd,alpha,beta,svd->V,v,u,u_1,svd->nconv,svd->n,swork,wv);
176: } else {
177: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,svd->n,swork,wv,wu);
178: }
180: /* compute SVD of bidiagonal matrix */
181: n = svd->n - svd->nconv;
182: PetscMemzero(PT,sizeof(PetscReal)*n*n);
183: PetscMemzero(Q,sizeof(PetscReal)*n*n);
184: for (i=0;i<n;i++)
185: PT[i*n+i] = Q[i*n+i] = 1.0;
186: PetscLogEventBegin(SVD_Dense,0,0,0,0);
187: LAPACKbdsdc_("U","I",&n,alpha,beta,Q,&n,PT,&n,PETSC_NULL,PETSC_NULL,work,iwork,&info,1,1);
188: PetscLogEventEnd(SVD_Dense,0,0,0,0);
190: /* compute error estimates */
191: k = 0;
192: conv = PETSC_TRUE;
193: for (i=svd->nconv;i<svd->n;i++) {
194: if (svd->which == SVD_SMALLEST) j = n-i+svd->nconv-1;
195: else j = i-svd->nconv;
196: svd->sigma[i] = alpha[j];
197: svd->errest[i] = PetscAbsScalar(Q[j*n+n-1])*beta[n-1];
198: if (alpha[j] > svd->tol) svd->errest[i] /= alpha[j];
199: if (conv) {
200: if (svd->errest[i] < svd->tol) k++;
201: else conv = PETSC_FALSE;
202: }
203: }
204:
205: /* check convergence */
206: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
207: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
208:
209: /* allocate work space for converged singular vectors */
210: if (nwork<k) {
211: for (i=nwork;i<k;i++)
212: if (lanczos->oneside) { SVDMatGetVecs(svd,workV+i,PETSC_NULL); }
213: else { SVDMatGetVecs(svd,workV+i,workU+i); }
214: nwork = k;
215: }
216:
217: /* compute converged singular vectors */
218: for (i=0;i<k;i++) {
219: if (svd->which == SVD_SMALLEST) j = n-i-1;
220: else j = i;
221: VecSet(workV[i],0.0);
222: for (m=0;m<n;m++) swork[m] = PT[m*n+j];
223: VecMAXPY(workV[i],n,swork,svd->V+svd->nconv);
224: if (!lanczos->oneside) {
225: VecSet(workU[i],0.0);
226: #if !defined(PETSC_USE_COMPLEX)
227: VecMAXPY(workU[i],n,Q+j*n,svd->U+svd->nconv);
228: #else
229: for (m=0;m<n;m++) swork[m] = Q[j*n+m];
230: VecMAXPY(workU[i],n,swork,svd->U+svd->nconv);
231: #endif
232: }
233: }
234:
235: /* compute restart vector */
236: if (svd->reason == SVD_CONVERGED_ITERATING) {
237: if (svd->which == SVD_SMALLEST) j = n-k-1;
238: else j = k;
239: VecSet(v,0.0);
240: for (m=0;m<n;m++) swork[m] = PT[m*n+j];
241: VecMAXPY(v,m,swork,svd->V+svd->nconv);
242: VecCopy(v,svd->V[svd->nconv+k]);
243: }
244:
245: /* copy converged singular vectors from temporary space */
246: for (i=0;i<k;i++) {
247: VecCopy(workV[i],svd->V[i+svd->nconv]);
248: if (!lanczos->oneside) {
249: VecCopy(workU[i],svd->U[i+svd->nconv]);
250: }
251: }
252:
253: svd->nconv += k;
254: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->n);
255: }
256:
257: /* sort singular triplets */
258: PetscMalloc(sizeof(PetscInt)*svd->nconv,&perm);
259: PetscMalloc(sizeof(Vec)*svd->nconv,&permV);
260: if (!lanczos->oneside) { PetscMalloc(sizeof(Vec)*svd->nconv,&permU); }
261: for (i=0;i<svd->nconv;i++) {
262: alpha[i] = svd->sigma[i];
263: beta[i] = svd->errest[i];
264: permV[i] = svd->V[i];
265: if (!lanczos->oneside) permU[i] = svd->U[i];
266: perm[i] = i;
267: }
268: PetscSortRealWithPermutation(svd->nconv,svd->sigma,perm);
269: for (i=0;i<svd->nconv;i++) {
270: if (svd->which == SVD_SMALLEST) j = perm[i];
271: else j = perm[svd->nconv-i-1];
272: svd->sigma[i] = alpha[j];
273: svd->errest[i] = beta[j];
274: svd->V[i] = permV[j];
275: if (!lanczos->oneside) svd->U[i] = permU[j];
276: }
277:
278: /* free working space */
279: VecDestroy(v);
280: VecDestroy(wv);
281: for (i=0;i<nwork;i++) { VecDestroy(workV[i]); }
282: PetscFree(workV);
283: if (lanczos->oneside) {
284: VecDestroy(u);
285: VecDestroy(u_1);
286: } else {
287: for (i=0;i<nwork;i++) { VecDestroy(workU[i]); }
288: PetscFree(workU);
289: PetscFree(permU);
290: VecDestroy(wu);
291: }
292: PetscFree(alpha);
293: PetscFree(beta);
294: PetscFree(Q);
295: PetscFree(PT);
296: PetscFree(work);
297: PetscFree(iwork);
298: PetscFree(swork);
299: PetscFree(perm);
300: PetscFree(permV);
301: return(0);
302: #endif
303: }
307: PetscErrorCode SVDSetFromOptions_LANCZOS(SVD svd)
308: {
310: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
313: PetscOptionsBegin(svd->comm,svd->prefix,"LANCZOS Singular Value Solver Options","SVD");
314: PetscOptionsTruth("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",PETSC_FALSE,&lanczos->oneside,PETSC_NULL);
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: PetscErrorCode SVDView_LANCZOS(SVD svd,PetscViewer viewer)
376: {
378: SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;
381: PetscViewerASCIIPrintf(viewer,"Lanczos reorthogonalization: %s\n",lanczos->oneside ? "one-side" : "two-side");
382: return(0);
383: }
388: PetscErrorCode SVDCreate_LANCZOS(SVD svd)
389: {
391: SVD_LANCZOS *lanczos;
394: PetscNew(SVD_LANCZOS,&lanczos);
395: PetscLogObjectMemory(svd,sizeof(SVD_LANCZOS));
396: svd->data = (void *)lanczos;
397: svd->ops->setup = SVDSetUp_LANCZOS;
398: svd->ops->solve = SVDSolve_LANCZOS;
399: svd->ops->destroy = SVDDestroy_Default;
400: svd->ops->setfromoptions = SVDSetFromOptions_LANCZOS;
401: svd->ops->view = SVDView_LANCZOS;
402: lanczos->oneside = PETSC_FALSE;
403: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","SVDLanczosSetOneSide_LANCZOS",SVDLanczosSetOneSide_LANCZOS);
404: return(0);
405: }