Actual source code: gklanczos.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*

  3:    SLEPc singular value solver: "lanczos"

  5:    Method: Explicitly restarted Lanczos

  7:    Algorithm:

  9:        Golub-Kahan-Lanczos bidiagonalization with explicit 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:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 24:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 26:    This file is part of SLEPc.

 28:    SLEPc is free software: you can redistribute it and/or modify it under  the
 29:    terms of version 3 of the GNU Lesser General Public License as published by
 30:    the Free Software Foundation.

 32:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 33:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 34:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 35:    more details.

 37:    You  should have received a copy of the GNU Lesser General  Public  License
 38:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 39:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 40: */

 42: #include <slepc/private/svdimpl.h>                /*I "slepcsvd.h" I*/

 44: typedef struct {
 45:   PetscBool oneside;
 46: } SVD_LANCZOS;

 50: PetscErrorCode SVDSetUp_Lanczos(SVD svd)
 51: {
 53:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
 54:   PetscInt       N;

 57:   SVDMatGetSize(svd,NULL,&N);
 58:   SVDSetDimensions_Default(svd);
 59:   if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
 60:   if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
 61:   svd->leftbasis = PetscNot(lanczos->oneside);
 62:   SVDAllocateSolution(svd,1);
 63:   DSSetType(svd->ds,DSSVD);
 64:   DSSetCompact(svd->ds,PETSC_TRUE);
 65:   DSAllocate(svd->ds,svd->ncv);
 66:   return(0);
 67: }

 71: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt n)
 72: {
 74:   PetscInt       i;
 75:   Vec            u,v;

 78:   BVGetColumn(svd->V,k,&v);
 79:   BVGetColumn(svd->U,k,&u);
 80:   SVDMatMult(svd,PETSC_FALSE,v,u);
 81:   BVRestoreColumn(svd->V,k,&v);
 82:   BVRestoreColumn(svd->U,k,&u);
 83:   BVOrthogonalizeColumn(svd->U,k,NULL,alpha+k,NULL);
 84:   BVScaleColumn(U,k,1.0/alpha[k]);

 86:   for (i=k+1;i<n;i++) {
 87:     BVGetColumn(svd->V,i,&v);
 88:     BVGetColumn(svd->U,i-1,&u);
 89:     SVDMatMult(svd,PETSC_TRUE,u,v);
 90:     BVRestoreColumn(svd->V,i,&v);
 91:     BVRestoreColumn(svd->U,i-1,&u);
 92:     BVOrthogonalizeColumn(svd->V,i,NULL,beta+i-1,NULL);
 93:     BVScaleColumn(V,i,1.0/beta[i-1]);

 95:     BVGetColumn(svd->V,i,&v);
 96:     BVGetColumn(svd->U,i,&u);
 97:     SVDMatMult(svd,PETSC_FALSE,v,u);
 98:     BVRestoreColumn(svd->V,i,&v);
 99:     BVRestoreColumn(svd->U,i,&u);
100:     BVOrthogonalizeColumn(svd->U,i,NULL,alpha+i,NULL);
101:     BVScaleColumn(U,i,1.0/alpha[i]);
102:   }

104:   BVGetColumn(svd->V,n,&v);
105:   BVGetColumn(svd->U,n-1,&u);
106:   SVDMatMult(svd,PETSC_TRUE,u,v);
107:   BVRestoreColumn(svd->V,n,&v);
108:   BVRestoreColumn(svd->U,n-1,&u);
109:   BVOrthogonalizeColumn(svd->V,n,NULL,beta+n-1,NULL);
110:   return(0);
111: }

