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-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_LANCZOS;

 39: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
 40: {
 42:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS *)svd->data;
 43:   PetscInt       N;

 46:   SVDMatGetSize(svd,PETSC_NULL,&N);
 47:   if (svd->ncv) { /* ncv set */
 48:     if (svd->ncv<svd->nsv) SETERRQ(((PetscObject)svd)->comm,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(((PetscObject)svd)->comm,1,"The value of ncv must not be larger than nev+mpd");
 59:   if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
 60:   if (!lanczos->oneside && svd->ncv != svd->n) {
 61:     if (svd->n) { VecDestroyVecs(svd->n,&svd->U); }
 62:     VecDuplicateVecs(svd->tl,svd->ncv,&svd->U);
 63:   }
 64:   DSSetType(svd->ds,DSSVD);
 65:   DSSetCompact(svd->ds,PETSC_TRUE);
 66:   DSAllocate(svd->ds,svd->ncv);
 67:   return(0);
 68: }

 72: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec *U,PetscInt k,PetscInt n,PetscScalar* work)
 73: {
 75:   PetscInt       i;
 76: 
 78:   SVDMatMult(svd,PETSC_FALSE,V[k],U[k]);
 79:   IPOrthogonalize(svd->ip,0,PETSC_NULL,k,PETSC_NULL,U,U[k],work,alpha+k,PETSC_NULL);
 80:   VecScale(U[k],1.0/alpha[k]);
 81:   for (i=k+1;i<n;i++) {
 82:     SVDMatMult(svd,PETSC_TRUE,U[i-1],V[i]);
 83:     IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,beta+i-1,PETSC_NULL);
 84:     VecScale(V[i],1.0/beta[i-1]);

 86:     SVDMatMult(svd,PETSC_FALSE,V[i],U[i]);
 87:     IPOrthogonalize(svd->ip,0,PETSC_NULL,i,PETSC_NULL,U,U[i],work,alpha+i,PETSC_NULL);
 88:     VecScale(U[i],1.0/alpha[i]);
 89:   }
 90:   SVDMatMult(svd,PETSC_TRUE,U[n-1],v);
 91:   IPOrthogonalize(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,beta+n-1,PETSC_NULL);
 92:   return(0);
 93: }

 97: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,Vec *V,Vec v,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
 98: {
100:   PetscInt       i;
101:   PetscReal      a,b;
102:   Vec            temp;
103: 
105:   SVDMatMult(svd,PETSC_FALSE,V[k],u);
106:   for (i=k+1;i<n;i++) {
107:     SVDMatMult(svd,PETSC_TRUE,u,V[i]);
108:     IPNormBegin(svd->ip,u,&a);
109:     IPMInnerProductBegin(svd->ip,V[i],i,V,work);
110:     IPNormEnd(svd->ip,u,&a);
111:     IPMInnerProductEnd(svd->ip,V[i],i,V,work);
112: 
113:     VecScale(u,1.0/a);
114:     SlepcVecMAXPBY(V[i],1.0/a,-1.0/a,i,work,V);

116:     IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,i,PETSC_NULL,V,V[i],work,PETSC_NULL,&b);
117:     VecScale(V[i],1.0/b);
118: 
119:     SVDMatMult(svd,PETSC_FALSE,V[i],u_1);
120:     VecAXPY(u_1,-b,u);

122:     alpha[i-1] = a;
123:     beta[i-1] = b;
124:     temp = u;
125:     u = u_1;
126:     u_1 = temp;
127:   }
128:   SVDMatMult(svd,PETSC_TRUE,u,v);
129:   IPNormBegin(svd->ip,u,&a);
130:   IPMInnerProductBegin(svd->ip,v,n,V,work);
131:   IPNormEnd(svd->ip,u,&a);
132:   IPMInnerProductEnd(svd->ip,v,n,V,work);
133: 
134:   VecScale(u,1.0/a);
135:   SlepcVecMAXPBY(v,1.0/a,-1.0/a,n,work,V);

