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-2009, 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 src/eps/impls/blzpack/blzpackp.h

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

 64: PetscErrorCode EPSSetUp_BLZPACK(EPS eps)
 65: {
 67:   PetscInt       N, n, listor, lrstor, ncuv, k1, k2, k3, k4;
 68:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
 69:   PetscTruth     flg;
 70:   KSP            ksp;
 71:   PC             pc;

 74:   VecGetSize(eps->vec_initial,&N);
 75:   VecGetLocalSize(eps->vec_initial,&n);
 76:   if (eps->ncv) {
 77:     if( eps->ncv < PetscMin(eps->nev+10,eps->nev*2) )
 78:       SETERRQ(0,"Warning: BLZpack recommends that ncv be larger than min(nev+10,nev*2)");
 79:   }
 80:   else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
 81:   if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
 82:   if (!eps->max_it) eps->max_it = PetscMax(1000,N);

 84:   if (!eps->ishermitian)
 85:     SETERRQ(PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 86:   if (blz->slice || eps->isgeneralized) {
 87:     PetscTypeCompare((PetscObject)eps->OP,STSINV,&flg);
 88:     if (!flg)
 89:       SETERRQ(PETSC_ERR_SUP,"Shift-and-invert ST is needed for generalized problems or spectrum slicing");
 90:     STGetKSP(eps->OP,&ksp);
 91:     PetscTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
 92:     if (!flg)
 93:       SETERRQ(PETSC_ERR_SUP,"Preonly KSP is needed for generalized problems or spectrum slicing");
 94:     KSPGetPC(ksp,&pc);
 95:     PetscTypeCompare((PetscObject)pc,PCCHOLESKY,&flg);
 96:     if (!flg)
 97:       SETERRQ(PETSC_ERR_SUP,"Cholesky PC is needed for generalized problems or spectrum slicing");
 98:   }
 99:   if (eps->which!=EPS_SMALLEST_REAL)
100:     SETERRQ(1,"Wrong value of eps->which");

102:   k1 = PetscMin(N,180);
103:   k2 = blz->block_size;
104:   k4 = PetscMin(eps->ncv,N);
105:   k3 = 484+k1*(13+k1*2+k2+PetscMax(18,k2+2))+k2*k2*3+k4*2;

107:   listor = 123+k1*12;
108:   PetscFree(blz->istor);
109:   PetscMalloc((17+listor)*sizeof(PetscBLASInt),&blz->istor);
110:   blz->istor[14] = PetscBLASIntCast(listor);

112:   if (blz->slice) lrstor = n*(k2*4+k1*2+k4)+k3;
113:   else lrstor = n*(k2*4+k1)+k3;
114:   PetscFree(blz->rstor);
115:   PetscMalloc((4+lrstor)*sizeof(PetscReal),&blz->rstor);
116:   blz->rstor[3] = lrstor;

118:   ncuv = PetscMax(3,blz->block_size);
119:   PetscFree(blz->u);
120:   PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->u);
121:   PetscFree(blz->v);
122:   PetscMalloc(ncuv*n*sizeof(PetscScalar),&blz->v);

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

127:   if (eps->extraction) {
128:      PetscInfo(eps,"Warning: extraction type ignored\n");
129:   }

131:   EPSAllocateSolution(eps);
132:   EPSDefaultGetWork(eps,1);
133:   return(0);
134: }

