Actual source code: trlanczos.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc singular value solver: "trlanczos"

 13:    Method: Thick-restart Lanczos

 15:    Algorithm:

 17:        Golub-Kahan-Lanczos bidiagonalization with thick-restart.

 19:    References:

 21:        [1] G.H. Golub and W. Kahan, "Calculating the singular values
 22:            and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
 23:            B 2:205-224, 1965.

 25:        [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
 26:            efficient parallel SVD solver based on restarted Lanczos
 27:            bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
 28:            2008.
 29: */

 31: #include <slepc/private/svdimpl.h>

 33: static PetscBool  cited = PETSC_FALSE;
 34: static const char citation[] =
 35:   "@Article{slepc-svd,\n"
 36:   "   author = \"V. Hern{\\'a}ndez and J. E. Rom{\\'a}n and A. Tom{\\'a}s\",\n"
 37:   "   title = \"A robust and efficient parallel {SVD} solver based on restarted {Lanczos} bidiagonalization\",\n"
 38:   "   journal = \"Electron. Trans. Numer. Anal.\",\n"
 39:   "   volume = \"31\",\n"
 40:   "   pages = \"68--85\",\n"
 41:   "   year = \"2008\"\n"
 42:   "}\n";

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

 48: PetscErrorCode SVDSetUp_TRLanczos(SVD svd)
 49: {
 51:   PetscInt       N;

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

 67: static PetscErrorCode SVDOneSideTRLanczosMGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
 68: {
 70:   PetscReal      a,b;
 71:   PetscInt       i,k=nconv+l;
 72:   Vec            ui,ui1,vi;

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

 88:   for (i=k+1;i<n;i++) {
 89:     BVGetColumn(V,i,&vi);
 90:     BVGetColumn(U,i-1,&ui1);
 91:     MatMult(svd->AT,ui1,vi);
 92:     BVRestoreColumn(V,i,&vi);
 93:     BVRestoreColumn(U,i-1,&ui1);
 94:     BVOrthogonalizeColumn(V,i,NULL,&b,NULL);
 95:     BVScaleColumn(V,i,1.0/b);
 96:     beta[i-1] = b;

 98:     BVGetColumn(V,i,&vi);
 99:     BVGetColumn(U,i,&ui);
100:     MatMult(svd->A,vi,ui);
101:     BVRestoreColumn(V,i,&vi);
102:     BVGetColumn(U,i-1,&ui1);
103:     VecAXPY(ui,-b,ui1);
104:     BVRestoreColumn(U,i-1,&ui1);
105:     BVRestoreColumn(U,i,&ui);
106:     BVNormColumn(U,i,NORM_2,&a);
107:     BVScaleColumn(U,i,1.0/a);
108:     alpha[i] = a;
109:   }

111:   BVGetColumn(V,n,&vi);
112:   BVGetColumn(U,n-1,&ui1);
113:   MatMult(svd->AT,ui1,vi);
114:   BVRestoreColumn(V,n,&vi);
115:   BVRestoreColumn(U,n-1,&ui1);
116:   BVOrthogonalizeColumn(V,n,NULL,&b,NULL);
117:   beta[n-1] = b;
118:   return(0);
119: }

121: /*
122:   Custom CGS orthogonalization, preprocess after first orthogonalization
123: */
124: static PetscErrorCode SVDOrthogonalizeCGS(BV V,PetscInt i,PetscScalar* h,PetscReal a,BVOrthogRefineType refine,PetscReal eta,PetscReal *norm)
125: {
127:   PetscReal      sum,onorm;
128:   PetscScalar    dot;
129:   PetscInt       j;

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

165: static PetscErrorCode SVDOneSideTRLanczosCGS(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt nconv,PetscInt l,PetscInt n,PetscScalar* work)
166: {
167:   PetscErrorCode     ierr;
168:   PetscReal          a,b,eta;
169:   PetscInt           i,j,k=nconv+l;
170:   Vec                ui,ui1,vi;
171:   BVOrthogRefineType refine;

174:   BVGetColumn(V,k,&vi);
175:   BVGetColumn(U,k,&ui);
176:   MatMult(svd->A,vi,ui);
177:   BVRestoreColumn(V,k,&vi);
178:   BVRestoreColumn(U,k,&ui);
179:   if (l>0) {
180:     for (i=0;i<l;i++) work[i]=beta[i+nconv];
181:     BVMultColumn(U,-1.0,1.0,k,work);
182:   }
183:   BVGetOrthogonalization(V,NULL,&refine,&eta,NULL);

185:   for (i=k+1;i<n;i++) {
186:     BVGetColumn(V,i,&vi);
187:     BVGetColumn(U,i-1,&ui1);
188:     MatMult(svd->AT,ui1,vi);
189:     BVRestoreColumn(V,i,&vi);
190:     BVRestoreColumn(U,i-1,&ui1);
191:     BVNormColumnBegin(U,i-1,NORM_2,&a);
192:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
193:       BVSetActiveColumns(V,0,i+1);
194:       BVGetColumn(V,i,&vi);
195:       BVDotVecBegin(V,vi,work);
196:     } else {
197:       BVSetActiveColumns(V,0,i);
198:       BVDotColumnBegin(V,i,work);
199:     }
200:     BVNormColumnEnd(U,i-1,NORM_2,&a);
201:     if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
202:       BVDotVecEnd(V,vi,work);
203:       BVRestoreColumn(V,i,&vi);
204:       BVSetActiveColumns(V,0,i);
205:     } else {
206:       BVDotColumnEnd(V,i,work);
207:     }

209:     BVScaleColumn(U,i-1,1.0/a);
210:     for (j=0;j<i;j++) work[j] = work[j] / a;
211:     BVMultColumn(V,-1.0,1.0/a,i,work);
212:     SVDOrthogonalizeCGS(V,i,work,a,refine,eta,&b);
213:     BVScaleColumn(V,i,1.0/b);
214:     if (PetscAbsReal(b)<10*PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)svd),1,"Recurrence generated a zero vector; use a two-sided variant");

216:     BVGetColumn(V,i,&vi);
217:     BVGetColumn(U,i,&ui);
218:     BVGetColumn(U,i-1,&ui1);
219:     MatMult(svd->A,vi,ui);
220:     VecAXPY(ui,-b,ui1);
221:     BVRestoreColumn(V,i,&vi);
222:     BVRestoreColumn(U,i,&ui);
223:     BVRestoreColumn(U,i-1,&ui1);

225:     alpha[i-1] = a;
226:     beta[i-1] = b;
227:   }

