Actual source code: arpack.c

  1: /*
  2:        This file implements a wrapper to the ARPACK package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2011, Universitat 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 <private/epsimpl.h>        /*I "slepceps.h" I*/
 25: #include <private/stimpl.h>         /*I "slepcst.h" I*/
 26: #include <../src/eps/impls/external/arpack/arpackp.h>

 28: PetscErrorCode EPSSolve_ARPACK(EPS);

 32: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
 33: {
 35:   PetscInt       ncv;
 36:   EPS_ARPACK     *ar = (EPS_ARPACK *)eps->data;

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

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

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

 72:   if (eps->balance!=EPS_BALANCE_NONE)
 73:     SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Balancing not supported in the Arpack interface");

 75:   EPSAllocateSolution(eps);
 76:   EPSDefaultGetWork(eps,2);

 78:   /* dispatch solve method */
 79:   if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
 80:   eps->ops->solve = EPSSolve_ARPACK;
 81:   return(0);
 82: }

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

102:   nev = PetscBLASIntCast(eps->nev);
103:   ncv = PetscBLASIntCast(eps->ncv);
104:   fcomm = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm));
105:   n = PetscBLASIntCast(eps->nloc);
106:   VecCreateMPIWithArray(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,PETSC_NULL,&x);
107:   VecCreateMPIWithArray(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,PETSC_NULL,&y);
108:   VecGetArray(eps->V[0],&pV);
109:   EPSGetStartVector(eps,0,eps->work[1],PETSC_NULL);
110:   VecGetArray(eps->work[1],&resid);
111: 
112:   ido  = 0;            /* first call to reverse communication interface */
113:   info = 1;            /* indicates a initial vector is provided */
114:   iparam[0] = 1;       /* use exact shifts */
115:   iparam[2] = PetscBLASIntCast(eps->max_it);  /* maximum number of Arnoldi update iterations */
116:   iparam[3] = 1;       /* blocksize */
117:   iparam[4] = 0;       /* number of converged Ritz values */
118: 
119:   /*
120:      Computational modes ([]=not supported):
121:             symmetric    non-symmetric    complex
122:         1     1  'I'        1  'I'         1  'I'
123:         2     3  'I'        3  'I'         3  'I'
124:         3     2  'G'        2  'G'         2  'G'
125:         4     3  'G'        3  'G'         3  'G'
126:         5   [ 4  'G' ]    [ 3  'G' ]
127:         6   [ 5  'G' ]    [ 4  'G' ]
128:    */
129:   PetscTypeCompare((PetscObject)eps->OP,STSINVERT,&isSinv);
130:   PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
131:   STGetShift(eps->OP,&sigmar);
132:   STGetOperators(eps->OP,&A,PETSC_NULL);

134:   if (isSinv) {
135:     /* shift-and-invert mode */
136:     iparam[6] = 3;
137:     if (eps->ispositive) bmat[0] = 'G';
138:     else bmat[0] = 'I';
139:   } else if (isShift && eps->ispositive) {
140:     /* generalized shift mode with B positive definite */
141:     iparam[6] = 2;
142:     bmat[0] = 'G';
143:   } else {
144:     /* regular mode */
145:     if (eps->ishermitian && eps->isgeneralized)
146:       SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
147:     iparam[6] = 1;
148:     bmat[0] = 'I';
149:   }
150: 
151: #if !defined(PETSC_USE_COMPLEX)
152:     if (eps->ishermitian) {
153:       switch(eps->which) {
154:         case EPS_TARGET_MAGNITUDE:
155:         case EPS_LARGEST_MAGNITUDE:  which = "LM"; break;
156:         case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
157:         case EPS_TARGET_REAL:
158:         case EPS_LARGEST_REAL:       which = "LA"; break;
159:         case EPS_SMALLEST_REAL:      which = "SA"; break;
160:         default: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
161:       }
162:     } else {
163: #endif
164:       switch(eps->which) {
165:         case EPS_TARGET_MAGNITUDE:
166:         case EPS_LARGEST_MAGNITUDE:  which = "LM"; break;
167:         case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
168:         case EPS_TARGET_REAL:
169:         case EPS_LARGEST_REAL:       which = "LR"; break;
170:         case EPS_SMALLEST_REAL:      which = "SR"; break;
171:         case EPS_TARGET_IMAGINARY:
172:         case EPS_LARGEST_IMAGINARY:  which = "LI"; break;
173:         case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
174:         default: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
175:       }
176: #if !defined(PETSC_USE_COMPLEX)
177:     }
178: #endif

