Actual source code: trlanczos.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  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:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 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: static PetscBool  cited = PETSC_FALSE;
 45: static const char citation[] =
 46:   "@Article{slepc-svd,\n"
 47:   "   author = \"V. Hern{\\'a}ndez and J. E. Rom{\\'a}n and A. Tom{\\'a}s\",\n"
 48:   "   title = \"A robust and efficient parallel {SVD} solver based on restarted {Lanczos} bidiagonalization\",\n"
 49:   "   journal = \"Electron. Trans. Numer. Anal.\",\n"
 50:   "   volume = \"31\",\n"
 51:   "   pages = \"68--85\",\n"
 52:   "   year = \"2008\"\n"
 53:   "}\n";

 55: typedef struct {
 56:   PetscBool oneside;
 57: } SVD_TRLANCZOS;

 61: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
 62: {
 64:   PetscInt       N;

 67:   SVDMatGetSize(svd,NULL,&N);
 68:   SVDSetDimensions_Default(svd);
 69:   if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),1,"The value of ncv must not be larger than nev+mpd");
 70:   if (!svd->max_it) svd->max_it = PetscMax(N/svd->ncv,100);
 71:   svd->leftbasis = PETSC_TRUE;
 72:   SVDAllocateSolution(svd,1);
 73:   DSSetType(svd->ds,DSSVD);
 74:   DSSetCompact(svd->ds,PETSC_TRUE);
 75:   DSAllocate(svd->ds,svd->ncv);
 76:   return(0);
 77: }

 81: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
 82: {
 84:   PetscReal      a,b;
 85:   PetscInt       i,k=nconv+l;
 86:   Vec            ui,ui1,vi;

 89:   BVGetColumn(V,k,&vi);
 90:   BVGetColumn(U,k,&ui);
 91:   SVDMatMult(svd,PETSC_FALSE,vi,ui);
 92:   BVRestoreColumn(V,k,&vi);
 93:   BVRestoreColumn(U,k,&ui);
 94:   if (l>0) {
 95:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
 96:     BVMultColumn(U,-1.0,1.0,k,work);
 97:   }
 98:   BVNormColumn(U,k,NORM_2,&a);
 99:   BVScaleColumn(U,k,1.0/a);
100:   alpha[k] = a;

102:   for (i=k+1;i<n;i++) {
103:     BVGetColumn(V,i,&vi);
104:     BVGetColumn(U,i-1,&ui1);
105:     SVDMatMult(svd,PETSC_TRUE,ui1,vi);
106:     BVRestoreColumn(V,i,&vi);
107:     BVRestoreColumn(U,i-1,&ui1);
108:     BVOrthogonalizeColumn(V,i,NULL,&b,NULL);
109:     BVScaleColumn(V,i,1.0/b);
110:     beta[i-1] = b;

112:     BVGetColumn(V,i,&vi);
113:     BVGetColumn(U,i,&ui);
114:     SVDMatMult(svd,PETSC_FALSE,vi,ui);
115:     BVRestoreColumn(V,i,&vi);
116:     BVGetColumn(U,i-1,&ui1);
117:     VecAXPY(ui,-b,ui1);
118:     BVRestoreColumn(U,i-1,&ui1);
119:     BVRestoreColumn(U,i,&ui);
120:     BVNormColumn(U,i,NORM_2,&a);
121:     BVScaleColumn(U,i,1.0/a);
122:     alpha[i] = a;
123:   }

125:   BVGetColumn(V,n,&vi);
126:   BVGetColumn(U,n-1,&ui1);
127:   SVDMatMult(svd,PETSC_TRUE,ui1,vi);
128:   BVRestoreColumn(V,n,&vi);
129:   BVRestoreColumn(U,n-1,&ui1);
130:   BVOrthogonalizeColumn(V,n,NULL,&b,NULL);
131:   beta[n-1] = b;
132:   return(0);
133: }

