Actual source code: arpack.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    This file implements a wrapper to the ARPACK package
 12: */

 14: #include <slepc/private/epsimpl.h>
 15: #include <../src/eps/impls/external/arpack/arpackp.h>

 17: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
 18: {
 20:   PetscInt       ncv;
 21:   PetscBool      flg,istrivial;
 22:   EPS_ARPACK     *ar = (EPS_ARPACK*)eps->data;

 25:   if (eps->ncv) {
 26:     if (eps->ncv<eps->nev+2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The value of ncv must be at least nev+2");
 27:   } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n); /* set default value of ncv */
 28:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
 29:   if (!eps->max_it) eps->max_it = PetscMax(300,(PetscInt)(2*eps->n/eps->ncv));
 30:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;

 32:   ncv = eps->ncv;
 33: #if defined(PETSC_USE_COMPLEX)
 34:   PetscFree(ar->rwork);
 35:   PetscMalloc1(ncv,&ar->rwork);
 36:   PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscReal));
 37:   PetscBLASIntCast(3*ncv*ncv+5*ncv,&ar->lworkl);
 38:   PetscFree(ar->workev);
 39:   PetscMalloc1(3*ncv,&ar->workev);
 40:   PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
 41: #else
 42:   if (eps->ishermitian) {
 43:     PetscBLASIntCast(ncv*(ncv+8),&ar->lworkl);
 44:   } else {
 45:     PetscBLASIntCast(3*ncv*ncv+6*ncv,&ar->lworkl);
 46:     PetscFree(ar->workev);
 47:     PetscMalloc1(3*ncv,&ar->workev);
 48:     PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
 49:   }
 50: #endif
 51:   PetscFree(ar->workl);
 52:   PetscMalloc1(ar->lworkl,&ar->workl);
 53:   PetscLogObjectMemory((PetscObject)eps,ar->lworkl*sizeof(PetscScalar));
 54:   PetscFree(ar->select);
 55:   PetscMalloc1(ncv,&ar->select);
 56:   PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscBool));
 57:   PetscFree(ar->workd);
 58:   PetscMalloc1(3*eps->nloc,&ar->workd);
 59:   PetscLogObjectMemory((PetscObject)eps,3*eps->nloc*sizeof(PetscScalar));

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

 63:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in the Arpack interface");
 64:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
 65:   if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");

 67:   EPSAllocateSolution(eps,0);
 68:   EPS_SetInnerProduct(eps);
 69:   EPSSetWorkVecs(eps,2);

 71:   PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
 72:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
 73:   RGIsTrivial(eps->rg,&istrivial);
 74:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
 75:   return(0);
 76: }

 78: PetscErrorCode EPSSolve_ARPACK(EPS eps)
 79: {
 81:   EPS_ARPACK     *ar = (EPS_ARPACK*)eps->data;
 82:   char           bmat[1],howmny[] = "A";
 83:   const char     *which;
 84:   PetscBLASInt   n,iparam[11],ipntr[14],ido,info,nev,ncv;
 85: #if !defined(PETSC_HAVE_MPIUNI)
 86:   PetscBLASInt   fcomm;
 87: #endif
 88:   PetscScalar    sigmar,*pV,*resid;
 89:   Vec            v0,x,y,w = eps->work[0];
 90:   Mat            A;
 91:   PetscBool      isSinv,isShift,rvec;
 92: #if !defined(PETSC_USE_COMPLEX)
 93:   PetscScalar    sigmai = 0.0;
 94: #endif

 97:   PetscBLASIntCast(eps->nev,&nev);
 98:   PetscBLASIntCast(eps->ncv,&ncv);
 99: #if !defined(PETSC_HAVE_MPIUNI)
100:   PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&fcomm);
101: #endif
102:   PetscBLASIntCast(eps->nloc,&n);
103:   EPSGetStartVector(eps,0,NULL);
104:   BVSetActiveColumns(eps->V,0,0);  /* just for deflation space */
105:   BVGetColumn(eps->V,0,&v0);
106:   VecCopy(v0,eps->work[1]);
107:   VecGetArray(v0,&pV);
108:   VecGetArray(eps->work[1],&resid);

