Actual source code: lobpcg.c

slepc-3.7.2 2016-07-19
Report Typos and Errors
  1: /*

  3:    SLEPc eigensolver: "lobpcg"

  5:    Method: Locally Optimal Block Preconditioned Conjugate Gradient

  7:    Algorithm:

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

 12:    References:

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

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

 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/epsimpl.h>                /*I "slepceps.h" I*/

 44: typedef struct {
 45:   PetscInt  bs;        /* block size */
 46:   PetscBool lock;      /* soft locking active/inactive */
 47:   PetscReal restart;   /* restart parameter */
 48: } EPS_LOBPCG;

 52: PetscErrorCode EPSSetDimensions_LOBPCG(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
 53: {
 54:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
 55:   PetscInt   k;

 58:   k = PetscMax(3*ctx->bs,((eps->nev-1)/ctx->bs+3)*ctx->bs);
 59:   if (*ncv) { /* ncv set */
 60:     if (*ncv<k) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv is not sufficiently large");
 61:   } else *ncv = k;
 62:   if (!*mpd) *mpd = 3*ctx->bs;
 63:   else if (*mpd!=3*ctx->bs) SETERRQ(PetscObjectComm((PetscObject)eps),1,"This solver does not allow a value of mpd different from 3*blocksize");
 64:   return(0);
 65: }

 69: PetscErrorCode EPSSetUp_LOBPCG(EPS eps)
 70: {
 72:   EPS_LOBPCG     *ctx = (EPS_LOBPCG*)eps->data;
 73:   PetscBool      precond,istrivial;

 76:   if (!eps->ishermitian || (eps->isgeneralized && !eps->ispositive)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"LOBPCG only works for Hermitian problems");
 77:   if (!ctx->bs) ctx->bs = PetscMin(16,eps->nev);
 78:   if (eps->n-eps->nds<5*ctx->bs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The problem size is too small relative to the block size");
 79:   EPSSetDimensions_LOBPCG(eps,eps->nev,&eps->ncv,&eps->mpd);
 80:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 81:   if (!eps->which) eps->which = EPS_SMALLEST_REAL;
 82:   if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
 83:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 84:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
 85:   RGIsTrivial(eps->rg,&istrivial);
 86:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");

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

 90:   STSetUp(eps->st);
 91:   PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond);
 92:   if (!precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"LOBPCG only works with precond ST");

 94:   EPSAllocateSolution(eps,0);
 95:   EPS_SetInnerProduct(eps);
 96:   DSSetType(eps->ds,DSGHEP);
 97:   DSAllocate(eps->ds,eps->mpd);
 98:   EPSSetWorkVecs(eps,1);
 99:   return(0);
100: }

104: PetscErrorCode EPSSolve_LOBPCG(EPS eps)
105: {
107:   EPS_LOBPCG     *ctx = (EPS_LOBPCG*)eps->data;
108:   PetscInt       i,j,k,ld,nv,ini,kini,nmat,nc,nconv,locked,guard,its;
109:   PetscReal      norm;
110:   PetscScalar    *eigr;
111:   PetscBool      breakdown,countc;
112:   Mat            A,B,M;
113:   Vec            v,w=eps->work[0];
114:   BV             X,Y,Z,R,P,AX,BX;

117:   DSGetLeadingDimension(eps->ds,&ld);
118:   STGetNumMatrices(eps->st,&nmat);
119:   STGetOperators(eps->st,0,&A);
120:   if (nmat>1) { STGetOperators(eps->st,1,&B); }
121:   else B = NULL;

123:   guard = (PetscInt)((1.0-ctx->restart)*ctx->bs);  /* number of guard vectors */

125:   /* 1. Allocate memory */
126:   PetscCalloc1(3*ctx->bs,&eigr);
127:   BVDuplicateResize(eps->V,3*ctx->bs,&Z);
128:   BVDuplicateResize(eps->V,ctx->bs,&X);
129:   BVDuplicateResize(eps->V,ctx->bs,&R);
130:   BVDuplicateResize(eps->V,ctx->bs,&P);
131:   BVDuplicateResize(eps->V,ctx->bs,&AX);
132:   if (B) {
133:     BVDuplicateResize(eps->V,ctx->bs,&BX);
134:   }
135:   nc = eps->nds;
136:   if (nc>0 || eps->nev>ctx->bs-guard) {
137:     BVDuplicateResize(eps->V,nc+eps->nev,&Y);
138:   }
139:   if (nc>0) {
140:     for (j=0;j<nc;j++) {
141:       BVGetColumn(eps->V,-nc+j,&v);
142:       BVInsertVec(Y,j,v);
143:       BVRestoreColumn(eps->V,-nc+j,&v);
144:     }
145:     BVSetActiveColumns(Y,0,nc);
146:   }

148:   /* 2. Apply the constraints to the initial vectors */
149:   kini = eps->nini;
150:   while (kini<eps->ncv-ctx->bs) { /* Generate more initial vectors if necessary */
151:     BVSetRandomColumn(eps->V,kini);
152:     BVOrthogonalizeColumn(eps->V,kini,NULL,&norm,&breakdown);
153:     if (norm>0.0 && !breakdown) {
154:       BVScaleColumn(eps->V,kini,1.0/norm);
155:       kini++;
156:     }
157:   }
158:   nv = ctx->bs;
159:   BVSetActiveColumns(eps->V,0,nv);
160:   BVSetActiveColumns(Z,0,nv);
161:   BVCopy(eps->V,Z);
162:   BVCopy(Z,X);

164:   /* 3. B-orthogonalize initial vectors */

166:   /* 4. Compute initial Ritz vectors */
167:   BVMatMult(X,A,AX);
168:   DSSetDimensions(eps->ds,nv,0,0,0);
169:   DSGetMat(eps->ds,DS_MAT_A,&M);
170:   BVMatProject(AX,NULL,X,M);
171:   DSRestoreMat(eps->ds,DS_MAT_A,&M);
172:   DSSetIdentity(eps->ds,DS_MAT_B);
173:   DSSetState(eps->ds,DS_STATE_RAW);
174:   DSSolve(eps->ds,eps->eigr,NULL);
175:   DSSort(eps->ds,eps->eigr,NULL,NULL,NULL,NULL);
176:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
177:   DSGetMat(eps->ds,DS_MAT_X,&M);
178:   BVMultInPlace(X,M,0,nv);
179:   BVMultInPlace(AX,M,0,nv);
180:   DSRestoreMat(eps->ds,DS_MAT_X,&M);

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

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

190:     if (ctx->lock) {
191:       BVSetActiveColumns(R,nconv,ctx->bs);
192:       BVSetActiveColumns(AX,nconv,ctx->bs);
193:       if (B) {
194:         BVSetActiveColumns(BX,nconv,ctx->bs);
195:       }
196:     }

198:     /* 7. Compute residuals */
199:     DSGetMat(eps->ds,DS_MAT_A,&M);
200:     BVCopy(AX,R);
201:     if (B) {
202:       BVMatMult(X,B,BX);
203:       BVMult(R,-1.0,1.0,BX,M);
204:     } else {
205:       BVMult(R,-1.0,1.0,X,M);
206:     }
207:     DSRestoreMat(eps->ds,DS_MAT_A,&M);

209:     /* 8. Compute residual norms and update index set of active iterates */
210:     ini = (ctx->lock)? nconv: 0;
211:     k = ini;
212:     countc = PETSC_TRUE;
213:     for (j=ini;j<ctx->bs;j++) {
214:       i = locked+j;
215:       BVGetColumn(R,j,&v);
216:       VecNorm(v,NORM_2,&norm);
217:       BVRestoreColumn(R,j,&v);
218:       (*eps->converged)(eps,eps->eigr[i],eps->eigi[i],norm,&eps->errest[i],eps->convergedctx);
219:       if (countc) {
220:         if (eps->errest[i] < eps->tol) k++;
221:         else countc = PETSC_FALSE;
222:       }
223:       if (!countc && !eps->trackall) break;
224:     }
225:     nconv = k;
226:     eps->nconv = locked + nconv;
227:     if (its) {
228:       EPSMonitor(eps,eps->its+its,eps->nconv,eps->eigr,eps->eigi,eps->errest,locked+ctx->bs);
229:     }
230:     (*eps->stopping)(eps,eps->its+its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
231:     if (eps->reason != EPS_CONVERGED_ITERATING || nconv >= ctx->bs-guard) {
232:       BVSetActiveColumns(eps->V,locked,eps->nconv);
233:       BVSetActiveColumns(X,0,nconv);
234:       BVCopy(X,eps->V);
235:     }
236:     if (eps->reason != EPS_CONVERGED_ITERATING) {
237:       eps->its += its;
238:       break;
239:     } else if (nconv >= ctx->bs-guard) {
240:       eps->its += its;
241:       its = 0;
242:     } else its++;

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

246:       /* extend constraints */
247:       BVSetActiveColumns(Y,nc+locked,nc+locked+nconv);
248:       BVCopy(X,Y);
249:       for (j=0;j<nconv;j++) {
250:         BVOrthogonalizeColumn(Y,nc+locked+j,NULL,&norm,&breakdown);
251:         if (norm>0.0 && !breakdown) {
252:           BVScaleColumn(Y,nc+locked+j,1.0/norm);
253:         } else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Orthogonalization of constraints failed");
254:       }
255:       BVSetActiveColumns(Y,0,nc+locked+nconv);

257:       /* shift work BV's */
258:       for (j=nconv;j<ctx->bs;j++) {
259:         BVCopyColumn(X,j,j-nconv);
260:         BVCopyColumn(R,j,j-nconv);
261:         BVCopyColumn(P,j,j-nconv);
262:         BVCopyColumn(AX,j,j-nconv);
263:         if (B) {
264:           BVCopyColumn(BX,j,j-nconv);
265:         }
266:       }

268:       /* set new initial vectors */
269:       BVSetActiveColumns(eps->V,locked+ctx->bs,locked+ctx->bs+nconv);
270:       BVSetActiveColumns(X,ctx->bs-nconv,ctx->bs);
271:       BVCopy(eps->V,X);
272:       for (j=ctx->bs-nconv;j<ctx->bs;j++) {
273:         BVGetColumn(X,j,&v);
274:         BVOrthogonalizeVec(Y,v,NULL,&norm,&breakdown);
275:         if (norm>0.0 && !breakdown) {
276:           VecScale(v,1.0/norm);
277:         } else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Orthogonalization of initial vector failed");
278:         BVRestoreColumn(X,j,&v);
279:       }
280:       locked += nconv;
281:       nconv = 0;
282:       BVSetActiveColumns(X,nconv,ctx->bs);

284:       /* B-orthogonalize initial vectors */
285:       BVOrthogonalize(X,NULL);
286:       BVSetActiveColumns(Z,nconv,ctx->bs);
287:       BVSetActiveColumns(AX,nconv,ctx->bs);
288:       BVCopy(X,Z);

290:       /* compute initial Ritz vectors */
291:       nv = ctx->bs;
292:       BVMatMult(X,A,AX);
293:       DSSetDimensions(eps->ds,nv,0,0,0);
294:       DSGetMat(eps->ds,DS_MAT_A,&M);
295:       BVMatProject(AX,NULL,X,M);
296:       DSRestoreMat(eps->ds,DS_MAT_A,&M);
297:       DSSetIdentity(eps->ds,DS_MAT_B);
298:       DSSetState(eps->ds,DS_STATE_RAW);
299:       DSSolve(eps->ds,eigr,NULL);
300:       DSSort(eps->ds,eigr,NULL,NULL,NULL,NULL);
301:       for (j=0;j<nv;j++) if (locked+j<eps->ncv) eps->eigr[locked+j] = eigr[j];
302:       DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
303:       DSGetMat(eps->ds,DS_MAT_X,&M);
304:       BVMultInPlace(X,M,0,nv);
305:       BVMultInPlace(AX,M,0,nv);
306:       DSRestoreMat(eps->ds,DS_MAT_X,&M);

308:       continue;   /* skip the rest of the iteration */
309:     }

311:     ini = (ctx->lock)? nconv: 0;
312:     if (ctx->lock) {
313:       BVSetActiveColumns(R,nconv,ctx->bs);
314:       BVSetActiveColumns(P,nconv,ctx->bs);
315:       BVSetActiveColumns(AX,nconv,ctx->bs);
316:       if (B) {
317:         BVSetActiveColumns(BX,nconv,ctx->bs);
318:       }
319:     }

321:     /* 9. Apply preconditioner to the residuals */
322:     for (j=ini;j<ctx->bs;j++) {
323:       BVGetColumn(R,j,&v);
324:       STMatSolve(eps->st,v,w);
325:       if (nc+locked>0) {
326:         BVOrthogonalizeVec(Y,w,NULL,&norm,&breakdown);
327:         if (norm>0.0 && !breakdown) {
328:           VecScale(w,1.0/norm);
329:         } else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Orthogonalization of preconditioned residual failed");
330:       }
331:       VecCopy(w,v);
332:       BVRestoreColumn(R,j,&v);
333:     }

335:     /* 11. B-orthonormalize preconditioned residuals */
336:     BVOrthogonalize(R,NULL);

338:     /* 13-16. B-orthonormalize conjugate directions */
339:     if (its>1) {
340:       BVOrthogonalize(P,NULL);
341:     }

343:     /* 17-23. Compute symmetric Gram matrices */
344:     BVSetActiveColumns(Z,0,ctx->bs);
345:     BVSetActiveColumns(X,0,ctx->bs);
346:     BVCopy(X,Z);
347:     BVSetActiveColumns(Z,ctx->bs,2*ctx->bs-ini);
348:     BVCopy(R,Z);
349:     if (its>1) {
350:       BVSetActiveColumns(Z,2*ctx->bs-ini,3*ctx->bs-2*ini);
351:       BVCopy(P,Z);
352:     }

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

357:     BVSetActiveColumns(Z,0,nv);
358:     DSSetDimensions(eps->ds,nv,0,0,0);
359:     DSGetMat(eps->ds,DS_MAT_A,&M);
360:     BVMatProject(Z,A,Z,M);
361:     DSRestoreMat(eps->ds,DS_MAT_A,&M);
362:     DSGetMat(eps->ds,DS_MAT_B,&M);
363:     if (B) {
364:       BVMatProject(Z,B,Z,M);
365:     } else {
366:       BVDot(Z,Z,M);
367:     }
368:     DSRestoreMat(eps->ds,DS_MAT_B,&M);
369:     
370:     /* 24. Solve the generalized eigenvalue problem */
371:     DSSetState(eps->ds,DS_STATE_RAW);
372:     DSSolve(eps->ds,eigr,NULL);
373:     DSSort(eps->ds,eigr,NULL,NULL,NULL,NULL);
374:     for (j=0;j<nv;j++) if (locked+j<eps->ncv) eps->eigr[locked+j] = eigr[j];
375:     DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
376:     
377:     /* 25-33. Compute Ritz vectors */
378:     DSGetMat(eps->ds,DS_MAT_X,&M);
379:     BVSetActiveColumns(Z,ctx->bs,nv);
380:     if (ctx->lock) {
381:       BVSetActiveColumns(P,0,ctx->bs);
382:     }
383:     BVMult(P,1.0,0.0,Z,M);
384:     BVCopy(P,X);
385:     if (ctx->lock) {
386:       BVSetActiveColumns(P,nconv,ctx->bs);
387:     }
388:     BVSetActiveColumns(Z,0,ctx->bs);
389:     BVMult(X,1.0,1.0,Z,M);
390:     if (ctx->lock) {
391:       BVSetActiveColumns(X,nconv,ctx->bs);
392:     }
393:     BVMatMult(X,A,AX);
394:     DSRestoreMat(eps->ds,DS_MAT_X,&M);
395:   }

397:   PetscFree(eigr);
398:   BVDestroy(&Z);
399:   BVDestroy(&X);
400:   BVDestroy(&R);
401:   BVDestroy(&P);
402:   BVDestroy(&AX);
403:   if (B) {
404:     BVDestroy(&BX);
405:   }
406:   if (nc>0 || eps->nev>ctx->bs-guard) {
407:     BVDestroy(&Y);
408:   }
409:   return(0);
410: }

