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-2010, 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/external/blzpack/blzpackp.h
 25:  #include private/stimpl.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:   PetscTruth     flg;
 73:   KSP            ksp;
 74:   PC             pc;

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

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

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

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

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

120:   ncuv = PetscMax(3,blz->block_size);
121:   PetscFree(blz->u);
122:   PetscMalloc(ncuv*eps->nloc*sizeof(PetscScalar),&blz->u);
123:   PetscFree(blz->v);
124:   PetscMalloc(ncuv*eps->nloc*sizeof(PetscScalar),&blz->v);

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

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

133:   EPSAllocateSolution(eps);

135:   /* dispatch solve method */
136:   if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
137:   eps->ops->solve = EPSSolve_BLZPACK;
138:   return(0);
139: }

143: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
144: {
146:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
147:   PetscInt       nn;
148:   PetscBLASInt   i, nneig, lflag, nvopu;
149:   Vec            x, y;
150:   PetscScalar    sigma,*pV;
151:   Mat            A;
152:   KSP            ksp;
153:   PC             pc;
154: 

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

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

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

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

191:   do {

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

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

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

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

256:   for (i=0;i<eps->nconv;i++) {
257:     eps->eigr[i]=blz->eig[i];
258:   }

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

270:   return(0);
271: }

275: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
276: {
278:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

281:   if (!blz->slice && !eps->isgeneralized) {
282:     EPSBackTransform_Default(eps);
283:   }
284:   return(0);
285: }

289: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
290: {
292:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;

296:   PetscFree(blz->istor);
297:   PetscFree(blz->rstor);
298:   PetscFree(blz->u);
299:   PetscFree(blz->v);
300:   PetscFree(blz->eig);
301:   PetscFree(eps->data);
302:   EPSFreeSolution(eps);
303:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetBlockSize_C","",PETSC_NULL);
304:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetInterval_C","",PETSC_NULL);
305:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetNSteps_C","",PETSC_NULL);
306:   return(0);
307: }

311: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
312: {
314:   EPS_BLZPACK    *blz = (EPS_BLZPACK *) eps->data;
315:   PetscTruth     isascii;

318:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
319:   if (!isascii) {
320:     SETERRQ1(1,"Viewer type %s not supported for EPSBLZPACK",((PetscObject)viewer)->type_name);
321:   }
322:   PetscViewerASCIIPrintf(viewer,"block size of the block-Lanczos algorithm: %d\n",blz->block_size);
323:   PetscViewerASCIIPrintf(viewer,"computational interval: [%f,%f]\n",blz->initial,blz->final);
324:   return(0);
325: }

329: PetscErrorCode EPSSetFromOptions_BLZPACK(EPS eps)
330: {
332:   EPS_BLZPACK    *blz = (EPS_BLZPACK *)eps->data;
333:   PetscInt       bs,n;
334:   PetscReal      interval[2];
335:   PetscTruth     flg;
336:   KSP            ksp;
337:   PC             pc;

340:   PetscOptionsBegin(((PetscObject)eps)->comm,((PetscObject)eps)->prefix,"BLZPACK Options","EPS");

342:   bs = blz->block_size;
343:   PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
344:   if (flg) {EPSBlzpackSetBlockSize(eps,bs);}

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

350:   interval[0] = blz->initial;
351:   interval[1] = blz->final;
352:   n = 2;
353:   PetscOptionsRealArray("-eps_blzpack_interval","Computational interval","EPSBlzpackSetInterval",interval,&n,&flg);
354:   if (flg) {
355:     if (n==1) interval[1]=interval[0];
356:     EPSBlzpackSetInterval(eps,interval[0],interval[1]);
357:   }

359:   if (blz->slice || eps->isgeneralized) {
360:     STSetType(eps->OP,STSINVERT);
361:     STGetKSP(eps->OP,&ksp);
362:     KSPSetType(ksp,KSPPREONLY);
363:     KSPGetPC(ksp,&pc);
364:     PCSetType(pc,PCCHOLESKY);
365:   }

367:   PetscOptionsEnd();
368:   return(0);
369: }

374: PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,PetscInt bs)
375: {
376:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

379:   if (bs == PETSC_DEFAULT) blz->block_size = 3;
380:   else if (bs <= 0) {
381:     SETERRQ(1, "Incorrect block size");
382:   } else blz->block_size = PetscBLASIntCast(bs);
383:   return(0);
384: }