137: /*
138:   Custom CGS orthogonalization, preprocess after first orthogonalization
139: */
140: static PetscErrorCode SVDOrthogonalizeCGS(BV V,PetscInt i,PetscScalar* h,PetscReal a,BVOrthogRefineType refine,PetscReal eta,PetscReal *norm)
141: {
143:   PetscReal      sum,onorm;
144:   PetscScalar    dot;
145:   PetscInt       j;

148:   switch (refine) {
149:   case BV_ORTHOG_REFINE_NEVER:
150:     BVNormColumn(V,i,NORM_2,norm);
151:     break;
152:   case BV_ORTHOG_REFINE_ALWAYS:
153:     BVSetActiveColumns(V,0,i);
154:     BVDotColumn(V,i,h);
155:     BVMultColumn(V,-1.0,1.0,i,h);
156:     BVNormColumn(V,i,NORM_2,norm);
157:     break;
158:   case BV_ORTHOG_REFINE_IFNEEDED:
159:     dot = h[i];
160:     onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
161:     sum = 0.0;
162:     for (j=0;j<i;j++) {
163:       sum += PetscRealPart(h[j] * PetscConj(h[j]));
164:     }
165:     *norm = PetscRealPart(dot)/(a*a) - sum;
166:     if (*norm>0.0) *norm = PetscSqrtReal(*norm);
167:     else {
168:       BVNormColumn(V,i,NORM_2,norm);
169:     }
170:     if (*norm < eta*onorm) {
171:       BVSetActiveColumns(V,0,i);
172:       BVDotColumn(V,i,h);
173:       BVMultColumn(V,-1.0,1.0,i,h);
174:       BVNormColumn(V,i,NORM_2,norm);
175:     }
176:     break;
177:   }
178:   return(0);
179: }

183: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
184: {
185:   PetscErrorCode     ierr;
186:   PetscReal          a,b,eta;
187:   PetscInt           i,j,k=nconv+l;
188:   Vec                ui,ui1,vi;
189:   BVOrthogRefineType refine;

192:   BVGetColumn(V,k,&vi);
193:   BVGetColumn(U,k,&ui);
194:   SVDMatMult(svd,PETSC_FALSE,vi,ui);
195:   BVRestoreColumn(V,k,&vi);
196:   BVRestoreColumn(U,k,&ui);
197:   if (l>0) {
198:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
199:     BVMultColumn(U,-1.0,1.0,k,work);
200:   }
201:   BVGetOrthogonalization(V,NULL,&refine,&eta,NULL);

203:   for (i=k+1;i<n;i++) {
204:     BVGetColumn(V,i,&vi);
205:     BVGetColumn(U,i-1,&ui1);
206:     SVDMatMult(svd,PETSC_TRUE,ui1,vi);
207:     BVRestoreColumn(V,i,&vi);
208:     BVRestoreColumn(U,i-1,&ui1);
209:     BVNormColumnBegin(U,i-1,NORM_2,&a);
210:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
211:       BVSetActiveColumns(V,0,i+1);
212:       BVGetColumn(V,i,&vi);
213:       BVDotVecBegin(V,vi,work);
214:     } else {
215:       BVSetActiveColumns(V,0,i);
216:       BVDotColumnBegin(V,i,work);
217:     }
218:     BVNormColumnEnd(U,i-1,NORM_2,&a);
219:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
220:       BVDotVecEnd(V,vi,work);
221:       BVRestoreColumn(V,i,&vi);
222:       BVSetActiveColumns(V,0,i);
223:     } else {
224:       BVDotColumnEnd(V,i,work);
225:     }

227:     BVScaleColumn(U,i-1,1.0/a);
228:     for (j=0;j<i;j++) work[j] = work[j] / a;
229:     BVMultColumn(V,-1.0,1.0/a,i,work);
230:     SVDOrthogonalizeCGS(V,i,work,a,refine,eta,&b);
231:     BVScaleColumn(V,i,1.0/b);

233:     BVGetColumn(V,i,&vi);
234:     BVGetColumn(U,i,&ui);
235:     BVGetColumn(U,i-1,&ui1);
236:     SVDMatMult(svd,PETSC_FALSE,vi,ui);
237:     VecAXPY(ui,-b,ui1);
238:     BVRestoreColumn(V,i,&vi);
239:     BVRestoreColumn(U,i,&ui);
240:     BVRestoreColumn(U,i-1,&ui1);

242:     alpha[i-1] = a;
243:     beta[i-1] = b;
244:   }

