Actual source code: trlanczos.c

slepc-3.6.1 2015-09-03
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-2015, 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_TRLANCZOS;

 50: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
 51: {
 53:   PetscInt       N;

 56:   SVDMatGetSize(svd,NULL,&N);
 57:   SVDSetDimensions_Default(svd);
 58:   if (svd->ncv>svd->nsv+svd->mpd) SETERRQ(PetscObjectComm((PetscObject)svd),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:   svd->leftbasis = PETSC_TRUE;
 61:   SVDAllocateSolution(svd,1);
 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,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
 71: {
 73:   PetscReal      a,b;
 74:   PetscInt       i,k=nconv+l;
 75:   Vec            ui,ui1,vi;

 78:   BVGetColumn(V,k,&vi);
 79:   BVGetColumn(U,k,&ui);
 80:   SVDMatMult(svd,PETSC_FALSE,vi,ui);
 81:   BVRestoreColumn(V,k,&vi);
 82:   BVRestoreColumn(U,k,&ui);
 83:   if (l>0) {
 84:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
 85:     BVMultColumn(U,-1.0,1.0,k,work);
 86:   }
 87:   BVNormColumn(U,k,NORM_2,&a);
 88:   BVScaleColumn(U,k,1.0/a);
 89:   alpha[k] = a;

 91:   for (i=k+1;i<n;i++) {
 92:     BVGetColumn(V,i,&vi);
 93:     BVGetColumn(U,i-1,&ui1);
 94:     SVDMatMult(svd,PETSC_TRUE,ui1,vi);
 95:     BVRestoreColumn(V,i,&vi);
 96:     BVRestoreColumn(U,i-1,&ui1);
 97:     BVOrthogonalizeColumn(V,i,NULL,&b,NULL);
 98:     BVScaleColumn(V,i,1.0/b);
 99:     beta[i-1] = b;

101:     BVGetColumn(V,i,&vi);
102:     BVGetColumn(U,i,&ui);
103:     SVDMatMult(svd,PETSC_FALSE,vi,ui);
104:     BVRestoreColumn(V,i,&vi);
105:     BVGetColumn(U,i-1,&ui1);
106:     VecAXPY(ui,-b,ui1);
107:     BVRestoreColumn(U,i-1,&ui1);
108:     BVRestoreColumn(U,i,&ui);
109:     BVNormColumn(U,i,NORM_2,&a);
110:     BVScaleColumn(U,i,1.0/a);
111:     alpha[i] = a;
112:   }

114:   BVGetColumn(V,n,&vi);
115:   BVGetColumn(U,n-1,&ui1);
116:   SVDMatMult(svd,PETSC_TRUE,ui1,vi);
117:   BVRestoreColumn(V,n,&vi);
118:   BVRestoreColumn(U,n-1,&ui1);
119:   BVOrthogonalizeColumn(V,n,NULL,&b,NULL);
120:   beta[n-1] = b;
121:   return(0);
122: }

126: /*
127:   Custom CGS orthogonalization, preprocess after first orthogonalization
128: */
129: static PetscErrorCode SVDOrthogonalizeCGS(BV V,PetscInt i,PetscScalar* h,PetscReal a,BVOrthogRefineType refine,PetscReal eta,PetscReal *norm)
130: {
132:   PetscReal      sum,onorm;
133:   PetscScalar    dot;
134:   PetscInt       j;

137:   switch (refine) {
138:   case BV_ORTHOG_REFINE_NEVER:
139:     BVNormColumn(V,i,NORM_2,norm);
140:     break;
141:   case BV_ORTHOG_REFINE_ALWAYS:
142:     BVSetActiveColumns(V,0,i);
143:     BVDotColumn(V,i,h);
144:     BVMultColumn(V,-1.0,1.0,i,h);
145:     BVNormColumn(V,i,NORM_2,norm);
146:     break;
147:   case BV_ORTHOG_REFINE_IFNEEDED:
148:     dot = h[i];
149:     onorm = PetscSqrtReal(PetscRealPart(dot)) / a;
150:     sum = 0.0;
151:     for (j=0;j<i;j++) {
152:       sum += PetscRealPart(h[j] * PetscConj(h[j]));
153:     }
154:     *norm = PetscRealPart(dot)/(a*a) - sum;
155:     if (*norm>0.0) *norm = PetscSqrtReal(*norm);
156:     else {
157:       BVNormColumn(V,i,NORM_2,norm);
158:     }
159:     if (*norm < eta*onorm) {
160:       BVSetActiveColumns(V,0,i);
161:       BVDotColumn(V,i,h);
162:       BVMultColumn(V,-1.0,1.0,i,h);
163:       BVNormColumn(V,i,NORM_2,norm);
164:     }
165:     break;
166:   }
167:   return(0);
168: }

172: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
173: {
174:   PetscErrorCode     ierr;
175:   PetscReal          a,b,eta;
176:   PetscInt           i,j,k=nconv+l;
177:   Vec                ui,ui1,vi;
178:   BVOrthogRefineType refine;

181:   BVGetColumn(V,k,&vi);
182:   BVGetColumn(U,k,&ui);
183:   SVDMatMult(svd,PETSC_FALSE,vi,ui);
184:   BVRestoreColumn(V,k,&vi);
185:   BVRestoreColumn(U,k,&ui);
186:   if (l>0) {
187:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
188:     BVMultColumn(U,-1.0,1.0,k,work);
189:   }
190:   BVGetOrthogonalization(V,NULL,&refine,&eta,NULL);

192:   for (i=k+1;i<n;i++) {
193:     BVGetColumn(V,i,&vi);
194:     BVGetColumn(U,i-1,&ui1);
195:     SVDMatMult(svd,PETSC_TRUE,ui1,vi);
196:     BVRestoreColumn(V,i,&vi);
197:     BVRestoreColumn(U,i-1,&ui1);
198:     BVNormColumnBegin(U,i-1,NORM_2,&a);
199:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
200:       BVSetActiveColumns(V,0,i+1);
201:       BVGetColumn(V,i,&vi);
202:       BVDotVecBegin(V,vi,work);
203:     } else {
204:       BVSetActiveColumns(V,0,i);
205:       BVDotColumnBegin(V,i,work);
206:     }
207:     BVNormColumnEnd(U,i-1,NORM_2,&a);
208:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
209:       BVDotVecEnd(V,vi,work);
210:       BVRestoreColumn(V,i,&vi);
211:       BVSetActiveColumns(V,0,i);
212:     } else {
213:       BVDotColumnEnd(V,i,work);
214:     }

216:     BVScaleColumn(U,i-1,1.0/a);
217:     for (j=0;j<i;j++) work[j] = work[j] / a;
218:     BVMultColumn(V,-1.0,1.0/a,i,work);
219:     SVDOrthogonalizeCGS(V,i,work,a,refine,eta,&b);
220:     BVScaleColumn(V,i,1.0/b);

222:     BVGetColumn(V,i,&vi);
223:     BVGetColumn(U,i,&ui);
224:     BVGetColumn(U,i-1,&ui1);
225:     SVDMatMult(svd,PETSC_FALSE,vi,ui);
226:     VecAXPY(ui,-b,ui1);
227:     BVRestoreColumn(V,i,&vi);
228:     BVRestoreColumn(U,i,&ui);
229:     BVRestoreColumn(U,i-1,&ui1);

231:     alpha[i-1] = a;
232:     beta[i-1] = b;
233:   }

235:   BVGetColumn(V,n,&vi);
236:   BVGetColumn(U,n-1,&ui1);
237:   SVDMatMult(svd,PETSC_TRUE,ui1,vi);
238:   BVRestoreColumn(V,n,&vi);
239:   BVRestoreColumn(U,n-1,&ui1);

241:   BVNormColumnBegin(svd->U,n-1,NORM_2,&a);
242:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
243:     BVSetActiveColumns(V,0,n+1);
244:     BVGetColumn(V,n,&vi);
245:     BVDotVecBegin(V,vi,work);
246:   } else {
247:     BVSetActiveColumns(V,0,n);
248:     BVDotColumnBegin(V,n,work);
249:   }
250:   BVNormColumnEnd(svd->U,n-1,NORM_2,&a);
251:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
252:     BVDotVecEnd(V,vi,work);
253:     BVRestoreColumn(V,n,&vi);
254:   } else {
255:     BVDotColumnEnd(V,n,work);
256:   }

