Actual source code: blzpack.c

  1: /*
  2:        This file implements a wrapper to the BLZPACK package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

  8:       This file is part of SLEPc. See the README file for conditions of use
  9:       and additional information.
 10:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 11: */

 13:  #include src/eps/impls/blzpack/blzpackp.h

 15: const char* blzpack_error[33] = {
 16:   "",
 17:   "illegal data, LFLAG ",
 18:   "illegal data, dimension of (U), (V), (X) ",
 19:   "illegal data, leading dimension of (U), (V), (X) ",
 20:   "illegal data, leading dimension of (EIG) ",
 21:   "illegal data, number of required eigenpairs ",
 22:   "illegal data, Lanczos algorithm block size ",
 23:   "illegal data, maximum number of steps ",
 24:   "illegal data, number of starting vectors ",
 25:   "illegal data, number of eigenpairs provided ",
 26:   "illegal data, problem type flag ",
 27:   "illegal data, spectrum slicing flag ",
 28:   "illegal data, eigenvectors purification flag ",
 29:   "illegal data, level of output ",
 30:   "illegal data, output file unit ",
 31:   "illegal data, LCOMM (MPI or PVM) ",
 32:   "illegal data, dimension of ISTOR ",
 33:   "illegal data, convergence threshold ",
 34:   "illegal data, dimension of RSTOR ",
 35:   "illegal data on at least one PE ",
 36:   "ISTOR(3:14) must be equal on all PEs ",
 37:   "RSTOR(1:3) must be equal on all PEs ",
 38:   "not enough space in ISTOR to start eigensolution ",
 39:   "not enough space in RSTOR to start eigensolution ",
 40:   "illegal data, number of negative eigenvalues ",
 41:   "illegal data, entries of V ",
 42:   "illegal data, entries of X ",
 43:   "failure in computational subinterval ",
 44:   "file I/O error, blzpack.__.BQ ",
 45:   "file I/O error, blzpack.__.BX ",
 46:   "file I/O error, blzpack.__.Q ",
 47:   "file I/O error, blzpack.__.X ",
 48:   "parallel interface error "
 49: };

 53: PetscErrorCode EPSSetUp_BLZPACK(EPS eps)
 54: {
 56:   PetscInt       N, n;
 57:   int            listor, lrstor, ncuv, k1, k2, k3, k4;
 58:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
 59:   PetscTruth     flg;
 60:   KSP            ksp;
 61:   PC             pc;

 64:   VecGetSize(eps->vec_initial,&N);
 65:   VecGetLocalSize(eps->vec_initial,&n);
 66:   if (eps->ncv) {
 67:     if( eps->ncv < PetscMin(eps->nev+10,eps->nev*2) )
 68:       SETERRQ(0,"Warning: BLZpack recommends that ncv be larger than min(nev+10,nev*2)");
 69:   }
 70:   else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
 71:   if (!eps->max_it) eps->max_it = PetscMax(1000,N);

 73:   if (!eps->ishermitian)
 74:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 75:   if (blz->slice || eps->isgeneralized) {
 76:     PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
 77:     if (!flg)
 78:       SETERRQ(PETSC_ERR_SUP,"Shift-and-invert ST is needed for generalized problems or spectrum slicing");
 79:     STGetKSP(eps->OP,&ksp);
 80:     PetscTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
 81:     if (!flg)
 82:       SETERRQ(PETSC_ERR_SUP,"Preonly KSP is needed for generalized problems or spectrum slicing");
 83:     KSPGetPC(ksp,&pc);
 84:     PetscTypeCompare((PetscObject)pc,PCCHOLESKY,&flg);
 85:     if (!flg)
 86:       SETERRQ(PETSC_ERR_SUP,"Cholesky PC is needed for generalized problems or spectrum slicing");
 87:   }
 88:   if (eps->which!=EPS_SMALLEST_REAL)
 89:     SETERRQ(1,"Wrong value of eps->which");

 91:   k1 = PetscMin(N,180);
 92:   k2 = blz->block_size;
 93:   k4 = PetscMin(eps->ncv,N);
 94:   k3 = 484+k1*(13+k1*2+k2+PetscMax(18,k2+2))+k2*k2*3+k4*2;

 96:   listor = 123+k1*12;
 97:   PetscFree(blz->istor);
 98:   PetscMalloc((17+listor)*sizeof(int),&blz->istor);
 99:   blz->istor[14] = listor;

101:   if (blz->slice) lrstor = n*(k2*4+k1*2+k4)+k3;
102:   else lrstor = n*(k2*4+k1)+k3;
103:   PetscFree(blz->rstor);
104:   PetscMalloc((4+lrstor)*sizeof(PetscReal),&blz->rstor);
105:   blz->rstor[3] = lrstor;

107:   ncuv = PetscMax(3,blz->block_size);
108:   PetscFree(blz->u);
109:   PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->u);
110:   PetscFree(blz->v);
111:   PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->v);

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

