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-2013, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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 <slepc-private/epsimpl.h>    /*I "slepceps.h" I*/
 25: #include <slepc-private/stimpl.h>     /*I "slepcst.h" I*/
 26: #include <../src/eps/impls/external/blzpack/blzpackp.h>

 28: PetscErrorCode EPSSolve_BLZPACK(EPS);

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

 68: PetscErrorCode EPSSetUp_BLZPACK(EPS eps)
 69: {
 71:   PetscInt       listor,lrstor,ncuv,k1,k2,k3,k4;
 72:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;
 73:   PetscBool      issinv;

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

 82:   if (!blz->block_size) blz->block_size = 3;
 83:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requested method is only available for Hermitian problems");
 84:   if (eps->which==EPS_ALL) {
 85:     if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Must define a computational interval when using EPS_ALL");
 86:     blz->slice = 1;
 87:   }
 88:   PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&issinv);
 89:   if (blz->slice || eps->isgeneralized) {
 90:     if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert ST is needed for generalized problems or spectrum slicing");
 91:   }
 92:   if (blz->slice) {
 93:     if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
 94:       if (eps->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The defined computational interval should have at least one of their sides bounded");
 95:       STSetDefaultShift(eps->st,eps->inta);
 96:     } else {
 97:       STSetDefaultShift(eps->st,eps->intb);
 98:     }
 99:   }
100:   if (!eps->which) {
101:     if (issinv) eps->which = EPS_TARGET_REAL;
102:     else eps->which = EPS_SMALLEST_REAL;
103:   }
104:   if ((issinv && eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_MAGNITUDE && eps->which!=EPS_ALL) || (!issinv && eps->which!=EPS_SMALLEST_REAL)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
105:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

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

112:   listor = 123+k1*12;
113:   PetscFree(blz->istor);
114:   PetscMalloc((17+listor)*sizeof(PetscBLASInt),&blz->istor);
115:   PetscLogObjectMemory(eps,(17+listor)*sizeof(PetscBLASInt));
116:   PetscBLASIntCast(listor,&blz->istor[14]);

118:   if (blz->slice) lrstor = eps->nloc*(k2*4+k1*2+k4)+k3;
119:   else lrstor = eps->nloc*(k2*4+k1)+k3;
120: lrstor*=10;
121:   PetscFree(blz->rstor);
122:   PetscMalloc((4+lrstor)*sizeof(PetscReal),&blz->rstor);
123:   PetscLogObjectMemory(eps,(4+lrstor)*sizeof(PetscReal));
124:   blz->rstor[3] = lrstor;

126:   ncuv = PetscMax(3,blz->block_size);
127:   PetscFree(blz->u);
128:   PetscMalloc(ncuv*eps->nloc*sizeof(PetscScalar),&blz->u);
129:   PetscLogObjectMemory(eps,ncuv*eps->nloc*sizeof(PetscScalar));
130:   PetscFree(blz->v);
131:   PetscMalloc(ncuv*eps->nloc*sizeof(PetscScalar),&blz->v);
132:   PetscLogObjectMemory(eps,ncuv*eps->nloc*sizeof(PetscScalar));

134:   PetscFree(blz->eig);
135:   PetscMalloc(2*eps->ncv*sizeof(PetscReal),&blz->eig);
136:   PetscLogObjectMemory(eps,2*eps->ncv*sizeof(PetscReal));

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

140:   EPSAllocateSolution(eps);

142:   /* dispatch solve method */
143:   if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
144:   eps->ops->solve = EPSSolve_BLZPACK;
145:   return(0);
146: }

150: PetscErrorCode EPSSolve_BLZPACK(EPS eps)
151: {
153:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;
154:   PetscInt       nn;
155:   PetscBLASInt   i,nneig,lflag,nvopu;
156:   Vec            x,y;
157:   PetscScalar    sigma,*pV;
158:   Mat            A;
159:   KSP            ksp;
160:   PC             pc;

163:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
164:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
165:   VecGetArray(eps->V[0],&pV);

167:   if (eps->isgeneralized && !blz->slice) {
168:     STGetShift(eps->st,&sigma); /* shift of origin */
169:     blz->rstor[0]  = sigma;        /* lower limit of eigenvalue interval */
170:     blz->rstor[1]  = sigma;        /* upper limit of eigenvalue interval */
171:   } else {
172:     sigma = 0.0;
173:     blz->rstor[0]  = eps->inta;    /* lower limit of eigenvalue interval */
174:     blz->rstor[1]  = eps->intb;    /* upper limit of eigenvalue interval */
175:   }
176:   nneig = 0;                       /* no. of eigs less than sigma */

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

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

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

197:   do {
198:     BLZpack_(blz->istor,blz->rstor,&sigma,&nneig,blz->u,blz->v,&lflag,&nvopu,blz->eig,pV);

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

253:   } while (lflag > 0);

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

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

260:   for (i=0;i<eps->nconv;i++) {
261:     eps->eigr[i]=blz->eig[i];
262:   }

264:   if (lflag!=0) {
265:     char msg[2048] = "";
266:     for (i = 0; i < 33; i++) {
267:       if (blz->istor[15] & (1 << i)) PetscStrcat(msg,blzpack_error[i]);
268:     }
269:     SETERRQ2(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15],msg);
270:   }
271:   VecDestroy(&x);
272:   VecDestroy(&y);
273:   return(0);
274: }

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

