Actual source code: lobpcg.c

slepc-3.20.2 2024-03-15
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 eigensolver: "lobpcg"

 13:    Method: Locally Optimal Block Preconditioned Conjugate Gradient

 15:    Algorithm:

 17:        LOBPCG with soft and hard locking. Follows the implementation
 18:        in BLOPEX [2].

 20:    References:

 22:        [1] A. V. Knyazev, "Toward the optimal preconditioned eigensolver:
 23:            locally optimal block preconditioned conjugate gradient method",
 24:            SIAM J. Sci. Comput. 23(2):517-541, 2001.

 26:        [2] A. V. Knyazev et al., "Block Locally Optimal Preconditioned
 27:            Eigenvalue Xolvers (BLOPEX) in Hypre and PETSc", SIAM J. Sci.
 28:            Comput. 29(5):2224-2239, 2007.
 29: */

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

 33: typedef struct {
 34:   PetscInt  bs;        /* block size */
 35:   PetscBool lock;      /* soft locking active/inactive */
 36:   PetscReal restart;   /* restart parameter */
 37:   PetscInt  guard;     /* number of guard vectors */
 38: } EPS_LOBPCG;

 40: static PetscErrorCode EPSSetDimensions_LOBPCG(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
 41: {
 42:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
 43:   PetscInt   k;

 45:   PetscFunctionBegin;
 46:   k = PetscMax(3*ctx->bs,((eps->nev-1)/ctx->bs+3)*ctx->bs);
 47:   if (*ncv!=PETSC_DEFAULT) { /* ncv set */
 48:     PetscCheck(*ncv>=k,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv is not sufficiently large");
 49:   } else *ncv = k;
 50:   if (*mpd==PETSC_DEFAULT) *mpd = 3*ctx->bs;
 51:   else PetscCheck(*mpd==3*ctx->bs,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"This solver does not allow a value of mpd different from 3*blocksize");
 52:   PetscFunctionReturn(PETSC_SUCCESS);
 53: }

 55: static PetscErrorCode EPSSetUp_LOBPCG(EPS eps)
 56: {
 57:   EPS_LOBPCG     *ctx = (EPS_LOBPCG*)eps->data;

 59:   PetscFunctionBegin;
 60:   EPSCheckHermitianDefinite(eps);
 61:   if (!ctx->bs) ctx->bs = PetscMin(16,eps->nev);
 62:   PetscCheck(eps->n-eps->nds>=5*ctx->bs,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The problem size is too small relative to the block size");
 63:   PetscCall(EPSSetDimensions_LOBPCG(eps,eps->nev,&eps->ncv,&eps->mpd));
 64:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 65:   if (!eps->which) eps->which = EPS_SMALLEST_REAL;
 66:   PetscCheck(eps->which==EPS_SMALLEST_REAL || eps->which==EPS_LARGEST_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only smallest real or largest real eigenvalues");
 67:   EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_EXTRACTION);
 68:   EPSCheckIgnored(eps,EPS_FEATURE_BALANCE);

 70:   if (!ctx->restart) ctx->restart = 0.9;

 72:   /* number of guard vectors */
 73:   if (ctx->bs==1) ctx->guard = 0;
 74:   else ctx->guard = PetscMin((PetscInt)((1.0-ctx->restart)*ctx->bs+0.45),ctx->bs-1);

 76:   PetscCall(EPSAllocateSolution(eps,0));
 77:   PetscCall(EPS_SetInnerProduct(eps));
 78:   PetscCall(DSSetType(eps->ds,DSGHEP));
 79:   PetscCall(DSAllocate(eps->ds,eps->mpd));
 80:   PetscCall(EPSSetWorkVecs(eps,1));
 81:   PetscFunctionReturn(PETSC_SUCCESS);
 82: }

 84: static PetscErrorCode EPSSolve_LOBPCG(EPS eps)
 85: {
 86:   EPS_LOBPCG     *ctx = (EPS_LOBPCG*)eps->data;
 87:   PetscInt       i,j,k,nv,ini,nmat,nc,nconv,locked,its,prev=0;
 88:   PetscReal      norm;
 89:   PetscScalar    *eigr,dot;
 90:   PetscBool      breakdown,countc,flip=PETSC_FALSE,checkprecond=PETSC_FALSE;
 91:   Mat            A,B,M,V=NULL,W=NULL;
 92:   Vec            v,z,w=eps->work[0];
 93:   BV             X,Y=NULL,Z,R,P,AX,BX;
 94:   SlepcSC        sc;

 96:   PetscFunctionBegin;
 97:   PetscCall(STGetNumMatrices(eps->st,&nmat));
 98:   PetscCall(STGetMatrix(eps->st,0,&A));
 99:   if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
100:   else B = NULL;

102:   if (eps->which==EPS_LARGEST_REAL) {  /* flip spectrum */
103:     flip = PETSC_TRUE;
104:     PetscCall(DSGetSlepcSC(eps->ds,&sc));
105:     sc->comparison = SlepcCompareSmallestReal;
106:   }

108:   /* undocumented option to check for a positive-definite preconditioner (turn-off by default) */
109:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-eps_lobpcg_checkprecond",&checkprecond,NULL));

111:   /* 1. Allocate memory */
112:   PetscCall(PetscCalloc1(3*ctx->bs,&eigr));
113:   PetscCall(BVDuplicateResize(eps->V,3*ctx->bs,&Z));
114:   PetscCall(BVDuplicateResize(eps->V,ctx->bs,&X));
115:   PetscCall(BVDuplicateResize(eps->V,ctx->bs,&R));
116:   PetscCall(BVDuplicateResize(eps->V,ctx->bs,&P));
117:   PetscCall(BVDuplicateResize(eps->V,ctx->bs,&AX));
118:   if (B) PetscCall(BVDuplicateResize(eps->V,ctx->bs,&BX));
119:   nc = eps->nds;
120:   if (nc>0 || eps->nev>ctx->bs-ctx->guard) PetscCall(BVDuplicateResize(eps->V,nc+eps->nev,&Y));
121:   if (nc>0) {
122:     for (j=0;j<nc;j++) {
123:       PetscCall(BVGetColumn(eps->V,-nc+j,&v));
124:       PetscCall(BVInsertVec(Y,j,v));
125:       PetscCall(BVRestoreColumn(eps->V,-nc+j,&v));
126:     }
127:     PetscCall(BVSetActiveColumns(Y,0,nc));
128:   }

130:   /* 2. Apply the constraints to the initial vectors */
131:   /* 3. B-orthogonalize initial vectors */
132:   for (k=eps->nini;k<eps->ncv-ctx->bs;k++) { /* Generate more initial vectors if necessary */
133:     PetscCall(BVSetRandomColumn(eps->V,k));
134:     PetscCall(BVOrthonormalizeColumn(eps->V,k,PETSC_TRUE,NULL,NULL));
135:   }
136:   nv = ctx->bs;
137:   PetscCall(BVSetActiveColumns(eps->V,0,nv));
138:   PetscCall(BVSetActiveColumns(Z,0,nv));
139:   PetscCall(BVCopy(eps->V,Z));
140:   PetscCall(BVCopy(Z,X));

142:   /* 4. Compute initial Ritz vectors */
143:   PetscCall(BVMatMult(X,A,AX));
144:   PetscCall(DSSetDimensions(eps->ds,nv,0,0));
145:   PetscCall(DSGetMat(eps->ds,DS_MAT_A,&M));
146:   PetscCall(BVMatProject(AX,NULL,X,M));
147:   if (flip) PetscCall(MatScale(M,-1.0));
148:   PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&M));
149:   PetscCall(DSSetIdentity(eps->ds,DS_MAT_B));
150:   PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
151:   PetscCall(DSSolve(eps->ds,eigr,NULL));
152:   PetscCall(DSSort(eps->ds,eigr,NULL,NULL,NULL,NULL));
153:   PetscCall(DSSynchronize(eps->ds,eigr,NULL));
154:   for (j=0;j<nv;j++) eps->eigr[j] = flip? -eigr[j]: eigr[j];
155:   PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
156:   PetscCall(DSGetMat(eps->ds,DS_MAT_X,&M));
157:   PetscCall(BVMultInPlace(X,M,0,nv));
158:   PetscCall(BVMultInPlace(AX,M,0,nv));
159:   PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&M));

