Actual source code: blzpack.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    This file implements a wrapper to the BLZPACK package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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 <../src/eps/impls/external/blzpack/blzpackp.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:   PetscBool      issinv,istrivial,flg;

 75:   if (eps->ncv) {
 76:     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)");
 77:   } else eps->ncv = PetscMin(eps->nev+10,eps->nev*2);
 78:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 79:   if (!eps->max_it) eps->max_it = PetscMax(1000,eps->n);

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

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

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

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

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

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

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

139:   EPSAllocateSolution(eps,0);
140:   EPS_SetInnerProduct(eps);
141:   PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
142:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
143:   RGIsTrivial(eps->rg,&istrivial);
144:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
145:   if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");

147:   /* dispatch solve method */
148:   eps->ops->solve = EPSSolve_BLZPACK;
149:   return(0);
150: }

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

167:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
168:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
169:   EPSGetStartVector(eps,0,NULL);
170:   BVSetActiveColumns(eps->V,0,0);  /* just for deflation space */
171:   BVGetColumn(eps->V,0,&v0);
172:   VecGetArray(v0,&pV);

174:   if (eps->isgeneralized && !blz->slice) {
175:     STGetShift(eps->st,&sigma); /* shift of origin */
176:     blz->rstor[0]  = sigma;        /* lower limit of eigenvalue interval */
177:     blz->rstor[1]  = sigma;        /* upper limit of eigenvalue interval */
178:   } else {
179:     sigma = 0.0;
180:     blz->rstor[0]  = eps->inta;    /* lower limit of eigenvalue interval */
181:     blz->rstor[1]  = eps->intb;    /* upper limit of eigenvalue interval */
182:   }
183:   nneig = 0;                       /* no. of eigs less than sigma */

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

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

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

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

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

260:   } while (lflag > 0);

262:   VecRestoreArray(v0,&pV);
263:   BVRestoreColumn(eps->V,0,&v0);

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

268:   for (i=0;i<eps->nconv;i++) {
269:     eps->eigr[i]=blz->eig[i];
270:   }

272:   if (lflag!=0) {
273:     char msg[2048] = "";
274:     for (i = 0; i < 33; i++) {
275:       if (blz->istor[15] & (1 << i)) PetscStrcat(msg,blzpack_error[i]);
276:     }
277:     SETERRQ2(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in BLZPACK (code=%d): '%s'",blz->istor[15],msg);
278:   }
279:   VecDestroy(&x);
280:   VecDestroy(&y);
281:   return(0);
282: }

286: PetscErrorCode EPSBackTransform_BLZPACK(EPS eps)
287: {
289:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;

292:   if (!blz->slice && !eps->isgeneralized) {
293:     EPSBackTransform_Default(eps);
294:   }
295:   return(0);
296: }

300: PetscErrorCode EPSReset_BLZPACK(EPS eps)
301: {
303:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;

306:   PetscFree(blz->istor);
307:   PetscFree(blz->rstor);
308:   PetscFree(blz->u);
309:   PetscFree(blz->v);
310:   PetscFree(blz->eig);
311:   return(0);
312: }

316: PetscErrorCode EPSDestroy_BLZPACK(EPS eps)
317: {

321:   PetscFree(eps->data);
322:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",NULL);
323:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",NULL);
324:   return(0);
325: }

329: PetscErrorCode EPSView_BLZPACK(EPS eps,PetscViewer viewer)
330: {
332:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;
333:   PetscBool      isascii;

336:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
337:   if (isascii) {
338:     PetscViewerASCIIPrintf(viewer,"  BLZPACK: block size=%d\n",blz->block_size);
339:     PetscViewerASCIIPrintf(viewer,"  BLZPACK: maximum number of steps per run=%d\n",blz->nsteps);
340:     if (blz->slice) {
341:       PetscViewerASCIIPrintf(viewer,"  BLZPACK: computational interval [%f,%f]\n",eps->inta,eps->intb);
342:     }
343:   }
344:   return(0);
345: }