137:   IPOrthogonalizeCGS1(svd->ip,0,PETSC_NULL,n,PETSC_NULL,V,v,work,PETSC_NULL,&b);
138: 
139:   alpha[n-1] = a;
140:   beta[n-1] = b;
141:   return(0);
142: }

146: PetscErrorCode SVDSolve_Lanczos(SVD svd)
147: {
149:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS *)svd->data;
150:   PetscReal      *alpha,*beta,lastbeta;
151:   PetscScalar    *swork,*w,*Q,*PT;
152:   PetscInt       i,k,j,nv,ld,off;
153:   Vec            v,u=0,u_1=0;
154:   PetscBool      conv;
155: 
157:   /* allocate working space */
158:   DSGetLeadingDimension(svd->ds,&ld);
159:   PetscMalloc(sizeof(PetscScalar)*ld,&w);
160:   PetscMalloc(sizeof(PetscScalar)*svd->n,&swork);

162:   VecDuplicate(svd->V[0],&v);
163:   if (lanczos->oneside) {
164:     SVDMatGetVecs(svd,PETSC_NULL,&u);
165:     SVDMatGetVecs(svd,PETSC_NULL,&u_1);
166:   }
167: 
168:   /* normalize start vector */
169:   if (svd->nini==0) {
170:     SlepcVecSetRandom(svd->V[0],svd->rand);
171:   }
172:   VecNormalize(svd->V[0],PETSC_NULL);
173: 
174:   while (svd->reason == SVD_CONVERGED_ITERATING) {
175:     svd->its++;

177:     /* inner loop */
178:     nv = PetscMin(svd->nconv+svd->mpd,svd->n);
179:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
180:     beta = alpha + ld;
181:     if (lanczos->oneside) {
182:       SVDOneSideLanczos(svd,alpha,beta,svd->V,v,u,u_1,svd->nconv,nv,swork);
183:     } else {
184:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,v,svd->U,svd->nconv,nv,swork);
185:     }
186:     lastbeta = beta[nv-1];
187:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);

189:     /* compute SVD of bidiagonal matrix */
190:     DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
191:     DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
192:     DSSolve(svd->ds,w,PETSC_NULL);
193:     DSSort(svd->ds,w,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);

195:     /* compute error estimates */
196:     k = 0;
197:     conv = PETSC_TRUE;
198:     DSGetArray(svd->ds,DS_MAT_U,&Q);
199:     for (i=svd->nconv;i<nv;i++) {
200:       svd->sigma[i] = PetscRealPart(w[i]);
201:       svd->errest[i] = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
202:       if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
203:       if (conv) {
204:         if (svd->errest[i] < svd->tol) k++;
205:         else conv = PETSC_FALSE;
206:       }
207:     }
208: 
209:     /* check convergence */
210:     if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
211:     if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
212: 
213:     /* compute restart vector */
214:     DSGetArray(svd->ds,DS_MAT_VT,&PT);
215:     if (svd->reason == SVD_CONVERGED_ITERATING) {
216:       for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
217:       SlepcVecMAXPBY(v,0.0,1.0,nv-svd->nconv,swork,svd->V+svd->nconv);
218:     }
219: 
220:     /* compute converged singular vectors */
221:     off = svd->nconv+svd->nconv*ld;
222:     SlepcUpdateVectors(nv-svd->nconv,svd->V+svd->nconv,0,k,PT+off,ld,PETSC_TRUE);
223:     if (!lanczos->oneside) {
224:       SlepcUpdateVectors(nv-svd->nconv,svd->U+svd->nconv,0,k,Q+off,ld,PETSC_FALSE);
225:     }
226:     DSRestoreArray(svd->ds,DS_MAT_VT,&PT);
227:     DSRestoreArray(svd->ds,DS_MAT_U,&Q);
228: 
229:     /* copy restart vector from temporary space */
230:     if (svd->reason == SVD_CONVERGED_ITERATING) {
231:       VecCopy(v,svd->V[svd->nconv+k]);
232:     }
233: 
234:     svd->nconv += k;
235:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
236:   }
237: 
238:   /* free working space */
239:   VecDestroy(&v);
240:   VecDestroy(&u);
241:   VecDestroy(&u_1);
242:   PetscFree(w);
243:   PetscFree(swork);
244:   return(0);
245: }

