Actual source code: blzpack.c
2: /*
3: This file implements a wrapper to the BLZPACK package
4: */
5: #include src/eps/impls/blzpack/blzpackp.h
7: const char* blzpack_error[33] = {
8: "",
9: "illegal data, LFLAG ",
10: "illegal data, dimension of (U), (V), (X) ",
11: "illegal data, leading dimension of (U), (V), (X) ",
12: "illegal data, leading dimension of (EIG) ",
13: "illegal data, number of required eigenpairs ",
14: "illegal data, Lanczos algorithm block size ",
15: "illegal data, maximum number of steps ",
16: "illegal data, number of starting vectors ",
17: "illegal data, number of eigenpairs provided ",
18: "illegal data, problem type flag ",
19: "illegal data, spectrum slicing flag ",
20: "illegal data, eigenvectors purification flag ",
21: "illegal data, level of output ",
22: "illegal data, output file unit ",
23: "illegal data, LCOMM (MPI or PVM) ",
24: "illegal data, dimension of ISTOR ",
25: "illegal data, convergence threshold ",
26: "illegal data, dimension of RSTOR ",
27: "illegal data on at least one PE ",
28: "ISTOR(3:14) must be equal on all PEs ",
29: "RSTOR(1:3) must be equal on all PEs ",
30: "not enough space in ISTOR to start eigensolution ",
31: "not enough space in RSTOR to start eigensolution ",
32: "illegal data, number of negative eigenvalues ",
33: "illegal data, entries of V ",
34: "illegal data, entries of X ",
35: "failure in computational subinterval ",
36: "file I/O error, blzpack.__.BQ ",
37: "file I/O error, blzpack.__.BX ",
38: "file I/O error, blzpack.__.Q ",
39: "file I/O error, blzpack.__.X ",
40: "parallel interface error "
41: };
45: static int EPSSetUp_BLZPACK(EPS eps)
46: {
47: int ierr, listor, lrstor, ncuv, N, n, k1, k2, k3, k4;
48: EPS_BLZPACK *blz = (EPS_BLZPACK *)eps->data;
51: #if defined(PETSC_USE_COMPLEX)
52: SETERRQ(PETSC_ERR_SUP,"Requested method is not available for complex problems");
53: #endif
54: if (!eps->ishermitian)
55: SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
57: if (eps->which!=EPS_SMALLEST_REAL)
58: SETERRQ(1,"Wrong value of eps->which");
60: VecGetSize(eps->vec_initial,&N);
61: VecGetLocalSize(eps->vec_initial,&n);
63: k1 = PetscMin(N,180);
64: k2 = blz->block_size;
65: k4 = PetscMin(eps->ncv,N);
66: k3 = 484+k1*(13+k1*2+k2+PetscMax(18,k2+2))+k2*k2*3+k4*2;
68: listor = 123+k1*12;
69: PetscMalloc((17+listor)*sizeof(int),&blz->istor);
70: blz->istor[14] = listor;
72: if( eps->isgeneralized ) lrstor = n*(k2*4+k1*2+k4)+k3;
73: else lrstor = n*(k2*4+k1)+k3;
74: PetscMalloc((4+lrstor)*sizeof(PetscReal),&blz->rstor);
75: blz->rstor[3] = lrstor;
77: ncuv = PetscMax(3,blz->block_size);
78: PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->u);
79: PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->v);
81: PetscMalloc(2*eps->ncv*sizeof(PetscReal),&blz->eig);
83: return(0);
84: }
88: static int EPSSetDefaults_BLZPACK(EPS eps)
89: {
90: int ierr, N;
93: VecGetSize(eps->vec_initial,&N);
94: if (eps->ncv) {
95: if( eps->ncv < PetscMin(eps->nev+10,eps->nev*2) )
96: SETERRQ(0,"Warning: BLZpack recommends that ncv be larger than min(nev+10,nev*2)");
97: }
98: else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
99: if (!eps->max_it) eps->max_it = PetscMax(100,N);
100: if (!eps->tol) eps->tol = 1.e-7;
101: return(0);
102: }
106: static int EPSSolve_BLZPACK(EPS eps)
107: {
108: EPS_BLZPACK *blz = (EPS_BLZPACK *)eps->data;
109: int i, n, nneig, lflag, nvopu, ierr;
110: Vec x, y;
111: PetscScalar sigma,*pV;
112: PetscTruth isSinv;
113: Mat A;
114:
117: VecGetLocalSize(eps->vec_initial,&n);
118: VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&x);
119: VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&y);
120: VecGetArray(eps->V[0],&pV);
121:
122: blz->istor[0] = n; /* number of rows of U, V, X*/
123: blz->istor[1] = n; /* leading dimension of U, V, X */
124: blz->istor[2] = eps->nev; /* number of required eigenpairs */
125: blz->istor[3] = eps->ncv; /* number of working eigenpairs */
126: blz->istor[4] = blz->block_size; /* number of vectors in a block */
127: blz->istor[5] = 0; /* maximun number of steps per run */
128: blz->istor[6] = 1; /* number of starting vectors as input */
129: blz->istor[7] = 0; /* number of eigenpairs given as input */
130: blz->istor[8] = eps->isgeneralized? 1: 0; /* problem type */
131: blz->istor[9] = (blz->initial == blz->final) ? 0 : 1; /* spectrum slicing */
132: blz->istor[10] = 0; /* solutions refinement */
133: blz->istor[11] = 3; /* level of printing */
134: blz->istor[12] = 6; /* file unit for output */
135: blz->istor[13] = MPI_Comm_c2f(eps->comm); /* communicator */
137: PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);
138: if(isSinv) { STGetShift(eps->OP,&sigma); }
139: else sigma = 0.0; /* shift from origin */
140: nneig = 0; /* no. of eigs less than sigma */
142: blz->rstor[0] = blz->initial; /* lower limit of eigenvalue interval */
143: blz->rstor[1] = blz->final; /* upper limit of eigenvalue interval */
144: blz->rstor[2] = eps->tol; /* threshold for convergence */
146: lflag = 0; /* reverse communication interface flag */
147: eps->its = 0;
149: for(;;) {
151: BLZpack_( blz->istor, blz->rstor, &sigma, &nneig, blz->u, blz->v,
152: &lflag, &nvopu, blz->eig, pV );
154: eps->its = eps->its + 1;
156: if( lflag == 1 ) {
157: /* compute v = OP u */
158: for (i=0;i<nvopu;i++) {
159: VecPlaceArray( x, blz->u+i*n );
160: VecPlaceArray( y, blz->v+i*n );
161: if (isSinv && eps->isgeneralized) {
162: STApplyNoB( eps->OP, x, y );
163: } else {
164: STApply( eps->OP, x, y );
165: }
166: }
167: }
168: else if( lflag == 2 ) {
169: /* compute v = B u */
170: for (i=0;i<nvopu;i++) {
171: VecPlaceArray( x, blz->u+i*n );
172: VecPlaceArray( y, blz->v+i*n );
173: STApplyB( eps->OP, x, y );
174: }
175: }
176: else if( lflag == 3 ) {
177: Mat B;
178: IS is;
179: int n;
180: MatFactorInfo info;
181: /* update shift */
182: STSetShift( eps->OP, sigma );
183: STGetOperators(eps->OP,&A,PETSC_NULL);
184: MatDuplicate(A, MAT_COPY_VALUES, &B);
185: MatGetSize(A, &n, PETSC_NULL);
186: ISCreateStride(eps->comm, n, 0, 1, &is);
187: PetscMemzero(&info,sizeof(MatFactorInfo));
188: info.fill = 1.0;
189: MatCholeskyFactor(B,is,&info);
190: MatGetInertia(B,&nneig,PETSC_NULL,PETSC_NULL);
191: MatDestroy(B);
192: ISDestroy(is);
193: }
194: else if( lflag == 4 ) {
195: /* copy the initial vector */
196: VecPlaceArray(x,blz->v);
197: VecCopy(eps->vec_initial,x);
198: }
199: else break;
200: }
202: VecRestoreArray( eps->V[0], &pV );
204: eps->nconv = BLZistorr_(blz->istor,"NTEIG",5);
205: eps->reason = EPS_CONVERGED_TOL;
207: for (i=0;i<eps->nconv;i++) {
208: eps->eigr[i]=blz->eig[i];
209: eps->eigi[i]=0.0;
210: }
212: if(isSinv) {
213: for( i=0; i<eps->nconv; i++ )
214: eps->eigr[i] = 1.0 / (eps->eigr[i] - sigma);
215: }
216: if (lflag!=0) {
217: char msg[2048] = "";
218: for (i = 0; i < 33; i++) {
219: if (blz->istor[15] & (1 << i)) PetscStrcat(msg, blzpack_error[i]);
220: }
221: SETERRQ2(PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15], msg);
222: }
223: VecDestroy(x);
224: VecDestroy(y);
226: return(0);
227: }
231: int EPSDestroy_BLZPACK(EPS eps)
232: {
233: EPS_BLZPACK *blz = (EPS_BLZPACK *)eps->data;
234: int ierr;
238: if (blz->istor) { PetscFree(blz->istor); }
239: if (blz->rstor) { PetscFree(blz->rstor); }
240: if (blz->u) { PetscFree(blz->u); }
241: if (blz->v) { PetscFree(blz->v); }
242: if (blz->eig) { PetscFree(blz->eig); }
243: EPSDefaultDestroy(eps);
244: return(0);
245: }
249: static int EPSView_BLZPACK(EPS eps,PetscViewer viewer)
250: {
251: EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;
252: int ierr;
253: PetscTruth isascii;
256: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
257: if (!isascii) {
258: SETERRQ1(1,"Viewer type %s not supported for EPSBLZPACK",((PetscObject)viewer)->type_name);
259: }
260: PetscViewerASCIIPrintf(viewer,"block size of the block-Lanczos algorithm: %d\n",blz->block_size);
261: PetscViewerASCIIPrintf(viewer,"computational interval: [%f,%f]\n",blz->initial,blz->final);
262: return(0);
263: }
267: static int EPSSetFromOptions_BLZPACK(EPS eps)
268: {
269: EPS_BLZPACK *blz = (EPS_BLZPACK *)eps->data;
270: int ierr,bs,n;
271: PetscReal interval[2];
272: PetscTruth flg;
275: PetscOptionsHead("BLZPACK options");
277: bs = blz->block_size;
278: PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
279: if (flg) {EPSBlzpackSetBlockSize(eps,bs);}
281: interval[0] = blz->initial;
282: interval[1] = blz->final;
283: n = 2;
284: PetscOptionsRealArray("-eps_blzpack_interval","Computational interval","EPSBlzpackSetInterval",interval,&n,&flg);
285: if (flg) {
286: if (n==1) interval[1]=interval[0];
287: EPSBlzpackSetInterval(eps,interval[0],interval[1]);
288: }
290: PetscOptionsTail();
291: return(0);
292: }
294: EXTERN_C_BEGIN
297: int EPSBlzpackSetBlockSize_BLZPACK(EPS eps,int bs)
298: {
299: EPS_BLZPACK *blz;
302: blz = (EPS_BLZPACK *) eps->data;
303: if (bs == PETSC_DEFAULT) blz->block_size = 3;
304: else if (bs <= 0) {
305: SETERRQ(1, "Incorrect block size");
306: } else blz->block_size = bs;
307: return(0);
308: }
309: EXTERN_C_END
313: /*@
314: EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.
316: Collective on EPS
318: Input Parameters:
319: + eps - the eigenproblem solver context
320: - bs - block size
322: Options Database Key:
323: . -eps_blzpack_block_size - Sets the value of the block size
325: Level: advanced
327: .seealso: EPSBlzpackSetInterval
328: @*/
329: int EPSBlzpackSetBlockSize(EPS eps,int bs)
330: {
331: int ierr, (*f)(EPS,int);
335: PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",(void (**)())&f);
336: if (f) {
337: (*f)(eps,bs);
338: }
339: return(0);
340: }
342: EXTERN_C_BEGIN
345: int EPSBlzpackSetInterval_BLZPACK(EPS eps,PetscReal initial,PetscReal final)
346: {
347: EPS_BLZPACK *blz;
350: blz = (EPS_BLZPACK *) eps->data;
351: blz->initial = initial;
352: blz->final = final;
353: return(0);
354: }
355: EXTERN_C_END
359: /*@
360: EPSBlzpackSetInterval - Sets the computational interval for the BLZPACK
361: package.
363: Collective on EPS
365: Input Parameters:
366: + eps - the eigenproblem solver context
367: . initial - lower bound of the interval
368: - final - upper bound of the interval
370: Options Database Key:
371: . -eps_blzpack_interval - Sets the bounds of the interval (two values
372: separated by commas)
374: Note:
375: The following possibilities are accepted (see Blzpack user's guide for
376: details).
377: initial>final: start seeking for eigenpairs in the upper bound
378: initial<final: start in the lower bound
379: initial=final: run around a single value (no interval)
380:
381: Level: advanced
383: .seealso: EPSBlzpackSetBlockSize()
384: @*/
385: int EPSBlzpackSetInterval(EPS eps,PetscReal initial,PetscReal final)
386: {
387: int ierr, (*f)(EPS,PetscReal,PetscReal);
391: PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetInterval_C",(void (**)())&f);
392: if (f) {
393: (*f)(eps,initial,final);
394: }
395: return(0);
396: }
398: EXTERN_C_BEGIN
401: int EPSCreate_BLZPACK(EPS eps)
402: {
403: EPS_BLZPACK *blzpack;
404: int ierr;
407: PetscNew(EPS_BLZPACK,&blzpack);
408: PetscMemzero(blzpack,sizeof(EPS_BLZPACK));
409: PetscLogObjectMemory(eps,sizeof(EPS_BLZPACK));
410: eps->data = (void *) blzpack;
411: eps->ops->setup = EPSSetUp_BLZPACK;
412: eps->ops->setdefaults = EPSSetDefaults_BLZPACK;
413: eps->ops->solve = EPSSolve_BLZPACK;
414: eps->ops->destroy = EPSDestroy_BLZPACK;
415: eps->ops->setfromoptions = EPSSetFromOptions_BLZPACK;
416: eps->ops->view = EPSView_BLZPACK;
417: eps->ops->backtransform = EPSBackTransform_Default;
419: blzpack->block_size = 3;
420: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetBlockSize_C","EPSBlzpackSetBlockSize_BLZPACK",EPSBlzpackSetBlockSize_BLZPACK);
421: blzpack->initial = 0.0;
422: blzpack->final = 0.0;
423: PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetInterval_C","EPSBlzpackSetInterval_BLZPACK",EPSBlzpackSetInterval_BLZPACK);
425: eps->which = EPS_SMALLEST_REAL;
427: return(0);
428: }
429: EXTERN_C_END