414: static PetscErrorCode EPSLOBPCGSetBlockSize_LOBPCG(EPS eps,PetscInt bs)
415: {
416:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

419:   ctx->bs = bs;
420:   return(0);
421: }

425: /*@
426:    EPSLOBPCGSetBlockSize - Sets the block size of the LOBPCG method.

428:    Logically Collective on EPS

430:    Input Parameters:
431: +  eps - the eigenproblem solver context
432: -  bs  - the block size

434:    Options Database Key:
435: .  -eps_lobpcg_blocksize - Sets the block size

437:    Level: advanced

439: .seealso: EPSLOBPCGGetBlockSize()
440: @*/
441: PetscErrorCode EPSLOBPCGSetBlockSize(EPS eps,PetscInt bs)
442: {

448:   PetscTryMethod(eps,"EPSLOBPCGSetBlockSize_C",(EPS,PetscInt),(eps,bs));
449:   return(0);
450: }

454: static PetscErrorCode EPSLOBPCGGetBlockSize_LOBPCG(EPS eps,PetscInt *bs)
455: {
456:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

459:   *bs = ctx->bs;
460:   return(0);
461: }

465: /*@
466:    EPSLOBPCGGetBlockSize - Gets the block size used in the LOBPCG method.

468:    Not Collective

470:    Input Parameter:
471: .  eps - the eigenproblem solver context

473:    Output Parameter:
474: .  bs - the block size

476:    Level: advanced

478: .seealso: EPSLOBPCGSetBlockSize()
479: @*/
480: PetscErrorCode EPSLOBPCGGetBlockSize(EPS eps,PetscInt *bs)
481: {

487:   PetscUseMethod(eps,"EPSLOBPCGGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
488:   return(0);
489: }

493: static PetscErrorCode EPSLOBPCGSetRestart_LOBPCG(EPS eps,PetscReal restart)
494: {
495:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

498:   if (restart==PETSC_DEFAULT) ctx->restart = 0.6;
499:   else {
500:     if (restart<0.1 || restart>1.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The restart argument must be in the range [0.1,1.0]");
501:     ctx->restart = restart;
502:   }
503:   return(0);
504: }

508: /*@
509:    EPSLOBPCGSetRestart - Sets the restart parameter for the LOBPCG method.
510:    The meaning of this parameter is the proportion of vectors within the
511:    current block iterate that must have converged in order to force a
512:    restart with hard locking.

514:    Logically Collective on EPS

516:    Input Parameters:
517: +  eps - the eigenproblem solver context
518: -  restart - the percentage of the block of vectors to force a restart

520:    Options Database Key:
521: .  -eps_lobpcg_restart - Sets the restart parameter

523:    Notes:
524:    Allowed values are in the range [0.1,1.0]. The default is 0.6.

526:    Level: advanced

528: .seealso: EPSLOBPCGGetRestart()
529: @*/
530: PetscErrorCode EPSLOBPCGSetRestart(EPS eps,PetscReal restart)
531: {

537:   PetscTryMethod(eps,"EPSLOBPCGSetRestart_C",(EPS,PetscReal),(eps,restart));
538:   return(0);
539: }

543: static PetscErrorCode EPSLOBPCGGetRestart_LOBPCG(EPS eps,PetscReal *restart)
544: {
545:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

548:   *restart = ctx->restart;
549:   return(0);
550: }

554: /*@
555:    EPSLOBPCGGetRestart - Gets the restart parameter used in the LOBPCG method.

557:    Not Collective

559:    Input Parameter:
560: .  eps - the eigenproblem solver context

562:    Output Parameter:
563: .  restart - the restart parameter

565:    Level: advanced

567: .seealso: EPSLOBPCGSetRestart()
568: @*/
569: PetscErrorCode EPSLOBPCGGetRestart(EPS eps,PetscReal *restart)
570: {

576:   PetscUseMethod(eps,"EPSLOBPCGGetRestart_C",(EPS,PetscReal*),(eps,restart));
577:   return(0);
578: }

582: static PetscErrorCode EPSLOBPCGSetLocking_LOBPCG(EPS eps,PetscBool lock)
583: {
584:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

587:   ctx->lock = lock;
588:   return(0);
589: }

593: /*@
594:    EPSLOBPCGSetLocking - Choose between locking and non-locking variants of
595:    the LOBPCG method.

597:    Logically Collective on EPS

599:    Input Parameters:
600: +  eps  - the eigenproblem solver context
601: -  lock - true if the locking variant must be selected

603:    Options Database Key:
604: .  -eps_lobpcg_locking - Sets the locking flag

606:    Notes:
607:    This flag refers to soft locking (converged vectors within the current
608:    block iterate), since hard locking is always used (when nev is larger
609:    than the block size).

611:    Level: advanced

613: .seealso: EPSLOBPCGGetLocking()
614: @*/
615: PetscErrorCode EPSLOBPCGSetLocking(EPS eps,PetscBool lock)
616: {

622:   PetscTryMethod(eps,"EPSLOBPCGSetLocking_C",(EPS,PetscBool),(eps,lock));
623:   return(0);
624: }

628: static PetscErrorCode EPSLOBPCGGetLocking_LOBPCG(EPS eps,PetscBool *lock)
629: {
630:   EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;

633:   *lock = ctx->lock;
634:   return(0);
635: }

639: /*@
640:    EPSLOBPCGGetLocking - Gets the locking flag used in the LOBPCG method.

642:    Not Collective

644:    Input Parameter:
645: .  eps - the eigenproblem solver context

647:    Output Parameter:
648: .  lock - the locking flag

650:    Level: advanced

652: .seealso: EPSLOBPCGSetLocking()
653: @*/
654: PetscErrorCode EPSLOBPCGGetLocking(EPS eps,PetscBool *lock)
655: {

661:   PetscUseMethod(eps,"EPSLOBPCGGetLocking_C",(EPS,PetscBool*),(eps,lock));
662:   return(0);
663: }

667: PetscErrorCode EPSView_LOBPCG(EPS eps,PetscViewer viewer)
668: {
670:   EPS_LOBPCG     *ctx = (EPS_LOBPCG*)eps->data;
671:   PetscBool      isascii;

674:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
675:   if (isascii) {
676:     PetscViewerASCIIPrintf(viewer,"  LOBPCG: block size %D\n",ctx->bs);
677:     PetscViewerASCIIPrintf(viewer,"  LOBPCG: restart parameter=%g (using %d guard vectors)\n",(double)ctx->restart,(int)((1.0-ctx->restart)*ctx->bs));
678:     PetscViewerASCIIPrintf(viewer,"  LOBPCG: soft locking %sactivated\n",ctx->lock?"":"de");
679:   }
680:   return(0);
681: }

685: PetscErrorCode EPSSetFromOptions_LOBPCG(PetscOptionItems *PetscOptionsObject,EPS eps)
686: {
688:   PetscBool      lock,flg;
689:   PetscInt       bs;
690:   PetscReal      restart;
691:   KSP            ksp;

694:   PetscOptionsHead(PetscOptionsObject,"EPS LOBPCG Options");
695:   PetscOptionsInt("-eps_lobpcg_blocksize","LOBPCG block size","EPSLOBPCGSetBlockSize",20,&bs,&flg);
696:   if (flg) {
697:     EPSLOBPCGSetBlockSize(eps,bs);
698:   }
699:   PetscOptionsReal("-eps_lobpcg_restart","Percentage of the block of vectors to force a restart","EPSLOBPCGSetRestart",0.5,&restart,&flg);
700:   if (flg) {
701:     EPSLOBPCGSetRestart(eps,restart);
702:   }
703:   PetscOptionsBool("-eps_lobpcg_locking","Choose between locking and non-locking variants","EPSLOBPCGSetLocking",PETSC_TRUE,&lock,&flg);
704:   if (flg) {
705:     EPSLOBPCGSetLocking(eps,lock);
706:   }

708:   /* Set STPrecond as the default ST */
709:   if (!((PetscObject)eps->st)->type_name) {
710:     STSetType(eps->st,STPRECOND);
711:   }

713:   /* Set the default options of the KSP */
714:   STGetKSP(eps->st,&ksp);
715:   if (!((PetscObject)ksp)->type_name) {
716:     KSPSetType(ksp,KSPPREONLY);
717:   }
718:   PetscOptionsTail();
719:   return(0);
720: }

724: PetscErrorCode EPSDestroy_LOBPCG(EPS eps)
725: {

729:   PetscFree(eps->data);
730:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetBlockSize_C",NULL);
731:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetBlockSize_C",NULL);
732:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetRestart_C",NULL);
733:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetRestart_C",NULL);
734:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetLocking_C",NULL);
735:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetLocking_C",NULL);
736:   return(0);
737: }

741: PETSC_EXTERN PetscErrorCode EPSCreate_LOBPCG(EPS eps)
742: {
743:   EPS_LOBPCG     *lobpcg;

747:   PetscNewLog(eps,&lobpcg);
748:   eps->data = (void*)lobpcg;
749:   lobpcg->lock = PETSC_TRUE;

751:   eps->ops->setup          = EPSSetUp_LOBPCG;
752:   eps->ops->solve          = EPSSolve_LOBPCG;
753:   eps->ops->setfromoptions = EPSSetFromOptions_LOBPCG;
754:   eps->ops->destroy        = EPSDestroy_LOBPCG;
755:   eps->ops->view           = EPSView_LOBPCG;
756:   eps->ops->backtransform  = EPSBackTransform_Default;
757:   STSetType(eps->st,STPRECOND);
758:   STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
759:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetBlockSize_C",EPSLOBPCGSetBlockSize_LOBPCG);
760:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetBlockSize_C",EPSLOBPCGGetBlockSize_LOBPCG);
761:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetRestart_C",EPSLOBPCGSetRestart_LOBPCG);
762:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetRestart_C",EPSLOBPCGGetRestart_LOBPCG);
763:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetLocking_C",EPSLOBPCGSetLocking_LOBPCG);
764:   PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetLocking_C",EPSLOBPCGGetLocking_LOBPCG);
765:   return(0);
766: }