Actual source code: primme.c
1: /*
2: This file implements a wrapper to the PRIMME library
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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 "petsc.h"
25: #include private/epsimpl.h
26: #include private/stimpl.h
28: PetscErrorCode EPSSolve_PRIMME(EPS);
31: #include "primme.h"
34: typedef struct {
35: primme_params primme; /* param struc */
36: primme_preset_method method; /* primme method */
37: Mat A; /* problem matrix */
38: EPS eps; /* EPS current context */
39: KSP ksp; /* preconditioner */
40: Vec x,y; /* auxiliar vectors */
41: } EPS_PRIMME;
43: const char *methodList[] = {
44: "dynamic",
45: "default_min_time",
46: "default_min_matvecs",
47: "arnoldi",
48: "gd",
49: "gd_plusk",
50: "gd_olsen_plusk",
51: "jd_olsen_plusk",
52: "rqi",
53: "jdqr",
54: "jdqmr",
55: "jdqmr_etol",
56: "subspace_iteration",
57: "lobpcg_orthobasis",
58: "lobpcg_orthobasis_window"
59: };
60: EPSPRIMMEMethod methodN[] = {
61: EPS_PRIMME_DYNAMIC,
62: EPS_PRIMME_DEFAULT_MIN_TIME,
63: EPS_PRIMME_DEFAULT_MIN_MATVECS,
64: EPS_PRIMME_ARNOLDI,
65: EPS_PRIMME_GD,
66: EPS_PRIMME_GD_PLUSK,
67: EPS_PRIMME_GD_OLSEN_PLUSK,
68: EPS_PRIMME_JD_OLSEN_PLUSK,
69: EPS_PRIMME_RQI,
70: EPS_PRIMME_JDQR,
71: EPS_PRIMME_JDQMR,
72: EPS_PRIMME_JDQMR_ETOL,
73: EPS_PRIMME_SUBSPACE_ITERATION,
74: EPS_PRIMME_LOBPCG_ORTHOBASIS,
75: EPS_PRIMME_LOBPCG_ORTHOBASISW
76: };
78: static void multMatvec_PRIMME(void *in, void *out, int *blockSize, primme_params *primme);
79: static void applyPreconditioner_PRIMME(void *in, void *out, int *blockSize, struct primme_params *primme);
81: void par_GlobalSumDouble(void *sendBuf, void *recvBuf, int *count, primme_params *primme) {
83: MPI_Allreduce((double*)sendBuf, (double*)recvBuf, *count, MPI_DOUBLE, MPI_SUM, ((PetscObject)(primme->commInfo))->comm);CHKERRABORT(((PetscObject)(primme->commInfo))->comm,ierr);
84: }
88: PetscErrorCode EPSSetUp_PRIMME(EPS eps)
89: {
91: PetscMPIInt numProcs, procID;
92: EPS_PRIMME *ops = (EPS_PRIMME *)eps->data;
93: primme_params *primme = &(((EPS_PRIMME *)eps->data)->primme);
94: PetscTruth t;
98: MPI_Comm_size(((PetscObject)eps)->comm,&numProcs);
99: MPI_Comm_rank(((PetscObject)eps)->comm,&procID);
100:
101: /* Check some constraints and set some default values */
102: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
103: STGetOperators(eps->OP, &ops->A, PETSC_NULL);
104: if (!ops->A) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"The problem matrix has to be specified first");
105: if (!eps->ishermitian)
106: SETERRQ(PETSC_ERR_SUP,"PRIMME is only available for Hermitian problems");
107: if (eps->isgeneralized)
108: SETERRQ(PETSC_ERR_SUP,"PRIMME is not available for generalized problems");
109: if (!eps->which) eps->which = EPS_LARGEST_REAL;
111: /* Change the default sigma to inf if necessary */
112: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
113: eps->which == EPS_LARGEST_IMAGINARY) {
114: STSetDefaultShift(eps->OP, 3e300);
115: }
117: STSetUp(eps->OP);
118: PetscTypeCompare((PetscObject)eps->OP, STPRECOND, &t);
119: if (t == PETSC_FALSE) SETERRQ(PETSC_ERR_SUP, "PRIMME only works with precond spectral transformation");
121: /* Transfer SLEPc options to PRIMME options */
122: primme->n = eps->n;
123: primme->nLocal = eps->nloc;
124: primme->numEvals = eps->nev;
125: primme->matrix = ops;
126: primme->commInfo = eps;
127: primme->maxMatvecs = eps->max_it;
128: primme->eps = eps->tol;
129: primme->numProcs = numProcs;
130: primme->procID = procID;
131: primme->printLevel = 0;
132: primme->correctionParams.precondition = 1;
134: if (!eps->which) eps->which = EPS_LARGEST_REAL;
135: switch(eps->which) {
136: case EPS_LARGEST_REAL:
137: primme->target = primme_largest;
138: break;
139: case EPS_SMALLEST_REAL:
140: primme->target = primme_smallest;
141: break;
142: default:
143: SETERRQ(PETSC_ERR_SUP,"PRIMME only allows EPS_LARGEST_REAL and EPS_SMALLEST_REAL for 'which' value");
144: break;
145: }
146:
147: if (primme_set_method(ops->method, primme) < 0)
148: SETERRQ(PETSC_ERR_SUP,"PRIMME method not valid");
149:
150: /* If user sets ncv, maxBasisSize is modified. If not, ncv is set as maxBasisSize */
151: if (eps->ncv) primme->maxBasisSize = eps->ncv;
152: else eps->ncv = primme->maxBasisSize;
153: if (eps->ncv < eps->nev+primme->maxBlockSize)
154: SETERRQ(PETSC_ERR_SUP,"PRIMME needs ncv >= nev+maxBlockSize");
155: if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
157: if (eps->extraction) {
158: PetscInfo(eps,"Warning: extraction type ignored\n");
159: }
161: /* Set workspace */
162: EPSAllocateSolution(eps);
164: /* Setup the preconditioner */
165: ops->eps = eps;
166: if (primme->correctionParams.precondition) {
167: STGetKSP(eps->OP, &ops->ksp);
168: PetscTypeCompare((PetscObject)ops->ksp, KSPPREONLY, &t);
169: if (t == PETSC_FALSE) SETERRQ(PETSC_ERR_SUP, "PRIMME only works with preonly ksp of the spectral transformation");
170: primme->preconditioner = PETSC_NULL;
171: primme->applyPreconditioner = applyPreconditioner_PRIMME;
172: }
174: /* Prepare auxiliary vectors */
175: VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->x);
176: VecCreateMPIWithArray(PETSC_COMM_WORLD,eps->nloc,eps->n,PETSC_NULL,&ops->y);
177:
178: /* dispatch solve method */
179: if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
180: eps->ops->solve = EPSSolve_PRIMME;
181: return(0);
182: }
186: PetscErrorCode EPSSolve_PRIMME(EPS eps)
187: {
189: EPS_PRIMME *ops = (EPS_PRIMME *)eps->data;
190: PetscScalar *a;
191: #ifdef PETSC_USE_COMPLEX
192: PetscInt i;
193: PetscReal *evals;
194: #endif
198: /* Reset some parameters left from previous runs */
199: ops->primme.aNorm = 0.0;
200: ops->primme.initSize = eps->nini;
201: ops->primme.iseed[0] = -1;
203: /* Call PRIMME solver */
204: VecGetArray(eps->V[0], &a);
205: #ifndef PETSC_USE_COMPLEX
206: dprimme(eps->eigr, a, eps->errest, &ops->primme);
207: #else
208: /* PRIMME returns real eigenvalues, but SLEPc works with complex ones */
209: PetscMalloc(eps->ncv*sizeof(PetscReal),&evals);
210: zprimme(evals, (Complex_Z*)a, eps->errest, &ops->primme);
211: for (i=0;i<eps->ncv;i++)
212: eps->eigr[i] = evals[i];
213: PetscFree(evals);
214: #endif
215: VecRestoreArray(eps->V[0], &a);
216:
217: switch(ierr) {
218: case 0: /* Successful */
219: break;
221: case -1:
222: SETERRQ(PETSC_ERR_SUP,"PRIMME: Failed to open output file");
223: break;
225: case -2:
226: SETERRQ(PETSC_ERR_SUP,"PRIMME: Insufficient integer or real workspace allocated");
227: break;
229: case -3:
230: SETERRQ(PETSC_ERR_SUP,"PRIMME: main_iter encountered a problem");
231: break;
233: default:
234: SETERRQ(PETSC_ERR_SUP,"PRIMME: some parameters wrong configured");
235: break;
236: }
238: eps->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
239: eps->reason = eps->ncv >= eps->nev ? EPS_CONVERGED_TOL : EPS_DIVERGED_ITS;
240: eps->its = ops->primme.stats.numOuterIterations;
241: eps->OP->applys = ops->primme.stats.numMatvecs;
243: return(0);
244: }
248: static void multMatvec_PRIMME(void *in, void *out, int *blockSize, primme_params *primme)
249: {
251: PetscInt i, N = primme->n;
252: EPS_PRIMME *ops = (EPS_PRIMME *)primme->matrix;
253: Vec x = ops->x;
254: Vec y = ops->y;
255: Mat A = ops->A;
258: for (i=0;i<*blockSize;i++) {
259: /* build vectors using 'in' an 'out' workspace */
260: VecPlaceArray(x, (PetscScalar*)in+N*i ); CHKERRABORT(PETSC_COMM_WORLD,ierr);
261: VecPlaceArray(y, (PetscScalar*)out+N*i ); CHKERRABORT(PETSC_COMM_WORLD,ierr);
263: MatMult(A, x, y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
264:
265: VecResetArray(x); CHKERRABORT(PETSC_COMM_WORLD,ierr);
266: VecResetArray(y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
267: }
268: PetscFunctionReturnVoid();
269: }
273: static void applyPreconditioner_PRIMME(void *in, void *out, int *blockSize, struct primme_params *primme)
274: {
276: PetscInt i, N = primme->n, lits;
277: EPS_PRIMME *ops = (EPS_PRIMME *)primme->matrix;
278: Vec x = ops->x;
279: Vec y = ops->y;
280:
282: for (i=0;i<*blockSize;i++) {
283: /* build vectors using 'in' an 'out' workspace */
284: VecPlaceArray(x, (PetscScalar*)in+N*i ); CHKERRABORT(PETSC_COMM_WORLD,ierr);
285: VecPlaceArray(y, (PetscScalar*)out+N*i ); CHKERRABORT(PETSC_COMM_WORLD,ierr);
287: KSPSolve(ops->ksp, x, y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
288: KSPGetIterationNumber(ops->ksp, &lits); CHKERRABORT(PETSC_COMM_WORLD,ierr);
289: ops->eps->OP->lineariterations+= lits;
290:
291: VecResetArray(x); CHKERRABORT(PETSC_COMM_WORLD,ierr);
292: VecResetArray(y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
293: }
294: PetscFunctionReturnVoid();
295: }
300: PetscErrorCode EPSDestroy_PRIMME(EPS eps)
301: {
303: EPS_PRIMME *ops = (EPS_PRIMME *)eps->data;
307:
308: primme_Free(&ops->primme);
309: VecDestroy(ops->x);
310: VecDestroy(ops->y);
311: PetscFree(eps->data);
312: EPSFreeSolution(eps);
313:
314: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","",PETSC_NULL);
315: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","",PETSC_NULL);
316: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","",PETSC_NULL);
317: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","",PETSC_NULL);
319: return(0);
320: }
324: PetscErrorCode EPSView_PRIMME(EPS eps,PetscViewer viewer)
325: {
327: PetscTruth isascii;
328: primme_params *primme = &((EPS_PRIMME *)eps->data)->primme;
329: EPSPRIMMEMethod methodn;
332: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
333: if (!isascii) {
334: SETERRQ1(1,"Viewer type %s not supported for EPSPRIMME",((PetscObject)viewer)->type_name);
335: }
336:
337: PetscViewerASCIIPrintf(viewer,"PRIMME solver block size: %d\n",primme->maxBlockSize);
338: EPSPRIMMEGetMethod(eps, &methodn);
339: PetscViewerASCIIPrintf(viewer,"PRIMME solver method: %s\n", methodList[methodn]);
341: /* Display PRIMME params */
342: primme_display_params(*primme);
344: return(0);
345: }
349: PetscErrorCode EPSSetFromOptions_PRIMME(EPS eps)
350: {
352: EPS_PRIMME *ops = (EPS_PRIMME *)eps->data;
353: PetscInt op;
354: PetscTruth flg;
357:
358: PetscOptionsBegin(((PetscObject)eps)->comm,((PetscObject)eps)->prefix,"PRIMME Options","EPS");
360: op = ops->primme.maxBlockSize;
361: PetscOptionsInt("-eps_primme_block_size"," maximum block size","EPSPRIMMESetBlockSize",op,&op,&flg);
362: if (flg) {EPSPRIMMESetBlockSize(eps,op);}
363: op = 0;
364: PetscOptionsEList("-eps_primme_method","set method for solving the eigenproblem",
365: "EPSPRIMMESetMethod",methodList,15,methodList[1],&op,&flg);
366: if (flg) {EPSPRIMMESetMethod(eps, methodN[op]);}
367:
368: PetscOptionsEnd();
369:
370: return(0);
371: }
376: PetscErrorCode EPSPRIMMESetBlockSize_PRIMME(EPS eps,PetscInt bs)
377: {
378: EPS_PRIMME *ops = (EPS_PRIMME *) eps->data;
382: if (bs == PETSC_DEFAULT) ops->primme.maxBlockSize = 1;
383: else if (bs <= 0) {
384: SETERRQ(1, "PRIMME: wrong block size");
385: } else ops->primme.maxBlockSize = bs;
386:
387: return(0);
388: }
393: /*@
394: EPSPRIMMESetBlockSize - The maximum block size the code will try to use.
395: The user should set
396: this based on the architecture specifics of the target computer,
397: as well as any a priori knowledge of multiplicities. The code does
398: NOT require BlockSize > 1 to find multiple eigenvalues. For some
399: methods, keeping BlockSize = 1 yields the best overall performance.
401: Collective on EPS
403: Input Parameters:
404: + eps - the eigenproblem solver context
405: - bs - block size
407: Options Database Key:
408: . -eps_primme_block_size - Sets the max allowed block size value
410: Notes:
411: If the block size is not set, the value established by primme_initialize
412: is used.
414: Level: advanced
415: .seealso: EPSPRIMMEGetBlockSize()
416: @*/
417: PetscErrorCode EPSPRIMMESetBlockSize(EPS eps,PetscInt bs)
418: {
419: PetscErrorCode ierr, (*f)(EPS,PetscInt);
422:
424: PetscObjectQueryFunction((PetscObject)eps,"EPSPRIMMESetBlockSize_C",(void (**)())&f);
425: if (f) {
426: (*f)(eps,bs);
427: }
428:
429: return(0);
430: }
435: PetscErrorCode EPSPRIMMEGetBlockSize_PRIMME(EPS eps,PetscInt *bs)
436: {
437: EPS_PRIMME *ops = (EPS_PRIMME *) eps->data;
441: if (bs) *bs = ops->primme.maxBlockSize;
442:
443: return(0);
444: }
449: /*@
450: EPSPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
451:
452: Collective on EPS
454: Input Parameters:
455: . eps - the eigenproblem solver context
456:
457: Output Parameters:
458: . bs - returned block size
460: Level: advanced
461: .seealso: EPSPRIMMESetBlockSize()
462: @*/
463: PetscErrorCode EPSPRIMMEGetBlockSize(EPS eps,PetscInt *bs)
464: {
465: PetscErrorCode ierr, (*f)(EPS,PetscInt*);
468:
470: PetscObjectQueryFunction((PetscObject)eps,"EPSPRIMMEGetBlockSize_C",(void (**)())&f);
471: if (f) {
472: (*f)(eps,bs);
473: }
474:
475: return(0);
476: }
481: PetscErrorCode EPSPRIMMESetMethod_PRIMME(EPS eps, EPSPRIMMEMethod method)
482: {
483: EPS_PRIMME *ops = (EPS_PRIMME *) eps->data;
487: if (method == PETSC_DEFAULT) ops->method = DEFAULT_MIN_TIME;
488: else ops->method = (primme_preset_method)method;
490: return(0);
491: }
496: /*@
497: EPSPRIMMESetMethod - Sets the method for the PRIMME library.
499: Collective on EPS
501: Input Parameters:
502: + eps - the eigenproblem solver context
503: - method - method that will be used by PRIMME. It must be one of:
504: EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
505: EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
506: EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK,
507: EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR,
508: EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
509: EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW
511: Options Database Key:
512: . -eps_primme_set_method - Sets the method for the PRIMME library.
514: Note:
515: If not set, the method defaults to EPS_PRIMME_DEFAULT_MIN_TIME.
517: Level: advanced
519: .seealso: EPSPRIMMEGetMethod(), EPSPRIMMEMethod
520: @*/
521: PetscErrorCode EPSPRIMMESetMethod(EPS eps, EPSPRIMMEMethod method)
522: {
523: PetscErrorCode ierr, (*f)(EPS,EPSPRIMMEMethod);
526:
528: PetscObjectQueryFunction((PetscObject)eps,"EPSPRIMMESetMethod_C",(void (**)())&f);
529: if (f) {
530: (*f)(eps,method);
531: }
532:
533: return(0);
534: }
539: PetscErrorCode EPSPRIMMEGetMethod_PRIMME(EPS eps, EPSPRIMMEMethod *method)
540: {
541: EPS_PRIMME *ops = (EPS_PRIMME *) eps->data;
545: if (method)
546: *method = (EPSPRIMMEMethod)ops->method;
548: return(0);
549: }
554: /*@C
555: EPSPRIMMEGetMethod - Gets the method for the PRIMME library.
557: Mon Collective on EPS
559: Input Parameters:
560: . eps - the eigenproblem solver context
561:
562: Output Parameters:
563: . method - method that will be used by PRIMME. It must be one of:
564: EPS_PRIMME_DYNAMIC, EPS_PRIMME_DEFAULT_MIN_TIME(EPS_PRIMME_JDQMR_ETOL),
565: EPS_PRIMME_DEFAULT_MIN_MATVECS(EPS_PRIMME_GD_OLSEN_PLUSK), EPS_PRIMME_ARNOLDI,
566: EPS_PRIMME_GD, EPS_PRIMME_GD_PLUSK, EPS_PRIMME_GD_OLSEN_PLUSK,
567: EPS_PRIMME_JD_OLSEN_PLUSK, EPS_PRIMME_RQI, EPS_PRIMME_JDQR, EPS_PRIMME_JDQMR,
568: EPS_PRIMME_JDQMR_ETOL, EPS_PRIMME_SUBSPACE_ITERATION,
569: EPS_PRIMME_LOBPCG_ORTHOBASIS, EPS_PRIMME_LOBPCG_ORTHOBASISW
571: Level: advanced
573: .seealso: EPSPRIMMESetMethod(), EPSPRIMMEMethod
574: @*/
575: PetscErrorCode EPSPRIMMEGetMethod(EPS eps, EPSPRIMMEMethod *method)
576: {
577: PetscErrorCode ierr, (*f)(EPS,EPSPRIMMEMethod*);
580:
582: PetscObjectQueryFunction((PetscObject)eps,"EPSPRIMMEGetMethod_C",(void (**)())&f);
583: if (f) {
584: (*f)(eps,method);
585: }
586:
587: return(0);
588: }
593: PetscErrorCode EPSCreate_PRIMME(EPS eps)
594: {
596: EPS_PRIMME *primme;
599:
600: STSetType(eps->OP, STPRECOND);
601: STPrecondSetKSPHasMat(eps->OP, PETSC_TRUE);
603: PetscNew(EPS_PRIMME,&primme);
604: PetscLogObjectMemory(eps,sizeof(EPS_PRIMME));
605: eps->data = (void *) primme;
606: eps->ops->setup = EPSSetUp_PRIMME;
607: eps->ops->setfromoptions = EPSSetFromOptions_PRIMME;
608: eps->ops->destroy = EPSDestroy_PRIMME;
609: eps->ops->view = EPSView_PRIMME;
610: eps->ops->backtransform = EPSBackTransform_Default;
611: eps->ops->computevectors = EPSComputeVectors_Default;
613: primme_initialize(&primme->primme);
614: primme->primme.matrixMatvec = multMatvec_PRIMME;
615: primme->primme.globalSumDouble = par_GlobalSumDouble;
616: primme->method = (primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME;
618: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetBlockSize_C","EPSPRIMMESetBlockSize_PRIMME",EPSPRIMMESetBlockSize_PRIMME);
619: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMESetMethod_C","EPSPRIMMESetMethod_PRIMME",EPSPRIMMESetMethod_PRIMME);
620: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetBlockSize_C","EPSPRIMMEGetBlockSize_PRIMME",EPSPRIMMEGetBlockSize_PRIMME);
621: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSPRIMMEGetMethod_C","EPSPRIMMEGetMethod_PRIMME",EPSPRIMMEGetMethod_PRIMME);
622: return(0);
623: }