229:   BVGetColumn(V,n,&vi);
230:   BVGetColumn(U,n-1,&ui1);
231:   MatMult(svd->AT,ui1,vi);
232:   BVRestoreColumn(V,n,&vi);
233:   BVRestoreColumn(U,n-1,&ui1);

235:   BVNormColumnBegin(svd->U,n-1,NORM_2,&a);
236:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
237:     BVSetActiveColumns(V,0,n+1);
238:     BVGetColumn(V,n,&vi);
239:     BVDotVecBegin(V,vi,work);
240:   } else {
241:     BVSetActiveColumns(V,0,n);
242:     BVDotColumnBegin(V,n,work);
243:   }
244:   BVNormColumnEnd(svd->U,n-1,NORM_2,&a);
245:   if (refine == BV_ORTHOG_REFINE_IFNEEDED) {
246:     BVDotVecEnd(V,vi,work);
247:     BVRestoreColumn(V,n,&vi);
248:   } else {
249:     BVDotColumnEnd(V,n,work);
250:   }

252:   BVScaleColumn(U,n-1,1.0/a);
253:   for (j=0;j<n;j++) work[j] = work[j] / a;
254:   BVMultColumn(V,-1.0,1.0/a,n,work);
255:   SVDOrthogonalizeCGS(V,n,work,a,refine,eta,&b);
256:   BVSetActiveColumns(V,nconv,n);
257:   alpha[n-1] = a;
258:   beta[n-1] = b;
259:   return(0);
260: }