246:   BVGetColumn(V,n,&vi);
247:   BVGetColumn(U,n-1,&ui1);
248:   SVDMatMult(svd,PETSC_TRUE,ui1,vi);
249:   BVRestoreColumn(V,n,&vi);
250:   BVRestoreColumn(U,n-1,&ui1);

252:   BVNormColumnBegin(svd->U,n-1,NORM_2,&a);
253:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
254:     BVSetActiveColumns(V,0,n+1);
255:     BVGetColumn(V,n,&vi);
256:     BVDotVecBegin(V,vi,work);
257:   } else {
258:     BVSetActiveColumns(V,0,n);
259:     BVDotColumnBegin(V,n,work);
260:   }
261:   BVNormColumnEnd(svd->U,n-1,NORM_2,&a);
262:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
263:     BVDotVecEnd(V,vi,work);
264:     BVRestoreColumn(V,n,&vi);
265:   } else {
266:     BVDotColumnEnd(V,n,work);
267:   }

269:   BVScaleColumn(U,n-1,1.0/a);
270:   for (j=0;j<n;j++) work[j] = work[j] / a;
271:   BVMultColumn(V,-1.0,1.0/a,n,work);
272:   SVDOrthogonalizeCGS(V,n,work,a,refine,eta,&b);
273:   BVSetActiveColumns(V,nconv,n);
274:   alpha[n-1] = a;
275:   beta[n-1] = b;
276:   return(0);
277: }

281: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
282: {
284:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
285:   PetscReal      *alpha,*beta,lastbeta,norm,resnorm;
286:   PetscScalar    *Q,*swork=NULL,*w;
287:   PetscInt       i,k,l,nv,ld;
288:   Mat            U,VT;
289:   PetscBool      conv;
290:   BVOrthogType   orthog;

293:   PetscCitationsRegister(citation,&cited);
294:   /* allocate working space */
295:   DSGetLeadingDimension(svd->ds,&ld);
296:   BVGetOrthogonalization(svd->V,&orthog,NULL,NULL,NULL);
297:   PetscMalloc1(ld,&w);
298:   if (lanczos->oneside) {
299:     PetscMalloc1(svd->ncv+1,&swork);
300:   }

302:   /* normalize start vector */
303:   if (!svd->nini) {
304:     BVSetRandomColumn(svd->V,0);
305:     BVNormColumn(svd->V,0,NORM_2,&norm);
306:     BVScaleColumn(svd->V,0,1.0/norm);
307:   }

309:   l = 0;
310:   while (svd->reason == SVD_CONVERGED_ITERATING) {
311:     svd->its++;

313:     /* inner loop */
314:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
315:     BVSetActiveColumns(svd->V,svd->nconv,nv);
316:     BVSetActiveColumns(svd->U,svd->nconv,nv);
317:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
318:     beta = alpha + ld;
319:     if (lanczos->oneside) {
320:       if (orthog == BV_ORTHOG_MGS) {
321:         SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
322:       } else {
323:         SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
324:       }
325:     } else {
326:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv+l,nv);
327:     }
328:     lastbeta = beta[nv-1];
329:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
330:     BVScaleColumn(svd->V,nv,1.0/lastbeta);

332:     /* compute SVD of general matrix */
333:     DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);
334:     if (l==0) {
335:       DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
336:     } else {
337:       DSSetState(svd->ds,DS_STATE_RAW);
338:     }
339:     DSSolve(svd->ds,w,NULL);
340:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);

342:     /* compute error estimates */
343:     k = 0;
344:     conv = PETSC_TRUE;
345:     DSGetArray(svd->ds,DS_MAT_U,&Q);
346:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
347:     beta = alpha + ld;
348:     for (i=svd->nconv;i<nv;i++) {
349:       svd->sigma[i] = PetscRealPart(w[i]);
350:       beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta;
351:       resnorm = PetscAbsReal(beta[i]);
352:       (*svd->converged)(svd,svd->sigma[i],resnorm,&svd->errest[i],svd->convergedctx);
353:       if (conv) {
354:         if (svd->errest[i] < svd->tol) k++;
355:         else conv = PETSC_FALSE;
356:       }
357:     }
358:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
359:     DSRestoreArray(svd->ds,DS_MAT_U,&Q);

