Actual source code: primme.c
slepc-3.7.3 2016-09-29
1: /*
2: This file implements a wrapper to the PRIMME package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
26: PetscErrorCode EPSSolve_PRIMME(EPS);
28: EXTERN_C_BEGIN
29: #include <primme.h>
30: EXTERN_C_END
32: typedef struct {
33: primme_params primme; /* param struc */
34: primme_preset_method method; /* primme method */
35: Mat A; /* problem matrix */
36: EPS eps; /* EPS current context */
37: KSP ksp; /* linear solver and preconditioner */
38: Vec x,y; /* auxiliary vectors */
39: PetscReal target; /* a copy of eps's target */
40: } EPS_PRIMME;
42: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme);
43: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme);
45: static void par_GlobalSumDouble(void *sendBuf,void *recvBuf,int *count,primme_params *primme)
46: {
48: MPI_Allreduce((double*)sendBuf,(double*)recvBuf,*count,MPI_DOUBLE,MPI_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),ierr);
49: }
53: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
54: {
56: PetscMPIInt numProcs,procID;
57: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
58: primme_params *primme = &ops->primme;
59: PetscBool istrivial,flg;
62: MPI_Comm_size(PetscObjectComm((PetscObject)eps),&numProcs);
63: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&procID);
65: /* Check some constraints and set some default values */
66: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
67: STGetOperators(eps->st,0,&ops->A);
68: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
69: if (eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
70: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
71: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
72: if (!eps->which) eps->which = EPS_LARGEST_REAL;
73: if (eps->converged != EPSConvergedAbsolute) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only supports absolute convergence test");
74: RGIsTrivial(eps->rg,&istrivial);
75: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
77: STSetUp(eps->st);
78: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&flg);
79: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with STPRECOND");
81: /* Transfer SLEPc options to PRIMME options */
82: primme->n = eps->n;
83: primme->nLocal = eps->nloc;
84: primme->numEvals = eps->nev;
85: primme->matrix = ops;
86: primme->commInfo = eps;
87: primme->maxMatvecs = eps->max_it;
88: primme->eps = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
89: primme->numProcs = numProcs;
90: primme->procID = procID;
91: primme->printLevel = 0;
92: primme->correctionParams.precondition = 1;
94: switch (eps->which) {
95: case EPS_LARGEST_REAL:
96: primme->target = primme_largest;
97: break;
98: case EPS_SMALLEST_REAL:
99: primme->target = primme_smallest;
100: break;
101: case EPS_TARGET_MAGNITUDE:
102: case EPS_TARGET_REAL:
103: primme->target = primme_closest_abs;
104: primme->numTargetShifts = 1;
105: ops->target = PetscRealPart(eps->target);
106: primme->targetShifts = &ops->target;
107: break;
108: default:
109: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
110: break;
111: }
113: if (primme_set_method(ops->method,primme) < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME method not valid");
115: /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
116: if (eps->ncv) primme->maxBasisSize = eps->ncv;
117: else eps->ncv = primme->maxBasisSize;
118: if (eps->ncv < eps->nev+primme->maxBlockSize) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
119: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
121: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
123: /* Set workspace */
124: EPSAllocateSolution(eps,0);
125: PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
126: if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
128: /* Setup the preconditioner */
129: ops->eps = eps;
130: if (primme->correctionParams.precondition) {
131: STGetKSP(eps->st,&ops->ksp);
132: PetscObjectTypeCompare((PetscObject)ops->ksp,KSPPREONLY,&flg);
133: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"PRIMME only works with KSPPREONLY");
134: primme->preconditioner = NULL;
135: primme->applyPreconditioner = applyPreconditioner_PRIMME;
136: }
138: /* Prepare auxiliary vectors */
139: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,eps->n,NULL,&ops->x);
140: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,eps->n,NULL,&ops->y);
141: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->x);
142: PetscLogObjectParent((PetscObject)eps,(PetscObject)ops->y);
144: /* dispatch solve method */
145: eps->ops->solve = EPSSolve_PRIMME;
146: return(0);
147: }
151: PetscErrorCode EPSSolve_PRIMME(EPS eps)
152: {
154: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
155: PetscScalar *a;
156: Vec v0;
157: #if defined(PETSC_USE_COMPLEX)
158: PetscInt i;
159: PetscReal *evals;
160: #endif
163: /* Reset some parameters left from previous runs */
164: ops->primme.aNorm = 1.0;
165: ops->primme.initSize = eps->nini;
166: ops->primme.iseed[0] = -1;
168: /* Call PRIMME solver */
169: BVGetColumn(eps->V,0,&v0);
170: VecGetArray(v0,&a);
171: #if !defined(PETSC_USE_COMPLEX)
172: dprimme(eps->eigr,a,eps->errest,&ops->primme);
173: if (ierr) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d",ierr);
174: #else
175: /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
176: PetscMalloc1(eps->ncv,&evals);
177: zprimme(evals,(Complex_Z*)a,eps->errest,&ops->primme);
178: if (ierr) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"PRIMME library failed with error code=%d",ierr);
179: for (i=0;i<eps->ncv;i++) eps->eigr[i] = evals[i];
180: PetscFree(evals);
181: #endif
182: VecRestoreArray(v0,&a);
183: BVRestoreColumn(eps->V,0,&v0);
185: eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
186: eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL: EPS_DIVERGED_ITS;
187: eps->its = ops->primme.stats.numOuterIterations;
188: return(0);
189: }
193: static void multMatvec_PRIMME(void *in,void *out,int *blockSize,primme_params *primme)
194: {
196: PetscInt i,N = primme->n;
197: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
198: Vec x = ops->x,y = ops->y;
199: Mat A = ops->A;
202: for (i=0;i<*blockSize;i++) {
203: /* build vectors using 'in' an 'out' workspace */
204: VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
205: VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
207: MatMult(A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
209: VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
210: VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)A),ierr);
211: }
212: PetscFunctionReturnVoid();
213: }
217: static void applyPreconditioner_PRIMME(void *in,void *out,int *blockSize,struct primme_params *primme)
218: {
220: PetscInt i,N = primme->n;
221: EPS_PRIMME *ops = (EPS_PRIMME*)primme->matrix;
222: Vec x = ops->x,y = ops->y;
225: for (i=0;i<*blockSize;i++) {
226: /* build vectors using 'in' an 'out' workspace */
227: VecPlaceArray(x,(PetscScalar*)in+N*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
228: VecPlaceArray(y,(PetscScalar*)out+N*i);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
229: KSPSolve(ops->ksp,x,y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
230: VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
231: VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)ops->ksp),ierr);
232: }
233: PetscFunctionReturnVoid();
234: }
238: PetscErrorCode EPSReset_PRIMME(EPS eps)
239: {
241: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
244: primme_Free(&ops->primme);
245: VecDestroy(&ops->x);
246: VecDestroy(&ops->y);
247: return(0);
248: }
252: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
253: {
257: PetscFree(eps->data);
258: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",NULL);
259: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",NULL);
260: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",NULL);
261: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",NULL);
262: return(0);
263: }
267: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
268: {
269: PetscErrorCode ierr;
270: PetscBool isascii;
271: primme_params *primme = &((EPS_PRIMME*)eps->data)->primme;
272: EPSPRIMMEMethod methodn;
273: PetscMPIInt rank;
276: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
277: if (isascii) {
278: PetscViewerASCIIPrintf(viewer," PRIMME: block size=%D\n",primme->maxBlockSize);
279: EPSPRIMMEGetMethod(eps,&methodn);
280: PetscViewerASCIIPrintf(viewer," PRIMME: solver method: %s\n",EPSPRIMMEMethods[methodn]);
282: /* Display PRIMME params */
283: MPI_Comm_rank(PetscObjectComm((PetscObject)eps),&rank);
284: if (!rank) primme_display_params(*primme);
285: }
286: return(0);
287: }
291: PetscErrorCode EPSSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,EPS eps)
292: {
293: PetscErrorCode ierr;
294: EPS_PRIMME *ctx = (EPS_PRIMME*)eps->data;
295: PetscInt bs;
296: EPSPRIMMEMethod meth;
297: PetscBool flg;
298: KSP ksp;
301: PetscOptionsHead(PetscOptionsObject,"EPS PRIMME Options");
302: PetscOptionsInt("-eps_primme_block_size","Maximum block size","EPSPRIMMESetBlockSize",ctx->primme.maxBlockSize,&bs,&flg);
303: if (flg) {
304: EPSPRIMMESetBlockSize(eps,bs);
305: }
306: PetscOptionsEnum("-eps_primme_method","Method for solving the eigenproblem","EPSPRIMMESetMethod",EPSPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
307: if (flg) {
308: EPSPRIMMESetMethod(eps,meth);
309: }
311: /* Set STPrecond as the default ST */
312: if (!((PetscObject)eps->st)->type_name) {
313: STSetType(eps->st,STPRECOND);
314: }
315: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
317: /* Set the default options of the KSP */
318: STGetKSP(eps->st,&ksp);
319: if (!((PetscObject)ksp)->type_name) {
320: KSPSetType(ksp,KSPPREONLY);
321: }
322: PetscOptionsTail();
323: return(0);
324: }
328: static PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
329: {
330: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
333: if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
334: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
335: else ops->primme.maxBlockSize = bs;
336: return(0);
337: }
341: /*@
342: EPSPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
344: Logically Collective on EPS
346: Input Parameters:
347: + eps - the eigenproblem solver context
348: - bs - block size
350: Options Database Key:
351: . -eps_primme_block_size - Sets the max allowed block size value
353: Notes:
354: If the block size is not set, the value established by primme_initialize
355: is used.
357: The user should set the block size based on the architecture specifics
358: of the target computer, as well as any a priori knowledge of multiplicities.
359: The code does NOT require bs > 1 to find multiple eigenvalues. For some
360: methods, keeping bs = 1 yields the best overall performance.
362: Level: advanced
364: .seealso: EPSPRIMMEGetBlockSize()
365: @*/
366: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
367: {
373: PetscTryMethod(eps,"EPSPRIMMESetBlockSize_C",(EPS,PetscInt),(eps,bs));
374: return(0);
375: }
379: static PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
380: {
381: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
384: *bs = ops->primme.maxBlockSize;
385: return(0);
386: }
390: /*@
391: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
393: Not Collective
395: Input Parameter:
396: . eps - the eigenproblem solver context
398: Output Parameter:
399: . bs - returned block size
401: Level: advanced
403: .seealso: EPSPRIMMESetBlockSize()
404: @*/
405: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
406: {
412: PetscUseMethod(eps,"EPSPRIMMEGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
413: return(0);
414: }
418: static PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps,EPSPRIMMEMethod method)
419: {
420: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
423: ops->method = (primme_preset_method)method;
424: return(0);
425: }
429: /*@
430: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
432: Logically Collective on EPS
434: Input Parameters:
435: + eps - the eigenproblem solver context
436: - method - method that will be used by PRIMME
438: Options Database Key:
439: . -eps_primme_method - Sets the method for the PRIMME library
441: Note:
442: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
444: Level: advanced
446: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
447: @*/
448: PetscErrorCode EPSPRIMMESetMethod(EPS eps,EPSPRIMMEMethod method)
449: {
455: PetscTryMethod(eps,"EPSPRIMMESetMethod_C",(EPS,EPSPRIMMEMethod),(eps,method));
456: return(0);
457: }
461: static PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps,EPSPRIMMEMethod *method)
462: {
463: EPS_PRIMME *ops = (EPS_PRIMME*)eps->data;
466: *method = (EPSPRIMMEMethod)ops->method;
467: return(0);
468: }
472: /*@
473: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
475: Not Collective
477: Input Parameter:
478: . eps - the eigenproblem solver context
480: Output Parameter:
481: . method - method that will be used by PRIMME
483: Level: advanced
485: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
486: @*/
487: PetscErrorCode EPSPRIMMEGetMethod(EPS eps,EPSPRIMMEMethod *method)
488: {
494: PetscUseMethod(eps,"EPSPRIMMEGetMethod_C",(EPS,EPSPRIMMEMethod*),(eps,method));
495: return(0);
496: }
500: PETSC_EXTERN PetscErrorCode EPSCreate_PRIMME(EPS eps)
501: {
503: EPS_PRIMME *primme;
506: PetscNewLog(eps,&primme);
507: eps->data = (void*)primme;
509: eps->ops->setup = EPSSetUp_PRIMME;
510: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
511: eps->ops->destroy = EPSDestroy_PRIMME;
512: eps->ops->reset = EPSReset_PRIMME;
513: eps->ops->view = EPSView_PRIMME;
514: eps->ops->backtransform = EPSBackTransform_Default;
516: primme_initialize(&primme->primme);
517: primme->primme.matrixMatvec = multMatvec_PRIMME;
518: primme->primme.globalSumDouble = par_GlobalSumDouble;
519: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
521: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",EPSPRIMMESetBlockSize_PRIMME);
522: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",EPSPRIMMESetMethod_PRIMME);
523: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",EPSPRIMMEGetBlockSize_PRIMME);
524: PetscObjectComposeFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",EPSPRIMMEGetMethod_PRIMME);
525: return(0);
526: }