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-2012, 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 <slepc-private/svdimpl.h> /*I "slepcsvd.h" I*/
30: #include <slepc-private/ipimpl.h>
31: #include <slepcblaslapack.h>
33: typedef struct {
34: PetscBool oneside;
35: } SVD_TRLANCZOS;
39: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
40: {
42: PetscInt N;
45: SVDMatGetSize(svd,PETSC_NULL,&N);
46: if (svd->ncv) { /* ncv set */
47: if (svd->ncv<svd->nsv) SETERRQ(((PetscObject)svd)->comm,1,"The value of ncv must be at least nsv");
48: }
49: else if (svd->mpd) { /* mpd set */
50: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
51: }
52: else { /* neither set: defaults depend on nsv being small or large */
53: if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
54: else { svd->mpd = 500; svd->ncv = PetscMin(N,svd->nsv+svd->mpd); }
55: }
56: if (!svd->mpd) svd->mpd = svd->ncv;
57: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(((PetscObject)svd)->comm,1,"The value of ncv must not be larger than nev+mpd");
58: if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
59: if (svd->ncv!=svd->n) {
60: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
61: }
62: DSSetType(svd->ds,DSSVD);
63: DSSetCompact(svd->ds,PETSC_TRUE);
64: DSAllocate(svd->ds,svd->ncv);
65: return(0);
66: }
70: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec* U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
71: {
73: PetscReal a,b;
74: PetscInt i,k=nconv+l;
77: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
78: if (l>0) {
79: for (i=0;i<l;i++) work[i]=beta[i+nconv];
80: SlepcVecMAXPBY(U[k],1.0,-1.0,l,work,U+nconv);
81: }
82: IPNorm(svd->ip,U[k],&a);
83: VecScale(U[k],1.0/a);
84: alpha[k] = a;
85: for (i=k+1;i<n;i++) {
86: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
87: IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,&b,PETSC_NULL);
88: VecScale(V[i],1.0/b);
89: beta[i-1] = b;
90:
91: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
92: VecAXPY(U[i],-b,U[i-1]);
93: IPNorm(svd->ip,U[i],&a);
94: VecScale(U[i],1.0/a);
95: alpha[i] = a;
96: }
97: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
98: IPOrthogonalize(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,&b,PETSC_NULL);
99: beta[n-1] = b;
100: return(0);
101: }
105: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec* U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
106: {
107: PetscErrorCode ierr;
108: PetscReal a,b,sum,onorm,eta;
109: PetscScalar dot;
110: PetscInt i,j,k=nconv+l;
111: IPOrthogRefineType refine;
114: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
115: if (l>0) {
116: for (i=0;i<l;i++) work[i]=beta[i+nconv];
117: SlepcVecMAXPBY(U[k],1.0,-1.0,l,work,U+nconv);
118: }
119: IPGetOrthogonalization(svd->ip,PETSC_NULL,&refine,&eta);
120: for (i=k+1;i<n;i++) {
121: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
122: IPNormBegin(svd->ip,U[i-1],&a);
123: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
124: IPInnerProductBegin(svd->ip,V[i],V[i],&dot);
125: }
126: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
127: IPNormEnd(svd->ip,U[i-1],&a);
128: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
129: IPInnerProductEnd(svd->ip,V[i],V[i],&dot);
130: }
131: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
132:
133: VecScale(U[i-1],1.0/a);
134: for (j=0;j<i;j++) work[j] = work[j] / a;
135: SlepcVecMAXPBY(V[i],1.0/a,-1.0,i,work,V);
137: switch (refine) {
138: case IP_ORTHOG_REFINE_NEVER:
139: IPNorm(svd->ip,V[i],&b);
140: break;
141: case IP_ORTHOG_REFINE_ALWAYS:
142: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
143: break;
144: case IP_ORTHOG_REFINE_IFNEEDED:
145: onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
146: sum = 0.0;
147: for (j=0;j<i;j++) {
148: sum += PetscRealPart(work[j] * PetscConj(work[j]));
149: }
150: b = PetscRealPart(dot)/(a*a) - sum;
151: if (b>0.0) b = PetscSqrtReal(b);
152: else {
153: IPNorm(svd->ip,V[i],&b);
154: }
155: if (b < eta*onorm) {
156: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
157: }
158: break;
159: }
160:
161: VecScale(V[i],1.0/b);
162:
163: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
164: VecAXPY(U[i],-b,U[i-1]);
166: alpha[i-1] = a;
167: beta[i-1] = b;
168: }
169: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
170: IPNormBegin(svd->ip,U[n-1],&a);
171: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
172: IPInnerProductBegin(svd->ip,v,v,&dot);
173: }
174: IPMInnerProductBegin(svd->ip,v,n,V,work);
175: IPNormEnd(svd->ip,U[n-1],&a);
176: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
177: IPInnerProductEnd(svd->ip,v,v,&dot);
178: }
179: IPMInnerProductEnd(svd->ip,v,n,V,work);
180:
181: VecScale(U[n-1],1.0/a);
182: for (j=0;j<n;j++) work[j] = work[j] / a;
183: SlepcVecMAXPBY(v,1.0/a,-1.0,n,work,V);
185: switch (refine) {
186: case IP_ORTHOG_REFINE_NEVER:
187: IPNorm(svd->ip,v,&b);
188: break;
189: case IP_ORTHOG_REFINE_ALWAYS:
190: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
191: break;
192: case IP_ORTHOG_REFINE_IFNEEDED:
193: onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
194: sum = 0.0;
195: for (j=0;j<i;j++) {
196: sum += PetscRealPart(work[j] * PetscConj(work[j]));
197: }
198: b = PetscRealPart(dot)/(a*a) - sum;
199: if (b>0.0) b = PetscSqrtReal(b);
200: else {
201: IPNorm(svd->ip,v,&b);
202: }
203: if (b < eta*onorm) {
204: IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
205: }
206: break;
207: }
208:
209: alpha[n-1] = a;
210: beta[n-1] = b;
211: return(0);
212: }
216: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
217: {
219: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
220: PetscReal *alpha,*beta,lastbeta,norm;
221: PetscScalar *Q,*PT,*swork,*w;
222: PetscInt i,k,l,nv,ld,off;
223: Vec v;
224: PetscBool conv;
225: IPOrthogType orthog;
226:
228: /* allocate working space */
229: DSGetLeadingDimension(svd->ds,&ld);
230: PetscMalloc(sizeof(PetscScalar)*ld,&w);
231: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
232: VecDuplicate(svd->V[0],&v);
233: IPGetOrthogonalization(svd->ip,&orthog,PETSC_NULL,PETSC_NULL);
234:
235: /* normalize start vector */
236: if (svd->nini==0) {
237: SlepcVecSetRandom(svd->V[0],svd->rand);
238: }
239: VecNormalize(svd->V[0],&norm);
240:
241: l = 0;
242: while (svd->reason == SVD_CONVERGED_ITERATING) {
243: svd->its++;
245: /* inner loop */
246: nv = PetscMin(svd->nconv+svd->mpd,svd->n);
247: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
248: beta = alpha + ld;
249: if (lanczos->oneside) {
250: if (orthog == IP_ORTHOG_MGS) {
251: SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,l,nv,swork);
252: } else {
253: SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,l,nv,swork);
254: }
255: } else {
256: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv+l,nv,swork);
257: }
258: lastbeta = beta[nv-1];
259: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
260: VecScale(v,1.0/lastbeta);
261:
262: /* compute SVD of general matrix */
263: DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);
264: if (l==0) {
265: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
266: } else {
267: DSSetState(svd->ds,DS_STATE_RAW);
268: }
269: DSSolve(svd->ds,w,PETSC_NULL);
270: DSSort(svd->ds,w,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
272: /* compute error estimates */
273: k = 0;
274: conv = PETSC_TRUE;
275: DSGetArray(svd->ds,DS_MAT_U,&Q);
276: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
277: beta = alpha + ld;
278: for (i=svd->nconv;i<nv;i++) {
279: svd->sigma[i] = PetscRealPart(w[i]);
280: beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta;
281: svd->errest[i] = PetscAbsScalar(beta[i]);
282: if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
283: if (conv) {
284: if (svd->errest[i] < svd->tol) k++;
285: else conv = PETSC_FALSE;
286: }
287: }
288: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
289:
290: /* check convergence and update l */
291: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
292: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
293: if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
294: else l = PetscMax((nv-svd->nconv-k)/2,0);
295:
296: /* compute converged singular vectors and restart vectors */
297: off = svd->nconv+svd->nconv*ld;
298: DSGetArray(svd->ds,DS_MAT_VT,&PT);
299: SlepcUpdateVectors(nv-svd->nconv,svd->V+svd->nconv,0,k+l,PT+off,ld,PETSC_TRUE);
300: SlepcUpdateVectors(nv-svd->nconv,svd->U+svd->nconv,0,k+l,Q+off,ld,PETSC_FALSE);
301: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
302: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
303:
304: /* copy the last vector to be the next initial vector */
305: if (svd->reason == SVD_CONVERGED_ITERATING) {
306: VecCopy(v,svd->V[svd->nconv+k+l]);
307: }
308:
309: svd->nconv += k;
310: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
311: }
312:
313: /* orthonormalize U columns in one side method */
314: if (lanczos->oneside) {
315: for (i=0;i<svd->nconv;i++) {
316: IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,svd->U,svd->U[i],PETSC_NULL,&norm,PETSC_NULL);
317: VecScale(svd->U[i],1.0/norm);
318: }
319: }
320:
321: /* free working space */
322: VecDestroy(&v);
323: PetscFree(w);
324: PetscFree(swork);
325: return(0);
326: }
330: PetscErrorCode SVDSetFromOptions_TRLanczos(SVD svd)
331: {
333: PetscBool set,val;
334: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
337: PetscOptionsHead("SVD TRLanczos Options");
338: PetscOptionsBool("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
339: if (set) {
340: SVDTRLanczosSetOneSide(svd,val);
341: }
342: PetscOptionsTail();
343: return(0);
344: }
346: EXTERN_C_BEGIN
349: PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
350: {
351: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
354: lanczos->oneside = oneside;
355: return(0);
356: }
357: EXTERN_C_END
361: /*@
362: SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
363: to be used is one-sided or two-sided.
365: Logically Collective on SVD
367: Input Parameters:
368: + svd - singular value solver
369: - oneside - boolean flag indicating if the method is one-sided or not
371: Options Database Key:
372: . -svd_trlanczos_oneside <boolean> - Indicates the boolean flag
374: Note:
375: By default, a two-sided variant is selected, which is sometimes slightly
376: more robust. However, the one-sided variant is faster because it avoids
377: the orthogonalization associated to left singular vectors.
379: Level: advanced
381: .seealso: SVDLanczosSetOneSide()
382: @*/
383: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
384: {
390: PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
391: return(0);
392: }
396: /*@
397: SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
398: to be used is one-sided or two-sided.
400: Not Collective
402: Input Parameters:
403: . svd - singular value solver
405: Output Parameters:
406: . oneside - boolean flag indicating if the method is one-sided or not
408: Level: advanced
410: .seealso: SVDTRLanczosSetOneSide()
411: @*/
412: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
413: {
419: PetscTryMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
420: return(0);
421: }
423: EXTERN_C_BEGIN
426: PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
427: {
428: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
431: *oneside = lanczos->oneside;
432: return(0);
433: }
434: EXTERN_C_END
438: PetscErrorCode SVDReset_TRLanczos(SVD svd)
439: {
443: VecDestroyVecs(svd->n,&svd->U);
444: return(0);
445: }
449: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
450: {
454: PetscFree(svd->data);
455: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","",PETSC_NULL);
456: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosGetOneSide_C","",PETSC_NULL);
457: return(0);
458: }
462: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
463: {
465: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS *)svd->data;
468: PetscViewerASCIIPrintf(viewer," TRLanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
469: return(0);
470: }
472: EXTERN_C_BEGIN
475: PetscErrorCode SVDCreate_TRLanczos(SVD svd)
476: {
480: PetscNewLog(svd,SVD_TRLANCZOS,&svd->data);
481: svd->ops->setup = SVDSetUp_TRLanczos;
482: svd->ops->solve = SVDSolve_TRLanczos;
483: svd->ops->destroy = SVDDestroy_TRLanczos;
484: svd->ops->reset = SVDReset_TRLanczos;
485: svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
486: svd->ops->view = SVDView_TRLanczos;
487: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosSetOneSide_C","SVDTRLanczosSetOneSide_TRLanczos",SVDTRLanczosSetOneSide_TRLanczos);
488: PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDTRLanczosGetOneSide_C","SVDTRLanczosGetOneSide_TRLanczos",SVDTRLanczosGetOneSide_TRLanczos);
489: return(0);
490: }
491: EXTERN_C_END