Actual source code: trlanczos.c
1: /*
3: SLEPc singular value solver: "trlanczos"
5: Method: Golub-Kahan-Lanczos bidiagonalization with thick-restart
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_TRLANCZOS;
27: PetscErrorCode SVDSetUp_TRLANCZOS(SVD svd)
28: {
29: PetscErrorCode ierr;
30: PetscInt N;
31: int i;
34: SVDMatGetSize(svd,PETSC_NULL,&N);
35: if (svd->ncv == PETSC_DECIDE)
36: svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
37: if (svd->max_it == PETSC_DECIDE)
38: svd->max_it = PetscMax(N/svd->ncv,100);
39: if (svd->ncv!=svd->n) {
40: if (svd->U) {
41: for (i=0;i<svd->n;i++) { VecDestroy(svd->U[i]); }
42: PetscFree(svd->U);
43: }
44: PetscMalloc(sizeof(Vec)*svd->ncv,&svd->U);
45: for (i=0;i<svd->ncv;i++) { SVDMatGetVecs(svd,PETSC_NULL,svd->U+i); }
46: }
47: return(0);
48: }
52: static PetscErrorCode SVDOneSideTRLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,PetscScalar* bb,Vec *V,Vec v,Vec* U,int nconv,int l,int n,PetscScalar* work,Vec wv,Vec wu)
53: {
55: PetscReal a,b;
56: int i,j,k=nconv+l;
59: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
60: if (l>0) {
61: VecSet(wu,0.0);
62: VecMAXPY(wu,l,bb,U+nconv);
63: VecAXPY(U[k],-1.0,wu);
64: }
65: for (i=k+1;i<n;i++) {
66: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
67: IPNormBegin(svd->ip,U[i-1],&a);
68: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
69: IPNormEnd(svd->ip,U[i-1],&a);
70: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
71:
72: VecScale(U[i-1],1.0/a);
73: VecScale(V[i],1.0/a);
74: for (j=0;j<i;j++) work[j] = - work[j] / a;
75: VecMAXPY(V[i],i,work,V);
77: IPOrthogonalizeCGS(svd->ip,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b,wv);
78: VecScale(V[i],1.0/b);
79:
80: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
81: VecAXPY(U[i],-b,U[i-1]);
83: alpha[i-k-1] = a;
84: beta[i-k-1] = b;
85: }
86: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
87: IPNormBegin(svd->ip,U[n-1],&a);
88: IPMInnerProductBegin(svd->ip,v,n,V,work);
89: IPNormEnd(svd->ip,U[n-1],&a);
90: IPMInnerProductEnd(svd->ip,v,n,V,work);
91:
92: VecScale(U[n-1],1.0/a);
93: VecScale(v,1.0/a);
94: for (j=0;j<n;j++) work[j] = - work[j] / a;
95: VecMAXPY(v,n,work,V);
97: IPOrthogonalizeCGS(svd->ip,n,PETSC_NULL,V,v,work,PETSC_NULL,&b,wv);
98:
99: alpha[n-k-1] = a;
100: beta[n-k-1] = b;
101: return(0);
102: }
106: PetscErrorCode SVDSolve_TRLANCZOS(SVD svd)
107: {
109: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
110: PetscReal *alpha,*beta,norm;
111: PetscScalar *b,*Q,*PT,*swork;
112: PetscInt *perm;
113: int i,j,k,l,m,n,nwork=0;
114: Vec v,wv,wu,*workV,*workU,*permV,*permU;
115: PetscTruth conv;
116:
118: /* allocate working space */
119: PetscMalloc(sizeof(PetscReal)*svd->n,&alpha);
120: PetscMalloc(sizeof(PetscReal)*svd->n,&beta);
121: PetscMalloc(sizeof(PetscScalar)*svd->n,&b);
122: PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&Q);
123: PetscMalloc(sizeof(PetscScalar)*svd->n*svd->n,&PT);
124: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
125: VecDuplicate(svd->V[0],&v);
126: VecDuplicate(svd->V[0],&wv);
127: VecDuplicate(svd->U[0],&wu);
128: PetscMalloc(sizeof(Vec)*svd->n,&workV);
129: PetscMalloc(sizeof(Vec)*svd->n,&workU);
130:
131: /* normalize start vector */
132: VecCopy(svd->vec_initial,svd->V[0]);
133: VecNormalize(svd->V[0],&norm);
134:
135: l = 0;
136: while (svd->reason == SVD_CONVERGED_ITERATING) {
137: svd->its++;
139: /* inner loop */
140: if (lanczos->oneside) {
141: SVDOneSideTRLanczos(svd,alpha,beta,b+svd->nconv,svd->V,v,svd->U,svd->nconv,l,svd->n,swork,wv,wu);
142: } else {
143: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv+l,svd->n,swork,wv,wu);
144: }
145: VecScale(v,1.0/beta[svd->n-svd->nconv-l-1]);
146:
147: /* compute SVD of general matrix */
148: n = svd->n - svd->nconv;
149: /* first l columns */
150: for (j=0;j<l;j++) {
151: for (i=0;i<j;i++) Q[j*n+i] = 0.0;
152: Q[j*n+j] = svd->sigma[svd->nconv+j];
153: for (i=j+1;i<n;i++) Q[j*n+i] = 0.0;
154: }
155: /* l+1 column */
156: for (i=0;i<l;i++) Q[l*n+i] = b[i+svd->nconv];
157: Q[l*n+l] = alpha[0];
158: for (i=l+1;i<n;i++) Q[l*n+i] = 0.0;
159: /* rest of matrix */
160: for (j=l+1;j<n;j++) {
161: for (i=0;i<j-1;i++) Q[j*n+i] = 0.0;
162: Q[j*n+j-1] = beta[j-l-1];
163: Q[j*n+j] = alpha[j-l];
164: for (i=j+1;i<n;i++) Q[j*n+i] = 0.0;
165: }
166: SVDDense(n,n,Q,alpha,PETSC_NULL,PT);
168: /* compute error estimates */
169: k = 0;
170: conv = PETSC_TRUE;
171: for (i=svd->nconv;i<svd->n;i++) {
172: if (svd->which == SVD_SMALLEST) j = n-i+svd->nconv-1;
173: else j = i-svd->nconv;
174: svd->sigma[i] = alpha[j];
175: b[i] = Q[j*n+n-1]*beta[n-l-1];
176: svd->errest[i] = PetscAbsScalar(b[i]);
177: if (alpha[j] > svd->tol) svd->errest[i] /= alpha[j];
178: if (conv) {
179: if (svd->errest[i] < svd->tol) k++;
180: else conv = PETSC_FALSE;
181: }
182: }
183:
184: /* check convergence and update l */
185: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
186: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
187: if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
188: else l = PetscMax((svd->n - svd->nconv - k) / 2,1);
189:
190: /* allocate work space for converged singular and restart vectors */
191: if (nwork<k+l) {
192: for (i=nwork;i<k+l;i++) {
193: SVDMatGetVecs(svd,workV+i,workU+i);
194: }
195: nwork = k+l;
196: }
197:
198: /* compute converged singular vectors and restart vectors*/
199: for (i=0;i<k+l;i++) {
200: if (svd->which == SVD_SMALLEST) j = n-i-1;
201: else j = i;
202: VecSet(workV[i],0.0);
203: for (m=0;m<n;m++) swork[m] = PT[m*n+j];
204: VecMAXPY(workV[i],n,swork,svd->V+svd->nconv);
205: VecSet(workU[i],0.0);
206: VecMAXPY(workU[i],n,Q+j*n,svd->U+svd->nconv);
207: }
208:
209: /* copy the last vector to be the next initial vector */
210: if (svd->reason == SVD_CONVERGED_ITERATING) {
211: VecCopy(v,svd->V[svd->nconv+k+l]);
212: }
213:
214: /* copy converged singular vectors and restart vectors from temporary space */
215: for (i=0;i<k+l;i++) {
216: VecCopy(workV[i],svd->V[i+svd->nconv]);
217: VecCopy(workU[i],svd->U[i+svd->nconv]);
218: }
219:
220: svd->nconv += k;
221: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,svd->n);
222: }
223:
224: /* sort singular triplets */
225: PetscMalloc(sizeof(PetscInt)*svd->nconv,&perm);
226: PetscMalloc(sizeof(Vec)*svd->nconv,&permV);
227: PetscMalloc(sizeof(Vec)*svd->nconv,&permU);
228: for (i=0;i<svd->nconv;i++) {
229: alpha[i] = svd->sigma[i];
230: beta[i] = svd->errest[i];
231: permV[i] = svd->V[i];
232: permU[i] = svd->U[i];
233: perm[i] = i;
234: }
235: PetscSortRealWithPermutation(svd->nconv,svd->sigma,perm);
236: for (i=0;i<svd->nconv;i++) {
237: if (svd->which == SVD_SMALLEST) j = perm[i];
238: else j = perm[svd->nconv-i-1];
239: svd->sigma[i] = alpha[j];
240: svd->errest[i] = beta[j];
241: svd->V[i] = permV[j];
242: svd->U[i] = permU[j];
243: }
244:
245: /* free working space */
246: VecDestroy(v);
247: VecDestroy(wv);
248: VecDestroy(wu);
249: for (i=0;i<nwork;i++) { VecDestroy(workV[i]); }
250: PetscFree(workV);
251: for (i=0;i<nwork;i++) { VecDestroy(workU[i]); }
252: PetscFree(workU);
254: PetscFree(alpha);
255: PetscFree(beta);
256: PetscFree(b);
257: PetscFree(Q);
258: PetscFree(PT);
259: PetscFree(swork);
260: PetscFree(perm);
261: PetscFree(permV);
262: PetscFree(permU);
263: return(0);
264: }
268: PetscErrorCode SVDSetFromOptions_TRLANCZOS(SVD svd)
269: {
271: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
274: PetscOptionsBegin(svd->comm,svd->prefix,"TRLANCZOS Singular Value Solver Options","SVD");
275: PetscOptionsTruth("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",PETSC_FALSE,&lanczos->oneside,PETSC_NULL);
276: PetscOptionsEnd();
277: return(0);
278: }
283: PetscErrorCode SVDTRLanczosSetOneSide_TRLANCZOS(SVD svd,PetscTruth oneside)
284: {
285: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
288: lanczos->oneside = oneside;
289: return(0);
290: }
295: /*@
296: SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
297: to be used is one-sided or two-sided.
299: Collective on SVD
301: Input Parameters:
302: + svd - singular value solver
303: - oneside - boolean flag indicating if the method is one-sided or not
305: Options Database Key:
306: . -svd_trlanczos_oneside <boolean> - Indicates the boolean flag
308: Note:
309: By default, a two-sided variant is selected, which is sometimes slightly
310: more robust. However, the one-sided variant is faster because it avoids
311: the orthogonalization associated to left singular vectors.
313: Level: advanced
315: .seealso: SVDLanczosSetOneSide()
316: @*/
317: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscTruth oneside)
318: {
319: PetscErrorCode ierr, (*f)(SVD,PetscTruth);
323: PetscObjectQueryFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",(void (**)())&f);
324: if (f) {
325: (*f)(svd,oneside);
326: }
327: return(0);
328: }
332: PetscErrorCode SVDView_TRLANCZOS(SVD svd,PetscViewer viewer)
333: {
335: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
338: PetscViewerASCIIPrintf(viewer,"Lanczos reorthogonalization: %s\n",lanczos->oneside ? "one-side" : "two-side");
339: return(0);
340: }
345: PetscErrorCode SVDCreate_TRLANCZOS(SVD svd)
346: {
348: SVD_TRLANCZOS *lanczos;
351: PetscNew(SVD_TRLANCZOS,&lanczos);
352: PetscLogObjectMemory(svd,sizeof(SVD_TRLANCZOS));
353: svd->data = (void *)lanczos;
354: svd->ops->setup = SVDSetUp_TRLANCZOS;
355: svd->ops->solve = SVDSolve_TRLANCZOS;
356: svd->ops->destroy = SVDDestroy_Default;
357: svd->ops->setfromoptions = SVDSetFromOptions_TRLANCZOS;
358: svd->ops->view = SVDView_TRLANCZOS;
359: lanczos->oneside = PETSC_FALSE;
360: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","SVDTRLanczosSetOneSide_TRLANCZOS",SVDTRLanczosSetOneSide_TRLANCZOS);
361: return(0);
362: }