Actual source code: trlanczos.c
1: /*
3: SLEPc singular value solver: "trlanczos"
5: Method: Thick-restart Lanczos
7: Algorithm:
9: Golub-Kahan-Lanczos bidiagonalization with thick-restart.
11: References:
13: [1] G.H. Golub and W. Kahan, "Calculating the singular values
14: and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
15: B 2:205-224, 1965.
17: [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
18: efficient parallel SVD solver based on restarted Lanczos
19: bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
20: 2008.
22: Last update: Jun 2007
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: SLEPc - Scalable Library for Eigenvalue Problem Computations
26: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
28: This file is part of SLEPc.
30: SLEPc is free software: you can redistribute it and/or modify it under the
31: terms of version 3 of the GNU Lesser General Public License as published by
32: the Free Software Foundation.
34: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
35: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
36: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
37: more details.
39: You should have received a copy of the GNU Lesser General Public License
40: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
41: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42: */
44: #include <slepc-private/svdimpl.h> /*I "slepcsvd.h" I*/
45: #include <slepc-private/ipimpl.h>
46: #include <slepcblaslapack.h>
48: typedef struct {
49: PetscBool oneside;
50: } SVD_TRLANCZOS;
54: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
55: {
57: PetscInt N;
60: SVDMatGetSize(svd,NULL,&N);
61: if (svd->ncv) { /* ncv set */
62: if (svd->ncv<svd->nsv) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must be at least nsv");
63: } else if (svd->mpd) { /* mpd set */
64: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
65: } else { /* neither set: defaults depend on nsv being small or large */
66: if (svd->nsv<500) svd->ncv = PetscMin(N,PetscMax(2*svd->nsv,10));
67: else {
68: svd->mpd = 500;
69: svd->ncv = PetscMin(N,svd->nsv+svd->mpd);
70: }
71: }
72: if (!svd->mpd) svd->mpd = svd->ncv;
73: if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
74: if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
75: if (svd->ncv!=svd->n) {
76: VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
77: PetscLogObjectParents(svd,svd->ncv,svd->U);
78: }
79: DSSetType(svd->ds,DSSVD);
80: DSSetCompact(svd->ds,PETSC_TRUE);
81: DSAllocate(svd->ds,svd->ncv);
82: return(0);
83: }
87: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec* U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
88: {
90: PetscReal a,b;
91: PetscInt i,k=nconv+l;
94: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
95: if (l>0) {
96: for (i=0;i<l;i++) work[i]=beta[i+nconv];
97: SlepcVecMAXPBY(U[k],1.0,-1.0,l,work,U+nconv);
98: }
99: IPNorm(svd->ip,U[k],&a);
100: VecScale(U[k],1.0/a);
101: alpha[k] = a;
102: for (i=k+1;i<n;i++) {
103: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
104: IPOrthogonalize(svd->ip,0,NULL,i,NULL,V,V[i],work,&b,NULL);
105: VecScale(V[i],1.0/b);
106: beta[i-1] = b;
108: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
109: VecAXPY(U[i],-b,U[i-1]);
110: IPNorm(svd->ip,U[i],&a);
111: VecScale(U[i],1.0/a);
112: alpha[i] = a;
113: }
114: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
115: IPOrthogonalize(svd->ip,0,NULL,n,NULL,V,v,work,&b,NULL);
116: beta[n-1] = b;
117: return(0);
118: }
122: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec* U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
123: {
124: PetscErrorCode ierr;
125: PetscReal a,b,sum,onorm,eta;
126: PetscScalar dot;
127: PetscInt i,j,k=nconv+l;
128: IPOrthogRefineType refine;
131: SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
132: if (l>0) {
133: for (i=0;i<l;i++) work[i]=beta[i+nconv];
134: SlepcVecMAXPBY(U[k],1.0,-1.0,l,work,U+nconv);
135: }
136: IPGetOrthogonalization(svd->ip,NULL,&refine,&eta);
137: for (i=k+1;i<n;i++) {
138: SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
139: IPNormBegin(svd->ip,U[i-1],&a);
140: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
141: IPInnerProductBegin(svd->ip,V[i],V[i],&dot);
142: }
143: IPMInnerProductBegin(svd->ip,V[i],i,V,work);
144: IPNormEnd(svd->ip,U[i-1],&a);
145: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
146: IPInnerProductEnd(svd->ip,V[i],V[i],&dot);
147: }
148: IPMInnerProductEnd(svd->ip,V[i],i,V,work);
150: VecScale(U[i-1],1.0/a);
151: for (j=0;j<i;j++) work[j] = work[j] / a;
152: SlepcVecMAXPBY(V[i],1.0/a,-1.0,i,work,V);
154: switch (refine) {
155: case IP_ORTHOG_REFINE_NEVER:
156: IPNorm(svd->ip,V[i],&b);
157: break;
158: case IP_ORTHOG_REFINE_ALWAYS:
159: IPOrthogonalizeCGS1(svd->ip,0,NULL,i,NULL,V,V[i],work,NULL,&b);
160: break;
161: case IP_ORTHOG_REFINE_IFNEEDED:
162: onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
163: sum = 0.0;
164: for (j=0;j<i;j++) {
165: sum += PetscRealPart(work[j] * PetscConj(work[j]));
166: }
167: b = PetscRealPart(dot)/(a*a) - sum;
168: if (b>0.0) b = PetscSqrtReal(b);
169: else {
170: IPNorm(svd->ip,V[i],&b);
171: }
172: if (b < eta*onorm) {
173: IPOrthogonalizeCGS1(svd->ip,0,NULL,i,NULL,V,V[i],work,NULL,&b);
174: }
175: break;
176: }
178: VecScale(V[i],1.0/b);
180: SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
181: VecAXPY(U[i],-b,U[i-1]);
183: alpha[i-1] = a;
184: beta[i-1] = b;
185: }
186: SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
187: IPNormBegin(svd->ip,U[n-1],&a);
188: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
189: IPInnerProductBegin(svd->ip,v,v,&dot);
190: }
191: IPMInnerProductBegin(svd->ip,v,n,V,work);
192: IPNormEnd(svd->ip,U[n-1],&a);
193: if (refine == IP_ORTHOG_REFINE_IFNEEDED) {
194: IPInnerProductEnd(svd->ip,v,v,&dot);
195: }
196: IPMInnerProductEnd(svd->ip,v,n,V,work);
198: VecScale(U[n-1],1.0/a);
199: for (j=0;j<n;j++) work[j] = work[j] / a;
200: SlepcVecMAXPBY(v,1.0/a,-1.0,n,work,V);
202: switch (refine) {
203: case IP_ORTHOG_REFINE_NEVER:
204: IPNorm(svd->ip,v,&b);
205: break;
206: case IP_ORTHOG_REFINE_ALWAYS:
207: IPOrthogonalizeCGS1(svd->ip,0,NULL,n,NULL,V,v,work,NULL,&b);
208: break;
209: case IP_ORTHOG_REFINE_IFNEEDED:
210: onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
211: sum = 0.0;
212: for (j=0;j<i;j++) {
213: sum += PetscRealPart(work[j] * PetscConj(work[j]));
214: }
215: b = PetscRealPart(dot)/(a*a) - sum;
216: if (b>0.0) b = PetscSqrtReal(b);
217: else {
218: IPNorm(svd->ip,v,&b);
219: }
220: if (b < eta*onorm) {
221: IPOrthogonalizeCGS1(svd->ip,0,NULL,n,NULL,V,v,work,NULL,&b);
222: }
223: break;
224: }
226: alpha[n-1] = a;
227: beta[n-1] = b;
228: return(0);
229: }
233: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
234: {
236: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
237: PetscReal *alpha,*beta,lastbeta,norm;
238: PetscScalar *Q,*PT,*swork,*w;
239: PetscInt i,k,l,nv,ld,off;
240: Vec v;
241: PetscBool conv;
242: IPOrthogType orthog;
245: /* allocate working space */
246: DSGetLeadingDimension(svd->ds,&ld);
247: PetscMalloc(sizeof(PetscScalar)*ld,&w);
248: PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);
249: VecDuplicate(svd->V[0],&v);
250: IPGetOrthogonalization(svd->ip,&orthog,NULL,NULL);
252: /* normalize start vector */
253: if (!svd->nini) {
254: SlepcVecSetRandom(svd->V[0],svd->rand);
255: }
256: VecNormalize(svd->V[0],&norm);
258: l = 0;
259: while (svd->reason == SVD_CONVERGED_ITERATING) {
260: svd->its++;
262: /* inner loop */
263: nv = PetscMin(svd->nconv+svd->mpd,svd->n);
264: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
265: beta = alpha + ld;
266: if (lanczos->oneside) {
267: if (orthog == IP_ORTHOG_MGS) {
268: SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,l,nv,swork);
269: } else {
270: SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,l,nv,swork);
271: }
272: } else {
273: SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv+l,nv,swork);
274: }
275: lastbeta = beta[nv-1];
276: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
277: VecScale(v,1.0/lastbeta);
279: /* compute SVD of general matrix */
280: DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);
281: if (l==0) {
282: DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
283: } else {
284: DSSetState(svd->ds,DS_STATE_RAW);
285: }
286: DSSolve(svd->ds,w,NULL);
287: DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
289: /* compute error estimates */
290: k = 0;
291: conv = PETSC_TRUE;
292: DSGetArray(svd->ds,DS_MAT_U,&Q);
293: DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
294: beta = alpha + ld;
295: for (i=svd->nconv;i<nv;i++) {
296: svd->sigma[i] = PetscRealPart(w[i]);
297: beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta;
298: svd->errest[i] = PetscAbsScalar(beta[i]);
299: if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
300: if (conv) {
301: if (svd->errest[i] < svd->tol) k++;
302: else conv = PETSC_FALSE;
303: }
304: }
305: DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
307: /* check convergence and update l */
308: if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
309: if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
310: if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
311: else l = PetscMax((nv-svd->nconv-k)/2,0);
313: /* compute converged singular vectors and restart vectors */
314: off = svd->nconv+svd->nconv*ld;
315: DSGetArray(svd->ds,DS_MAT_VT,&PT);
316: SlepcUpdateVectors(nv-svd->nconv,svd->V+svd->nconv,0,k+l,PT+off,ld,PETSC_TRUE);
317: SlepcUpdateVectors(nv-svd->nconv,svd->U+svd->nconv,0,k+l,Q+off,ld,PETSC_FALSE);
318: DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
319: DSRestoreArray(svd->ds,DS_MAT_U,&Q);
321: /* copy the last vector to be the next initial vector */
322: if (svd->reason == SVD_CONVERGED_ITERATING) {
323: VecCopy(v,svd->V[svd->nconv+k+l]);
324: }
326: svd->nconv += k;
327: SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
328: }
330: /* orthonormalize U columns in one side method */
331: if (lanczos->oneside) {
332: for (i=0;i<svd->nconv;i++) {
333: IPOrthogonalize(svd->ip,0,NULL,i,NULL,svd->U,svd->U[i],NULL,&norm,NULL);
334: VecScale(svd->U[i],1.0/norm);
335: }
336: }
338: /* free working space */
339: VecDestroy(&v);
340: PetscFree(w);
341: PetscFree(swork);
342: return(0);
343: }
347: PetscErrorCode SVDSetFromOptions_TRLanczos(SVD svd)
348: {
350: PetscBool set,val;
351: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
354: PetscOptionsHead("SVD TRLanczos Options");
355: PetscOptionsBool("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
356: if (set) {
357: SVDTRLanczosSetOneSide(svd,val);
358: }
359: PetscOptionsTail();
360: return(0);
361: }
365: static PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
366: {
367: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
370: lanczos->oneside = oneside;
371: return(0);
372: }
376: /*@
377: SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
378: to be used is one-sided or two-sided.
380: Logically Collective on SVD
382: Input Parameters:
383: + svd - singular value solver
384: - oneside - boolean flag indicating if the method is one-sided or not
386: Options Database Key:
387: . -svd_trlanczos_oneside <boolean> - Indicates the boolean flag
389: Note:
390: By default, a two-sided variant is selected, which is sometimes slightly
391: more robust. However, the one-sided variant is faster because it avoids
392: the orthogonalization associated to left singular vectors.
394: Level: advanced
396: .seealso: SVDLanczosSetOneSide()
397: @*/
398: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
399: {
405: PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
406: return(0);
407: }
411: /*@
412: SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
413: to be used is one-sided or two-sided.
415: Not Collective
417: Input Parameters:
418: . svd - singular value solver
420: Output Parameters:
421: . oneside - boolean flag indicating if the method is one-sided or not
423: Level: advanced
425: .seealso: SVDTRLanczosSetOneSide()
426: @*/
427: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
428: {
434: PetscTryMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
435: return(0);
436: }
440: static PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
441: {
442: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
445: *oneside = lanczos->oneside;
446: return(0);
447: }
451: PetscErrorCode SVDReset_TRLanczos(SVD svd)
452: {
456: VecDestroyVecs(svd->n,&svd->U);
457: return(0);
458: }
462: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
463: {
467: PetscFree(svd->data);
468: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",NULL);
469: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",NULL);
470: return(0);
471: }
475: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
476: {
478: SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;
481: PetscViewerASCIIPrintf(viewer," TRLanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
482: return(0);
483: }
487: PETSC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD svd)
488: {
492: PetscNewLog(svd,SVD_TRLANCZOS,&svd->data);
493: svd->ops->setup = SVDSetUp_TRLanczos;
494: svd->ops->solve = SVDSolve_TRLanczos;
495: svd->ops->destroy = SVDDestroy_TRLanczos;
496: svd->ops->reset = SVDReset_TRLanczos;
497: svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
498: svd->ops->view = SVDView_TRLanczos;
499: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",SVDTRLanczosSetOneSide_TRLanczos);
500: PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",SVDTRLanczosGetOneSide_TRLanczos);
501: return(0);
502: }