389: /*@
390:    EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.

392:    Collective on EPS

394:    Input Parameters:
395: +  eps - the eigenproblem solver context
396: -  bs - block size

398:    Options Database Key:
399: .  -eps_blzpack_block_size - Sets the value of the block size

401:    Level: advanced

403: .seealso: EPSBlzpackSetInterval()
404: @*/
405: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,PetscInt bs)
406: {
407:   PetscErrorCode ierr, (*f)(EPS,PetscInt);

411:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",(void (**)())&f);
412:   if (f) {
413:     (*f)(eps,bs);
414:   }
415:   return(0);
416: }

421: PetscErrorCode EPSBlzpackSetInterval_BLZPACK(EPS eps,PetscReal initial,PetscReal final)
422: {
424:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;;

427:   blz->initial    = initial;
428:   blz->final      = final;
429:   blz->slice      = 1;
430:   STSetShift(eps->OP,initial);
431:   return(0);
432: }

437: /*@
438:    EPSBlzpackSetInterval - Sets the computational interval for the BLZPACK
439:    package.

441:    Collective on EPS

443:    Input Parameters:
444: +  eps     - the eigenproblem solver context
445: .  initial - lower bound of the interval
446: -  final   - upper bound of the interval

448:    Options Database Key:
449: .  -eps_blzpack_interval - Sets the bounds of the interval (two values
450:    separated by commas)

452:    Note:
453:    The following possibilities are accepted (see Blzpack user's guide for
454:    details).
455:      initial>final: start seeking for eigenpairs in the upper bound
456:      initial<final: start in the lower bound
457:      initial=final: run around a single value (no interval)
458:    
459:    Level: advanced

461: .seealso: EPSBlzpackSetBlockSize()
462: @*/
463: PetscErrorCode EPSBlzpackSetInterval(EPS eps,PetscReal initial,PetscReal final)
464: {
465:   PetscErrorCode ierr, (*f)(EPS,PetscReal,PetscReal);

469:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetInterval_C",(void (**)())&f);
470:   if (f) {
471:     (*f)(eps,initial,final);
472:   }
473:   return(0);
474: }

479: PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,PetscInt nsteps)
480: {
481:   EPS_BLZPACK *blz = (EPS_BLZPACK *) eps->data;

484:   if (nsteps == PETSC_DEFAULT) blz->nsteps = 0;
485:   else { blz->nsteps = PetscBLASIntCast(nsteps); }
486:   return(0);
487: }

492: /*@
493:    EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
494:    package.

496:    Collective on EPS

498:    Input Parameters:
499: +  eps     - the eigenproblem solver context
500: -  nsteps  - maximum number of steps

502:    Options Database Key:
503: .  -eps_blzpack_nsteps - Sets the maximum number of steps per run

505:    Level: advanced

507: @*/
508: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,PetscInt nsteps)
509: {
510:   PetscErrorCode ierr, (*f)(EPS,PetscInt);

514:   PetscObjectQueryFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",(void (**)())&f);
515:   if (f) {
516:     (*f)(eps,nsteps);
517:   }
518:   return(0);
519: }

524: PetscErrorCode EPSCreate_BLZPACK(EPS eps)
525: {
527:   EPS_BLZPACK    *blzpack;

530:   PetscNew(EPS_BLZPACK,&blzpack);
531:   PetscLogObjectMemory(eps,sizeof(EPS_BLZPACK));
532:   eps->data                      = (void *) blzpack;
533:   eps->ops->setup                = EPSSetUp_BLZPACK;
534:   eps->ops->setfromoptions       = EPSSetFromOptions_BLZPACK;
535:   eps->ops->destroy              = EPSDestroy_BLZPACK;
536:   eps->ops->view                 = EPSView_BLZPACK;
537:   eps->ops->backtransform        = EPSBackTransform_BLZPACK;
538:   eps->ops->computevectors       = EPSComputeVectors_Default;

540:   blzpack->block_size = 3;
541:   blzpack->initial = 0.0;
542:   blzpack->final = 0.0;
543:   blzpack->slice = 0;
544:   blzpack->nsteps = 0;

546:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetBlockSize_C","EPSBlzpackSetBlockSize_BLZPACK",EPSBlzpackSetBlockSize_BLZPACK);
547:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetInterval_C","EPSBlzpackSetInterval_BLZPACK",EPSBlzpackSetInterval_BLZPACK);
548:   PetscObjectComposeFunctionDynamic((PetscObject)eps,"EPSBlzpackSetNSteps_C","EPSBlzpackSetNSteps_BLZPACK",EPSBlzpackSetNSteps_BLZPACK);
549:   return(0);
550: }