Actual source code: gklanczos.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: "lanczos"

 13:    Method: Explicitly restarted Lanczos

 15:    Algorithm:

 17:        Golub-Kahan-Lanczos bidiagonalization with explicit 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: typedef struct {
 34:   PetscBool oneside;
 35: } SVD_LANCZOS;

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

 42:   SVDCheckStandard(svd);
 43:   SVDCheckDefinite(svd);
 44:   MatGetSize(svd->A,NULL,&N);
 45:   SVDSetDimensions_Default(svd);
 47:   if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
 48:   svd->leftbasis = PetscNot(lanczos->oneside);
 49:   SVDAllocateSolution(svd,1);
 50:   DSSetType(svd->ds,DSSVD);
 51:   DSSetCompact(svd->ds,PETSC_TRUE);
 52:   DSSetExtraRow(svd->ds,PETSC_TRUE);
 53:   DSAllocate(svd->ds,svd->ncv+1);
 54:   return 0;
 55: }

 57: PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
 58: {
 59:   PetscInt       i;
 60:   Vec            u,v;
 61:   PetscBool      lindep=PETSC_FALSE;

 63:   BVGetColumn(svd->V,k,&v);
 64:   BVGetColumn(svd->U,k,&u);
 65:   MatMult(svd->A,v,u);
 66:   BVRestoreColumn(svd->V,k,&v);
 67:   BVRestoreColumn(svd->U,k,&u);
 68:   BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep);
 69:   if (PetscUnlikely(lindep)) {
 70:     *n = k;
 71:     if (breakdown) *breakdown = lindep;
 72:     return 0;
 73:   }

 75:   for (i=k+1;i<*n;i++) {
 76:     BVGetColumn(svd->V,i,&v);
 77:     BVGetColumn(svd->U,i-1,&u);
 78:     MatMult(svd->AT,u,v);
 79:     BVRestoreColumn(svd->V,i,&v);
 80:     BVRestoreColumn(svd->U,i-1,&u);
 81:     BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep);
 82:     if (PetscUnlikely(lindep)) {
 83:       *n = i;
 84:       break;
 85:     }
 86:     BVGetColumn(svd->V,i,&v);
 87:     BVGetColumn(svd->U,i,&u);
 88:     MatMult(svd->A,v,u);
 89:     BVRestoreColumn(svd->V,i,&v);
 90:     BVRestoreColumn(svd->U,i,&u);
 91:     BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep);
 92:     if (PetscUnlikely(lindep)) {
 93:       *n = i;
 94:       break;
 95:     }
 96:   }

 98:   if (!lindep) {
 99:     BVGetColumn(svd->V,*n,&v);
100:     BVGetColumn(svd->U,*n-1,&u);
101:     MatMult(svd->AT,u,v);
102:     BVRestoreColumn(svd->V,*n,&v);
103:     BVRestoreColumn(svd->U,*n-1,&u);
104:     BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep);
105:   }
106:   if (PetscUnlikely(breakdown)) *breakdown = lindep;
107:   return 0;
108: }

110: static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
111: {
112:   PetscInt       i,bvl,bvk;
113:   PetscReal      a,b;
114:   Vec            z,temp;

116:   BVGetActiveColumns(V,&bvl,&bvk);
117:   BVGetColumn(V,k,&z);
118:   MatMult(svd->A,z,u);
119:   BVRestoreColumn(V,k,&z);

121:   for (i=k+1;i<n;i++) {
122:     BVGetColumn(V,i,&z);
123:     MatMult(svd->AT,u,z);
124:     BVRestoreColumn(V,i,&z);
125:     VecNormBegin(u,NORM_2,&a);
126:     BVSetActiveColumns(V,0,i);
127:     BVDotColumnBegin(V,i,work);
128:     VecNormEnd(u,NORM_2,&a);
129:     BVDotColumnEnd(V,i,work);
130:     VecScale(u,1.0/a);
131:     BVMultColumn(V,-1.0/a,1.0/a,i,work);

133:     /* h = V^* z, z = z - V h  */
134:     BVDotColumn(V,i,work);
135:     BVMultColumn(V,-1.0,1.0,i,work);
136:     BVNormColumn(V,i,NORM_2,&b);
138:     BVScaleColumn(V,i,1.0/b);

140:     BVGetColumn(V,i,&z);
141:     MatMult(svd->A,z,u_1);
142:     BVRestoreColumn(V,i,&z);
143:     VecAXPY(u_1,-b,u);
144:     alpha[i-1] = a;
145:     beta[i-1] = b;
146:     temp = u;
147:     u = u_1;
148:     u_1 = temp;
149:   }

151:   BVGetColumn(V,n,&z);
152:   MatMult(svd->AT,u,z);
153:   BVRestoreColumn(V,n,&z);
154:   VecNormBegin(u,NORM_2,&a);
155:   BVDotColumnBegin(V,n,work);
156:   VecNormEnd(u,NORM_2,&a);
157:   BVDotColumnEnd(V,n,work);
158:   VecScale(u,1.0/a);
159:   BVMultColumn(V,-1.0/a,1.0/a,n,work);

161:   /* h = V^* z, z = z - V h  */
162:   BVDotColumn(V,n,work);
163:   BVMultColumn(V,-1.0,1.0,n,work);
164:   BVNormColumn(V,i,NORM_2,&b);

166:   alpha[n-1] = a;
167:   beta[n-1] = b;
168:   BVSetActiveColumns(V,bvl,bvk);
169:   return 0;
170: }