258:   BVScaleColumn(U,n-1,1.0/a);
259:   for (j=0;j<n;j++) work[j] = work[j] / a;
260:   BVMultColumn(V,-1.0,1.0/a,n,work);
261:   SVDOrthogonalizeCGS(V,n,work,a,refine,eta,&b);
262:   BVSetActiveColumns(V,nconv,n);
263:   alpha[n-1] = a;
264:   beta[n-1] = b;
265:   return(0);
266: }

270: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
271: {
273:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
274:   PetscReal      *alpha,*beta,lastbeta,norm;
275:   PetscScalar    *Q,*swork=NULL,*w;
276:   PetscInt       i,k,l,nv,ld;
277:   Mat            U,VT;
278:   PetscBool      conv;
279:   BVOrthogType   orthog;

282:   /* allocate working space */
283:   DSGetLeadingDimension(svd->ds,&ld);
284:   BVGetOrthogonalization(svd->V,&orthog,NULL,NULL,NULL);
285:   PetscMalloc1(ld,&w);
286:   if (lanczos->oneside) {
287:     PetscMalloc1(svd->ncv+1,&swork);
288:   }

290:   /* normalize start vector */
291:   if (!svd->nini) {
292:     BVSetRandomColumn(svd->V,0,svd->rand);
293:     BVNormColumn(svd->V,0,NORM_2,&norm);
294:     BVScaleColumn(svd->V,0,1.0/norm);
295:   }

297:   l = 0;
298:   while (svd->reason == SVD_CONVERGED_ITERATING) {
299:     svd->its++;

301:     /* inner loop */
302:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
303:     BVSetActiveColumns(svd->V,svd->nconv,nv);
304:     BVSetActiveColumns(svd->U,svd->nconv,nv);
305:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
306:     beta = alpha + ld;
307:     if (lanczos->oneside) {
308:       if (orthog == BV_ORTHOG_MGS) {
309:         SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
310:       } else {
311:         SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
312:       }
313:     } else {
314:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv+l,nv);
315:     }
316:     lastbeta = beta[nv-1];
317:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
318:     BVScaleColumn(svd->V,nv,1.0/lastbeta);

320:     /* compute SVD of general matrix */
321:     DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);
322:     if (l==0) {
323:       DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
324:     } else {
325:       DSSetState(svd->ds,DS_STATE_RAW);
326:     }
327:     DSSolve(svd->ds,w,NULL);
328:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);