161:   /* 5. Initialize range of active iterates */
162:   locked = 0;  /* hard-locked vectors, the leading locked columns of V are eigenvectors */
163:   nconv  = 0;  /* number of converged eigenvalues in the current block */
164:   its    = 0;  /* iterations for the current block */

166:   /* 6. Main loop */
167:   while (eps->reason == EPS_CONVERGED_ITERATING) {

169:     if (ctx->lock) {
170:       PetscCall(BVSetActiveColumns(R,nconv,ctx->bs));
171:       PetscCall(BVSetActiveColumns(AX,nconv,ctx->bs));
172:       if (B) PetscCall(BVSetActiveColumns(BX,nconv,ctx->bs));
173:     }

175:     /* 7. Compute residuals */
176:     ini = (ctx->lock)? nconv: 0;
177:     PetscCall(BVCopy(AX,R));
178:     if (B) PetscCall(BVMatMult(X,B,BX));
179:     for (j=ini;j<ctx->bs;j++) {
180:       PetscCall(BVGetColumn(R,j,&v));
181:       PetscCall(BVGetColumn(B?BX:X,j,&z));
182:       PetscCall(VecAXPY(v,-eps->eigr[locked+j],z));
183:       PetscCall(BVRestoreColumn(R,j,&v));
184:       PetscCall(BVRestoreColumn(B?BX:X,j,&z));
185:     }

187:     /* 8. Compute residual norms and update index set of active iterates */
188:     k = ini;
189:     countc = PETSC_TRUE;
190:     for (j=ini;j<ctx->bs;j++) {
191:       i = locked+j;
192:       PetscCall(BVGetColumn(R,j,&v));
193:       PetscCall(VecNorm(v,NORM_2,&norm));
194:       PetscCall(BVRestoreColumn(R,j,&v));
195:       PetscCall((*eps->converged)(eps,eps->eigr[i],eps->eigi[i],norm,&eps->errest[i],eps->convergedctx));
196:       if (countc) {
197:         if (eps->errest[i] < eps->tol) k++;
198:         else countc = PETSC_FALSE;
199:       }
200:       if (!countc && !eps->trackall) break;
201:     }
202:     nconv = k;
203:     eps->nconv = locked + nconv;
204:     if (its) PetscCall(EPSMonitor(eps,eps->its+its,eps->nconv,eps->eigr,eps->eigi,eps->errest,locked+ctx->bs));
205:     PetscCall((*eps->stopping)(eps,eps->its+its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
206:     if (eps->reason != EPS_CONVERGED_ITERATING || nconv >= ctx->bs-ctx->guard) {
207:       PetscCall(BVSetActiveColumns(eps->V,locked,eps->nconv));
208:       PetscCall(BVSetActiveColumns(X,0,nconv));
209:       PetscCall(BVCopy(X,eps->V));
210:     }
211:     if (eps->reason != EPS_CONVERGED_ITERATING) {
212:       break;
213:     } else if (nconv >= ctx->bs-ctx->guard) {
214:       eps->its += its-1;
215:       its = 0;
216:     } else its++;

218:     if (nconv >= ctx->bs-ctx->guard) {  /* force hard locking of vectors and compute new R */

220:       /* extend constraints */
221:       PetscCall(BVSetActiveColumns(Y,nc+locked,nc+locked+nconv));
222:       PetscCall(BVCopy(X,Y));
223:       PetscCall(BVSetActiveColumns(Y,0,nc+locked+nconv));

225:       /* shift work BV's */
226:       for (j=nconv;j<ctx->bs;j++) {
227:         PetscCall(BVCopyColumn(X,j,j-nconv));
228:         PetscCall(BVCopyColumn(R,j,j-nconv));
229:         PetscCall(BVCopyColumn(P,j,j-nconv));
230:         PetscCall(BVCopyColumn(AX,j,j-nconv));
231:         if (B) PetscCall(BVCopyColumn(BX,j,j-nconv));
232:       }

234:       /* set new initial vectors */
235:       PetscCall(BVSetActiveColumns(eps->V,locked+ctx->bs,locked+ctx->bs+nconv));
236:       PetscCall(BVSetActiveColumns(X,ctx->bs-nconv,ctx->bs));
237:       PetscCall(BVCopy(eps->V,X));
238:       for (j=ctx->bs-nconv;j<ctx->bs;j++) {
239:         PetscCall(BVGetColumn(X,j,&v));
240:         PetscCall(BVOrthogonalizeVec(Y,v,NULL,&norm,&breakdown));
241:         if (norm>0.0 && !breakdown) PetscCall(VecScale(v,1.0/norm));
242:         else {
243:           PetscCall(PetscInfo(eps,"Orthogonalization of initial vector failed\n"));
244:           eps->reason = EPS_DIVERGED_BREAKDOWN;
245:           goto diverged;
246:         }
247:         PetscCall(BVRestoreColumn(X,j,&v));
248:       }
249:       locked += nconv;
250:       nconv = 0;
251:       PetscCall(BVSetActiveColumns(X,nconv,ctx->bs));

253:       /* B-orthogonalize initial vectors */
254:       PetscCall(BVOrthogonalize(X,NULL));
255:       PetscCall(BVSetActiveColumns(Z,nconv,ctx->bs));
256:       PetscCall(BVSetActiveColumns(AX,nconv,ctx->bs));
257:       PetscCall(BVCopy(X,Z));

259:       /* compute initial Ritz vectors */
260:       nv = ctx->bs;
261:       PetscCall(BVMatMult(X,A,AX));
262:       PetscCall(DSSetDimensions(eps->ds,nv,0,0));
263:       PetscCall(DSGetMat(eps->ds,DS_MAT_A,&M));
264:       PetscCall(BVMatProject(AX,NULL,X,M));
265:       if (flip) PetscCall(MatScale(M,-1.0));
266:       PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&M));
267:       PetscCall(DSSetIdentity(eps->ds,DS_MAT_B));
268:       PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
269:       PetscCall(DSSolve(eps->ds,eigr,NULL));
270:       PetscCall(DSSort(eps->ds,eigr,NULL,NULL,NULL,NULL));
271:       PetscCall(DSSynchronize(eps->ds,eigr,NULL));
272:       for (j=0;j<nv;j++) if (locked+j<eps->ncv) eps->eigr[locked+j] = flip? -eigr[j]: eigr[j];
273:       PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
274:       PetscCall(DSGetMat(eps->ds,DS_MAT_X,&M));
275:       PetscCall(BVMultInPlace(X,M,0,nv));
276:       PetscCall(BVMultInPlace(AX,M,0,nv));
277:       PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&M));

