Actual source code: blopex.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: This file implements a wrapper to the BLOPEX package
12: */
14: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
15: #include "slepc-interface.h"
16: #include <blopex_lobpcg.h>
17: #include <blopex_interpreter.h>
18: #include <blopex_multivector.h>
19: #include <blopex_temp_multivector.h>
21: PetscInt slepc_blopex_useconstr = -1;
23: typedef struct {
24: lobpcg_Tolerance tol;
25: lobpcg_BLASLAPACKFunctions blap_fn;
26: mv_InterfaceInterpreter ii;
27: ST st;
28: Vec w;
29: PetscInt bs; /* block size */
30: } EPS_BLOPEX;
32: static void Precond_FnSingleVector(void *data,void *x,void *y)
33: {
35: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
36: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
37: KSP ksp;
40: STGetKSP(blopex->st,&ksp);CHKERRABORT(comm,ierr);
41: KSPSolve(ksp,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
42: PetscFunctionReturnVoid();
43: }
45: static void Precond_FnMultiVector(void *data,void *x,void *y)
46: {
47: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
50: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
51: PetscFunctionReturnVoid();
52: }
54: static void OperatorASingleVector(void *data,void *x,void *y)
55: {
57: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
58: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
59: Mat A,B;
60: PetscScalar sigma;
61: PetscInt nmat;
64: STGetNumMatrices(blopex->st,&nmat);CHKERRABORT(comm,ierr);
65: STGetMatrix(blopex->st,0,&A);CHKERRABORT(comm,ierr);
66: if (nmat>1) { STGetMatrix(blopex->st,1,&B);CHKERRABORT(comm,ierr); }
67: MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
68: STGetShift(blopex->st,&sigma);CHKERRABORT(comm,ierr);
69: if (sigma != 0.0) {
70: if (nmat>1) {
71: MatMult(B,(Vec)x,blopex->w);CHKERRABORT(comm,ierr);
72: } else {
73: VecCopy((Vec)x,blopex->w);CHKERRABORT(comm,ierr);
74: }
75: VecAXPY((Vec)y,-sigma,blopex->w);CHKERRABORT(comm,ierr);
76: }
77: PetscFunctionReturnVoid();
78: }
80: static void OperatorAMultiVector(void *data,void *x,void *y)
81: {
82: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
85: blopex->ii.Eval(OperatorASingleVector,data,x,y);
86: PetscFunctionReturnVoid();
87: }
89: static void OperatorBSingleVector(void *data,void *x,void *y)
90: {
92: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
93: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
94: Mat B;
97: STGetMatrix(blopex->st,1,&B);CHKERRABORT(comm,ierr);
98: MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
99: PetscFunctionReturnVoid();
100: }
102: static void OperatorBMultiVector(void *data,void *x,void *y)
103: {
104: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
107: blopex->ii.Eval(OperatorBSingleVector,data,x,y);
108: PetscFunctionReturnVoid();
109: }
111: PetscErrorCode EPSSetDimensions_BLOPEX(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
112: {
114: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
115: PetscInt k;
118: k = ((eps->nev-1)/ctx->bs+1)*ctx->bs;
119: if (*ncv) { /* ncv set */
120: if (*ncv<k) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv is not sufficiently large");
121: } else *ncv = k;
122: if (!*mpd) *mpd = *ncv;
123: else { PetscInfo(eps,"Warning: given value of mpd ignored\n"); }
124: return(0);
125: }
127: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
128: {
129: #if defined(PETSC_MISSING_LAPACK_POTRF) || defined(PETSC_MISSING_LAPACK_SYGV)
131: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF/SYGV - Lapack routines are unavailable");
132: #else
134: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
135: PetscBool istrivial,flg;
138: if (!eps->ishermitian || (eps->isgeneralized && !eps->ispositive)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"blopex only works for Hermitian problems");
139: if (!blopex->bs) blopex->bs = PetscMin(16,eps->nev);
140: EPSSetDimensions_BLOPEX(eps,eps->nev,&eps->ncv,&eps->mpd);
141: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
142: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
143: if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
144: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
145: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
146: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
147: RGIsTrivial(eps->rg,&istrivial);
148: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
150: blopex->st = eps->st;
152: if (eps->converged == EPSConvergedRelative) {
153: blopex->tol.absolute = 0.0;
154: blopex->tol.relative = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
155: } else if (eps->converged == EPSConvergedAbsolute) {
156: blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
157: blopex->tol.relative = 0.0;
158: } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Convergence test not supported in this solver");
160: SLEPCSetupInterpreter(&blopex->ii);
162: /* allocate memory */
163: if (!eps->V) { EPSGetBV(eps,&eps->V); }
164: PetscObjectTypeCompareAny((PetscObject)eps->V,&flg,BVVECS,BVCONTIGUOUS,"");
165: if (!flg) { /* blopex only works with BVVECS or BVCONTIGUOUS */
166: BVSetType(eps->V,BVCONTIGUOUS);
167: }
168: EPSAllocateSolution(eps,0);
169: BVCreateVec(eps->V,&blopex->w);
170: PetscLogObjectParent((PetscObject)eps,(PetscObject)blopex->w);
172: #if defined(PETSC_USE_COMPLEX)
173: blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
174: blopex->blap_fn.zhegv = PETSC_zsygv_interface;
175: #else
176: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
177: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
178: #endif
179: return(0);
180: #endif
181: }
183: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
184: {
185: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
186: PetscScalar sigma,*eigr=NULL;
187: PetscReal *errest=NULL;
188: int i,j,info,its,nconv;
189: double *residhist=NULL;
190: PetscErrorCode ierr;
191: mv_MultiVectorPtr eigenvectors,constraints;
192: #if defined(PETSC_USE_COMPLEX)
193: komplex *lambda=NULL,*lambdahist=NULL;
194: #else
195: double *lambda=NULL,*lambdahist=NULL;
196: #endif
199: STGetShift(eps->st,&sigma);
200: PetscMalloc1(blopex->bs,&lambda);
201: if (eps->numbermonitors>0) {
202: PetscMalloc4(blopex->bs*(eps->max_it+1),&lambdahist,eps->ncv,&eigr,blopex->bs*(eps->max_it+1),&residhist,eps->ncv,&errest);
203: }
205: /* Complete the initial basis with random vectors */
206: for (i=eps->nini;i<eps->ncv;i++) {
207: BVSetRandomColumn(eps->V,i);
208: }
210: while (eps->reason == EPS_CONVERGED_ITERATING) {
212: /* Create multivector of constraints from leading columns of V */
213: PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,1);
214: BVSetActiveColumns(eps->V,0,eps->nconv);
215: constraints = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds+eps->nconv,eps->V);
217: /* Create multivector where eigenvectors of this run will be stored */
218: PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,0);
219: BVSetActiveColumns(eps->V,eps->nconv,eps->nconv+blopex->bs);
220: eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,blopex->bs,eps->V);
222: #if defined(PETSC_USE_COMPLEX)
223: info = lobpcg_solve_complex(eigenvectors,blopex,OperatorAMultiVector,
224: eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
225: blopex,Precond_FnMultiVector,constraints,
226: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
227: lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs);
228: #else
229: info = lobpcg_solve_double(eigenvectors,blopex,OperatorAMultiVector,
230: eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
231: blopex,Precond_FnMultiVector,constraints,
232: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
233: lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs);
234: #endif
235: if (info>0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"BLOPEX failed with exit code=%d",info);
236: mv_MultiVectorDestroy(constraints);
237: mv_MultiVectorDestroy(eigenvectors);
239: for (j=0;j<blopex->bs;j++) {
240: #if defined(PETSC_USE_COMPLEX)
241: eps->eigr[eps->nconv+j] = lambda[j].real+PETSC_i*lambda[j].imag;
242: #else
243: eps->eigr[eps->nconv+j] = lambda[j];
244: #endif
245: }
247: if (eps->numbermonitors>0) {
248: for (i=0;i<its;i++) {
249: nconv = 0;
250: for (j=0;j<blopex->bs;j++) {
251: #if defined(PETSC_USE_COMPLEX)
252: eigr[eps->nconv+j] = lambdahist[j+i*blopex->bs].real+PETSC_i*lambdahist[j+i*blopex->bs].imag;
253: #else
254: eigr[eps->nconv+j] = lambdahist[j+i*blopex->bs];
255: #endif
256: errest[eps->nconv+j] = residhist[j+i*blopex->bs];
257: if (residhist[j+i*blopex->bs]<=eps->tol) nconv++;
258: }
259: EPSMonitor(eps,eps->its+i,eps->nconv+nconv,eigr,eps->eigi,errest,eps->nconv+blopex->bs);
260: }
261: }
263: eps->its += its;
264: if (info==-1) {
265: eps->reason = EPS_DIVERGED_ITS;
266: break;
267: } else {
268: for (i=0;i<blopex->bs;i++) {
269: if (sigma != 0.0) eps->eigr[eps->nconv+i] += sigma;
270: }
271: eps->nconv += blopex->bs;
272: if (eps->nconv>=eps->nev) eps->reason = EPS_CONVERGED_TOL;
273: }
274: }
276: PetscFree(lambda);
277: if (eps->numbermonitors>0) {
278: PetscFree4(lambdahist,eigr,residhist,errest);
279: }
280: return(0);
281: }
283: static PetscErrorCode EPSBLOPEXSetBlockSize_BLOPEX(EPS eps,PetscInt bs)
284: {
285: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
288: ctx->bs = bs;
289: return(0);
290: }
292: /*@
293: EPSBLOPEXSetBlockSize - Sets the block size of the BLOPEX solver.
295: Logically Collective on EPS
297: Input Parameters:
298: + eps - the eigenproblem solver context
299: - bs - the block size
301: Options Database Key:
302: . -eps_blopex_blocksize - Sets the block size
304: Level: advanced
306: .seealso: EPSBLOPEXGetBlockSize()
307: @*/
308: PetscErrorCode EPSBLOPEXSetBlockSize(EPS eps,PetscInt bs)
309: {
315: PetscTryMethod(eps,"EPSBLOPEXSetBlockSize_C",(EPS,PetscInt),(eps,bs));
316: return(0);
317: }
319: static PetscErrorCode EPSBLOPEXGetBlockSize_BLOPEX(EPS eps,PetscInt *bs)
320: {
321: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
324: *bs = ctx->bs;
325: return(0);
326: }
328: /*@
329: EPSBLOPEXGetBlockSize - Gets the block size used in the BLOPEX solver.
331: Not Collective
333: Input Parameter:
334: . eps - the eigenproblem solver context
336: Output Parameter:
337: . bs - the block size
339: Level: advanced
341: .seealso: EPSBLOPEXSetBlockSize()
342: @*/
343: PetscErrorCode EPSBLOPEXGetBlockSize(EPS eps,PetscInt *bs)
344: {
350: PetscUseMethod(eps,"EPSBLOPEXGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
351: return(0);
352: }
354: PetscErrorCode EPSReset_BLOPEX(EPS eps)
355: {
357: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
360: VecDestroy(&blopex->w);
361: return(0);
362: }
364: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
365: {
369: LOBPCG_DestroyRandomContext();
370: PetscFree(eps->data);
371: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",NULL);
372: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",NULL);
373: return(0);
374: }
376: PetscErrorCode EPSView_BLOPEX(EPS eps,PetscViewer viewer)
377: {
379: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
380: PetscBool isascii;
383: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
384: if (isascii) {
385: PetscViewerASCIIPrintf(viewer," block size %D\n",ctx->bs);
386: }
387: return(0);
388: }
390: PetscErrorCode EPSSetFromOptions_BLOPEX(PetscOptionItems *PetscOptionsObject,EPS eps)
391: {
393: PetscBool flg;
394: PetscInt bs;
397: PetscOptionsHead(PetscOptionsObject,"EPS BLOPEX Options");
399: PetscOptionsInt("-eps_blopex_blocksize","Block size","EPSBLOPEXSetBlockSize",20,&bs,&flg);
400: if (flg) { EPSBLOPEXSetBlockSize(eps,bs); }
402: PetscOptionsTail();
404: LOBPCG_SetFromOptionsRandomContext();
405: return(0);
406: }
408: SLEPC_EXTERN PetscErrorCode EPSCreate_BLOPEX(EPS eps)
409: {
410: EPS_BLOPEX *ctx;
414: PetscNewLog(eps,&ctx);
415: eps->data = (void*)ctx;
417: eps->categ = EPS_CATEGORY_PRECOND;
419: eps->ops->solve = EPSSolve_BLOPEX;
420: eps->ops->setup = EPSSetUp_BLOPEX;
421: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
422: eps->ops->destroy = EPSDestroy_BLOPEX;
423: eps->ops->reset = EPSReset_BLOPEX;
424: eps->ops->view = EPSView_BLOPEX;
425: eps->ops->backtransform = EPSBackTransform_Default;
426: eps->ops->setdefaultst = EPSSetDefaultST_Precond;
428: LOBPCG_InitRandomContext(PetscObjectComm((PetscObject)eps),NULL);
429: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",EPSBLOPEXSetBlockSize_BLOPEX);
430: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",EPSBLOPEXGetBlockSize_BLOPEX);
431: if (slepc_blopex_useconstr < 0) { PetscObjectComposedDataRegister(&slepc_blopex_useconstr); }
432: return(0);
433: }