Actual source code: arpack.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    This file implements a wrapper to the ARPACK package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, 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>
 25: #include <../src/eps/impls/external/arpack/arpackp.h>

 27: PetscErrorCode EPSSolve_ARPACK(EPS);

 31: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
 32: {
 34:   PetscInt       ncv;
 35:   PetscBool      flg,istrivial;
 36:   EPS_ARPACK     *ar = (EPS_ARPACK*)eps->data;

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

 46:   ncv = eps->ncv;
 47: #if defined(PETSC_USE_COMPLEX)
 48:   PetscFree(ar->rwork);
 49:   PetscMalloc1(ncv,&ar->rwork);
 50:   PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscReal));
 51:   PetscBLASIntCast(3*ncv*ncv+5*ncv,&ar->lworkl);
 52:   PetscFree(ar->workev);
 53:   PetscMalloc1(3*ncv,&ar->workev);
 54:   PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
 55: #else
 56:   if (eps->ishermitian) {
 57:     PetscBLASIntCast(ncv*(ncv+8),&ar->lworkl);
 58:   } else {
 59:     PetscBLASIntCast(3*ncv*ncv+6*ncv,&ar->lworkl);
 60:     PetscFree(ar->workev);
 61:     PetscMalloc1(3*ncv,&ar->workev);
 62:     PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
 63:   }
 64: #endif
 65:   PetscFree(ar->workl);
 66:   PetscMalloc1(ar->lworkl,&ar->workl);
 67:   PetscLogObjectMemory((PetscObject)eps,ar->lworkl*sizeof(PetscScalar));
 68:   PetscFree(ar->select);
 69:   PetscMalloc1(ncv,&ar->select);
 70:   PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscBool));
 71:   PetscFree(ar->workd);
 72:   PetscMalloc1(3*eps->nloc,&ar->workd);
 73:   PetscLogObjectMemory((PetscObject)eps,3*eps->nloc*sizeof(PetscScalar));

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

 77:   if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in the Arpack interface");
 78:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

 80:   EPSAllocateSolution(eps,0);
 81:   EPSSetWorkVecs(eps,2);

 83:   PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
 84:   if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
 85:   RGIsTrivial(eps->rg,&istrivial);
 86:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");

 88:   /* dispatch solve method */
 89:   eps->ops->solve = EPSSolve_ARPACK;
 90:   return(0);
 91: }

 95: PetscErrorCode EPSSolve_ARPACK(EPS eps)
 96: {
 98:   EPS_ARPACK     *ar = (EPS_ARPACK*)eps->data;
 99:   char           bmat[1],howmny[] = "A";
100:   const char     *which;
101:   PetscBLASInt   n,iparam[11],ipntr[14],ido,info,nev,ncv;
102: #if !defined(PETSC_HAVE_MPIUNI)
103:   PetscBLASInt   fcomm;
104: #endif
105:   PetscScalar    sigmar,*pV,*resid;
106:   Vec            v0,x,y,w = eps->work[0];
107:   Mat            A;
108:   PetscBool      isSinv,isShift,rvec;
109: #if !defined(PETSC_USE_COMPLEX)
110:   PetscScalar    sigmai = 0.0;
111: #endif

114:   PetscBLASIntCast(eps->nev,&nev);
115:   PetscBLASIntCast(eps->ncv,&ncv);
116: #if !defined(PETSC_HAVE_MPIUNI)
117:   PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&fcomm);
118: #endif
119:   PetscBLASIntCast(eps->nloc,&n);
120:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
121:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
122:   EPSGetStartVector(eps,0,NULL);
123:   BVSetActiveColumns(eps->V,0,0);  /* just for deflation space */
124:   BVGetColumn(eps->V,0,&v0);
125:   VecCopy(v0,eps->work[1]);
126:   VecGetArray(v0,&pV);
127:   VecGetArray(eps->work[1],&resid);

129:   ido  = 0;            /* first call to reverse communication interface */
130:   info = 1;            /* indicates a initial vector is provided */
131:   iparam[0] = 1;       /* use exact shifts */
132:   PetscBLASIntCast(eps->max_it,&iparam[2]);  /* max Arnoldi iterations */
133:   iparam[3] = 1;       /* blocksize */
134:   iparam[4] = 0;       /* number of converged Ritz values */

136:   /*
137:      Computational modes ([]=not supported):
138:             symmetric    non-symmetric    complex
139:         1     1  'I'        1  'I'         1  'I'
140:         2     3  'I'        3  'I'         3  'I'
141:         3     2  'G'        2  'G'         2  'G'
142:         4     3  'G'        3  'G'         3  'G'
143:         5   [ 4  'G' ]    [ 3  'G' ]
144:         6   [ 5  'G' ]    [ 4  'G' ]
145:    */
146:   PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
147:   PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift);
148:   STGetShift(eps->st,&sigmar);
149:   STGetOperators(eps->st,0,&A);

151:   if (isSinv) {
152:     /* shift-and-invert mode */
153:     iparam[6] = 3;
154:     if (eps->ispositive) bmat[0] = 'G';
155:     else bmat[0] = 'I';
156:   } else if (isShift && eps->ispositive) {
157:     /* generalized shift mode with B positive definite */
158:     iparam[6] = 2;
159:     bmat[0] = 'G';
160:   } else {
161:     /* regular mode */
162:     if (eps->ishermitian && eps->isgeneralized)
163:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
164:     iparam[6] = 1;
165:     bmat[0] = 'I';
166:   }

