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-2012, 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 <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
25: #include <slepc-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,PETSC_ERR_ARG_OUTOFRANGE,"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) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Balancing not supported in the Arpack interface");
73: if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
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;
93: #if !defined(PETSC_HAVE_MPIUNI)
94: PetscBLASInt fcomm;
95: #endif
96: PetscScalar sigmar,*pV,*resid;
97: Vec x,y,w = eps->work[0];
98: Mat A;
99: PetscBool isSinv,isShift,rvec;
100: #if !defined(PETSC_USE_COMPLEX)
101: PetscScalar sigmai = 0.0;
102: #endif
105: nev = PetscBLASIntCast(eps->nev);
106: ncv = PetscBLASIntCast(eps->ncv);
107: #if !defined(PETSC_HAVE_MPIUNI)
108: fcomm = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm));
109: #endif
110: n = PetscBLASIntCast(eps->nloc);
111: VecCreateMPIWithArray(((PetscObject)eps)->comm,1,eps->nloc,PETSC_DECIDE,PETSC_NULL,&x);
112: VecCreateMPIWithArray(((PetscObject)eps)->comm,1,eps->nloc,PETSC_DECIDE,PETSC_NULL,&y);
113: VecGetArray(eps->V[0],&pV);
114: EPSGetStartVector(eps,0,eps->work[1],PETSC_NULL);
115: VecGetArray(eps->work[1],&resid);
116:
117: ido = 0; /* first call to reverse communication interface */
118: info = 1; /* indicates a initial vector is provided */
119: iparam[0] = 1; /* use exact shifts */
120: iparam[2] = PetscBLASIntCast(eps->max_it); /* maximum number of Arnoldi update iterations */
121: iparam[3] = 1; /* blocksize */
122: iparam[4] = 0; /* number of converged Ritz values */
123:
124: /*
125: Computational modes ([]=not supported):
126: symmetric non-symmetric complex
127: 1 1 'I' 1 'I' 1 'I'
128: 2 3 'I' 3 'I' 3 'I'
129: 3 2 'G' 2 'G' 2 'G'
130: 4 3 'G' 3 'G' 3 'G'
131: 5 [ 4 'G' ] [ 3 'G' ]
132: 6 [ 5 'G' ] [ 4 'G' ]
133: */
134: PetscObjectTypeCompare((PetscObject)eps->OP,STSINVERT,&isSinv);
135: PetscObjectTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
136: STGetShift(eps->OP,&sigmar);
137: STGetOperators(eps->OP,&A,PETSC_NULL);
139: if (isSinv) {
140: /* shift-and-invert mode */
141: iparam[6] = 3;
142: if (eps->ispositive) bmat[0] = 'G';
143: else bmat[0] = 'I';
144: } else if (isShift && eps->ispositive) {
145: /* generalized shift mode with B positive definite */
146: iparam[6] = 2;
147: bmat[0] = 'G';
148: } else {
149: /* regular mode */
150: if (eps->ishermitian && eps->isgeneralized)
151: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
152: iparam[6] = 1;
153: bmat[0] = 'I';
154: }
155:
156: #if !defined(PETSC_USE_COMPLEX)
157: if (eps->ishermitian) {
158: switch(eps->which) {
159: case EPS_TARGET_MAGNITUDE:
160: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
161: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
162: case EPS_TARGET_REAL:
163: case EPS_LARGEST_REAL: which = "LA"; break;
164: case EPS_SMALLEST_REAL: which = "SA"; break;
165: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
166: }
167: } else {
168: #endif
169: switch(eps->which) {
170: case EPS_TARGET_MAGNITUDE:
171: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
172: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
173: case EPS_TARGET_REAL:
174: case EPS_LARGEST_REAL: which = "LR"; break;
175: case EPS_SMALLEST_REAL: which = "SR"; break;
176: case EPS_TARGET_IMAGINARY:
177: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
178: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
179: default: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
180: }
181: #if !defined(PETSC_USE_COMPLEX)
182: }
183: #endif
185: do {
187: #if !defined(PETSC_USE_COMPLEX)
188: if (eps->ishermitian) {
189: ARsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,
190: resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
191: ar->workl,&ar->lworkl,&info);
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,&info);
197: }
198: #else
199: ARnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,
200: resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
201: ar->workl,&ar->lworkl,ar->rwork,&info);
202: #endif
203:
204: if (ido == -1 || ido == 1 || ido == 2) {
205: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
206: /* special case for shift-and-invert with B semi-positive definite*/
207: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
208: } else {
209: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
210: }
211: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
212:
213: if (ido == -1) {
214: /* Y = OP * X for for the initialization phase to
215: force the starting vector into the range of OP */
216: STApply(eps->OP,x,y);
217: } else if (ido == 2) {
218: /* Y = B * X */
219: IPApplyMatrix(eps->ip,x,y);
220: } else { /* ido == 1 */
221: if (iparam[6] == 3 && bmat[0] == 'G') {
222: /* Y = OP * X for shift-and-invert with B semi-positive definite */
223: STAssociatedKSPSolve(eps->OP,x,y);
224: } else if (iparam[6] == 2) {
225: /* X=A*X Y=B^-1*X for shift with B positive definite */
226: MatMult(A,x,y);
227: if (sigmar != 0.0) {
228: IPApplyMatrix(eps->ip,x,w);
229: VecAXPY(y,sigmar,w);
230: }
231: VecCopy(y,x);
232: STAssociatedKSPSolve(eps->OP,x,y);
233: } else {
234: /* Y = OP * X */
235: STApply(eps->OP,x,y);
236: }
237: IPOrthogonalize(eps->ip,0,PETSC_NULL,eps->nds,PETSC_NULL,eps->defl,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
238: }
239:
240: VecResetArray(x);
241: VecResetArray(y);
242: } else if (ido != 99) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Internal error in ARPACK reverse comunication interface (ido=%d)",ido);
243:
244: } while (ido != 99);
246: eps->nconv = iparam[4];
247: eps->its = iparam[2];
248:
249: if (info==3) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_LIB,"No shift could be applied in xxAUPD.\nTry increasing the size of NCV relative to NEV");
250: else if (info!=0 && info!=1) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);
252: rvec = PETSC_TRUE;
254: if (eps->nconv > 0) {
255: #if !defined(PETSC_USE_COMPLEX)
256: if (eps->ishermitian) {
257: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
258: ARseupd_ (&fcomm,&rvec,howmny,ar->select,eps->eigr,
259: pV,&n,&sigmar,
260: bmat,&n,which,&nev,&eps->tol,
261: resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
262: ar->workl,&ar->lworkl,&info);
263: }
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: ARneupd_ (&fcomm,&rvec,howmny,ar->select,eps->eigr,eps->eigi,
267: pV,&n,&sigmar,&sigmai,ar->workev,
268: bmat,&n,which,&nev,&eps->tol,
269: resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
270: ar->workl,&ar->lworkl,&info);
271: }
272: #else
273: EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
274: ARneupd_ (&fcomm,&rvec,howmny,ar->select,eps->eigr,
275: pV,&n,&sigmar,ar->workev,
276: bmat,&n,which,&nev,&eps->tol,
277: resid,&ncv,pV,&n,iparam,ipntr,ar->workd,
278: ar->workl,&ar->lworkl,ar->rwork,&info);
279: #endif
280: if (info!=0) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info);
281: }
283: VecRestoreArray(eps->V[0],&pV);
284: VecRestoreArray(eps->work[1],&resid);
285: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
286: else eps->reason = EPS_DIVERGED_ITS;
288: if (eps->ishermitian) {
289: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
290: } else {
291: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
292: }
294: VecDestroy(&x);
295: VecDestroy(&y);
296: return(0);
297: }
301: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
302: {
304: PetscBool isSinv;
307: PetscObjectTypeCompare((PetscObject)eps->OP,STSINVERT,&isSinv);
308: if (!isSinv) {
309: EPSBackTransform_Default(eps);
310: }
311: return(0);
312: }
316: PetscErrorCode EPSReset_ARPACK(EPS eps)
317: {
319: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
322: PetscFree(ar->workev);
323: PetscFree(ar->workl);
324: PetscFree(ar->select);
325: PetscFree(ar->workd);
326: #if defined(PETSC_USE_COMPLEX)
327: PetscFree(ar->rwork);
328: #endif
329: EPSDefaultFreeWork(eps);
330: EPSFreeSolution(eps);
331: return(0);
332: }
336: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
337: {
341: PetscFree(eps->data);
342: return(0);
343: }
345: EXTERN_C_BEGIN
348: PetscErrorCode EPSCreate_ARPACK(EPS eps)
349: {
353: PetscNewLog(eps,EPS_ARPACK,&eps->data);
354: eps->ops->setup = EPSSetUp_ARPACK;
355: eps->ops->destroy = EPSDestroy_ARPACK;
356: eps->ops->reset = EPSReset_ARPACK;
357: eps->ops->backtransform = EPSBackTransform_ARPACK;
358: eps->ops->computevectors = EPSComputeVectors_Default;
359: return(0);
360: }
361: EXTERN_C_END