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: PetscErrorCode EPSSetUp_BLZPACK(EPS eps)
 46: {
 48:   PetscInt       N, n;
 49:   int            listor, lrstor, ncuv, k1, k2, k3, k4;
 50:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
 51:   PetscTruth     flg;
 52:   KSP            ksp;
 53:   PC             pc;

 56:   VecGetSize(eps->vec_initial,&N);
 57:   VecGetLocalSize(eps->vec_initial,&n);
 58:   if (eps->ncv) {
 59:     if( eps->ncv < PetscMin(eps->nev+10,eps->nev*2) )
 60:       SETERRQ(0,"Warning: BLZpack recommends that ncv be larger than min(nev+10,nev*2)");
 61:   }
 62:   else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
 63:   if (!eps->max_it) eps->max_it = PetscMax(100,N);
 64:   if (!eps->tol) eps->tol = 1.e-7;

 66:   if (!eps->ishermitian)
 67:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 68:   if (blz->slice) {
 69:     PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
 70:     if (!flg)
 71:       SETERRQ(PETSC_ERR_SUP,"Shift-and-invert ST is needed for spectrum slicing");
 72:     STGetKSP(eps->OP,&ksp);
 73:     PetscTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
 74:     if (!flg)
 75:       SETERRQ(PETSC_ERR_SUP,"Preonly KSP is needed for spectrum slicing");
 76:     KSPGetPC(ksp,&pc);
 77:     PetscTypeCompare((PetscObject)pc,PCCHOLESKY,&flg);
 78:     if (!flg)
 79:       SETERRQ(PETSC_ERR_SUP,"Cholesky PC is needed for spectrum slicing");
 80:   }
 81:   if (eps->which!=EPS_SMALLEST_REAL)
 82:     SETERRQ(1,"Wrong value of eps->which");

 84:   k1 = PetscMin(N,180);
 85:   k2 = blz->block_size;
 86:   k4 = PetscMin(eps->ncv,N);
 87:   k3 = 484+k1*(13+k1*2+k2+PetscMax(18,k2+2))+k2*k2*3+k4*2;

 89:   listor = 123+k1*12;
 90:   PetscFree(blz->istor);
 91:   PetscMalloc((17+listor)*sizeof(int),&blz->istor);
 92:   blz->istor[14] = listor;

 94:   if (blz->slice) lrstor = n*(k2*4+k1*2+k4)+k3;
 95:   else lrstor = n*(k2*4+k1)+k3;
 96:   PetscFree(blz->rstor);
 97:   PetscMalloc((4+lrstor)*sizeof(PetscReal),&blz->rstor);
 98:   blz->rstor[3] = lrstor;

100:   ncuv = PetscMax(3,blz->block_size);
101:   PetscFree(blz->u);
102:   PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->u);
103:   PetscFree(blz->v);
104:   PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->v);

106:   PetscFree(blz->eig);
107:   PetscMalloc(2*eps->ncv*sizeof(PetscReal),&blz->eig);

109:   EPSAllocateSolutionContiguous(eps);
110:   return(0);
111: }

115: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
116: {
118:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
119:   PetscInt       n, nn;
120:   int            i, nneig, lflag, nvopu;
121:   Vec            x, y;
122:   PetscScalar    sigma,*pV;
123:   Mat            A;
124:   KSP            ksp;
125:   PC             pc;
126: 

129:   VecGetLocalSize(eps->vec_initial,&n);
130:   VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&x);
131:   VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&y);
132:   VecGetArray(eps->V[0],&pV);
133: 
134:   if (blz->slice) { STGetShift(eps->OP,&sigma); }
135:   else sigma = 0.0;              /* shift of origin */
136:   nneig = 0;                     /* no. of eigs less than sigma */