284:   if (!blz->slice && !eps->isgeneralized) {
285:     EPSBackTransform_Default(eps);
286:   }
287:   return(0);
288: }

292: PetscErrorCode EPSReset_BLZPACK(EPS eps)
293: {
295:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;

298:   PetscFree(blz->istor);
299:   PetscFree(blz->rstor);
300:   PetscFree(blz->u);
301:   PetscFree(blz->v);
302:   PetscFree(blz->eig);
303:   EPSFreeSolution(eps);
304:   return(0);
305: }

309: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
310: {

314:   PetscFree(eps->data);
315:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",NULL);
316:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",NULL);
317:   return(0);
318: }

322: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
323: {
325:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;
326:   PetscBool      isascii;

329:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
330:   if (isascii) {
331:     PetscViewerASCIIPrintf(viewer,"  BLZPACK: block size=%d\n",blz->block_size);
332:     PetscViewerASCIIPrintf(viewer,"  BLZPACK: maximum number of steps per run=%d\n",blz->nsteps);
333:     if (blz->slice) {
334:       PetscViewerASCIIPrintf(viewer,"  BLZPACK: computational interval [%f,%f]\n",eps->inta,eps->intb);
335:     }
336:   }
337:   return(0);
338: }

342: PetscErrorCode EPSSetFromOptions_BLZPACK(EPS eps)
343: {
345:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;
346:   PetscInt       bs,n;
347:   PetscBool      flg;

350:   PetscOptionsHead("EPS BLZPACK Options");

352:   bs = blz->block_size;
353:   PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
354:   if (flg) {
355:     EPSBlzpackSetBlockSize(eps,bs);
356:   }

358:   n = blz->nsteps;
359:   PetscOptionsInt("-eps_blzpack_nsteps","Number of steps","EPSBlzpackSetNSteps",n,&n,&flg);
360:   if (flg) {
361:     EPSBlzpackSetNSteps(eps,n);
362:   }

364:   PetscOptionsTail();
365:   return(0);
366: }

370: static PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,PetscInt bs)
371: {
372:   EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;;

375:   if (bs == PETSC_DEFAULT) blz->block_size = 3;
376:   else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Block size must be positive");
377:   else {
378:     PetscBLASIntCast(bs,&blz->block_size);
379:   }
380:   return(0);
381: }

385: /*@
386:    EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.

388:    Collective on EPS

390:    Input Parameters:
391: +  eps - the eigenproblem solver context
392: -  bs - block size

394:    Options Database Key:
395: .  -eps_blzpack_block_size - Sets the value of the block size

397:    Level: advanced
398: @*/
399: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,PetscInt bs)
400: {

406:   PetscTryMethod(eps,"EPSBlzpackSetBlockSize_C",(EPS,PetscInt),(eps,bs));
407:   return(0);
408: }

412: static PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,PetscInt nsteps)
413: {
414:   EPS_BLZPACK *blz = (EPS_BLZPACK*)eps->data;

417:   if (nsteps == PETSC_DEFAULT) blz->nsteps = 0;
418:   else {
419:     PetscBLASIntCast(nsteps,&blz->nsteps);
420:   }
421:   return(0);
422: }

426: /*@
427:    EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
428:    package.

430:    Collective on EPS

432:    Input Parameters:
433: +  eps     - the eigenproblem solver context
434: -  nsteps  - maximum number of steps

436:    Options Database Key:
437: .  -eps_blzpack_nsteps - Sets the maximum number of steps per run

439:    Level: advanced

441: @*/
442: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,PetscInt nsteps)
443: {

449:   PetscTryMethod(eps,"EPSBlzpackSetNSteps_C",(EPS,PetscInt),(eps,nsteps));
450:   return(0);
451: }

455: PETSC_EXTERN PetscErrorCode EPSCreate_BLZPACK(EPS eps)
456: {
458:   EPS_BLZPACK    *blzpack;

461:   PetscNewLog(eps,EPS_BLZPACK,&blzpack);
462:   eps->data                      = (void*)blzpack;
463:   eps->ops->setup                = EPSSetUp_BLZPACK;
464:   eps->ops->setfromoptions       = EPSSetFromOptions_BLZPACK;
465:   eps->ops->destroy              = EPSDestroy_BLZPACK;
466:   eps->ops->reset                = EPSReset_BLZPACK;
467:   eps->ops->view                 = EPSView_BLZPACK;
468:   eps->ops->backtransform        = EPSBackTransform_BLZPACK;
469:   eps->ops->computevectors       = EPSComputeVectors_Default;
470:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",EPSBlzpackSetBlockSize_BLZPACK);
471:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",EPSBlzpackSetNSteps_BLZPACK);
472:   return(0);
473: }