330:     /* compute error estimates */
331:     k = 0;
332:     conv = PETSC_TRUE;
333:     DSGetArray(svd->ds,DS_MAT_U,&Q);
334:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
335:     beta = alpha + ld;
336:     for (i=svd->nconv;i<nv;i++) {
337:       svd->sigma[i] = PetscRealPart(w[i]);
338:       beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta;
339:       svd->errest[i] = PetscAbsReal(beta[i]);
340:       if (svd->sigma[i] > svd->tol) svd->errest[i] /= svd->sigma[i];
341:       if (conv) {
342:         if (svd->errest[i] < svd->tol) k++;
343:         else conv = PETSC_FALSE;
344:       }
345:     }
346:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
347:     DSRestoreArray(svd->ds,DS_MAT_U,&Q);

349:     /* check convergence and update l */
350:     if (svd->its >= svd->max_it) svd->reason = SVD_DIVERGED_ITS;
351:     if (svd->nconv+k >= svd->nsv) svd->reason = SVD_CONVERGED_TOL;
352:     if (svd->reason != SVD_CONVERGED_ITERATING) l = 0;
353:     else l = PetscMax((nv-svd->nconv-k)/2,0);

355:     /* compute converged singular vectors and restart vectors */
356:     DSGetMat(svd->ds,DS_MAT_VT,&VT);
357:     BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k+l);
358:     MatDestroy(&VT);
359:     DSGetMat(svd->ds,DS_MAT_U,&U);
360:     BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k+l);
361:     MatDestroy(&U);