279:       continue;   /* skip the rest of the iteration */
280:     }

282:     ini = (ctx->lock)? nconv: 0;
283:     if (ctx->lock) {
284:       PetscCall(BVSetActiveColumns(R,nconv,ctx->bs));
285:       PetscCall(BVSetActiveColumns(P,nconv,ctx->bs));
286:       PetscCall(BVSetActiveColumns(AX,nconv,ctx->bs));
287:       if (B) PetscCall(BVSetActiveColumns(BX,nconv,ctx->bs));
288:     }

290:     /* 9. Apply preconditioner to the residuals */
291:     PetscCall(BVGetMat(R,&V));
292:     if (prev != ctx->bs-ini) {
293:       prev = ctx->bs-ini;
294:       PetscCall(MatDestroy(&W));
295:       PetscCall(MatDuplicate(V,MAT_SHARE_NONZERO_PATTERN,&W));
296:     }
297:     PetscCall(STApplyMat(eps->st,V,W));
298:     if (checkprecond) {
299:       for (j=ini;j<ctx->bs;j++) {
300:         PetscCall(MatDenseGetColumnVecRead(V,j-ini,&v));
301:         PetscCall(MatDenseGetColumnVecRead(W,j-ini,&w));
302:         PetscCall(VecDot(v,w,&dot));
303:         PetscCall(MatDenseRestoreColumnVecRead(W,j-ini,&w));
304:         PetscCall(MatDenseRestoreColumnVecRead(V,j-ini,&v));
305:         if (PetscRealPart(dot)<0.0) {
306:           PetscCall(PetscInfo(eps,"The preconditioner is not positive-definite\n"));
307:           eps->reason = EPS_DIVERGED_BREAKDOWN;
308:           goto diverged;
309:         }
310:       }
311:     }
312:     if (nc+locked>0) {
313:       for (j=ini;j<ctx->bs;j++) {
314:         PetscCall(MatDenseGetColumnVecWrite(W,j-ini,&w));
315:         PetscCall(BVOrthogonalizeVec(Y,w,NULL,&norm,&breakdown));
316:         if (norm>0.0 && !breakdown) PetscCall(VecScale(w,1.0/norm));
317:         PetscCall(MatDenseRestoreColumnVecWrite(W,j-ini,&w));
318:         if (norm<=0.0 || breakdown) {
319:           PetscCall(PetscInfo(eps,"Orthogonalization of preconditioned residual failed\n"));
320:           eps->reason = EPS_DIVERGED_BREAKDOWN;
321:           goto diverged;
322:         }
323:       }
324:     }
325:     PetscCall(MatCopy(W,V,SAME_NONZERO_PATTERN));
326:     PetscCall(BVRestoreMat(R,&V));