116:   EPSAllocateSolutionContiguous(eps);
117:   EPSDefaultGetWork(eps,1);
118:   return(0);
119: }

123: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
124: {
126:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
127:   PetscInt       n, nn;
128:   int            i, nneig, lflag, nvopu;
129:   Vec            x, y;
130:   PetscScalar    sigma,*pV;
131:   Mat            A;
132:   KSP            ksp;
133:   PC             pc;
134: 

137:   VecGetLocalSize(eps->vec_initial,&n);
138:   VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&x);
139:   VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&y);
140:   VecGetArray(eps->V[0],&pV);
141: 
142:   if (eps->isgeneralized && !blz->slice) {
143:     STGetShift(eps->OP,&sigma); /* shift of origin */
144:     blz->rstor[0]  = sigma;        /* lower limit of eigenvalue interval */
145:     blz->rstor[1]  = sigma;        /* upper limit of eigenvalue interval */
146:   } else {
147:     sigma = 0.0;
148:     blz->rstor[0]  = blz->initial; /* lower limit of eigenvalue interval */
149:     blz->rstor[1]  = blz->final;   /* upper limit of eigenvalue interval */
150:   }
151:   nneig = 0;                     /* no. of eigs less than sigma */

153:   blz->istor[0]  = n;            /* number of rows of U, V, X*/
154:   blz->istor[1]  = n;            /* leading dimension of U, V, X */
155:   blz->istor[2]  = eps->nev;     /* number of required eigenpairs */
156:   blz->istor[3]  = eps->ncv;     /* number of working eigenpairs */
157:   blz->istor[4]  = blz->block_size;    /* number of vectors in a block */
158:   blz->istor[5]  = blz->nsteps;  /* maximun number of steps per run */
159:   blz->istor[6]  = 1;            /* number of starting vectors as input */
160:   blz->istor[7]  = 0;            /* number of eigenpairs given as input */
161:   blz->istor[8]  = (blz->slice || eps->isgeneralized) ? 1 : 0;   /* problem type */
162:   blz->istor[9]  = blz->slice;   /* spectrum slicing */
163:   blz->istor[10] = eps->isgeneralized ? 1 : 0;   /* solutions refinement (purify) */
164:   blz->istor[11] = 0;            /* level of printing */
165:   blz->istor[12] = 6;            /* file unit for output */
166:   blz->istor[13] = MPI_Comm_c2f(eps->comm);    /* communicator */

168:   blz->rstor[2]  = eps->tol;     /* threshold for convergence */

170:   lflag = 0;           /* reverse communication interface flag */

172:   do {

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

177:     switch (lflag) {
178:     case 1:
179:       /* compute v = OP u */
180:       for (i=0;i<nvopu;i++) {
181:         VecPlaceArray( x, blz->u+i*n );
182:         VecPlaceArray( y, blz->v+i*n );
183:         if (blz->slice || eps->isgeneralized) {
184:           STAssociatedKSPSolve( eps->OP, x, y );
185:         } else {
186:           STApply( eps->OP, x, y );
187:         }
188:         IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0]);
189:         VecResetArray(x);
190:         VecResetArray(y);
191:       }
192:       /* monitor */
193:       eps->nconv  = BLZistorr_(blz->istor,"NTEIG",5);
194:       EPSMonitor(eps,eps->its,eps->nconv,
195:         blz->rstor+BLZistorr_(blz->istor,"IRITZ",5),
196:         eps->eigi,
197:         blz->rstor+BLZistorr_(blz->istor,"IRITZ",5)+BLZistorr_(blz->istor,"JT",2),
198:         BLZistorr_(blz->istor,"NRITZ",5));
199:       eps->its = eps->its + 1;
200:       if (eps->its >= eps->max_it || eps->nconv >= eps->nev) lflag = 5;
201:       break;
202:     case 2:
203:       /* compute v = B u */
204:       for (i=0;i<nvopu;i++) {
205:         VecPlaceArray( x, blz->u+i*n );
206:         VecPlaceArray( y, blz->v+i*n );
207:         IPApplyMatrix(eps->ip, x, y );
208:         VecResetArray(x);
209:         VecResetArray(y);
210:       }
211:       break;
212:     case 3:
213:       /* update shift */
214:       PetscInfo1(eps,"Factorization update (sigma=%g)\n",sigma);
215:       STSetShift(eps->OP,sigma);
216:       STGetKSP(eps->OP,&ksp);
217:       KSPGetPC(ksp,&pc);
218:       PCGetFactoredMatrix(pc,&A);
219:       MatGetInertia(A,&nn,PETSC_NULL,PETSC_NULL);
220:       nneig = nn;
221:       break;
222:     case 4:
223:       /* copy the initial vector */
224:       VecPlaceArray(x,blz->v);
225:       EPSGetStartVector(eps,0,x,PETSC_NULL);
226:       VecResetArray(x);
227:       break;
228:     }
229: 
230:   } while (lflag > 0);

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

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