172: /*
173:    SVDKrylovConvergence - Implements the loop that checks for convergence
174:    in Krylov methods.

176:    Input Parameters:
177:      svd     - the solver; some error estimates are updated in svd->errest
178:      getall  - whether all residuals must be computed
179:      kini    - initial value of k (the loop variable)
180:      nits    - number of iterations of the loop
181:      normr   - norm of triangular factor of qr([A;B]), used only in GSVD

183:    Output Parameter:
184:      kout  - the first index where the convergence test failed
185: */
186: PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal normr,PetscInt *kout)
187: {
188:   PetscInt       k,marker,ld;
189:   PetscReal      *alpha,*beta,*betah,resnorm;
190:   PetscBool      extra;

192:   if (PetscUnlikely(svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it)) *kout = svd->nsv;
193:   else {
194:     DSGetLeadingDimension(svd->ds,&ld);
195:     DSGetExtraRow(svd->ds,&extra);
197:     marker = -1;
198:     if (svd->trackall) getall = PETSC_TRUE;
199:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
200:     beta = alpha + ld;
201:     betah = alpha + 2*ld;
202:     for (k=kini;k<kini+nits;k++) {
203:       if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k])*normr;
204:       else resnorm = PetscAbsReal(beta[k]);
205:       (*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx);
206:       if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
207:       if (marker!=-1 && !getall) break;
208:     }
209:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);
210:     if (marker!=-1) k = marker;
211:     *kout = k;
212:   }
213:   return 0;
214: }

216: PetscErrorCode SVDSolve_Lanczos(SVD svd)
217: {
218:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
219:   PetscReal      *alpha,*beta;
220:   PetscScalar    *swork,*w,*P;
221:   PetscInt       i,k,j,nv,ld;
222:   Vec            u=NULL,u_1=NULL;
223:   Mat            U,V;

225:   /* allocate working space */
226:   DSGetLeadingDimension(svd->ds,&ld);
227:   PetscMalloc2(ld,&w,svd->ncv,&swork);

229:   if (lanczos->oneside) {
230:     MatCreateVecs(svd->A,NULL,&u);
231:     MatCreateVecs(svd->A,NULL,&u_1);
232:   }

234:   /* normalize start vector */
235:   if (!svd->nini) {
236:     BVSetRandomColumn(svd->V,0);
237:     BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL);
238:   }

240:   while (svd->reason == SVD_CONVERGED_ITERATING) {
241:     svd->its++;

243:     /* inner loop */
244:     nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
245:     BVSetActiveColumns(svd->V,svd->nconv,nv);
246:     DSGetArrayReal(svd->ds,DS_MAT_T,&alpha);
247:     beta = alpha + ld;
248:     if (lanczos->oneside) SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork);
249:     else {
250:       BVSetActiveColumns(svd->U,svd->nconv,nv);
251:       SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL);
252:     }
253:     DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha);

255:     /* compute SVD of bidiagonal matrix */
256:     DSSetDimensions(svd->ds,nv,svd->nconv,0);
257:     DSSVDSetDimensions(svd->ds,nv);
258:     DSSetState(svd->ds,DS_STATE_INTERMEDIATE);
259:     DSSolve(svd->ds,w,NULL);
260:     DSSort(svd->ds,w,NULL,NULL,NULL,NULL);
261:     DSUpdateExtraRow(svd->ds);
262:     DSSynchronize(svd->ds,w,NULL);
263:     for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);

