Actual source code: lobpcg.c
slepc-3.6.3 2016-03-29
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-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/epsimpl.h> /*I "slepceps.h" I*/
43: #include <slepc/private/dsimpl.h> /*I "slepcds.h" I*/
45: typedef struct {
46: PetscInt bs; /* block size */
47: PetscBool lock; /* soft locking active/inactive */
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: STSetUp(eps->st);
89: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond);
90: if (!precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"LOBPCG only works with precond ST");
92: EPSAllocateSolution(eps,0);
93: EPS_SetInnerProduct(eps);
94: DSSetType(eps->ds,DSGHEP);
95: DSAllocate(eps->ds,eps->mpd);
96: EPSSetWorkVecs(eps,1);
97: return(0);
98: }
102: PetscErrorCode EPSSolve_LOBPCG(EPS eps)
103: {
105: EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
106: PetscInt i,j,k,ld,nv,ini,kini,nmat,nc,nconv,bdone,its;
107: PetscReal norm;
108: PetscBool breakdown,countc;
109: Mat A,B,M;
110: Vec v,w=eps->work[0];
111: BV X,Y,Z,R,P,AX,BX;
114: DSGetLeadingDimension(eps->ds,&ld);
115: STGetNumMatrices(eps->st,&nmat);
116: STGetOperators(eps->st,0,&A);
117: if (nmat>1) { STGetOperators(eps->st,1,&B); }
118: else B = NULL;
120: /* 1. Allocate memory */
121: BVDuplicateResize(eps->V,3*ctx->bs,&Z);
122: BVDuplicateResize(eps->V,ctx->bs,&X);
123: BVDuplicateResize(eps->V,ctx->bs,&R);
124: BVDuplicateResize(eps->V,ctx->bs,&P);
125: BVDuplicateResize(eps->V,ctx->bs,&AX);
126: if (B) {
127: BVDuplicateResize(eps->V,ctx->bs,&BX);
128: }
129: nc = eps->nds;
130: if (nc>0 || eps->nev>ctx->bs) {
131: BVDuplicateResize(eps->V,nc+eps->nev,&Y);
132: }
133: if (nc>0) {
134: for (j=0;j<nc;j++) {
135: BVGetColumn(eps->V,-nc+j,&v);
136: BVInsertVec(Y,j,v);
137: BVRestoreColumn(eps->V,-nc+j,&v);
138: }
139: BVSetActiveColumns(Y,0,nc);
140: }
142: /* 2. Apply the constraints to the initial vectors */
143: kini = eps->nini;
144: while (kini<eps->ncv-2*ctx->bs) { /* Generate more initial vectors if necessary */
145: BVSetRandomColumn(eps->V,kini,eps->rand);
146: BVOrthogonalizeColumn(eps->V,kini,NULL,&norm,&breakdown);
147: if (norm>0.0 && !breakdown) {
148: BVScaleColumn(eps->V,kini,1.0/norm);
149: kini++;
150: }
151: }
152: nv = ctx->bs;
153: BVSetActiveColumns(eps->V,0,nv);
154: BVSetActiveColumns(Z,0,nv);
155: BVCopy(eps->V,Z);
156: BVCopy(Z,X);
158: /* 3. B-orthogonalize initial vectors */
159: if (B) {
160: BVMatMult(X,B,BX);
161: }
163: /* 4. Compute initial Ritz vectors */
164: BVMatMult(X,A,AX);
165: DSSetDimensions(eps->ds,nv,0,0,0);
166: DSGetMat(eps->ds,DS_MAT_A,&M);
167: BVMatProject(AX,NULL,X,M);
168: DSRestoreMat(eps->ds,DS_MAT_A,&M);
169: DSSetIdentity(eps->ds,DS_MAT_B);
170: DSSetState(eps->ds,DS_STATE_RAW);
171: DSSolve(eps->ds,eps->eigr,eps->eigi);
172: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
173: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
174: DSGetMat(eps->ds,DS_MAT_X,&M);
175: BVMultInPlace(X,M,0,nv);
176: BVMultInPlace(AX,M,0,nv);
177: DSRestoreMat(eps->ds,DS_MAT_X,&M);
179: /* 5. Initialize range of active iterates */
180: bdone = 0; /* finished blocks, the leading bdone*bs columns of V are eigenvectors */
181: nconv = 0; /* number of converged eigenvalues in the current block */
182: its = 0; /* iterations for the current block */
184: /* 6. Main loop */
185: while (eps->reason == EPS_CONVERGED_ITERATING) {
187: if (ctx->lock) {
188: BVSetActiveColumns(R,nconv,ctx->bs);
189: BVSetActiveColumns(AX,nconv,ctx->bs);
190: if (B) {
191: BVSetActiveColumns(BX,nconv,ctx->bs);
192: }
193: }
195: /* 7. Compute residuals */
196: DSGetMat(eps->ds,DS_MAT_A,&M);
197: BVCopy(AX,R);
198: if (B) {
199: BVMult(R,-1.0,1.0,BX,M);
200: } else {
201: BVMult(R,-1.0,1.0,X,M);
202: }
203: DSRestoreMat(eps->ds,DS_MAT_A,&M);
205: /* 8. Compute residual norms and update index set of active iterates */
206: ini = (ctx->lock)? nconv: 0;
207: k = ini;
208: countc = PETSC_TRUE;
209: for (j=ini;j<ctx->bs;j++) {
210: i = bdone*ctx->bs+j;
211: BVGetColumn(R,j,&v);
212: VecNorm(v,NORM_2,&norm);
213: BVRestoreColumn(R,j,&v);
214: (*eps->converged)(eps,eps->eigr[i],eps->eigi[i],norm,&eps->errest[i],eps->convergedctx);
215: if (countc) {
216: if (eps->errest[i] < eps->tol) k++;
217: else countc = PETSC_FALSE;
218: }
219: if (!countc && !eps->trackall) break;
220: }
221: nconv = k;
222: eps->nconv = bdone*ctx->bs + nconv;
223: if (eps->its+its) {
224: EPSMonitor(eps,eps->its+its,eps->nconv,eps->eigr,eps->eigi,eps->errest,(bdone+1)*ctx->bs);
225: }
226: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
227: else if (eps->its+its >= eps->max_it) {
228: eps->its += its;
229: eps->reason = EPS_DIVERGED_ITS;
230: }
231: if (eps->reason != EPS_CONVERGED_ITERATING || nconv == ctx->bs) {
232: BVSetActiveColumns(eps->V,bdone*ctx->bs,eps->nconv);
233: BVSetActiveColumns(Z,0,nconv);
234: BVSetActiveColumns(X,0,nconv);
235: BVCopy(X,eps->V);
236: BVCopy(X,Z);
237: }
238: if (eps->reason != EPS_CONVERGED_ITERATING) {
239: eps->its += its;
240: break;
241: } else if (nconv == ctx->bs) {
242: eps->its += its;
243: its = 0;
244: }
245: its++;
247: if (nconv == ctx->bs) { /* block finished: lock eigenvectors and compute new R */
249: /* extend constraints */
250: BVSetActiveColumns(Y,nc+bdone*ctx->bs,nc+(bdone+1)*ctx->bs);
251: BVCopy(Z,Y);
252: for (j=0;j<ctx->bs;j++) {
253: BVOrthogonalizeColumn(Y,nc+bdone*ctx->bs+j,NULL,&norm,&breakdown);
254: if (norm>0.0 && !breakdown) {
255: BVScaleColumn(Y,nc+bdone*ctx->bs+j,1.0/norm);
256: } else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Orthogonalization of constraints failed");
257: }
258: BVSetActiveColumns(Y,0,nc+(bdone+1)*ctx->bs);
260: bdone++;
261: nconv = 0;
263: /* set new initial vectors */
264: BVSetActiveColumns(eps->V,bdone*ctx->bs,(bdone+1)*ctx->bs);
265: BVCopy(eps->V,Z);
266: for (j=0;j<ctx->bs;j++) {
267: BVGetColumn(Z,j,&v);
268: BVOrthogonalizeVec(Y,v,NULL,NULL,NULL);
269: VecNormalize(v,NULL);
270: BVRestoreColumn(Z,j,&v);
271: }
272: BVSetActiveColumns(X,nconv,ctx->bs);
273: BVSetActiveColumns(AX,nconv,ctx->bs);
274: BVCopy(Z,X);
276: /* B-orthogonalize initial vectors */
277: if (B) {
278: BVOrthogonalize(X,NULL);
279: }
281: /* Compute initial Ritz vectors */
282: nv = ctx->bs;
283: BVMatMult(X,A,AX);
284: DSSetDimensions(eps->ds,nv,0,0,0);
285: DSGetMat(eps->ds,DS_MAT_A,&M);
286: BVMatProject(AX,NULL,X,M);
287: DSRestoreMat(eps->ds,DS_MAT_A,&M);
288: DSSetIdentity(eps->ds,DS_MAT_B);
289: DSSetState(eps->ds,DS_STATE_RAW);
290: DSSolve(eps->ds,eps->eigr+bdone*ctx->bs,eps->eigi);
291: DSSort(eps->ds,eps->eigr+bdone*ctx->bs,eps->eigi,NULL,NULL,NULL);
292: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
293: DSGetMat(eps->ds,DS_MAT_X,&M);
294: BVMultInPlace(X,M,0,nv);
295: BVMultInPlace(AX,M,0,nv);
296: DSRestoreMat(eps->ds,DS_MAT_X,&M);
298: EPSMonitor(eps,eps->its+its-1,eps->nconv,eps->eigr,eps->eigi,eps->errest,(bdone+1)*ctx->bs);
300: if (ctx->lock) {
301: BVSetActiveColumns(R,nconv,ctx->bs);
302: BVSetActiveColumns(AX,nconv,ctx->bs);
303: if (B) {
304: BVSetActiveColumns(BX,nconv,ctx->bs);
305: }
306: }
308: /* compute residuals */
309: DSGetMat(eps->ds,DS_MAT_A,&M);
310: BVCopy(AX,R);
311: if (B) {
312: BVMatMult(X,B,BX);
313: BVMult(R,-1.0,1.0,BX,M);
314: } else {
315: BVMult(R,-1.0,1.0,X,M);
316: }
317: DSRestoreMat(eps->ds,DS_MAT_A,&M);
318: }
320: ini = (ctx->lock)? nconv: 0;
321: if (ctx->lock) {
322: BVSetActiveColumns(R,nconv,ctx->bs);
323: BVSetActiveColumns(P,nconv,ctx->bs);
324: BVSetActiveColumns(AX,nconv,ctx->bs);
325: if (B) {
326: BVSetActiveColumns(BX,nconv,ctx->bs);
327: }
328: }
330: /* 9. Apply preconditioner to the residuals */
331: for (j=ini;j<ctx->bs;j++) {
332: BVGetColumn(R,j,&v);
333: STMatSolve(eps->st,v,w);
334: if (nc+bdone*ctx->bs>0) {
335: BVOrthogonalizeVec(Y,w,NULL,NULL,NULL);
336: VecNormalize(w,NULL);
337: }
338: VecCopy(w,v);
339: BVRestoreColumn(R,j,&v);
340: }
342: /* 11. B-orthonormalize preconditioned residuals */
343: BVOrthogonalize(R,NULL);
345: /* 13-16. B-orthonormalize conjugate directions */
346: if (its>1) {
347: BVOrthogonalize(P,NULL);
348: }
350: /* 17-23. Compute symmetric Gram matrices */
351: BVSetActiveColumns(Z,0,ctx->bs);
352: BVSetActiveColumns(X,0,ctx->bs);
353: BVCopy(X,Z);
354: BVSetActiveColumns(Z,ctx->bs,2*ctx->bs-ini);
355: BVCopy(R,Z);
356: if (its>1) {
357: BVSetActiveColumns(Z,2*ctx->bs-ini,3*ctx->bs-2*ini);
358: BVCopy(P,Z);
359: }
361: if (its>1) nv = 3*ctx->bs-2*ini;
362: else nv = 2*ctx->bs-ini;
364: BVSetActiveColumns(Z,0,nv);
365: DSSetDimensions(eps->ds,nv,0,0,0);
366: DSGetMat(eps->ds,DS_MAT_A,&M);
367: BVMatProject(Z,A,Z,M);
368: DSRestoreMat(eps->ds,DS_MAT_A,&M);
369: DSGetMat(eps->ds,DS_MAT_B,&M);
370: if (B) {
371: BVMatProject(Z,B,Z,M);
372: } else {
373: BVDot(Z,Z,M);
374: }
375: DSRestoreMat(eps->ds,DS_MAT_B,&M);
376:
377: /* 24. Solve the generalized eigenvalue problem */
378: DSSetState(eps->ds,DS_STATE_RAW);
379: DSSolve(eps->ds,eps->eigr+bdone*ctx->bs,eps->eigi);
380: DSSort(eps->ds,eps->eigr+bdone*ctx->bs,eps->eigi,NULL,NULL,NULL);
381: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
382:
383: /* 25-33. Compute Ritz vectors */
384: DSGetMat(eps->ds,DS_MAT_X,&M);
385: BVSetActiveColumns(Z,ctx->bs,nv);
386: if (ctx->lock) {
387: BVSetActiveColumns(P,0,ctx->bs);
388: }
389: BVMult(P,1.0,0.0,Z,M);
390: BVCopy(P,X);
391: if (ctx->lock) {
392: BVSetActiveColumns(P,nconv,ctx->bs);
393: }
394: BVSetActiveColumns(Z,0,ctx->bs);
395: BVMult(X,1.0,1.0,Z,M);
396: if (ctx->lock) {
397: BVSetActiveColumns(X,nconv,ctx->bs);
398: }
399: BVMatMult(X,A,AX);
400: if (B) {
401: BVMatMult(X,B,BX);
402: }
403: DSRestoreMat(eps->ds,DS_MAT_X,&M);
404: }
406: BVDestroy(&Z);
407: BVDestroy(&X);
408: BVDestroy(&R);
409: BVDestroy(&P);
410: BVDestroy(&AX);
411: if (B) {
412: BVDestroy(&BX);
413: }
414: if (nc>0 || eps->nev>ctx->bs) {
415: BVDestroy(&Y);
416: }
417: return(0);
418: }
422: static PetscErrorCode EPSLOBPCGSetBlockSize_LOBPCG(EPS eps,PetscInt bs)
423: {
424: EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
427: ctx->bs = bs;
428: return(0);
429: }
433: /*@
434: EPSLOBPCGSetBlockSize - Sets the block size of the LOBPCG method.
436: Logically Collective on EPS
438: Input Parameters:
439: + eps - the eigenproblem solver context
440: - bs - the block size
442: Options Database Key:
443: . -eps_lobpcg_blocksize - Sets the block size
445: Level: advanced
447: .seealso: EPSLOBPCGGetBlockSize()
448: @*/
449: PetscErrorCode EPSLOBPCGSetBlockSize(EPS eps,PetscInt bs)
450: {
456: PetscTryMethod(eps,"EPSLOBPCGSetBlockSize_C",(EPS,PetscInt),(eps,bs));
457: return(0);
458: }
462: static PetscErrorCode EPSLOBPCGGetBlockSize_LOBPCG(EPS eps,PetscInt *bs)
463: {
464: EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
467: *bs = ctx->bs;
468: return(0);
469: }
473: /*@
474: EPSLOBPCGGetBlockSize - Gets the block size used in the LOBPCG method.
476: Not Collective
478: Input Parameter:
479: . eps - the eigenproblem solver context
481: Output Parameter:
482: . bs - the block size
484: Level: advanced
486: .seealso: EPSLOBPCGSetBlockSize()
487: @*/
488: PetscErrorCode EPSLOBPCGGetBlockSize(EPS eps,PetscInt *bs)
489: {
495: PetscTryMethod(eps,"EPSLOBPCGGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
496: return(0);
497: }
501: static PetscErrorCode EPSLOBPCGSetLocking_LOBPCG(EPS eps,PetscBool lock)
502: {
503: EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
506: ctx->lock = lock;
507: return(0);
508: }
512: /*@
513: EPSLOBPCGSetLocking - Choose between locking and non-locking variants of
514: the LOBPCG method.
516: Logically Collective on EPS
518: Input Parameters:
519: + eps - the eigenproblem solver context
520: - lock - true if the locking variant must be selected
522: Options Database Key:
523: . -eps_lobpcg_locking - Sets the locking flag
525: Notes:
526: This flag refers to soft locking (converged vectors within the current
527: block iterate), since hard locking is always used (when nev is larger
528: than the block size).
530: Level: advanced
532: .seealso: EPSLOBPCGGetLocking()
533: @*/
534: PetscErrorCode EPSLOBPCGSetLocking(EPS eps,PetscBool lock)
535: {
541: PetscTryMethod(eps,"EPSLOBPCGSetLocking_C",(EPS,PetscBool),(eps,lock));
542: return(0);
543: }
547: static PetscErrorCode EPSLOBPCGGetLocking_LOBPCG(EPS eps,PetscBool *lock)
548: {
549: EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
552: *lock = ctx->lock;
553: return(0);
554: }
558: /*@
559: EPSLOBPCGGetLocking - Gets the locking flag used in the LOBPCG method.
561: Not Collective
563: Input Parameter:
564: . eps - the eigenproblem solver context
566: Output Parameter:
567: . lock - the locking flag
569: Level: advanced
571: .seealso: EPSLOBPCGSetLocking()
572: @*/
573: PetscErrorCode EPSLOBPCGGetLocking(EPS eps,PetscBool *lock)
574: {
580: PetscTryMethod(eps,"EPSLOBPCGGetLocking_C",(EPS,PetscBool*),(eps,lock));
581: return(0);
582: }
586: PetscErrorCode EPSView_LOBPCG(EPS eps,PetscViewer viewer)
587: {
589: EPS_LOBPCG *ctx = (EPS_LOBPCG*)eps->data;
590: PetscBool isascii;
593: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
594: if (isascii) {
595: PetscViewerASCIIPrintf(viewer," LOBPCG: block size %D\n",ctx->bs);
596: PetscViewerASCIIPrintf(viewer," LOBPCG: soft locking %sactivated\n",ctx->lock?"":"de");
597: }
598: return(0);
599: }
603: PetscErrorCode EPSSetFromOptions_LOBPCG(PetscOptions *PetscOptionsObject,EPS eps)
604: {
606: PetscBool lock,flg;
607: PetscInt bs;
608: KSP ksp;
611: PetscOptionsHead(PetscOptionsObject,"EPS LOBPCG Options");
612: PetscOptionsInt("-eps_lobpcg_blocksize","LOBPCG block size","EPSLOBPCGSetBlockSize",20,&bs,&flg);
613: if (flg) {
614: EPSLOBPCGSetBlockSize(eps,bs);
615: }
616: PetscOptionsBool("-eps_lobpcg_locking","Choose between locking and non-locking variants","EPSLOBPCGSetLocking",PETSC_TRUE,&lock,&flg);
617: if (flg) {
618: EPSLOBPCGSetLocking(eps,lock);
619: }
621: /* Set STPrecond as the default ST */
622: if (!((PetscObject)eps->st)->type_name) {
623: STSetType(eps->st,STPRECOND);
624: }
626: /* Set the default options of the KSP */
627: STGetKSP(eps->st,&ksp);
628: if (!((PetscObject)ksp)->type_name) {
629: KSPSetType(ksp,KSPPREONLY);
630: }
631: PetscOptionsTail();
632: return(0);
633: }
637: PetscErrorCode EPSDestroy_LOBPCG(EPS eps)
638: {
642: PetscFree(eps->data);
643: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetBlockSize_C",NULL);
644: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetBlockSize_C",NULL);
645: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetLocking_C",NULL);
646: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetLocking_C",NULL);
647: return(0);
648: }
652: PETSC_EXTERN PetscErrorCode EPSCreate_LOBPCG(EPS eps)
653: {
654: EPS_LOBPCG *lobpcg;
658: PetscNewLog(eps,&lobpcg);
659: eps->data = (void*)lobpcg;
660: lobpcg->lock = PETSC_TRUE;
662: eps->ops->setup = EPSSetUp_LOBPCG;
663: eps->ops->solve = EPSSolve_LOBPCG;
664: eps->ops->setfromoptions = EPSSetFromOptions_LOBPCG;
665: eps->ops->destroy = EPSDestroy_LOBPCG;
666: eps->ops->view = EPSView_LOBPCG;
667: eps->ops->backtransform = EPSBackTransform_Default;
668: STSetType(eps->st,STPRECOND);
669: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
670: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetBlockSize_C",EPSLOBPCGSetBlockSize_LOBPCG);
671: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetBlockSize_C",EPSLOBPCGGetBlockSize_LOBPCG);
672: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGSetLocking_C",EPSLOBPCGSetLocking_LOBPCG);
673: PetscObjectComposeFunction((PetscObject)eps,"EPSLOBPCGGetLocking_C",EPSLOBPCGGetLocking_LOBPCG);
674: return(0);
675: }