328:     /* 11. B-orthonormalize preconditioned residuals */
329:     PetscCall(BVOrthogonalize(R,NULL));

331:     /* 13-16. B-orthonormalize conjugate directions */
332:     if (its>1) PetscCall(BVOrthogonalize(P,NULL));

334:     /* 17-23. Compute symmetric Gram matrices */
335:     PetscCall(BVSetActiveColumns(Z,0,ctx->bs));
336:     PetscCall(BVSetActiveColumns(X,0,ctx->bs));
337:     PetscCall(BVCopy(X,Z));
338:     PetscCall(BVSetActiveColumns(Z,ctx->bs,2*ctx->bs-ini));
339:     PetscCall(BVCopy(R,Z));
340:     if (its>1) {
341:       PetscCall(BVSetActiveColumns(Z,2*ctx->bs-ini,3*ctx->bs-2*ini));
342:       PetscCall(BVCopy(P,Z));
343:     }

345:     if (its>1) nv = 3*ctx->bs-2*ini;
346:     else nv = 2*ctx->bs-ini;

348:     PetscCall(BVSetActiveColumns(Z,0,nv));
349:     PetscCall(DSSetDimensions(eps->ds,nv,0,0));
350:     PetscCall(DSGetMat(eps->ds,DS_MAT_A,&M));
351:     PetscCall(BVMatProject(Z,A,Z,M));
352:     if (flip) PetscCall(MatScale(M,-1.0));
353:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_A,&M));
354:     PetscCall(DSGetMat(eps->ds,DS_MAT_B,&M));
355:     PetscCall(BVMatProject(Z,B,Z,M)); /* covers also the case B=NULL */
356:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_B,&M));