110:   ido  = 0;            /* first call to reverse communication interface */
111:   info = 1;            /* indicates an initial vector is provided */
112:   iparam[0] = 1;       /* use exact shifts */
113:   PetscBLASIntCast(eps->max_it,&iparam[2]);  /* max Arnoldi iterations */
114:   iparam[3] = 1;       /* blocksize */
115:   iparam[4] = 0;       /* number of converged Ritz values */

117:   /*
118:      Computational modes ([]=not supported):
119:             symmetric    non-symmetric    complex
120:         1     1  'I'        1  'I'         1  'I'
121:         2     3  'I'        3  'I'         3  'I'
122:         3     2  'G'        2  'G'         2  'G'
123:         4     3  'G'        3  'G'         3  'G'
124:         5   [ 4  'G' ]    [ 3  'G' ]
125:         6   [ 5  'G' ]    [ 4  'G' ]
126:    */
127:   PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
128:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift);
129:   STGetShift(eps->st,&sigmar);
130:   STGetMatrix(eps->st,0,&A);
131:   MatCreateVecsEmpty(A,&x,&y);

133:   if (isSinv) {
134:     /* shift-and-invert mode */
135:     iparam[6] = 3;
136:     if (eps->ispositive) bmat[0] = 'G';
137:     else bmat[0] = 'I';
138:   } else if (isShift && eps->ispositive) {
139:     /* generalized shift mode with B positive definite */
140:     iparam[6] = 2;
141:     bmat[0] = 'G';
142:   } else {
143:     /* regular mode */
144:     if (eps->ishermitian && eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
145:     iparam[6] = 1;
146:     bmat[0] = 'I';
147:   }

149: #if !defined(PETSC_USE_COMPLEX)
150:     if (eps->ishermitian) {
151:       switch (eps->which) {
152:         case EPS_TARGET_MAGNITUDE:
153:         case EPS_LARGEST_MAGNITUDE:  which = "LM"; break;
154:         case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
155:         case EPS_TARGET_REAL:
156:         case EPS_LARGEST_REAL:       which = "LA"; break;
157:         case EPS_SMALLEST_REAL:      which = "SA"; break;
158:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
159:       }
160:     } else {
161: #endif
162:       switch (eps->which) {
163:         case EPS_TARGET_MAGNITUDE:
164:         case EPS_LARGEST_MAGNITUDE:  which = "LM"; break;
165:         case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
166:         case EPS_TARGET_REAL:
167:         case EPS_LARGEST_REAL:       which = "LR"; break;
168:         case EPS_SMALLEST_REAL:      which = "SR"; break;
169:         case EPS_TARGET_IMAGINARY:
170:         case EPS_LARGEST_IMAGINARY:  which = "LI"; break;
171:         case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
172:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
173:       }
174: #if !defined(PETSC_USE_COMPLEX)
175:     }
176: #endif

178:   do {

180: #if !defined(PETSC_USE_COMPLEX)
181:     if (eps->ishermitian) {
182:       PetscStackCall("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
183:     } else {
184:       PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
185:     }
186: #else
187:     PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,ar->rwork,&info));
188: #endif

190:     if (ido == -1 || ido == 1 || ido == 2) {
191:       if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
192:         /* special case for shift-and-invert with B semi-positive definite*/
193:         VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
194:       } else {
195:         VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
196:       }
197:       VecPlaceArray(y,&ar->workd[ipntr[1]-1]);

199:       if (ido == -1) {
200:         /* Y = OP * X for for the initialization phase to
201:            force the starting vector into the range of OP */
202:         STApply(eps->st,x,y);
203:       } else if (ido == 2) {
204:         /* Y = B * X */
205:         BVApplyMatrix(eps->V,x,y);
206:       } else { /* ido == 1 */
207:         if (iparam[6] == 3 && bmat[0] == 'G') {
208:           /* Y = OP * X for shift-and-invert with B semi-positive definite */
209:           STMatSolve(eps->st,x,y);
210:         } else if (iparam[6] == 2) {
211:           /* X=A*X Y=B^-1*X for shift with B positive definite */
212:           MatMult(A,x,y);
213:           if (sigmar != 0.0) {
214:             BVApplyMatrix(eps->V,x,w);
215:             VecAXPY(y,sigmar,w);
216:           }
217:           VecCopy(y,x);
218:           STMatSolve(eps->st,x,y);
219:         } else {
220:           /* Y = OP * X */
221:           STApply(eps->st,x,y);
222:         }
223:         BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
224:       }

226:       VecResetArray(x);
227:       VecResetArray(y);
228:     } else if (ido != 99) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in ARPACK reverse comunication interface (ido=%d)",ido);