262: PetscErrorCode SVDSolve_TRLanczos(SVD svd)
263: {
265:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
266:   PetscReal      *alpha,*beta,lastbeta,resnorm;
267:   PetscScalar    *Q,*swork=NULL,*w;
268:   PetscInt       i,k,l,nv,ld;
269:   Mat            U,VT;
270:   PetscBool      conv;
271:   BVOrthogType   orthog;

274:   PetscCitationsRegister(citation,&cited);
275:   /* allocate working space */
276:   DSGetLeadingDimension(svd->ds,&ld);
277:   BVGetOrthogonalization(svd->V,&orthog,NULL,NULL,NULL);
278:   PetscMalloc1(ld,&w);
279:   if (lanczos->oneside) {
280:     PetscMalloc1(svd->ncv+1,&swork);
281:   }

283:   /* normalize start vector */
284:   if (!svd->nini) {
285:     BVSetRandomColumn(svd->V,0);
286:     BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL);
287:   }

289:   l = 0;
290:   while (svd->reason == SVD_CONVERGED_ITERATING) {
291:     svd->its++;

293:     /* inner loop */
294:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
295:     BVSetActiveColumns(svd->V,svd->nconv,nv);
296:     BVSetActiveColumns(svd->U,svd->nconv,nv);
297:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
298:     beta = alpha + ld;
299:     if (lanczos->oneside) {
300:       if (orthog == BV_ORTHOG_MGS) {
301:         SVDOneSideTRLanczosMGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
302:       } else {
303:         SVDOneSideTRLanczosCGS(svd,alpha,beta,svd->V,svd->U,svd->nconv,l,nv,swork);
304:       }
305:     } else {
306:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv+l,nv);
307:     }
308:     lastbeta = beta[nv-1];
309:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
310:     BVScaleColumn(svd->V,nv,1.0/lastbeta);

312:     /* compute SVD of general matrix */
313:     DSSetDimensions(svd->ds,nv,nv,svd->nconv,svd->nconv+l);
314:     if (l==0) {
315:       DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
316:     } else {
317:       DSSetState(svd->ds,DS_STATE_RAW);
318:     }
319:     DSSolve(svd->ds,w,NULL);
320:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
321:     DSSynchronize(svd->ds,w,NULL);

323:     /* compute error estimates */
324:     k = 0;
325:     conv = PETSC_TRUE;
326:     DSGetArray(svd->ds,DS_MAT_U,&Q);
327:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
328:     beta = alpha + ld;
329:     for (i=svd->nconv;i<nv;i++) {
330:       svd->sigma[i] = PetscRealPart(w[i]);
331:       beta[i] = PetscRealPart(Q[nv-1+i*ld])*lastbeta;
332:       resnorm = PetscAbsReal(beta[i]);
333:       (*svd->converged)(svd,svd->sigma[i],resnorm,&svd->errest[i],svd->convergedctx);
334:       if (conv) {
335:         if (svd->errest[i] < svd->tol) k++;
336:         else conv = PETSC_FALSE;
337:       }
338:     }
339:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
340:     DSRestoreArray(svd->ds,DS_MAT_U,&Q);
341:     if (svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it) k = svd->nsv;

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

348:     /* compute converged singular vectors and restart vectors */
349:     DSGetMat(svd->ds,DS_MAT_VT,&VT);
350:     BVMultInPlaceTranspose(svd->V,VT,svd->nconv,svd->nconv+k+l);
351:     MatDestroy(&VT);
352:     DSGetMat(svd->ds,DS_MAT_U,&U);
353:     BVMultInPlace(svd->U,U,svd->nconv,svd->nconv+k+l);
354:     MatDestroy(&U);

356:     /* copy the last vector to be the next initial vector */
357:     if (svd->reason == SVD_CONVERGED_ITERATING) {
358:       BVCopyColumn(svd->V,nv,svd->nconv+k+l);
359:     }

361:     svd->nconv += k;
362:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
363:   }

365:   /* orthonormalize U columns in one side method */
366:   if (lanczos->oneside) {
367:     for (i=0;i<svd->nconv;i++) {
368:       BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,NULL,NULL);
369:     }
370:   }

372:   /* free working space */
373:   PetscFree(w);
374:   if (swork) { PetscFree(swork); }
375:   return(0);
376: }