265:     /* check convergence */
266:     SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k);
267:     (*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx);

269:     /* compute restart vector */
270:     if (svd->reason == SVD_CONVERGED_ITERATING) {
271:       if (k<nv) {
272:         DSGetArray(svd->ds,DS_MAT_V,&P);
273:         for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
274:         DSRestoreArray(svd->ds,DS_MAT_V,&P);
275:         BVMultColumn(svd->V,1.0,0.0,nv,swork);
276:       } else {
277:         /* all approximations have converged, generate a new initial vector */
278:         BVSetRandomColumn(svd->V,nv);
279:         BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL);
280:       }
281:     }

283:     /* compute converged singular vectors */
284:     DSGetMat(svd->ds,DS_MAT_V,&V);
285:     BVMultInPlace(svd->V,V,svd->nconv,k);
286:     DSRestoreMat(svd->ds,DS_MAT_V,&V);
287:     if (!lanczos->oneside) {
288:       DSGetMat(svd->ds,DS_MAT_U,&U);
289:       BVMultInPlace(svd->U,U,svd->nconv,k);
290:       DSRestoreMat(svd->ds,DS_MAT_U,&U);
291:     }

293:     /* copy restart vector from the last column */
294:     if (svd->reason == SVD_CONVERGED_ITERATING) BVCopyColumn(svd->V,nv,k);

296:     svd->nconv = k;
297:     SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv);
298:   }

300:   /* free working space */
301:   VecDestroy(&u);
302:   VecDestroy(&u_1);
303:   PetscFree2(w,swork);
304:   return 0;
305: }

307: PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd,PetscOptionItems *PetscOptionsObject)
308: {
309:   PetscBool      set,val;
310:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;

312:   PetscOptionsHeadBegin(PetscOptionsObject,"SVD Lanczos Options");

314:     PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set);
315:     if (set) SVDLanczosSetOneSide(svd,val);

317:   PetscOptionsHeadEnd();
318:   return 0;
319: }

321: static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
322: {
323:   SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;

325:   if (lanczos->oneside != oneside) {
326:     lanczos->oneside = oneside;
327:     svd->state = SVD_STATE_INITIAL;
328:   }
329:   return 0;
330: }

332: /*@
333:    SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
334:    to be used is one-sided or two-sided.

336:    Logically Collective on svd

338:    Input Parameters:
339: +  svd     - singular value solver
340: -  oneside - boolean flag indicating if the method is one-sided or not

342:    Options Database Key:
343: .  -svd_lanczos_oneside <boolean> - Indicates the boolean flag

345:    Note:
346:    By default, a two-sided variant is selected, which is sometimes slightly
347:    more robust. However, the one-sided variant is faster because it avoids
348:    the orthogonalization associated to left singular vectors. It also saves
349:    the memory required for storing such vectors.

351:    Level: advanced

353: .seealso: SVDTRLanczosSetOneSide()
354: @*/
355: PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
356: {
359:   PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
360:   return 0;
361: }

363: static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
364: {
365:   SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;

367:   *oneside = lanczos->oneside;
368:   return 0;
369: }

371: /*@
372:    SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
373:    to be used is one-sided or two-sided.

375:    Not Collective

377:    Input Parameters:
378: .  svd     - singular value solver

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

383:    Level: advanced

385: .seealso: SVDLanczosSetOneSide()
386: @*/
387: PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
388: {
391:   PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
392:   return 0;
393: }

395: PetscErrorCode SVDDestroy_Lanczos(SVD svd)
396: {
397:   PetscFree(svd->data);
398:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL);
399:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL);
400:   return 0;
401: }

403: PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
404: {
405:   SVD_LANCZOS    *lanczos = (SVD_LANCZOS*)svd->data;
406:   PetscBool      isascii;

408:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
409:   if (isascii) PetscViewerASCIIPrintf(viewer,"  %s-sided reorthogonalization\n",lanczos->oneside? "one": "two");
410:   return 0;
411: }

413: SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
414: {
415:   SVD_LANCZOS    *ctx;

417:   PetscNew(&ctx);
418:   svd->data = (void*)ctx;

420:   svd->ops->setup          = SVDSetUp_Lanczos;
421:   svd->ops->solve          = SVDSolve_Lanczos;
422:   svd->ops->destroy        = SVDDestroy_Lanczos;
423:   svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
424:   svd->ops->view           = SVDView_Lanczos;
425:   svd->ops->computevectors = SVDComputeVectors_Left;
426:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos);
427:   PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos);
428:   return 0;
429: }