230:   } while (ido != 99);

232:   eps->nconv = iparam[4];
233:   eps->its = iparam[2];

235:   if (info==3) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"No shift could be applied in xxAUPD.\nTry increasing the size of NCV relative to NEV");
236:   else if (info!=0 && info!=1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",(int)info);

238:   rvec = PETSC_TRUE;

240:   if (eps->nconv > 0) {
241: #if !defined(PETSC_USE_COMPLEX)
242:     if (eps->ishermitian) {
243:       EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
244:       PetscStackCall("ARPACKseupd",ARPACKseupd_(&fcomm,&rvec,howmny,ar->select,eps->eigr,pV,&n,&sigmar,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
245:     } else {
246:       EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
247:       PetscStackCall("ARPACKneupd",ARPACKneupd_(&fcomm,&rvec,howmny,ar->select,eps->eigr,eps->eigi,pV,&n,&sigmar,&sigmai,ar->workev,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
248:     }
249: #else
250:     EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
251:     PetscStackCall("ARPACKneupd",ARPACKneupd_(&fcomm,&rvec,howmny,ar->select,eps->eigr,pV,&n,&sigmar,ar->workev,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,ar->rwork,&info));
252: #endif
253:     if (info!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",(int)info);
254:   }

256:   VecRestoreArray(v0,&pV);
257:   BVRestoreColumn(eps->V,0,&v0);
258:   VecRestoreArray(eps->work[1],&resid);
259:   if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
260:   else eps->reason = EPS_DIVERGED_ITS;

262:   if (eps->ishermitian) {
263:     PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
264:   } else {
265:     PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
266:   }

268:   VecDestroy(&x);
269:   VecDestroy(&y);
270:   return(0);
271: }

273: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
274: {
276:   PetscBool      isSinv;

279:   PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
280:   if (!isSinv) {
281:     EPSBackTransform_Default(eps);
282:   }
283:   return(0);
284: }

286: PetscErrorCode EPSReset_ARPACK(EPS eps)
287: {
289:   EPS_ARPACK     *ar = (EPS_ARPACK*)eps->data;

292:   PetscFree(ar->workev);
293:   PetscFree(ar->workl);
294:   PetscFree(ar->select);
295:   PetscFree(ar->workd);
296: #if defined(PETSC_USE_COMPLEX)
297:   PetscFree(ar->rwork);
298: #endif
299:   return(0);
300: }

302: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
303: {

307:   PetscFree(eps->data);
308:   return(0);
309: }

311: PETSC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
312: {
313:   EPS_ARPACK     *ctx;

317:   PetscNewLog(eps,&ctx);
318:   eps->data = (void*)ctx;

320:   eps->ops->solve          = EPSSolve_ARPACK;
321:   eps->ops->setup          = EPSSetUp_ARPACK;
322:   eps->ops->destroy        = EPSDestroy_ARPACK;
323:   eps->ops->reset          = EPSReset_ARPACK;
324:   eps->ops->backtransform  = EPSBackTransform_ARPACK;
325:   return(0);
326: }