378: PetscErrorCode SVDSetFromOptions_TRLanczos(PetscOptionItems *PetscOptionsObject,SVD svd)
379: {
381:   PetscBool      set,val;
382:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;

385:   PetscOptionsHead(PetscOptionsObject,"SVD TRLanczos Options");

387:     PetscOptionsBool("-svd_trlanczos_oneside","Use one-side reorthogonalization","SVDTRLanczosSetOneSide",lanczos->oneside,&val,&set);
388:     if (set) { SVDTRLanczosSetOneSide(svd,val); }

390:   PetscOptionsTail();
391:   return(0);
392: }

394: static PetscErrorCode SVDTRLanczosSetOneSide_TRLanczos(SVD svd,PetscBool oneside)
395: {
396:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

399:   if (lanczos->oneside != oneside) {
400:     lanczos->oneside = oneside;
401:     svd->state = SVD_STATE_INITIAL;
402:   }
403:   return(0);
404: }

406: /*@
407:    SVDTRLanczosSetOneSide - Indicate if the variant of the Lanczos method
408:    to be used is one-sided or two-sided.

410:    Logically Collective on svd

412:    Input Parameters:
413: +  svd     - singular value solver
414: -  oneside - boolean flag indicating if the method is one-sided or not

416:    Options Database Key:
417: .  -svd_trlanczos_oneside <boolean> - Indicates the boolean flag

419:    Note:
420:    By default, a two-sided variant is selected, which is sometimes slightly
421:    more robust. However, the one-sided variant is faster because it avoids
422:    the orthogonalization associated to left singular vectors.

424:    Level: advanced

426: .seealso: SVDLanczosSetOneSide()
427: @*/
428: PetscErrorCode SVDTRLanczosSetOneSide(SVD svd,PetscBool oneside)
429: {

435:   PetscTryMethod(svd,"SVDTRLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
436:   return(0);
437: }

439: static PetscErrorCode SVDTRLanczosGetOneSide_TRLanczos(SVD svd,PetscBool *oneside)
440: {
441:   SVD_TRLANCZOS *lanczos = (SVD_TRLANCZOS*)svd->data;

444:   *oneside = lanczos->oneside;
445:   return(0);
446: }

448: /*@
449:    SVDTRLanczosGetOneSide - Gets if the variant of the Lanczos method
450:    to be used is one-sided or two-sided.

452:    Not Collective

454:    Input Parameters:
455: .  svd     - singular value solver

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

460:    Level: advanced

462: .seealso: SVDTRLanczosSetOneSide()
463: @*/
464: PetscErrorCode SVDTRLanczosGetOneSide(SVD svd,PetscBool *oneside)
465: {

471:   PetscUseMethod(svd,"SVDTRLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
472:   return(0);
473: }

475: PetscErrorCode SVDDestroy_TRLanczos(SVD svd)
476: {

480:   PetscFree(svd->data);
481:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",NULL);
482:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",NULL);
483:   return(0);
484: }

486: PetscErrorCode SVDView_TRLanczos(SVD svd,PetscViewer viewer)
487: {
489:   SVD_TRLANCZOS  *lanczos = (SVD_TRLANCZOS*)svd->data;
490:   PetscBool      isascii;

493:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
494:   if (isascii) {
495:     PetscViewerASCIIPrintf(viewer,"  %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
496:   }
497:   return(0);
498: }

500: SLEPC_EXTERN PetscErrorCode SVDCreate_TRLanczos(SVD svd)
501: {
503:   SVD_TRLANCZOS  *ctx;

506:   PetscNewLog(svd,&ctx);
507:   svd->data = (void*)ctx;

509:   svd->ops->setup          = SVDSetUp_TRLanczos;
510:   svd->ops->solve          = SVDSolve_TRLanczos;
511:   svd->ops->destroy        = SVDDestroy_TRLanczos;
512:   svd->ops->setfromoptions = SVDSetFromOptions_TRLanczos;
513:   svd->ops->view           = SVDView_TRLanczos;
514:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosSetOneSide_C",SVDTRLanczosSetOneSide_TRLanczos);
515:   PetscObjectComposeFunction((PetscObject)svd,"SVDTRLanczosGetOneSide_C",SVDTRLanczosGetOneSide_TRLanczos);
516:   return(0);
517: }