138:   blz->istor[0]  = n;            /* number of rows of U, V, X*/
139:   blz->istor[1]  = n;            /* leading dimension of U, V, X */
140:   blz->istor[2]  = eps->nev;     /* number of required eigenpairs */
141:   blz->istor[3]  = eps->ncv;     /* number of working eigenpairs */
142:   blz->istor[4]  = blz->block_size;    /* number of vectors in a block */
143:   blz->istor[5]  = blz->nsteps;  /* maximun number of steps per run */
144:   blz->istor[6]  = 1;            /* number of starting vectors as input */
145:   blz->istor[7]  = 0;            /* number of eigenpairs given as input */
146:   blz->istor[8]  = blz->slice;   /* problem type */
147:   blz->istor[9]  = blz->slice;   /* spectrum slicing */
148:   blz->istor[10] = blz->slice;   /* solutions refinement (purify) */
149:   blz->istor[11] = 0;            /* level of printing */
150:   blz->istor[12] = 6;            /* file unit for output */
151:   blz->istor[13] = MPI_Comm_c2f(eps->comm);    /* communicator */

153:   blz->rstor[0]  = blz->initial; /* lower limit of eigenvalue interval */
154:   blz->rstor[1]  = blz->final;   /* upper limit of eigenvalue interval */
155:   blz->rstor[2]  = eps->tol;     /* threshold for convergence */

157:   lflag = 0;           /* reverse communication interface flag */
158:   eps->its  = 0;

160:   do {

162:     BLZpack_( blz->istor, blz->rstor, &sigma, &nneig, blz->u, blz->v,
163:               &lflag, &nvopu, blz->eig, pV );

165:     switch (lflag) {
166:     case 1:
167:       /* compute v = OP u */
168:       for (i=0;i<nvopu;i++) {
169:         VecPlaceArray( x, blz->u+i*n );
170:         VecPlaceArray( y, blz->v+i*n );
171:         if (blz->slice) {
172:           STApplyNoB( eps->OP, x, y );
173:         } else {
174:           STApply( eps->OP, x, y );
175:         }
176:         EPSOrthogonalize(eps,eps->nds,eps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
177:         VecResetArray(x);
178:         VecResetArray(y);
179:       }
180:       /* monitor */
181:       eps->nconv  = BLZistorr_(blz->istor,"NTEIG",5);
182:       EPSMonitor(eps,eps->its,eps->nconv,
183:         blz->rstor+BLZistorr_(blz->istor,"IRITZ",5),
184:         eps->eigi,
185:         blz->rstor+BLZistorr_(blz->istor,"IRITZ",5)+BLZistorr_(blz->istor,"JT",2),
186:         BLZistorr_(blz->istor,"NRITZ",5));
187:       eps->its = eps->its + 1;
188:       if (eps->its >= eps->max_it || eps->nconv >= eps->nev) lflag = 5;
189:       break;
190:     case 2:
191:       /* compute v = B u */
192:       for (i=0;i<nvopu;i++) {
193:         VecPlaceArray( x, blz->u+i*n );
194:         VecPlaceArray( y, blz->v+i*n );
195:         STApplyB( eps->OP, x, y );
196:         VecResetArray(x);
197:         VecResetArray(y);
198:       }
199:       break;
200:     case 3:
201:       /* update shift */
202:       STSetShift(eps->OP,sigma);
203:       STGetKSP(eps->OP,&ksp);
204:       KSPGetPC(ksp,&pc);
205:       PCGetFactoredMatrix(pc,&A);
206:       MatGetInertia(A,&nn,PETSC_NULL,PETSC_NULL);
207:       nneig = nn;
208:       break;
209:     case 4:
210:       /* copy the initial vector */
211:       VecPlaceArray(x,blz->v);
212:       VecCopy(eps->vec_initial,x);
213:       VecResetArray(x);
214:       break;
215:     }
216: 
217:   } while (lflag > 0);

219:   VecRestoreArray( eps->V[0], &pV );

221:   eps->nconv  = BLZistorr_(blz->istor,"NTEIG",5);
222:   eps->reason = EPS_CONVERGED_TOL;

224:   for (i=0;i<eps->nconv;i++) {
225:     eps->eigr[i]=blz->eig[i];
226:     eps->eigi[i]=0.0;
227:   }

229:   if (lflag!=0) {
230:     char msg[2048] = "";
231:     for (i = 0; i < 33; i++) {
232:       if (blz->istor[15] & (1 << i)) PetscStrcat(msg, blzpack_error[i]);
233:     }
234:     SETERRQ2(PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15], msg);
235:   }
236:   VecDestroy(x);
237:   VecDestroy(y);