138: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
139: {
141:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
142:   PetscInt       n, nn;
143:   PetscBLASInt   i, nneig, lflag, nvopu;
144:   Vec            x, y;
145:   PetscScalar    sigma,*pV;
146:   Mat            A;
147:   KSP            ksp;
148:   PC             pc;
149: 

152:   VecGetLocalSize(eps->vec_initial,&n);
153:   VecCreateMPIWithArray(((PetscObject)eps)->comm,n,PETSC_DECIDE,PETSC_NULL,&x);
154:   VecCreateMPIWithArray(((PetscObject)eps)->comm,n,PETSC_DECIDE,PETSC_NULL,&y);
155:   VecGetArray(eps->V[0],&pV);
156: 
157:   if (eps->isgeneralized && !blz->slice) {
158:     STGetShift(eps->OP,&sigma); /* shift of origin */
159:     blz->rstor[0]  = sigma;        /* lower limit of eigenvalue interval */
160:     blz->rstor[1]  = sigma;        /* upper limit of eigenvalue interval */
161:   } else {
162:     sigma = 0.0;
163:     blz->rstor[0]  = blz->initial; /* lower limit of eigenvalue interval */
164:     blz->rstor[1]  = blz->final;   /* upper limit of eigenvalue interval */
165:   }
166:   nneig = 0;                     /* no. of eigs less than sigma */

168:   blz->istor[0]  = PetscBLASIntCast(n); /* number of rows of U, V, X*/
169:   blz->istor[1]  = PetscBLASIntCast(n); /* leading dimension of U, V, X */
170:   blz->istor[2]  = PetscBLASIntCast(eps->nev); /* number of required eigenpairs */
171:   blz->istor[3]  = PetscBLASIntCast(eps->ncv); /* number of working eigenpairs */
172:   blz->istor[4]  = blz->block_size;    /* number of vectors in a block */
173:   blz->istor[5]  = blz->nsteps;  /* maximun number of steps per run */
174:   blz->istor[6]  = 1;            /* number of starting vectors as input */
175:   blz->istor[7]  = 0;            /* number of eigenpairs given as input */
176:   blz->istor[8]  = (blz->slice || eps->isgeneralized) ? 1 : 0;   /* problem type */
177:   blz->istor[9]  = blz->slice;   /* spectrum slicing */
178:   blz->istor[10] = eps->isgeneralized ? 1 : 0;   /* solutions refinement (purify) */
179:   blz->istor[11] = 0;            /* level of printing */
180:   blz->istor[12] = 6;            /* file unit for output */
181:   blz->istor[13] = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm)); /* communicator */

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

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

187:   do {

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

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

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

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

252:   for (i=0;i<eps->nconv;i++) {
253:     eps->eigr[i]=blz->eig[i];
254:   }

256:   if (lflag!=0) {
257:     char msg[2048] = "";
258:     for (i = 0; i < 33; i++) {
259:       if (blz->istor[15] & (1 << i)) PetscStrcat(msg, blzpack_error[i]);
260:     }
261:     SETERRQ2(PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15], msg);
262:   }
263:   VecDestroy(x);
264:   VecDestroy(y);

266:   return(0);
267: }

271: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
272: {
274:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

277:   if (!blz->slice && !eps->isgeneralized) {
278:     EPSBackTransform_Default(eps);
279:   }
280:   return(0);
281: }

285: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
286: {
288:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

292:   PetscFree(blz->istor);
293:   PetscFree(blz->rstor);
294:   PetscFree(blz->u);
295:   PetscFree(blz->v);
296:   PetscFree(blz->eig);
297:   PetscFree(eps->data);
298:   EPSFreeSolution(eps);
299:   return(0);
300: }

304: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
305: {
307:   EPS_BLZPACK    *blz = (EPS_BLZPACK *) eps->data;
308:   PetscTruth     isascii;

311:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
312:   if (!isascii) {
313:     SETERRQ1(1,"Viewer type %s not supported for EPSBLZPACK",((PetscObject)viewer)->type_name);
314:   }
315:   PetscViewerASCIIPrintf(viewer,"block size of the block-Lanczos algorithm: %d\n",blz->block_size);
316:   PetscViewerASCIIPrintf(viewer,"computational interval: [%f,%f]\n",blz->initial,blz->final);
317:   return(0);
318: }

322: PetscErrorCode EPSSetFromOptions_BLZPACK(EPS eps)
323: {
325:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
326:   PetscInt       bs,n;
327:   PetscReal      interval[2];
328:   PetscTruth     flg;
329:   KSP            ksp;
330:   PC             pc;

333:   PetscOptionsHead("BLZPACK options");

335:   bs = blz->block_size;
336:   PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
337:   if (flg) {EPSBlzpackSetBlockSize(eps,bs);}

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

343:   interval[0] = blz->initial;
344:   interval[1] = blz->final;
345:   n = 2;
346:   PetscOptionsRealArray("-eps_blzpack_interval","Computational interval","EPSBlzpackSetInterval",interval,&n,&flg);
347:   if (flg) {
348:     if (n==1) interval[1]=interval[0];
349:     EPSBlzpackSetInterval(eps,interval[0],interval[1]);
350:   }

352:   if (blz->slice || eps->isgeneralized) {
353:     STSetType(eps->OP,STSINV);
354:     STGetKSP(eps->OP,&ksp);
355:     KSPSetType(ksp,KSPPREONLY);
356:     KSPGetPC(ksp,&pc);
357:     PCSetType(pc,PCCHOLESKY);
358:   }

360:   PetscOptionsTail();
361:   return(0);
362: }

367: PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,PetscInt bs)
368: {
369:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

372:   if (bs == PETSC_DEFAULT) blz->block_size = 3;
373:   else if (bs <= 0) {
374:     SETERRQ(1, "Incorrect block size");
375:   } else blz->block_size = PetscBLASIntCast(bs);
376:   return(0);
377: }

