Actual source code: blzpack.c
slepc-3.6.1 2015-09-03
1: /*
2: This file implements a wrapper to the BLZPACK package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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*/
25: #include <../src/eps/impls/external/blzpack/blzpackp.h>
27: PetscErrorCode EPSSolve_BLZPACK(EPS);
29: const char* blzpack_error[33] = {
30: "",
31: "illegal data, LFLAG ",
32: "illegal data, dimension of (U), (V), (X) ",
33: "illegal data, leading dimension of (U), (V), (X) ",
34: "illegal data, leading dimension of (EIG) ",
35: "illegal data, number of required eigenpairs ",
36: "illegal data, Lanczos algorithm block size ",
37: "illegal data, maximum number of steps ",
38: "illegal data, number of starting vectors ",
39: "illegal data, number of eigenpairs provided ",
40: "illegal data, problem type flag ",
41: "illegal data, spectrum slicing flag ",
42: "illegal data, eigenvectors purification flag ",
43: "illegal data, level of output ",
44: "illegal data, output file unit ",
45: "illegal data, LCOMM (MPI or PVM) ",
46: "illegal data, dimension of ISTOR ",
47: "illegal data, convergence threshold ",
48: "illegal data, dimension of RSTOR ",
49: "illegal data on at least one PE ",
50: "ISTOR(3:14) must be equal on all PEs ",
51: "RSTOR(1:3) must be equal on all PEs ",
52: "not enough space in ISTOR to start eigensolution ",
53: "not enough space in RSTOR to start eigensolution ",
54: "illegal data, number of negative eigenvalues ",
55: "illegal data, entries of V ",
56: "illegal data, entries of X ",
57: "failure in computational subinterval ",
58: "file I/O error, blzpack.__.BQ ",
59: "file I/O error, blzpack.__.BX ",
60: "file I/O error, blzpack.__.Q ",
61: "file I/O error, blzpack.__.X ",
62: "parallel interface error "
63: };
67: PetscErrorCode EPSSetUp_BLZPACK(EPS eps)
68: {
70: PetscInt listor,lrstor,ncuv,k1,k2,k3,k4;
71: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
72: PetscBool issinv,istrivial,flg;
75: if (eps->ncv) {
76: if (eps->ncv < PetscMin(eps->nev+10,eps->nev*2)) SETERRQ(PetscObjectComm((PetscObject)eps),0,"Warning: BLZpack recommends that ncv be larger than min(nev+10,nev*2)");
77: } else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
78: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
79: if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);
81: if (!blz->block_size) blz->block_size = 3;
82: if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
83: if (eps->which==EPS_ALL) {
84: if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Must define a computational interval when using EPS_ALL");
85: blz->slice = 1;
86: }
87: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&issinv);
88: if (blz->slice || eps->isgeneralized) {
89: if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert ST is needed for generalized problems or spectrum slicing");
90: }
91: if (blz->slice) {
92: if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
93: if (eps->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The defined computational interval should have at least one of their sides bounded");
94: STSetDefaultShift(eps->st,eps->inta);
95: } else {
96: STSetDefaultShift(eps->st,eps->intb);
97: }
98: }
99: if (!eps->which) {
100: if (issinv) eps->which = EPS_TARGET_REAL;
101: else eps->which = EPS_SMALLEST_REAL;
102: }
103: if ((issinv && eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_MAGNITUDE && eps->which!=EPS_ALL) || (!issinv && eps->which!=EPS_SMALLEST_REAL)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
104: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
106: k1 = PetscMin(eps->n,180);
107: k2 = blz->block_size;
108: k4 = PetscMin(eps->ncv,eps->n);
109: k3 = 484+k1*(13+k1*2+k2+PetscMax(18,k2+2))+k2*k2*3+k4*2;
111: listor = 123+k1*12;
112: PetscFree(blz->istor);
113: PetscMalloc1((17+listor),&blz->istor);
114: PetscLogObjectMemory((PetscObject)eps,(17+listor)*sizeof(PetscBLASInt));
115: PetscBLASIntCast(listor,&blz->istor[14]);
117: if (blz->slice) lrstor = eps->nloc*(k2*4+k1*2+k4)+k3;
118: else lrstor = eps->nloc*(k2*4+k1)+k3;
119: lrstor*=10;
120: PetscFree(blz->rstor);
121: PetscMalloc1((4+lrstor),&blz->rstor);
122: PetscLogObjectMemory((PetscObject)eps,(4+lrstor)*sizeof(PetscReal));
123: blz->rstor[3] = lrstor;
125: ncuv = PetscMax(3,blz->block_size);
126: PetscFree(blz->u);
127: PetscMalloc1(ncuv*eps->nloc,&blz->u);
128: PetscLogObjectMemory((PetscObject)eps,ncuv*eps->nloc*sizeof(PetscScalar));
129: PetscFree(blz->v);
130: PetscMalloc1(ncuv*eps->nloc,&blz->v);
131: PetscLogObjectMemory((PetscObject)eps,ncuv*eps->nloc*sizeof(PetscScalar));
133: PetscFree(blz->eig);
134: PetscMalloc1(2*eps->ncv,&blz->eig);
135: PetscLogObjectMemory((PetscObject)eps,2*eps->ncv*sizeof(PetscReal));
137: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
139: EPSAllocateSolution(eps,0);
140: PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
141: if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
142: RGIsTrivial(eps->rg,&istrivial);
143: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
145: /* dispatch solve method */
146: eps->ops->solve = EPSSolve_BLZPACK;
147: return(0);
148: }
152: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
153: {
155: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
156: PetscInt nn;
157: PetscBLASInt i,nneig,lflag,nvopu;
158: Vec x,y,v0;
159: PetscScalar sigma,*pV;
160: Mat A;
161: KSP ksp;
162: PC pc;
165: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
166: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
167: EPSGetStartVector(eps,0,NULL);
168: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
169: BVGetColumn(eps->V,0,&v0);
170: VecGetArray(v0,&pV);
172: if (eps->isgeneralized && !blz->slice) {
173: STGetShift(eps->st,&sigma); /* shift of origin */
174: blz->rstor[0] = sigma; /* lower limit of eigenvalue interval */
175: blz->rstor[1] = sigma; /* upper limit of eigenvalue interval */
176: } else {
177: sigma = 0.0;
178: blz->rstor[0] = eps->inta; /* lower limit of eigenvalue interval */
179: blz->rstor[1] = eps->intb; /* upper limit of eigenvalue interval */
180: }
181: nneig = 0; /* no. of eigs less than sigma */
183: PetscBLASIntCast(eps->nloc,&blz->istor[0]); /* no. of rows of U, V, X */
184: PetscBLASIntCast(eps->nloc,&blz->istor[1]); /* leading dim of U, V, X */
185: PetscBLASIntCast(eps->nev,&blz->istor[2]); /* required eigenpairs */
186: PetscBLASIntCast(eps->ncv,&blz->istor[3]); /* working eigenpairs */
187: blz->istor[4] = blz->block_size; /* number of vectors in a block */
188: blz->istor[5] = blz->nsteps; /* maximun number of steps per run */
189: blz->istor[6] = 1; /* number of starting vectors as input */
190: blz->istor[7] = 0; /* number of eigenpairs given as input */
191: blz->istor[8] = (blz->slice || eps->isgeneralized) ? 1 : 0; /* problem type */
192: blz->istor[9] = blz->slice; /* spectrum slicing */
193: blz->istor[10] = eps->isgeneralized ? 1 : 0; /* solutions refinement (purify) */
194: blz->istor[11] = 0; /* level of printing */
195: blz->istor[12] = 6; /* file unit for output */
196: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&blz->istor[13]);
198: blz->rstor[2] = eps->tol; /* threshold for convergence */
200: lflag = 0; /* reverse communication interface flag */
202: do {
203: BLZpack_(blz->istor,blz->rstor,&sigma,&nneig,blz->u,blz->v,&lflag,&nvopu,blz->eig,pV);
205: switch (lflag) {
206: case 1:
207: /* compute v = OP u */
208: for (i=0;i<nvopu;i++) {
209: VecPlaceArray(x,blz->u+i*eps->nloc);
210: VecPlaceArray(y,blz->v+i*eps->nloc);
211: if (blz->slice || eps->isgeneralized) {
212: STMatSolve(eps->st,x,y);
213: } else {
214: STApply(eps->st,x,y);
215: }
216: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
217: VecResetArray(x);
218: VecResetArray(y);
219: }
220: /* monitor */
221: eps->nconv = BLZistorr_(blz->istor,"NTEIG",5);
222: EPSMonitor(eps,eps->its,eps->nconv,
223: blz->rstor+BLZistorr_(blz->istor,"IRITZ",5),
224: eps->eigi,
225: blz->rstor+BLZistorr_(blz->istor,"IRITZ",5)+BLZistorr_(blz->istor,"JT",2),
226: BLZistorr_(blz->istor,"NRITZ",5));
227: eps->its = eps->its + 1;
228: if (eps->its >= eps->max_it || eps->nconv >= eps->nev) lflag = 5;
229: break;
230: case 2:
231: /* compute v = B u */
232: for (i=0;i<nvopu;i++) {
233: VecPlaceArray(x,blz->u+i*eps->nloc);
234: VecPlaceArray(y,blz->v+i*eps->nloc);
235: BVApplyMatrix(eps->V,x,y);
236: VecResetArray(x);
237: VecResetArray(y);
238: }
239: break;
240: case 3:
241: /* update shift */
242: PetscInfo1(eps,"Factorization update (sigma=%g)\n",(double)sigma);
243: STSetShift(eps->st,sigma);
244: STGetKSP(eps->st,&ksp);
245: KSPGetPC(ksp,&pc);
246: PCFactorGetMatrix(pc,&A);
247: MatGetInertia(A,&nn,NULL,NULL);
248: PetscBLASIntCast(nn,&nneig);
249: break;
250: case 4:
251: /* copy the initial vector */
252: VecPlaceArray(x,blz->v);
253: VecCopy(v0,x);
254: VecResetArray(x);
255: break;
256: }
258: } while (lflag > 0);
260: VecRestoreArray(v0,&pV);
261: BVRestoreColumn(eps->V,0,&v0);
263: eps->nconv = BLZistorr_(blz->istor,"NTEIG",5);
264: eps->reason = EPS_CONVERGED_TOL;
266: for (i=0;i<eps->nconv;i++) {
267: eps->eigr[i]=blz->eig[i];
268: }
270: if (lflag!=0) {
271: char msg[2048] = "";
272: for (i = 0; i < 33; i++) {
273: if (blz->istor[15] & (1 << i)) PetscStrcat(msg,blzpack_error[i]);
274: }
275: SETERRQ2(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15],msg);
276: }
277: VecDestroy(&x);
278: VecDestroy(&y);
279: return(0);
280: }
284: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
285: {
287: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
290: if (!blz->slice && !eps->isgeneralized) {
291: EPSBackTransform_Default(eps);
292: }
293: return(0);
294: }
298: PetscErrorCode EPSReset_BLZPACK(EPS eps)
299: {
301: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
304: PetscFree(blz->istor);
305: PetscFree(blz->rstor);
306: PetscFree(blz->u);
307: PetscFree(blz->v);
308: PetscFree(blz->eig);
309: return(0);
310: }
314: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
315: {
319: PetscFree(eps->data);
320: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",NULL);
321: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",NULL);
322: return(0);
323: }
327: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
328: {
330: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
331: PetscBool isascii;
334: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
335: if (isascii) {
336: PetscViewerASCIIPrintf(viewer," BLZPACK: block size=%d\n",blz->block_size);
337: PetscViewerASCIIPrintf(viewer," BLZPACK: maximum number of steps per run=%d\n",blz->nsteps);
338: if (blz->slice) {
339: PetscViewerASCIIPrintf(viewer," BLZPACK: computational interval [%f,%f]\n",eps->inta,eps->intb);
340: }
341: }
342: return(0);
343: }
347: PetscErrorCode EPSSetFromOptions_BLZPACK(PetscOptions *PetscOptionsObject,EPS eps)
348: {
350: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
351: PetscInt bs,n;
352: PetscBool flg;
355: PetscOptionsHead(PetscOptionsObject,"EPS BLZPACK Options");
357: bs = blz->block_size;
358: PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
359: if (flg) {
360: EPSBlzpackSetBlockSize(eps,bs);
361: }
363: n = blz->nsteps;
364: PetscOptionsInt("-eps_blzpack_nsteps","Number of steps","EPSBlzpackSetNSteps",n,&n,&flg);
365: if (flg) {
366: EPSBlzpackSetNSteps(eps,n);
367: }
369: PetscOptionsTail();
370: return(0);
371: }
375: static PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,PetscInt bs)
376: {
378: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
381: if (bs == PETSC_DEFAULT) blz->block_size = 3;
382: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Block size must be positive");
383: else {
384: PetscBLASIntCast(bs,&blz->block_size);
385: }
386: return(0);
387: }
391: /*@
392: EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.
394: Collective on EPS
396: Input Parameters:
397: + eps - the eigenproblem solver context
398: - bs - block size
400: Options Database Key:
401: . -eps_blzpack_block_size - Sets the value of the block size
403: Level: advanced
404: @*/
405: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,PetscInt bs)
406: {
412: PetscTryMethod(eps,"EPSBlzpackSetBlockSize_C",(EPS,PetscInt),(eps,bs));
413: return(0);
414: }
418: static PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,PetscInt nsteps)
419: {
421: EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;
424: if (nsteps == PETSC_DEFAULT) blz->nsteps = 0;
425: else {
426: PetscBLASIntCast(nsteps,&blz->nsteps);
427: }
428: return(0);
429: }
433: /*@
434: EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
435: package.
437: Collective on EPS
439: Input Parameters:
440: + eps - the eigenproblem solver context
441: - nsteps - maximum number of steps
443: Options Database Key:
444: . -eps_blzpack_nsteps - Sets the maximum number of steps per run
446: Level: advanced
448: @*/
449: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,PetscInt nsteps)
450: {
456: PetscTryMethod(eps,"EPSBlzpackSetNSteps_C",(EPS,PetscInt),(eps,nsteps));
457: return(0);
458: }
462: PETSC_EXTERN PetscErrorCode EPSCreate_BLZPACK(EPS eps)
463: {
465: EPS_BLZPACK *blzpack;
468: PetscNewLog(eps,&blzpack);
469: eps->data = (void*)blzpack;
471: eps->ops->setup = EPSSetUp_BLZPACK;
472: eps->ops->setfromoptions = EPSSetFromOptions_BLZPACK;
473: eps->ops->destroy = EPSDestroy_BLZPACK;
474: eps->ops->reset = EPSReset_BLZPACK;
475: eps->ops->view = EPSView_BLZPACK;
476: eps->ops->backtransform = EPSBackTransform_BLZPACK;
477: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",EPSBlzpackSetBlockSize_BLZPACK);
478: PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",EPSBlzpackSetNSteps_BLZPACK);
479: return(0);
480: }