239:   return(0);
240: }

244: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
245: {
247:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

250:   if (!blz->slice) {
251:     EPSBackTransform_Default(eps);
252:   }
253:   return(0);
254: }

258: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
259: {
261:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

265:   PetscFree(blz->istor);
266:   PetscFree(blz->rstor);
267:   PetscFree(blz->u);
268:   PetscFree(blz->v);
269:   PetscFree(blz->eig);
270:   PetscFree(eps->data);
271:   EPSFreeSolutionContiguous(eps);
272:   return(0);
273: }

277: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
278: {
280:   EPS_BLZPACK    *blz = (EPS_BLZPACK *) eps->data;
281:   PetscTruth     isascii;

284:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
285:   if (!isascii) {
286:     SETERRQ1(1,"Viewer type %s not supported for EPSBLZPACK",((PetscObject)viewer)->type_name);
287:   }
288:   PetscViewerASCIIPrintf(viewer,"block size of the block-Lanczos algorithm: %d\n",blz->block_size);
289:   PetscViewerASCIIPrintf(viewer,"computational interval: [%f,%f]\n",blz->initial,blz->final);
290:   return(0);
291: }

295: PetscErrorCode EPSSetFromOptions_BLZPACK(EPS eps)
296: {
298:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
299:   PetscInt       bs,n;
300:   PetscReal      interval[2];
301:   PetscTruth     flg;
302:   KSP            ksp;
303:   PC             pc;

306:   PetscOptionsHead("BLZPACK options");

308:   bs = blz->block_size;
309:   PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
310:   if (flg) {EPSBlzpackSetBlockSize(eps,bs);}

312:   n = blz->nsteps;
313:   PetscOptionsInt("-eps_blzpack_nsteps","Number of steps","EPSBlzpackSetNSteps",n,&n,&flg);
314:   if (flg) {EPSBlzpackSetNSteps(eps,n);}

316:   interval[0] = blz->initial;
317:   interval[1] = blz->final;
318:   n = 2;
319:   PetscOptionsRealArray("-eps_blzpack_interval","Computational interval","EPSBlzpackSetInterval",interval,&n,&flg);
320:   if (flg) {
321:     if (n==1) interval[1]=interval[0];
322:     EPSBlzpackSetInterval(eps,interval[0],interval[1]);
323:   }

325:   if (blz->slice) {
326:     STSetType(eps->OP,STSINV);
327:     STGetKSP(eps->OP,&ksp);
328:     KSPSetType(ksp,KSPPREONLY);
329:     KSPGetPC(ksp,&pc);
330:     PCSetType(pc,PCCHOLESKY);
331:   }

333:   PetscOptionsTail();
334:   return(0);
335: }

340: PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,int bs)
341: {
342:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

345:   if (bs == PETSC_DEFAULT) blz->block_size = 3;
346:   else if (bs <= 0) {
347:     SETERRQ(1, "Incorrect block size");
348:   } else blz->block_size = bs;
349:   return(0);
350: }

355: /*@
356:    EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.

358:    Collective on EPS

360:    Input Parameters:
361: +  eps - the eigenproblem solver context
362: -  bs - block size

364:    Options Database Key:
365: .  -eps_blzpack_block_size - Sets the value of the block size

367:    Level: advanced

369: .seealso: EPSBlzpackSetInterval()
370: @*/
371: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,int bs)
372: {
373:   PetscErrorCode ierr, (*f)(EPS,int);

377:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",(void (**)())&f);
378:   if (f) {
379:     (*f)(eps,bs);
380:   }
381:   return(0);
382: }