249: PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd)
250: {
252:   PetscBool      set,val;
253:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS *)svd->data;

256:   PetscOptionsHead("SVD Lanczos Options");
257:   PetscOptionsBool("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
258:   if (set) {
259:     SVDLanczosSetOneSide(svd,val);
260:   }
261:   PetscOptionsTail();
262:   return(0);
263: }

265: EXTERN_C_BEGIN
268: PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
269: {
270:   SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;

273:   if (lanczos->oneside != oneside) {
274:     lanczos->oneside = oneside;
275:     svd->setupcalled = 0;
276:   }
277:   return(0);
278: }
279: EXTERN_C_END

283: /*@
284:    SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method 
285:    to be used is one-sided or two-sided.

287:    Logically Collective on SVD

289:    Input Parameters:
290: +  svd     - singular value solver
291: -  oneside - boolean flag indicating if the method is one-sided or not

293:    Options Database Key:
294: .  -svd_lanczos_oneside <boolean> - Indicates the boolean flag

296:    Note:
297:    By default, a two-sided variant is selected, which is sometimes slightly
298:    more robust. However, the one-sided variant is faster because it avoids 
299:    the orthogonalization associated to left singular vectors. It also saves
300:    the memory required for storing such vectors.

302:    Level: advanced

304: .seealso: SVDTRLanczosSetOneSide()
305: @*/
306: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
307: {

313:   PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
314:   return(0);
315: }

319: /*@
320:    SVDLanczosGetOneSide - Gets if the variant of the Lanczos method 
321:    to be used is one-sided or two-sided.

323:    Not Collective

325:    Input Parameters:
326: .  svd     - singular value solver

328:    Output Parameters:
329: .  oneside - boolean flag indicating if the method is one-sided or not

331:    Level: advanced

333: .seealso: SVDLanczosSetOneSide()
334: @*/
335: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
336: {

342:   PetscTryMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
343:   return(0);
344: }

346: EXTERN_C_BEGIN
349: PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
350: {
351:   SVD_LANCZOS *lanczos = (SVD_LANCZOS *)svd->data;

354:   *oneside = lanczos->oneside;
355:   return(0);
356: }
357: EXTERN_C_END

361: PetscErrorCode SVDReset_Lanczos(SVD svd)
362: {

366:   VecDestroyVecs(svd->n,&svd->U);
367:   return(0);
368: }

372: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
373: {

377:   PetscFree(svd->data);
378:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","",PETSC_NULL);
379:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosGetOneSide_C","",PETSC_NULL);
380:   return(0);
381: }

385: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
386: {
388:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS *)svd->data;

391:   PetscViewerASCIIPrintf(viewer,"  Lanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
392:   return(0);
393: }

395: EXTERN_C_BEGIN
398: PetscErrorCode SVDCreate_Lanczos(SVD svd)
399: {

403:   PetscNewLog(svd,SVD_LANCZOS,&svd->data);
404:   svd->ops->setup          = SVDSetUp_Lanczos;
405:   svd->ops->solve          = SVDSolve_Lanczos;
406:   svd->ops->destroy        = SVDDestroy_Lanczos;
407:   svd->ops->reset          = SVDReset_Lanczos;
408:   svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
409:   svd->ops->view           = SVDView_Lanczos;
410:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosSetOneSide_C","SVDLanczosSetOneSide_Lanczos",SVDLanczosSetOneSide_Lanczos);
411:   PetscObjectComposeFunctionDynamic((PetscObject)svd,"SVDLanczosGetOneSide_C","SVDLanczosGetOneSide_Lanczos",SVDLanczosGetOneSide_Lanczos);
412:   return(0);
413: }
414: EXTERN_C_END