168: #if !defined(PETSC_USE_COMPLEX)
169:     if (eps->ishermitian) {
170:       switch (eps->which) {
171:         case EPS_TARGET_MAGNITUDE:
172:         case EPS_LARGEST_MAGNITUDE:  which = "LM"; break;
173:         case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
174:         case EPS_TARGET_REAL:
175:         case EPS_LARGEST_REAL:       which = "LA"; break;
176:         case EPS_SMALLEST_REAL:      which = "SA"; break;
177:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
178:       }
179:     } else {
180: #endif
181:       switch (eps->which) {
182:         case EPS_TARGET_MAGNITUDE:
183:         case EPS_LARGEST_MAGNITUDE:  which = "LM"; break;
184:         case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
185:         case EPS_TARGET_REAL:
186:         case EPS_LARGEST_REAL:       which = "LR"; break;
187:         case EPS_SMALLEST_REAL:      which = "SR"; break;
188:         case EPS_TARGET_IMAGINARY:
189:         case EPS_LARGEST_IMAGINARY:  which = "LI"; break;
190:         case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
191:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
192:       }
193: #if !defined(PETSC_USE_COMPLEX)
194:     }
195: #endif

197:   do {

199: #if !defined(PETSC_USE_COMPLEX)
200:     if (eps->ishermitian) {
201:       PetscStackCall("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
202:     } else {
203:       PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
204:     }
205: #else
206:     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));
207: #endif

209:     if (ido == -1 || ido == 1 || ido == 2) {
210:       if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
211:         /* special case for shift-and-invert with B semi-positive definite*/
212:         VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
213:       } else {
214:         VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
215:       }
216:       VecPlaceArray(y,&ar->workd[ipntr[1]-1]);

218:       if (ido == -1) {
219:         /* Y = OP * X for for the initialization phase to
220:            force the starting vector into the range of OP */
221:         STApply(eps->st,x,y);
222:       } else if (ido == 2) {
223:         /* Y = B * X */
224:         BVApplyMatrix(eps->V,x,y);
225:       } else { /* ido == 1 */
226:         if (iparam[6] == 3 && bmat[0] == 'G') {
227:           /* Y = OP * X for shift-and-invert with B semi-positive definite */
228:           STMatSolve(eps->st,x,y);
229:         } else if (iparam[6] == 2) {
230:           /* X=A*X Y=B^-1*X for shift with B positive definite */
231:           MatMult(A,x,y);
232:           if (sigmar != 0.0) {
233:             BVApplyMatrix(eps->V,x,w);
234:             VecAXPY(y,sigmar,w);
235:           }
236:           VecCopy(y,x);
237:           STMatSolve(eps->st,x,y);
238:         } else {
239:           /* Y = OP * X */
240:           STApply(eps->st,x,y);
241:         }
242:         BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
243:       }

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

249:   } while (ido != 99);

251:   eps->nconv = iparam[4];
252:   eps->its = iparam[2];

254:   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");
255:   else if (info!=0 && info!=1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);

257:   rvec = PETSC_TRUE;

259:   if (eps->nconv > 0) {
260: #if !defined(PETSC_USE_COMPLEX)
261:     if (eps->ishermitian) {
262:       EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
263:       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));
264:     } else {
265:       EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
266:       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));
267:     }
268: #else
269:     EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
270:     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));
271: #endif
272:     if (info!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info);
273:   }

275:   VecRestoreArray(v0,&pV);
276:   BVRestoreColumn(eps->V,0,&v0);
277:   VecRestoreArray(eps->work[1],&resid);
278:   if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
279:   else eps->reason = EPS_DIVERGED_ITS;

281:   if (eps->ishermitian) {
282:     PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
283:   } else {
284:     PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
285:   }

287:   VecDestroy(&x);
288:   VecDestroy(&y);
289:   return(0);
290: }

294: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
295: {
297:   PetscBool      isSinv;

300:   PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
301:   if (!isSinv) {
302:     EPSBackTransform_Default(eps);
303:   }
304:   return(0);
305: }

309: PetscErrorCode EPSReset_ARPACK(EPS eps)
310: {
312:   EPS_ARPACK     *ar = (EPS_ARPACK*)eps->data;

315:   PetscFree(ar->workev);
316:   PetscFree(ar->workl);
317:   PetscFree(ar->select);
318:   PetscFree(ar->workd);
319: #if defined(PETSC_USE_COMPLEX)
320:   PetscFree(ar->rwork);
321: #endif
322:   return(0);
323: }

327: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
328: {

332:   PetscFree(eps->data);
333:   return(0);
334: }

338: PETSC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
339: {
340:   EPS_ARPACK     *ctx;

344:   PetscNewLog(eps,&ctx);
345:   eps->data = (void*)ctx;

347:   eps->ops->setup                = EPSSetUp_ARPACK;
348:   eps->ops->destroy              = EPSDestroy_ARPACK;
349:   eps->ops->reset                = EPSReset_ARPACK;
350:   eps->ops->backtransform        = EPSBackTransform_ARPACK;
351:   return(0);
352: }