115: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
116: {
118:   PetscInt       i,bvl,bvk;
119:   PetscReal      a,b;
120:   Vec            z,temp;

123:   BVGetActiveColumns(V,&bvl,&bvk);
124:   BVGetColumn(V,k,&z);
125:   SVDMatMult(svd,PETSC_FALSE,z,u);
126:   BVRestoreColumn(V,k,&z);

128:   for (i=k+1;i<n;i++) {
129:     BVGetColumn(V,i,&z);
130:     SVDMatMult(svd,PETSC_TRUE,u,z);
131:     BVRestoreColumn(V,i,&z);
132:     VecNormBegin(u,NORM_2,&a);
133:     BVSetActiveColumns(V,0,i);
134:     BVDotColumnBegin(V,i,work);
135:     VecNormEnd(u,NORM_2,&a);
136:     BVDotColumnEnd(V,i,work);
137:     VecScale(u,1.0/a);
138:     BVMultColumn(V,-1.0/a,1.0/a,i,work);

140:     /* h = V^* z, z = z - V h  */
141:     BVDotColumn(V,i,work);
142:     BVMultColumn(V,-1.0,1.0,i,work);
143:     BVNormColumn(V,i,NORM_2,&b);
144:     BVScaleColumn(V,i,1.0/b);

146:     BVGetColumn(V,i,&z);
147:     SVDMatMult(svd,PETSC_FALSE,z,u_1);
148:     BVRestoreColumn(V,i,&z);
149:     VecAXPY(u_1,-b,u);
150:     alpha[i-1] = a;
151:     beta[i-1] = b;
152:     temp = u;
153:     u = u_1;
154:     u_1 = temp;
155:   }

157:   BVGetColumn(V,n,&z);
158:   SVDMatMult(svd,PETSC_TRUE,u,z);
159:   BVRestoreColumn(V,n,&z);
160:   VecNormBegin(u,NORM_2,&a);
161:   BVDotColumnBegin(V,n,work);
162:   VecNormEnd(u,NORM_2,&a);
163:   BVDotColumnEnd(V,n,work);
164:   VecScale(u,1.0/a);
165:   BVMultColumn(V,-1.0/a,1.0/a,n,work);

167:   /* h = V^* z, z = z - V h  */
168:   BVDotColumn(V,n,work);
169:   BVMultColumn(V,-1.0,1.0,n,work);
170:   BVNormColumn(V,i,NORM_2,&b);

172:   alpha[n-1] = a;
173:   beta[n-1] = b;
174:   BVSetActiveColumns(V,bvl,bvk);
175:   return(0);
176: }

180: PetscErrorCode SVDSolve_Lanczos(SVD svd)
181: {
183:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
184:   PetscReal      *alpha,*beta,lastbeta,norm,resnorm;
185:   PetscScalar    *swork,*w,*Q,*PT;
186:   PetscInt       i,k,j,nv,ld;
187:   Vec            u=0,u_1=0;
188:   Mat            U,VT;
189:   PetscBool      conv;

192:   /* allocate working space */
193:   DSGetLeadingDimension(svd->ds,&ld);
194:   PetscMalloc2(ld,&w,svd->ncv,&swork);

196:   if (lanczos->oneside) {
197:     SVDMatCreateVecs(svd,NULL,&u);
198:     SVDMatCreateVecs(svd,NULL,&u_1);
199:   }

201:   /* normalize start vector */
202:   if (!svd->nini) {
203:     BVSetRandomColumn(svd->V,0);
204:     BVNormColumn(svd->V,0,NORM_2,&norm);
205:     BVScaleColumn(svd->V,0,1.0/norm);
206:   }

208:   while (svd->reason == SVD_CONVERGED_ITERATING) {
209:     svd->its++;

211:     /* inner loop */
212:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
213:     BVSetActiveColumns(svd->V,svd->nconv,nv);
214:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
215:     beta = alpha + ld;
216:     if (lanczos->oneside) {
217:       SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
218:     } else {
219:       BVSetActiveColumns(svd->U,svd->nconv,nv);
220:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,nv);
221:     }
222:     lastbeta = beta[nv-1];
223:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);

225:     /* compute SVD of bidiagonal matrix */
226:     DSSetDimensions(svd->ds,nv,nv,svd->nconv,0);
227:     DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
228:     DSSolve(svd->ds,w,NULL);
229:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);

231:     /* compute error estimates */
232:     k = 0;
233:     conv = PETSC_TRUE;
234:     DSGetArray(svd->ds,DS_MAT_U,&Q);
235:     for (i=svd->nconv;i<nv;i++) {
236:       svd->sigma[i] = PetscRealPart(w[i]);
237:       resnorm = PetscAbsScalar(Q[nv-1+i*ld])*lastbeta;
238:       (*svd->converged)(svd,svd->sigma[i],resnorm,&svd->errest[i],svd->convergedctx);
239:       if (conv) {
240:         if (svd->errest[i] < svd->tol) k++;
241:         else conv = PETSC_FALSE;
242:       }
243:     }
244:     DSRestoreArray(svd->ds,DS_MAT_U,&Q);