180:   do {

182: #if !defined(PETSC_USE_COMPLEX)
183:     if (eps->ishermitian) {
184:       ARsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,
185:                resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
186:                ar->workl,&ar->lworkl,&info,1,2);
187:     }
188:     else {
189:       ARnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,
190:                resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
191:                ar->workl,&ar->lworkl,&info,1,2);
192:     }
193: #else
194:     ARnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,
195:              resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
196:              ar->workl,&ar->lworkl,ar->rwork,&info,1,2);
197: #endif
198: 
199:     if (ido == -1 || ido == 1 || ido == 2) {
200:       if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
201:         /* special case for shift-and-invert with B semi-positive definite*/
202:         VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
203:       } else {
204:         VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
205:       }
206:       VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
207: 
208:       if (ido == -1) {
209:         /* Y = OP * X for for the initialization phase to 
210:            force the starting vector into the range of OP */
211:         STApply(eps->OP,x,y);
212:       } else if (ido == 2) {
213:         /* Y = B * X */
214:         IPApplyMatrix(eps->ip,x,y);
215:       } else { /* ido == 1 */
216:         if (iparam[6] == 3 && bmat[0] == 'G') {
217:           /* Y = OP * X for shift-and-invert with B semi-positive definite */
218:           STAssociatedKSPSolve(eps->OP,x,y);
219:         } else if (iparam[6] == 2) {
220:           /* X=A*X Y=B^-1*X for shift with B positive definite */
221:           MatMult(A,x,y);
222:           if (sigmar != 0.0) {
223:             IPApplyMatrix(eps->ip,x,w);
224:             VecAXPY(y,sigmar,w);
225:           }
226:           VecCopy(y,x);
227:           STAssociatedKSPSolve(eps->OP,x,y);
228:         } else  {
229:           /* Y = OP * X */
230:           STApply(eps->OP,x,y);
231:         }
232:         IPOrthogonalize(eps->ip,0,PETSC_NULL,eps->nds,PETSC_NULL,eps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
233:       }
234: 
235:       VecResetArray(x);
236:       VecResetArray(y);
237:     } else if (ido != 99) {
238:       SETERRQ1(((PetscObject)eps)->comm,1,"Internal error in ARPACK reverse comunication interface (ido=%d)\n",ido);
239:     }
240: 
241:   } while (ido != 99);

243:   eps->nconv = iparam[4];
244:   eps->its = iparam[2];
245: 
246:   if (info==3) { SETERRQ(((PetscObject)eps)->comm,1,"No shift could be applied in xxAUPD.\n"
247:                            "Try increasing the size of NCV relative to NEV."); }
248:   else if (info!=0 && info!=1) { SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);}

250:   rvec = PETSC_TRUE;

252:   if (eps->nconv > 0) {
253: #if !defined(PETSC_USE_COMPLEX)
254:     if (eps->ishermitian) {
255:       EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
256:       ARseupd_ (&fcomm,&rvec,howmny,ar->select,eps->eigr,
257:                 pV,&n,&sigmar,
258:                 bmat,&n,which,&nev,&eps->tol,
259:                 resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
260:                 ar->workl,&ar->lworkl,&info,1,1,2);
261:     }
262:     else {
263:       EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
264:       ARneupd_ (&fcomm,&rvec,howmny,ar->select,eps->eigr,eps->eigi,
265:                 pV,&n,&sigmar,&sigmai,ar->workev,
266:                 bmat,&n,which,&nev,&eps->tol,
267:                 resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
268:                 ar->workl,&ar->lworkl,&info,1,1,2);
269:     }
270: #else
271:     EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
272:     ARneupd_ (&fcomm,&rvec,howmny,ar->select,eps->eigr,
273:               pV,&n,&sigmar,ar->workev,
274:               bmat,&n,which,&nev,&eps->tol,
275:               resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
276:               ar->workl,&ar->lworkl,ar->rwork,&info,1,1,2);
277: #endif
278:     if (info!=0) { SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info); }
279:   }

281:   VecRestoreArray(eps->V[0],&pV);
282:   VecRestoreArray(eps->work[1],&resid);
283:   if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
284:   else eps->reason = EPS_DIVERGED_ITS;

286:   if (eps->ishermitian) {
287:     PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
288:   } else {
289:     PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
290:   }

292:   VecDestroy(&x);
293:   VecDestroy(&y);
294:   return(0);
295: }

299: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
300: {
302:   PetscBool      isSinv;

305:   PetscTypeCompare((PetscObject)eps->OP,STSINVERT,&isSinv);
306:   if (!isSinv) {
307:     EPSBackTransform_Default(eps);
308:   }
309:   return(0);
310: }

314: PetscErrorCode EPSReset_ARPACK(EPS eps)
315: {
317:   EPS_ARPACK     *ar = (EPS_ARPACK *)eps->data;

320:   PetscFree(ar->workev);
321:   PetscFree(ar->workl);
322:   PetscFree(ar->select);
323:   PetscFree(ar->workd);
324: #if defined(PETSC_USE_COMPLEX)
325:   PetscFree(ar->rwork);
326: #endif
327:   EPSDefaultFreeWork(eps);
328:   EPSFreeSolution(eps);
329:   return(0);
330: }

334: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
335: {

339:   PetscFree(eps->data);
340:   return(0);
341: }

343: EXTERN_C_BEGIN
346: PetscErrorCode EPSCreate_ARPACK(EPS eps)
347: {

351:   PetscNewLog(eps,EPS_ARPACK,&eps->data);
352:   eps->ops->setup                = EPSSetUp_ARPACK;
353:   eps->ops->destroy              = EPSDestroy_ARPACK;
354:   eps->ops->reset                = EPSReset_ARPACK;
355:   eps->ops->backtransform        = EPSBackTransform_ARPACK;
356:   eps->ops->computevectors       = EPSComputeVectors_Default;
357:   return(0);
358: }
359: EXTERN_C_END