237:   for (i=0;i<eps->nconv;i++) {
238:     eps->eigr[i]=blz->eig[i];
239:   }

241:   if (lflag!=0) {
242:     char msg[2048] = "";
243:     for (i = 0; i < 33; i++) {
244:       if (blz->istor[15] & (1 << i)) PetscStrcat(msg, blzpack_error[i]);
245:     }
246:     SETERRQ2(PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15], msg);
247:   }
248:   VecDestroy(x);
249:   VecDestroy(y);

251:   return(0);
252: }

256: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
257: {
259:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

262:   if (!blz->slice && !eps->isgeneralized) {
263:     EPSBackTransform_Default(eps);
264:   }
265:   return(0);
266: }

270: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
271: {
273:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

277:   PetscFree(blz->istor);
278:   PetscFree(blz->rstor);
279:   PetscFree(blz->u);
280:   PetscFree(blz->v);
281:   PetscFree(blz->eig);
282:   PetscFree(eps->data);
283:   EPSFreeSolutionContiguous(eps);
284:   return(0);
285: }

289: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
290: {
292:   EPS_BLZPACK    *blz = (EPS_BLZPACK *) eps->data;
293:   PetscTruth     isascii;

296:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
297:   if (!isascii) {
298:     SETERRQ1(1,"Viewer type %s not supported for EPSBLZPACK",((PetscObject)viewer)->type_name);
299:   }
300:   PetscViewerASCIIPrintf(viewer,"block size of the block-Lanczos algorithm: %d\n",blz->block_size);
301:   PetscViewerASCIIPrintf(viewer,"computational interval: [%f,%f]\n",blz->initial,blz->final);
302:   return(0);
303: }

307: PetscErrorCode EPSSetFromOptions_BLZPACK(EPS eps)
308: {
310:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
311:   PetscInt       bs,n;
312:   PetscReal      interval[2];
313:   PetscTruth     flg;
314:   KSP            ksp;
315:   PC             pc;

318:   PetscOptionsHead("BLZPACK options");

320:   bs = blz->block_size;
321:   PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
322:   if (flg) {EPSBlzpackSetBlockSize(eps,bs);}

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

328:   interval[0] = blz->initial;
329:   interval[1] = blz->final;
330:   n = 2;
331:   PetscOptionsRealArray("-eps_blzpack_interval","Computational interval","EPSBlzpackSetInterval",interval,&n,&flg);
332:   if (flg) {
333:     if (n==1) interval[1]=interval[0];
334:     EPSBlzpackSetInterval(eps,interval[0],interval[1]);
335:   }

337:   if (blz->slice || eps->isgeneralized) {
338:     STSetType(eps->OP,STSINV);
339:     STGetKSP(eps->OP,&ksp);
340:     KSPSetType(ksp,KSPPREONLY);
341:     KSPGetPC(ksp,&pc);
342:     PCSetType(pc,PCCHOLESKY);
343:   }

345:   PetscOptionsTail();
346:   return(0);
347: }

352: PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,int bs)
353: {
354:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

357:   if (bs == PETSC_DEFAULT) blz->block_size = 3;
358:   else if (bs <= 0) {
359:     SETERRQ(1, "Incorrect block size");
360:   } else blz->block_size = bs;
361:   return(0);
362: }