246:     /* check convergence */
247:     (*svd->stopping)(svd,svd->its,svd->max_it,svd->nconv+k,svd->nsv,&svd->reason,svd->stoppingctx);

249:     /* compute restart vector */
250:     DSGetArray(svd->ds,DS_MAT_VT,&PT);
251:     if (svd->reason == SVD_CONVERGED_ITERATING) {
252:       for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PT[k+svd->nconv+j*ld];
253:       BVMultColumn(svd->V,1.0,0.0,nv,swork);
254:     }
255:     DSRestoreArray(svd->ds,DS_MAT_VT,&PT);

257:     /* compute converged singular vectors */
258:     DSGetMat(svd->ds,DS_MAT_VT,&VT);
259:     BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k);
260:     MatDestroy(&VT);
261:     if (!lanczos->oneside) {
262:       DSGetMat(svd->ds,DS_MAT_U,&U);
263:       BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k);
264:       MatDestroy(&U);
265:     }

267:     /* copy restart vector from the last column */
268:     if (svd->reason == SVD_CONVERGED_ITERATING) {
269:       BVCopyColumn(svd->V,nv,svd->nconv+k);
270:     }

272:     svd->nconv += k;
273:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
274:   }

276:   /* free working space */
277:   VecDestroy(&u);
278:   VecDestroy(&u_1);
279:   PetscFree2(w,swork);
280:   return(0);
281: }

285: PetscErrorCode SVDSetFromOptions_Lanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
286: {
288:   PetscBool      set,val;
289:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;

292:   PetscOptionsHead(PetscOptionsObject,"SVD Lanczos Options");
293:   PetscOptionsBool("-svd_lanczos_oneside","Lanczos one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
294:   if (set) {
295:     SVDLanczosSetOneSide(svd,val);
296:   }
297:   PetscOptionsTail();
298:   return(0);
299: }

303: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
304: {
305:   SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;

308:   if (lanczos->oneside != oneside) {
309:     lanczos->oneside = oneside;
310:     svd->state = SVD_STATE_INITIAL;
311:   }
312:   return(0);
313: }

317: /*@
318:    SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
319:    to be used is one-sided or two-sided.

321:    Logically Collective on SVD

323:    Input Parameters:
324: +  svd     - singular value solver
325: -  oneside - boolean flag indicating if the method is one-sided or not

327:    Options Database Key:
328: .  -svd_lanczos_oneside <boolean> - Indicates the boolean flag

330:    Note:
331:    By default, a two-sided variant is selected, which is sometimes slightly
332:    more robust. However, the one-sided variant is faster because it avoids
333:    the orthogonalization associated to left singular vectors. It also saves
334:    the memory required for storing such vectors.

336:    Level: advanced

338: .seealso: SVDTRLanczosSetOneSide()
339: @*/
340: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
341: {

347:   PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
348:   return(0);
349: }

353: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
354: {
355:   SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;

358:   *oneside = lanczos->oneside;
359:   return(0);
360: }

364: /*@
365:    SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
366:    to be used is one-sided or two-sided.

368:    Not Collective

370:    Input Parameters:
371: .  svd     - singular value solver

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

376:    Level: advanced

378: .seealso: SVDLanczosSetOneSide()
379: @*/
380: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
381: {

387:   PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
388:   return(0);
389: }

393: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
394: {

398:   PetscFree(svd->data);
399:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
400:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
401:   return(0);
402: }

406: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
407: {
409:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
410:   PetscBool      isascii;

413:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
414:   if (isascii) {
415:     PetscViewerASCIIPrintf(viewer,"  Lanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
416:   }
417:   return(0);
418: }

422: PETSC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
423: {
425:   SVD_LANCZOS    *ctx;

428:   PetscNewLog(svd,&ctx);
429:   svd->data = (void*)ctx;

431:   svd->ops->setup          = SVDSetUp_Lanczos;
432:   svd->ops->solve          = SVDSolve_Lanczos;
433:   svd->ops->destroy        = SVDDestroy_Lanczos;
434:   svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
435:   svd->ops->view           = SVDView_Lanczos;
436:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
437:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
438:   return(0);
439: }