349: PetscErrorCode EPSSetFromOptions_BLZPACK(PetscOptionItems *PetscOptionsObject,EPS eps)
350: {
352:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;
353:   PetscInt       bs,n;
354:   PetscBool      flg;

357:   PetscOptionsHead(PetscOptionsObject,"EPS BLZPACK Options");

359:   bs = blz->block_size;
360:   PetscOptionsInt("-eps_blzpack_block_size","Block size","EPSBlzpackSetBlockSize",bs,&bs,&flg);
361:   if (flg) {
362:     EPSBlzpackSetBlockSize(eps,bs);
363:   }

365:   n = blz->nsteps;
366:   PetscOptionsInt("-eps_blzpack_nsteps","Number of steps","EPSBlzpackSetNSteps",n,&n,&flg);
367:   if (flg) {
368:     EPSBlzpackSetNSteps(eps,n);
369:   }

371:   PetscOptionsTail();
372:   return(0);
373: }

377: static PetscErrorCode EPSBlzpackSetBlockSize_BLZPACK(EPS eps,PetscInt bs)
378: {
380:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;

383:   if (bs == PETSC_DEFAULT) blz->block_size = 3;
384:   else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Block size must be positive");
385:   else {
386:     PetscBLASIntCast(bs,&blz->block_size);
387:   }
388:   return(0);
389: }

393: /*@
394:    EPSBlzpackSetBlockSize - Sets the block size for the BLZPACK package.

396:    Collective on EPS

398:    Input Parameters:
399: +  eps - the eigenproblem solver context
400: -  bs - block size

402:    Options Database Key:
403: .  -eps_blzpack_block_size - Sets the value of the block size

405:    Level: advanced
406: @*/
407: PetscErrorCode EPSBlzpackSetBlockSize(EPS eps,PetscInt bs)
408: {

414:   PetscTryMethod(eps,"EPSBlzpackSetBlockSize_C",(EPS,PetscInt),(eps,bs));
415:   return(0);
416: }

420: static PetscErrorCode EPSBlzpackSetNSteps_BLZPACK(EPS eps,PetscInt nsteps)
421: {
423:   EPS_BLZPACK    *blz = (EPS_BLZPACK*)eps->data;

426:   if (nsteps == PETSC_DEFAULT) blz->nsteps = 0;
427:   else {
428:     PetscBLASIntCast(nsteps,&blz->nsteps);
429:   }
430:   return(0);
431: }

435: /*@
436:    EPSBlzpackSetNSteps - Sets the maximum number of steps per run for the BLZPACK
437:    package.

439:    Collective on EPS

441:    Input Parameters:
442: +  eps     - the eigenproblem solver context
443: -  nsteps  - maximum number of steps

445:    Options Database Key:
446: .  -eps_blzpack_nsteps - Sets the maximum number of steps per run

448:    Level: advanced

450: @*/
451: PetscErrorCode EPSBlzpackSetNSteps(EPS eps,PetscInt nsteps)
452: {

458:   PetscTryMethod(eps,"EPSBlzpackSetNSteps_C",(EPS,PetscInt),(eps,nsteps));
459:   return(0);
460: }

464: PETSC_EXTERN PetscErrorCode EPSCreate_BLZPACK(EPS eps)
465: {
467:   EPS_BLZPACK    *blzpack;

470:   PetscNewLog(eps,&blzpack);
471:   eps->data = (void*)blzpack;

473:   eps->ops->setup                = EPSSetUp_BLZPACK;
474:   eps->ops->setfromoptions       = EPSSetFromOptions_BLZPACK;
475:   eps->ops->destroy              = EPSDestroy_BLZPACK;
476:   eps->ops->reset                = EPSReset_BLZPACK;
477:   eps->ops->view                 = EPSView_BLZPACK;
478:   eps->ops->backtransform        = EPSBackTransform_BLZPACK;
479:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetBlockSize_C",EPSBlzpackSetBlockSize_BLZPACK);
480:   PetscObjectComposeFunction((PetscObject)eps,"EPSBlzpackSetNSteps_C",EPSBlzpackSetNSteps_BLZPACK);
481:   return(0);
482: }