358:     /* 24. Solve the generalized eigenvalue problem */
359:     PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
360:     PetscCall(DSSolve(eps->ds,eigr,NULL));
361:     PetscCall(DSSort(eps->ds,eigr,NULL,NULL,NULL,NULL));
362:     PetscCall(DSSynchronize(eps->ds,eigr,NULL));
363:     for (j=0;j<nv;j++) if (locked+j<eps->ncv) eps->eigr[locked+j] = flip? -eigr[j]: eigr[j];
364:     PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));

366:     /* 25-33. Compute Ritz vectors */
367:     PetscCall(DSGetMat(eps->ds,DS_MAT_X,&M));
368:     PetscCall(BVSetActiveColumns(Z,ctx->bs,nv));
369:     if (ctx->lock) PetscCall(BVSetActiveColumns(P,0,ctx->bs));
370:     PetscCall(BVMult(P,1.0,0.0,Z,M));
371:     PetscCall(BVCopy(P,X));
372:     if (ctx->lock) PetscCall(BVSetActiveColumns(P,nconv,ctx->bs));
373:     PetscCall(BVSetActiveColumns(Z,0,ctx->bs));
374:     PetscCall(BVMult(X,1.0,1.0,Z,M));
375:     if (ctx->lock) PetscCall(BVSetActiveColumns(X,nconv,ctx->bs));
376:     PetscCall(BVMatMult(X,A,AX));
377:     PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&M));
378:   }

380: diverged:
381:   eps->its += its;