361:     /* check convergence and update l */
362:     (*svd->stopping)(svd,svd->its,svd->max_it,svd->nconv+k,svd->nsv,&svd->reason,svd->stoppingctx);
363:     if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
364:     else l = PetscMax((nv-svd->nconv-k)/2,0);

366:     /* compute converged singular vectors and restart vectors */
367:     DSGetMat(svd->ds,DS_MAT_VT,&VT);
368:     BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k+l);
369:     MatDestroy(&VT);
370:     DSGetMat(svd->ds,DS_MAT_U,&U);
371:     BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k+l);
372:     MatDestroy(&U);

374:     /* copy the last vector to be the next initial vector */
375:     if (svd->reason == SVD_CONVERGED_ITERATING) {
376:       BVCopyColumn(svd->V,nv,svd->nconv+k+l);
377:     }

379:     svd->nconv += k;
380:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
381:   }

383:   /* orthonormalize U columns in one side method */
384:   if (lanczos->oneside) {
385:     for (i=0;i<svd->nconv;i++) {
386:       BVOrthogonalizeColumn(svd->U,i,NULL,&norm,NULL);
387:       BVScaleColumn(svd->U,i,1.0/norm);
388:     }
389:   }

391:   /* free working space */
392:   PetscFree(w);
393:   if (swork) { PetscFree(swork); }
394:   return(0);
395: }

399: PetscErrorCode SVDSetFromOptions_TRLanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
400: {
402:   PetscBool      set,val;
403:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

406:   PetscOptionsHead(PetscOptionsObject,"SVD TRLanczos Options");
407:   PetscOptionsBool("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
408:   if (set) {
409:     SVDTRLanczosSetOneSide(svd,val);
410:   }
411:   PetscOptionsTail();
412:   return(0);
413: }

417: static PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
418: {
419:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

422:   lanczos->oneside = oneside;
423:   return(0);
424: }

428: /*@
429:    SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
430:    to be used is one-sided or two-sided.

432:    Logically Collective on SVD

434:    Input Parameters:
435: +  svd     - singular value solver
436: -  oneside - boolean flag indicating if the method is one-sided or not

438:    Options Database Key:
439: .  -svd_trlanczos_oneside <boolean> - Indicates the boolean flag

441:    Note:
442:    By default, a two-sided variant is selected, which is sometimes slightly
443:    more robust. However, the one-sided variant is faster because it avoids
444:    the orthogonalization associated to left singular vectors.

446:    Level: advanced

448: .seealso: SVDLanczosSetOneSide()
449: @*/
450: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
451: {

457:   PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
458:   return(0);
459: }

463: static PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
464: {
465:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

468:   *oneside = lanczos->oneside;
469:   return(0);
470: }

474: /*@
475:    SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
476:    to be used is one-sided or two-sided.

478:    Not Collective

480:    Input Parameters:
481: .  svd     - singular value solver

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

486:    Level: advanced

488: .seealso: SVDTRLanczosSetOneSide()
489: @*/
490: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
491: {

497:   PetscUseMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
498:   return(0);
499: }

503: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
504: {

508:   PetscFree(svd->data);
509:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",NULL);
510:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",NULL);
511:   return(0);
512: }

516: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
517: {
519:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
520:   PetscBool      isascii;

523:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
524:   if (isascii) {
525:     PetscViewerASCIIPrintf(viewer,"  TRLanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
526:   }
527:   return(0);
528: }

532: PETSC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD svd)
533: {
535:   SVD_TRLANCZOS  *ctx;

538:   PetscNewLog(svd,&ctx);
539:   svd->data = (void*)ctx;

541:   svd->ops->setup          = SVDSetUp_TRLanczos;
542:   svd->ops->solve          = SVDSolve_TRLanczos;
543:   svd->ops->destroy        = SVDDestroy_TRLanczos;
544:   svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
545:   svd->ops->view           = SVDView_TRLanczos;
546:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",SVDTRLanczosSetOneSide_TRLanczos);
547:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",SVDTRLanczosGetOneSide_TRLanczos);
548:   return(0);
549: }