382: /*@
383:    EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.

385:    Collective on EPS

387:    Input Parameters:
388: +  eps - the eigenproblem solver context
389: -  bs - block size

391:    Options Database Key:
392: .  -eps_blzpack_block_size - Sets the value of the block size

394:    Level: advanced

396: .seealso: EPSBlzpackSetInterval()
397: @*/
398: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,PetscInt bs)
399: {
400:   PetscErrorCode ierr, (*f)(EPS,PetscInt);

404:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",(void (**)())&f);
405:   if (f) {
406:     (*f)(eps,bs);
407:   }
408:   return(0);
409: }

414: PetscErrorCode EPSBlzpackSetInterval_BLZPACK(EPS eps,PetscReal initial,PetscReal final)
415: {
417:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

420:   blz->initial    = initial;
421:   blz->final      = final;
422:   blz->slice      = 1;
423:   STSetShift(eps->OP,initial);
424:   return(0);
425: }

430: /*@
431:    EPSBlzpackSetInterval - Sets the computational interval for the BLZPACK
432:    package.

434:    Collective on EPS

436:    Input Parameters:
437: +  eps     - the eigenproblem solver context
438: .  initial - lower bound of the interval
439: -  final   - upper bound of the interval

441:    Options Database Key:
442: .  -eps_blzpack_interval - Sets the bounds of the interval (two values
443:    separated by commas)

445:    Note:
446:    The following possibilities are accepted (see Blzpack user's guide for
447:    details).
448:      initial>final: start seeking for eigenpairs in the upper bound
449:      initial<final: start in the lower bound
450:      initial=final: run around a single value (no interval)
451:    
452:    Level: advanced

454: .seealso: EPSBlzpackSetBlockSize()
455: @*/
456: PetscErrorCode EPSBlzpackSetInterval(EPS eps,PetscReal initial,PetscReal final)
457: {
458:   PetscErrorCode ierr, (*f)(EPS,PetscReal,PetscReal);

462:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetInterval_C",(void (**)())&f);
463:   if (f) {
464:     (*f)(eps,initial,final);
465:   }
466:   return(0);
467: }

472: PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,PetscInt nsteps)
473: {
474:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;

477:   if (nsteps == PETSC_DEFAULT) blz->nsteps = 0;
478:   else { blz->nsteps = PetscBLASIntCast(nsteps); }
479:   return(0);
480: }

485: /*@
486:    EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
487:    package.

489:    Collective on EPS

491:    Input Parameters:
492: +  eps     - the eigenproblem solver context
493: -  nsteps  - maximum number of steps

495:    Options Database Key:
496: .  -eps_blzpack_nsteps - Sets the maximum number of steps per run

498:    Level: advanced

500: @*/
501: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,PetscInt nsteps)
502: {
503:   PetscErrorCode ierr, (*f)(EPS,PetscInt);

507:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",(void (**)())&f);
508:   if (f) {
509:     (*f)(eps,nsteps);
510:   }
511:   return(0);
512: }

517: PetscErrorCode EPSCreate_BLZPACK(EPS eps)
518: {
520:   EPS_BLZPACK    *blzpack;

523:   PetscNew(EPS_BLZPACK,&blzpack);
524:   PetscLogObjectMemory(eps,sizeof(EPS_BLZPACK));
525:   eps->data                      = (void *) blzpack;
526:   eps->ops->solve                = EPSSolve_BLZPACK;
527:   eps->ops->setup                = EPSSetUp_BLZPACK;
528:   eps->ops->setfromoptions       = EPSSetFromOptions_BLZPACK;
529:   eps->ops->destroy              = EPSDestroy_BLZPACK;
530:   eps->ops->view                 = EPSView_BLZPACK;
531:   eps->ops->backtransform        = EPSBackTransform_BLZPACK;
532:   eps->ops->computevectors       = EPSComputeVectors_Default;

534:   blzpack->block_size = 3;
535:   blzpack->initial = 0.0;
536:   blzpack->final = 0.0;
537:   blzpack->slice = 0;
538:   blzpack->nsteps = 0;

540:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetBlockSize_C","EPSBlzpackSetBlockSize_BLZPACK",EPSBlzpackSetBlockSize_BLZPACK);
541:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetInterval_C","EPSBlzpackSetInterval_BLZPACK",EPSBlzpackSetInterval_BLZPACK);
542:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetNSteps_C","EPSBlzpackSetNSteps_BLZPACK",EPSBlzpackSetNSteps_BLZPACK);

544:   eps->which = EPS_SMALLEST_REAL;

546:   return(0);
547: }