Actual source code: feast.c
slepc-main 2024-11-15
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: This file implements a wrapper to the FEAST solver in MKL
12: */
14: #include <petscsys.h>
15: #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
16: #define MKL_ILP64
17: #endif
18: #include <mkl.h>
19: #include <slepc/private/epsimpl.h>
21: #if defined(PETSC_USE_COMPLEX)
22: # if defined(PETSC_USE_REAL_SINGLE)
23: # define FEAST_RCI cfeast_hrci
24: # define SCALAR_CAST (MKL_Complex8*)
25: # else
26: # define FEAST_RCI zfeast_hrci
27: # define SCALAR_CAST (MKL_Complex16*)
28: # endif
29: #else
30: # if defined(PETSC_USE_REAL_SINGLE)
31: # define FEAST_RCI sfeast_srci
32: # else
33: # define FEAST_RCI dfeast_srci
34: # endif
35: # define SCALAR_CAST
36: #endif
38: typedef struct {
39: PetscInt npoints; /* number of contour points */
40: PetscScalar *work1,*Aq,*Bq; /* workspace */
41: #if defined(PETSC_USE_REAL_SINGLE)
42: MKL_Complex8 *work2;
43: #else
44: MKL_Complex16 *work2;
45: #endif
46: } EPS_FEAST;
48: static PetscErrorCode EPSSetUp_FEAST(EPS eps)
49: {
50: PetscInt ncv;
51: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
52: PetscMPIInt size;
54: PetscFunctionBegin;
55: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)eps),&size));
56: PetscCheck(size==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The FEAST interface is supported for sequential runs only");
57: EPSCheckHermitianDefinite(eps);
58: EPSCheckSinvertCayley(eps);
59: if (eps->ncv!=PETSC_DETERMINE) {
60: PetscCheck(eps->ncv>=eps->nev+2,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The value of ncv must be at least nev+2");
61: } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n); /* set default value of ncv */
62: if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
63: if (eps->max_it==PETSC_DETERMINE) eps->max_it = 20;
64: if (!eps->which) eps->which = EPS_ALL;
65: PetscCheck(eps->which==EPS_ALL && eps->inta!=eps->intb,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver must be used with a computational interval");
66: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);
67: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
69: if (!ctx->npoints) ctx->npoints = 8;
71: ncv = eps->ncv;
72: PetscCall(PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq));
73: PetscCall(PetscMalloc4(eps->nloc*ncv,&ctx->work1,eps->nloc*ncv,&ctx->work2,ncv*ncv,&ctx->Aq,ncv*ncv,&ctx->Bq));
75: PetscCall(EPSAllocateSolution(eps,0));
76: PetscCall(EPSSetWorkVecs(eps,2));
77: PetscFunctionReturn(PETSC_SUCCESS);
78: }
80: static PetscErrorCode EPSSolve_FEAST(EPS eps)
81: {
82: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
83: MKL_INT fpm[128],ijob,n,ncv,nconv,loop,info;
84: PetscReal *evals,epsout=0.0;
85: PetscInt i,k,nmat,ld;
86: PetscScalar *pV,*pz,*X=NULL;
87: Vec x,y,w=eps->work[0],z=eps->work[1];
88: Mat A,B;
89: #if defined(PETSC_USE_REAL_SINGLE)
90: MKL_Complex8 Ze;
91: #else
92: MKL_Complex16 Ze;
93: #endif
95: PetscFunctionBegin;
96: #if !defined(PETSC_HAVE_64BIT_BLAS_INDICES)
97: PetscCall(PetscBLASIntCast(eps->ncv,&ncv));
98: PetscCall(PetscBLASIntCast(eps->nloc,&n));
99: #else
100: ncv = eps->ncv;
101: n = eps->nloc;
102: #endif
104: /* parameters */
105: feastinit(fpm);
106: fpm[0] = (eps->numbermonitors>0)? 1: 0; /* runtime comments */
107: fpm[1] = (MKL_INT)ctx->npoints; /* contour points */
108: #if !defined(PETSC_USE_REAL_SINGLE)
109: fpm[2] = -PetscLog10Real(eps->tol); /* tolerance for trace */
110: #endif
111: fpm[3] = (MKL_INT)eps->max_it; /* refinement loops */
112: fpm[5] = 1; /* second stopping criterion */
113: #if defined(PETSC_USE_REAL_SINGLE)
114: fpm[6] = -PetscLog10Real(eps->tol); /* tolerance for trace */
115: #endif
117: PetscCall(PetscMalloc1(eps->ncv,&evals));
118: PetscCall(BVGetLeadingDimension(eps->V,&ld));
119: PetscCall(BVGetArray(eps->V,&pV));
120: if (ld==n) X = pV;
121: else PetscCall(PetscMalloc1(eps->ncv*n,&X));
123: ijob = -1; /* first call to reverse communication interface */
124: PetscCall(STGetNumMatrices(eps->st,&nmat));
125: PetscCall(STGetMatrix(eps->st,0,&A));
126: if (nmat>1) PetscCall(STGetMatrix(eps->st,1,&B));
127: else B = NULL;
128: PetscCall(MatCreateVecsEmpty(A,&x,&y));
130: do {
132: FEAST_RCI(&ijob,&n,&Ze,SCALAR_CAST ctx->work1,ctx->work2,SCALAR_CAST ctx->Aq,SCALAR_CAST ctx->Bq,fpm,&epsout,&loop,&eps->inta,&eps->intb,&ncv,evals,SCALAR_CAST X,&nconv,eps->errest,&info);
134: PetscCheck(ncv==eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"FEAST changed value of ncv to %d",(int)ncv);
135: if (ijob == 10) {
136: /* set new quadrature point */
137: PetscCall(STSetShift(eps->st,Ze.real));
138: } else if (ijob == 20) {
139: /* use same quadrature point and factorization for transpose solve */
140: } else if (ijob == 11 || ijob == 21) {
141: /* linear solve (A-sigma*B)\work2, overwrite work2 */
142: for (k=0;k<ncv;k++) {
143: PetscCall(VecGetArray(z,&pz));
144: #if defined(PETSC_USE_COMPLEX)
145: for (i=0;i<eps->nloc;i++) pz[i] = PetscCMPLX(ctx->work2[eps->nloc*k+i].real,ctx->work2[eps->nloc*k+i].imag);
146: #else
147: for (i=0;i<eps->nloc;i++) pz[i] = ctx->work2[eps->nloc*k+i].real;
148: #endif
149: PetscCall(VecRestoreArray(z,&pz));
150: if (ijob == 11) PetscCall(STMatSolve(eps->st,z,w));
151: else PetscCall(STMatSolveHermitianTranspose(eps->st,z,w));
152: PetscCall(VecGetArray(w,&pz));
153: #if defined(PETSC_USE_COMPLEX)
154: for (i=0;i<eps->nloc;i++) {
155: ctx->work2[eps->nloc*k+i].real = PetscRealPart(pz[i]);
156: ctx->work2[eps->nloc*k+i].imag = PetscImaginaryPart(pz[i]);
157: }
158: #else
159: for (i=0;i<eps->nloc;i++) ctx->work2[eps->nloc*k+i].real = pz[i];
160: #endif
161: PetscCall(VecRestoreArray(w,&pz));
162: }
163: } else if (ijob == 30 || ijob == 40) {
164: /* multiplication A*V or B*V, result in work1 */
165: for (k=fpm[23]-1;k<fpm[23]+fpm[24]-1;k++) {
166: PetscCall(VecPlaceArray(x,&X[k*eps->nloc]));
167: PetscCall(VecPlaceArray(y,&ctx->work1[k*eps->nloc]));
168: if (ijob == 30) PetscCall(MatMult(A,x,y));
169: else if (nmat>1) PetscCall(MatMult(B,x,y));
170: else PetscCall(VecCopy(x,y));
171: PetscCall(VecResetArray(x));
172: PetscCall(VecResetArray(y));
173: }
174: } else PetscCheck(ijob==0 || ijob==-2,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in FEAST reverse communication interface (ijob=%d)",(int)ijob);
176: } while (ijob);
178: eps->reason = EPS_CONVERGED_TOL;
179: eps->its = loop;
180: eps->nconv = nconv;
181: if (info) {
182: switch (info) {
183: case 1: /* No eigenvalue has been found in the proposed search interval */
184: eps->nconv = 0;
185: break;
186: case 2: /* FEAST did not converge "yet" */
187: eps->reason = EPS_DIVERGED_ITS;
188: break;
189: default:
190: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by FEAST (%d)",(int)info);
191: }
192: }
194: for (i=0;i<eps->nconv;i++) eps->eigr[i] = evals[i];
195: if (ld!=n) {
196: for (i=0;i<eps->nconv;i++) PetscCall(PetscArraycpy(pV+i*ld,X+i*n,n));
197: PetscCall(PetscFree(X));
198: }
199: PetscCall(BVRestoreArray(eps->V,&pV));
200: PetscCall(VecDestroy(&x));
201: PetscCall(VecDestroy(&y));
202: PetscCall(PetscFree(evals));
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: static PetscErrorCode EPSReset_FEAST(EPS eps)
207: {
208: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
210: PetscFunctionBegin;
211: PetscCall(PetscFree4(ctx->work1,ctx->work2,ctx->Aq,ctx->Bq));
212: PetscFunctionReturn(PETSC_SUCCESS);
213: }
215: static PetscErrorCode EPSDestroy_FEAST(EPS eps)
216: {
217: PetscFunctionBegin;
218: PetscCall(PetscFree(eps->data));
219: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",NULL));
220: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",NULL));
221: PetscFunctionReturn(PETSC_SUCCESS);
222: }
224: static PetscErrorCode EPSSetFromOptions_FEAST(EPS eps,PetscOptionItems *PetscOptionsObject)
225: {
226: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
227: PetscInt n;
228: PetscBool flg;
230: PetscFunctionBegin;
231: PetscOptionsHeadBegin(PetscOptionsObject,"EPS FEAST Options");
233: n = ctx->npoints;
234: PetscCall(PetscOptionsInt("-eps_feast_num_points","Number of contour integration points","EPSFEASTSetNumPoints",n,&n,&flg));
235: if (flg) PetscCall(EPSFEASTSetNumPoints(eps,n));
237: PetscOptionsHeadEnd();
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: static PetscErrorCode EPSView_FEAST(EPS eps,PetscViewer viewer)
242: {
243: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
244: PetscBool isascii;
246: PetscFunctionBegin;
247: PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
248: if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer," number of contour integration points=%" PetscInt_FMT "\n",ctx->npoints));
249: PetscFunctionReturn(PETSC_SUCCESS);
250: }
252: static PetscErrorCode EPSSetDefaultST_FEAST(EPS eps)
253: {
254: PetscFunctionBegin;
255: if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSINVERT));
256: PetscFunctionReturn(PETSC_SUCCESS);
257: }
259: static PetscErrorCode EPSFEASTSetNumPoints_FEAST(EPS eps,PetscInt npoints)
260: {
261: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
263: PetscFunctionBegin;
264: if (npoints == PETSC_DEFAULT || npoints == PETSC_DECIDE) ctx->npoints = 8;
265: else ctx->npoints = npoints;
266: PetscFunctionReturn(PETSC_SUCCESS);
267: }
269: /*@
270: EPSFEASTSetNumPoints - Sets the number of contour integration points for
271: the FEAST package.
273: Logically Collective
275: Input Parameters:
276: + eps - the eigenproblem solver context
277: - npoints - number of contour integration points
279: Options Database Key:
280: . -eps_feast_num_points - Sets the number of points
282: Level: advanced
284: .seealso: EPSFEASTGetNumPoints()
285: @*/
286: PetscErrorCode EPSFEASTSetNumPoints(EPS eps,PetscInt npoints)
287: {
288: PetscFunctionBegin;
291: PetscTryMethod(eps,"EPSFEASTSetNumPoints_C",(EPS,PetscInt),(eps,npoints));
292: PetscFunctionReturn(PETSC_SUCCESS);
293: }
295: static PetscErrorCode EPSFEASTGetNumPoints_FEAST(EPS eps,PetscInt *npoints)
296: {
297: EPS_FEAST *ctx = (EPS_FEAST*)eps->data;
299: PetscFunctionBegin;
300: *npoints = ctx->npoints;
301: PetscFunctionReturn(PETSC_SUCCESS);
302: }
304: /*@
305: EPSFEASTGetNumPoints - Gets the number of contour integration points for
306: the FEAST package.
308: Not Collective
310: Input Parameter:
311: . eps - the eigenproblem solver context
313: Output Parameter:
314: . npoints - number of contour integration points
316: Level: advanced
318: .seealso: EPSFEASTSetNumPoints()
319: @*/
320: PetscErrorCode EPSFEASTGetNumPoints(EPS eps,PetscInt *npoints)
321: {
322: PetscFunctionBegin;
324: PetscAssertPointer(npoints,2);
325: PetscUseMethod(eps,"EPSFEASTGetNumPoints_C",(EPS,PetscInt*),(eps,npoints));
326: PetscFunctionReturn(PETSC_SUCCESS);
327: }
329: SLEPC_EXTERN PetscErrorCode EPSCreate_FEAST(EPS eps)
330: {
331: EPS_FEAST *ctx;
333: PetscFunctionBegin;
334: PetscCall(PetscNew(&ctx));
335: eps->data = (void*)ctx;
337: eps->categ = EPS_CATEGORY_CONTOUR;
339: eps->ops->solve = EPSSolve_FEAST;
340: eps->ops->setup = EPSSetUp_FEAST;
341: eps->ops->setupsort = EPSSetUpSort_Basic;
342: eps->ops->setfromoptions = EPSSetFromOptions_FEAST;
343: eps->ops->destroy = EPSDestroy_FEAST;
344: eps->ops->reset = EPSReset_FEAST;
345: eps->ops->view = EPSView_FEAST;
346: eps->ops->setdefaultst = EPSSetDefaultST_FEAST;
348: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTSetNumPoints_C",EPSFEASTSetNumPoints_FEAST));
349: PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSFEASTGetNumPoints_C",EPSFEASTGetNumPoints_FEAST));
350: PetscFunctionReturn(PETSC_SUCCESS);
351: }