383:   if (flip) sc->comparison = SlepcCompareLargestReal;
384:   PetscCall(PetscFree(eigr));
385:   PetscCall(MatDestroy(&W));
386:   if (V) PetscCall(BVRestoreMat(R,&V)); /* only needed when goto diverged is reached */
387:   PetscCall(BVDestroy(&Z));
388:   PetscCall(BVDestroy(&X));
389:   PetscCall(BVDestroy(&R));
390:   PetscCall(BVDestroy(&P));
391:   PetscCall(BVDestroy(&AX));
392:   if (B) PetscCall(BVDestroy(&BX));
393:   if (nc>0 || eps->nev>ctx->bs-ctx->guard) PetscCall(BVDestroy(&Y));
394:   PetscFunctionReturn(PETSC_SUCCESS);
395: }

397: static PetscErrorCode EPSLOBPCGSetBlockSize_LOBPCG(EPS eps,PetscInt bs)
398: {
399:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

401:   PetscFunctionBegin;
402:   if (bs == PETSC_DEFAULT || bs == PETSC_DECIDE) bs = 0;
403:   else PetscCheck(bs>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size %" PetscInt_FMT,bs);
404:   if (ctx->bs != bs) {
405:     ctx->bs = bs;
406:     eps->state = EPS_STATE_INITIAL;
407:   }
408:   PetscFunctionReturn(PETSC_SUCCESS);
409: }

411: /*@
412:    EPSLOBPCGSetBlockSize - Sets the block size of the LOBPCG method.

414:    Logically Collective

416:    Input Parameters:
417: +  eps - the eigenproblem solver context
418: -  bs  - the block size

420:    Options Database Key:
421: .  -eps_lobpcg_blocksize - Sets the block size

423:    Level: advanced

425: .seealso: EPSLOBPCGGetBlockSize()
426: @*/
427: PetscErrorCode EPSLOBPCGSetBlockSize(EPS eps,PetscInt bs)
428: {
429:   PetscFunctionBegin;
432:   PetscTryMethod(eps,"EPSLOBPCGSetBlockSize_C",(EPS,PetscInt),(eps,bs));
433:   PetscFunctionReturn(PETSC_SUCCESS);
434: }

436: static PetscErrorCode EPSLOBPCGGetBlockSize_LOBPCG(EPS eps,PetscInt *bs)
437: {
438:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

440:   PetscFunctionBegin;
441:   *bs = ctx->bs;
442:   PetscFunctionReturn(PETSC_SUCCESS);
443: }

445: /*@
446:    EPSLOBPCGGetBlockSize - Gets the block size used in the LOBPCG method.

448:    Not Collective

450:    Input Parameter:
451: .  eps - the eigenproblem solver context

453:    Output Parameter:
454: .  bs - the block size

456:    Level: advanced

458: .seealso: EPSLOBPCGSetBlockSize()
459: @*/
460: PetscErrorCode EPSLOBPCGGetBlockSize(EPS eps,PetscInt *bs)
461: {
462:   PetscFunctionBegin;
464:   PetscAssertPointer(bs,2);
465:   PetscUseMethod(eps,"EPSLOBPCGGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
466:   PetscFunctionReturn(PETSC_SUCCESS);
467: }

469: static PetscErrorCode EPSLOBPCGSetRestart_LOBPCG(EPS eps,PetscReal restart)
470: {
471:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

473:   PetscFunctionBegin;
474:   if (restart==(PetscReal)PETSC_DEFAULT) restart = 0.9;
475:   PetscCheck(restart>=0.1 && restart<=1.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The restart argument %g must be in the range [0.1,1.0]",(double)restart);
476:   if (restart != ctx->restart) {
477:     ctx->restart = restart;
478:     eps->state = EPS_STATE_INITIAL;
479:   }
480:   PetscFunctionReturn(PETSC_SUCCESS);
481: }

483: /*@
484:    EPSLOBPCGSetRestart - Sets the restart parameter for the LOBPCG method.

486:    Logically Collective

488:    Input Parameters:
489: +  eps - the eigenproblem solver context
490: -  restart - the percentage of the block of vectors to force a restart

492:    Options Database Key:
493: .  -eps_lobpcg_restart - Sets the restart parameter

495:    Notes:
496:    The meaning of this parameter is the proportion of vectors within the
497:    current block iterate that must have converged in order to force a
498:    restart with hard locking.
499:    Allowed values are in the range [0.1,1.0]. The default is 0.9.

501:    Level: advanced

503: .seealso: EPSLOBPCGGetRestart()
504: @*/
505: PetscErrorCode EPSLOBPCGSetRestart(EPS eps,PetscReal restart)
506: {
507:   PetscFunctionBegin;
510:   PetscTryMethod(eps,"EPSLOBPCGSetRestart_C",(EPS,PetscReal),(eps,restart));
511:   PetscFunctionReturn(PETSC_SUCCESS);
512: }

514: static PetscErrorCode EPSLOBPCGGetRestart_LOBPCG(EPS eps,PetscReal *restart)
515: {
516:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

518:   PetscFunctionBegin;
519:   *restart = ctx->restart;
520:   PetscFunctionReturn(PETSC_SUCCESS);
521: }

523: /*@
524:    EPSLOBPCGGetRestart - Gets the restart parameter used in the LOBPCG method.

526:    Not Collective

528:    Input Parameter:
529: .  eps - the eigenproblem solver context

531:    Output Parameter:
532: .  restart - the restart parameter

534:    Level: advanced

536: .seealso: EPSLOBPCGSetRestart()
537: @*/
538: PetscErrorCode EPSLOBPCGGetRestart(EPS eps,PetscReal *restart)
539: {
540:   PetscFunctionBegin;
542:   PetscAssertPointer(restart,2);
543:   PetscUseMethod(eps,"EPSLOBPCGGetRestart_C",(EPS,PetscReal*),(eps,restart));
544:   PetscFunctionReturn(PETSC_SUCCESS);
545: }

547: static PetscErrorCode EPSLOBPCGSetLocking_LOBPCG(EPS eps,PetscBool lock)
548: {
549:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

551:   PetscFunctionBegin;
552:   ctx->lock = lock;
553:   PetscFunctionReturn(PETSC_SUCCESS);
554: }

556: /*@
557:    EPSLOBPCGSetLocking - Choose between locking and non-locking variants of
558:    the LOBPCG method.

560:    Logically Collective

562:    Input Parameters:
563: +  eps  - the eigenproblem solver context
564: -  lock - true if the locking variant must be selected

566:    Options Database Key:
567: .  -eps_lobpcg_locking - Sets the locking flag

569:    Notes:
570:    This flag refers to soft locking (converged vectors within the current
571:    block iterate), since hard locking is always used (when nev is larger
572:    than the block size).

574:    Level: advanced

576: .seealso: EPSLOBPCGGetLocking()
577: @*/
578: PetscErrorCode EPSLOBPCGSetLocking(EPS eps,PetscBool lock)
579: {
580:   PetscFunctionBegin;
583:   PetscTryMethod(eps,"EPSLOBPCGSetLocking_C",(EPS,PetscBool),(eps,lock));
584:   PetscFunctionReturn(PETSC_SUCCESS);
585: }

587: static PetscErrorCode EPSLOBPCGGetLocking_LOBPCG(EPS eps,PetscBool *lock)
588: {
589:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

591:   PetscFunctionBegin;
592:   *lock = ctx->lock;
593:   PetscFunctionReturn(PETSC_SUCCESS);
594: }

596: /*@
597:    EPSLOBPCGGetLocking - Gets the locking flag used in the LOBPCG method.

599:    Not Collective

601:    Input Parameter:
602: .  eps - the eigenproblem solver context

604:    Output Parameter:
605: .  lock - the locking flag

607:    Level: advanced

609: .seealso: EPSLOBPCGSetLocking()
610: @*/
611: PetscErrorCode EPSLOBPCGGetLocking(EPS eps,PetscBool *lock)
612: {
613:   PetscFunctionBegin;
615:   PetscAssertPointer(lock,2);
616:   PetscUseMethod(eps,"EPSLOBPCGGetLocking_C",(EPS,PetscBool*),(eps,lock));
617:   PetscFunctionReturn(PETSC_SUCCESS);
618: }

620: static PetscErrorCode EPSView_LOBPCG(EPS eps,PetscViewer viewer)
621: {
622:   EPS_LOBPCG     *ctx = (EPS_LOBPCG*)eps->data;
623:   PetscBool      isascii;

625:   PetscFunctionBegin;
626:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
627:   if (isascii) {
628:     PetscCall(PetscViewerASCIIPrintf(viewer,"  block size %" PetscInt_FMT "\n",ctx->bs));
629:     PetscCall(PetscViewerASCIIPrintf(viewer,"  restart parameter=%g (using %" PetscInt_FMT " guard vectors)\n",(double)ctx->restart,ctx->guard));
630:     PetscCall(PetscViewerASCIIPrintf(viewer,"  soft locking %sactivated\n",ctx->lock?"":"de"));
631:   }
632:   PetscFunctionReturn(PETSC_SUCCESS);
633: }

635: static PetscErrorCode EPSSetFromOptions_LOBPCG(EPS eps,PetscOptionItems *PetscOptionsObject)
636: {
637:   PetscBool      lock,flg;
638:   PetscInt       bs;
639:   PetscReal      restart;

641:   PetscFunctionBegin;
642:   PetscOptionsHeadBegin(PetscOptionsObject,"EPS LOBPCG Options");

644:     PetscCall(PetscOptionsInt("-eps_lobpcg_blocksize","Block size","EPSLOBPCGSetBlockSize",20,&bs,&flg));
645:     if (flg) PetscCall(EPSLOBPCGSetBlockSize(eps,bs));

647:     PetscCall(PetscOptionsReal("-eps_lobpcg_restart","Percentage of the block of vectors to force a restart","EPSLOBPCGSetRestart",0.5,&restart,&flg));
648:     if (flg) PetscCall(EPSLOBPCGSetRestart(eps,restart));

650:     PetscCall(PetscOptionsBool("-eps_lobpcg_locking","Choose between locking and non-locking variants","EPSLOBPCGSetLocking",PETSC_TRUE,&lock,&flg));
651:     if (flg) PetscCall(EPSLOBPCGSetLocking(eps,lock));

653:   PetscOptionsHeadEnd();
654:   PetscFunctionReturn(PETSC_SUCCESS);
655: }

657: static PetscErrorCode EPSDestroy_LOBPCG(EPS eps)
658: {
659:   PetscFunctionBegin;
660:   PetscCall(PetscFree(eps->data));
661:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetBlockSize_C",NULL));
662:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetBlockSize_C",NULL));
663:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetRestart_C",NULL));
664:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetRestart_C",NULL));
665:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetLocking_C",NULL));
666:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetLocking_C",NULL));
667:   PetscFunctionReturn(PETSC_SUCCESS);
668: }

670: SLEPC_EXTERN PetscErrorCode EPSCreate_LOBPCG(EPS eps)
671: {
672:   EPS_LOBPCG     *lobpcg;

674:   PetscFunctionBegin;
675:   PetscCall(PetscNew(&lobpcg));
676:   eps->data = (void*)lobpcg;
677:   lobpcg->lock = PETSC_TRUE;

679:   eps->useds = PETSC_TRUE;
680:   eps->categ = EPS_CATEGORY_PRECOND;

682:   eps->ops->solve          = EPSSolve_LOBPCG;
683:   eps->ops->setup          = EPSSetUp_LOBPCG;
684:   eps->ops->setupsort      = EPSSetUpSort_Default;
685:   eps->ops->setfromoptions = EPSSetFromOptions_LOBPCG;
686:   eps->ops->destroy        = EPSDestroy_LOBPCG;
687:   eps->ops->view           = EPSView_LOBPCG;
688:   eps->ops->backtransform  = EPSBackTransform_Default;
689:   eps->ops->setdefaultst   = EPSSetDefaultST_GMRES;

691:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetBlockSize_C",EPSLOBPCGSetBlockSize_LOBPCG));
692:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetBlockSize_C",EPSLOBPCGGetBlockSize_LOBPCG));
693:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetRestart_C",EPSLOBPCGSetRestart_LOBPCG));
694:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetRestart_C",EPSLOBPCGGetRestart_LOBPCG));
695:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetLocking_C",EPSLOBPCGSetLocking_LOBPCG));
696:   PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetLocking_C",EPSLOBPCGGetLocking_LOBPCG));
697:   PetscFunctionReturn(PETSC_SUCCESS);
698: }