363:     /* copy the last vector to be the next initial vector */
364:     if (svd->reason == SVD_CONVERGED_ITERATING) {
365:       BVCopyColumn(svd->V,nv,svd->nconv+k+l);
366:     }

368:     svd->nconv += k;
369:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
370:   }

372:   /* orthonormalize U columns in one side method */
373:   if (lanczos->oneside) {
374:     for (i=0;i<svd->nconv;i++) {
375:       BVOrthogonalizeColumn(svd->U,i,NULL,&norm,NULL);
376:       BVScaleColumn(svd->U,i,1.0/norm);
377:     }
378:   }

380:   /* free working space */
381:   PetscFree(w);
382:   if (swork) { PetscFree(swork); }
383:   return(0);
384: }

388: PetscErrorCode SVDSetFromOptions_TRLanczos(PetscOptions *PetscOptionsObject,SVD svd)
389: {
391:   PetscBool      set,val;
392:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

395:   PetscOptionsHead(PetscOptionsObject,"SVD TRLanczos Options");
396:   PetscOptionsBool("-svd_trlanczos_oneside","Lanczos one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
397:   if (set) {
398:     SVDTRLanczosSetOneSide(svd,val);
399:   }
400:   PetscOptionsTail();
401:   return(0);
402: }

406: static PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
407: {
408:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

411:   lanczos->oneside = oneside;
412:   return(0);
413: }

417: /*@
418:    SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
419:    to be used is one-sided or two-sided.

421:    Logically Collective on SVD

423:    Input Parameters:
424: +  svd     - singular value solver
425: -  oneside - boolean flag indicating if the method is one-sided or not

427:    Options Database Key:
428: .  -svd_trlanczos_oneside <boolean> - Indicates the boolean flag

430:    Note:
431:    By default, a two-sided variant is selected, which is sometimes slightly
432:    more robust. However, the one-sided variant is faster because it avoids
433:    the orthogonalization associated to left singular vectors.

435:    Level: advanced

437: .seealso: SVDLanczosSetOneSide()
438: @*/
439: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
440: {

446:   PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
447:   return(0);
448: }

452: /*@
453:    SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
454:    to be used is one-sided or two-sided.

456:    Not Collective

458:    Input Parameters:
459: .  svd     - singular value solver

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

464:    Level: advanced

466: .seealso: SVDTRLanczosSetOneSide()
467: @*/
468: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
469: {

475:   PetscTryMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
476:   return(0);
477: }

481: static PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
482: {
483:   SVD_TRLANCZOS    *lanczos = (SVD_TRLANCZOS*)svd->data;

486:   *oneside = lanczos->oneside;
487:   return(0);
488: }

492: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
493: {

497:   PetscFree(svd->data);
498:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",NULL);
499:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",NULL);
500:   return(0);
501: }

505: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
506: {
508:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

511:   PetscViewerASCIIPrintf(viewer,"  TRLanczos: %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
512:   return(0);
513: }

517: PETSC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD svd)
518: {
520:   SVD_TRLANCZOS  *ctx;

523:   PetscNewLog(svd,&ctx);
524:   svd->data = (void*)ctx;

526:   svd->ops->setup          = SVDSetUp_TRLanczos;
527:   svd->ops->solve          = SVDSolve_TRLanczos;
528:   svd->ops->destroy        = SVDDestroy_TRLanczos;
529:   svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
530:   svd->ops->view           = SVDView_TRLanczos;
531:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",SVDTRLanczosSetOneSide_TRLanczos);
532:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",SVDTRLanczosGetOneSide_TRLanczos);
533:   return(0);
534: }