387: PetscErrorCode EPSBlzpackSetInterval_BLZPACK(EPS eps,PetscReal initial,PetscReal final)
388: {
389:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

392:   blz->initial    = initial;
393:   blz->final      = final;
394:   blz->slice      = 1;
395:   return(0);
396: }

401: /*@
402:    EPSBlzpackSetInterval - Sets the computational interval for the BLZPACK
403:    package.

405:    Collective on EPS

407:    Input Parameters:
408: +  eps     - the eigenproblem solver context
409: .  initial - lower bound of the interval
410: -  final   - upper bound of the interval

412:    Options Database Key:
413: .  -eps_blzpack_interval - Sets the bounds of the interval (two values
414:    separated by commas)

416:    Note:
417:    The following possibilities are accepted (see Blzpack user's guide for
418:    details).
419:      initial>final: start seeking for eigenpairs in the upper bound
420:      initial<final: start in the lower bound
421:      initial=final: run around a single value (no interval)
422:    
423:    Level: advanced

425: .seealso: EPSBlzpackSetBlockSize()
426: @*/
427: PetscErrorCode EPSBlzpackSetInterval(EPS eps,PetscReal initial,PetscReal final)
428: {
429:   PetscErrorCode ierr, (*f)(EPS,PetscReal,PetscReal);

433:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetInterval_C",(void (**)())&f);
434:   if (f) {
435:     (*f)(eps,initial,final);
436:   }
437:   return(0);
438: }

443: PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,int nsteps)
444: {
445:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;

448:   blz->nsteps = nsteps == PETSC_DEFAULT ? 0 : nsteps;
449:   return(0);
450: }

455: /*@
456:    EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
457:    package.

459:    Collective on EPS

461:    Input Parameters:
462: +  eps     - the eigenproblem solver context
463: -  nsteps  - maximum number of steps

465:    Options Database Key:
466: .  -eps_blzpack_nsteps - Sets the maximum number of steps per run

468:    Level: advanced

470: @*/
471: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,int nsteps)
472: {
473:   PetscErrorCode ierr, (*f)(EPS,int);

477:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",(void (**)())&f);
478:   if (f) {
479:     (*f)(eps,nsteps);
480:   }
481:   return(0);
482: }

487: PetscErrorCode EPSCreate_BLZPACK(EPS eps)
488: {
490:   EPS_BLZPACK    *blzpack;

493:   PetscNew(EPS_BLZPACK,&blzpack);
494:   PetscLogObjectMemory(eps,sizeof(EPS_BLZPACK));
495:   eps->data                      = (void *) blzpack;
496:   eps->ops->solve                = EPSSolve_BLZPACK;
497:   eps->ops->setup                = EPSSetUp_BLZPACK;
498:   eps->ops->setfromoptions       = EPSSetFromOptions_BLZPACK;
499:   eps->ops->destroy              = EPSDestroy_BLZPACK;
500:   eps->ops->view                 = EPSView_BLZPACK;
501:   eps->ops->backtransform        = EPSBackTransform_BLZPACK;
502:   eps->ops->computevectors       = EPSComputeVectors_Default;

504:   blzpack->block_size = 3;
505:   blzpack->initial = 0.0;
506:   blzpack->final = 0.0;
507:   blzpack->slice = 0;
508:   blzpack->nsteps = 0;

510:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetBlockSize_C","EPSBlzpackSetBlockSize_BLZPACK",EPSBlzpackSetBlockSize_BLZPACK);
511:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetInterval_C","EPSBlzpackSetInterval_BLZPACK",EPSBlzpackSetInterval_BLZPACK);
512:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetNSteps_C","EPSBlzpackSetNSteps_BLZPACK",EPSBlzpackSetNSteps_BLZPACK);

514:   eps->which = EPS_SMALLEST_REAL;

516:   return(0);
517: }