367: /*@
368:    EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.

370:    Collective on EPS

372:    Input Parameters:
373: +  eps - the eigenproblem solver context
374: -  bs - block size

376:    Options Database Key:
377: .  -eps_blzpack_block_size - Sets the value of the block size

379:    Level: advanced

381: .seealso: EPSBlzpackSetInterval()
382: @*/
383: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,int bs)
384: {
385:   PetscErrorCode ierr, (*f)(EPS,int);

389:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",(void (**)())&f);
390:   if (f) {
391:     (*f)(eps,bs);
392:   }
393:   return(0);
394: }

399: PetscErrorCode EPSBlzpackSetInterval_BLZPACK(EPS eps,PetscReal initial,PetscReal final)
400: {
402:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

405:   blz->initial    = initial;
406:   blz->final      = final;
407:   blz->slice      = 1;
408:   STSetShift(eps->OP,initial);
409:   return(0);
410: }

415: /*@
416:    EPSBlzpackSetInterval - Sets the computational interval for the BLZPACK
417:    package.

419:    Collective on EPS

421:    Input Parameters:
422: +  eps     - the eigenproblem solver context
423: .  initial - lower bound of the interval
424: -  final   - upper bound of the interval

426:    Options Database Key:
427: .  -eps_blzpack_interval - Sets the bounds of the interval (two values
428:    separated by commas)

430:    Note:
431:    The following possibilities are accepted (see Blzpack user's guide for
432:    details).
433:      initial>final: start seeking for eigenpairs in the upper bound
434:      initial<final: start in the lower bound
435:      initial=final: run around a single value (no interval)
436:    
437:    Level: advanced

439: .seealso: EPSBlzpackSetBlockSize()
440: @*/
441: PetscErrorCode EPSBlzpackSetInterval(EPS eps,PetscReal initial,PetscReal final)
442: {
443:   PetscErrorCode ierr, (*f)(EPS,PetscReal,PetscReal);

447:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetInterval_C",(void (**)())&f);
448:   if (f) {
449:     (*f)(eps,initial,final);
450:   }
451:   return(0);
452: }

457: PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,int nsteps)
458: {
459:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;

462:   blz->nsteps = nsteps == PETSC_DEFAULT ? 0 : nsteps;
463:   return(0);
464: }

469: /*@
470:    EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
471:    package.

473:    Collective on EPS

475:    Input Parameters:
476: +  eps     - the eigenproblem solver context
477: -  nsteps  - maximum number of steps

479:    Options Database Key:
480: .  -eps_blzpack_nsteps - Sets the maximum number of steps per run

482:    Level: advanced

484: @*/
485: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,int nsteps)
486: {
487:   PetscErrorCode ierr, (*f)(EPS,int);

491:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",(void (**)())&f);
492:   if (f) {
493:     (*f)(eps,nsteps);
494:   }
495:   return(0);
496: }

501: PetscErrorCode EPSCreate_BLZPACK(EPS eps)
502: {
504:   EPS_BLZPACK    *blzpack;

507:   PetscNew(EPS_BLZPACK,&blzpack);
508:   PetscLogObjectMemory(eps,sizeof(EPS_BLZPACK));
509:   eps->data                      = (void *) blzpack;
510:   eps->ops->solve                = EPSSolve_BLZPACK;
511:   eps->ops->setup                = EPSSetUp_BLZPACK;
512:   eps->ops->setfromoptions       = EPSSetFromOptions_BLZPACK;
513:   eps->ops->destroy              = EPSDestroy_BLZPACK;
514:   eps->ops->view                 = EPSView_BLZPACK;
515:   eps->ops->backtransform        = EPSBackTransform_BLZPACK;
516:   eps->ops->computevectors       = EPSComputeVectors_Default;

518:   blzpack->block_size = 3;
519:   blzpack->initial = 0.0;
520:   blzpack->final = 0.0;
521:   blzpack->slice = 0;
522:   blzpack->nsteps = 0;

524:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetBlockSize_C","EPSBlzpackSetBlockSize_BLZPACK",EPSBlzpackSetBlockSize_BLZPACK);
525:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetInterval_C","EPSBlzpackSetInterval_BLZPACK",EPSBlzpackSetInterval_BLZPACK);
526:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetNSteps_C","EPSBlzpackSetNSteps_BLZPACK",EPSBlzpackSetNSteps_BLZPACK);

528:   eps->which = EPS_SMALLEST_REAL;

530:   return(0);
531: }