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-2013, 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 <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(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: PetscMalloc(ncv*sizeof(PetscReal),&ar->rwork);
50: PetscLogObjectMemory(eps,ncv*sizeof(PetscReal));
51: PetscBLASIntCast(3*ncv*ncv+5*ncv,&ar->lworkl);
52: PetscFree(ar->workev);
53: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
54: PetscLogObjectMemory(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: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
62: PetscLogObjectMemory(eps,3*ncv*sizeof(PetscScalar));
63: }
64: #endif
65: PetscFree(ar->workl);
66: PetscMalloc(ar->lworkl*sizeof(PetscScalar),&ar->workl);
67: PetscLogObjectMemory(eps,ar->lworkl*sizeof(PetscScalar));
68: PetscFree(ar->select);
69: PetscMalloc(ncv*sizeof(PetscBool),&ar->select);
70: PetscLogObjectMemory(eps,ncv*sizeof(PetscBool));
71: PetscFree(ar->workd);
72: PetscMalloc(3*eps->nloc*sizeof(PetscScalar),&ar->workd);
73: PetscLogObjectMemory(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);
81: EPSSetWorkVecs(eps,2);
83: /* dispatch solve method */
84: if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
85: eps->ops->solve = EPSSolve_ARPACK;
86: return(0);
87: }
91: PetscErrorCode EPSSolve_ARPACK(EPS eps)
92: {
94: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
95: char bmat[1],howmny[] = "A";
96: const char *which;
97: PetscBLASInt n,iparam[11],ipntr[14],ido,info,nev,ncv;
98: #if !defined(PETSC_HAVE_MPIUNI)
99: PetscBLASInt fcomm;
100: #endif
101: PetscScalar sigmar,*pV,*resid;
102: Vec x,y,w = eps->work[0];
103: Mat A;
104: PetscBool isSinv,isShift,rvec;
105: #if !defined(PETSC_USE_COMPLEX)
106: PetscScalar sigmai = 0.0;
107: #endif
110: PetscBLASIntCast(eps->nev,&nev);
111: PetscBLASIntCast(eps->ncv,&ncv);
112: #if !defined(PETSC_HAVE_MPIUNI)
113: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&fcomm);
114: #endif
115: PetscBLASIntCast(eps->nloc,&n);
116: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
117: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
118: VecGetArray(eps->V[0],&pV);
119: EPSGetStartVector(eps,0,eps->work[1],NULL);
120: VecGetArray(eps->work[1],&resid);
122: ido = 0; /* first call to reverse communication interface */
123: info = 1; /* indicates a initial vector is provided */
124: iparam[0] = 1; /* use exact shifts */
125: PetscBLASIntCast(eps->max_it,&iparam[2]); /* max Arnoldi iterations */
126: iparam[3] = 1; /* blocksize */
127: iparam[4] = 0; /* number of converged Ritz values */
129: /*
130: Computational modes ([]=not supported):
131: symmetric non-symmetric complex
132: 1 1 'I' 1 'I' 1 'I'
133: 2 3 'I' 3 'I' 3 'I'
134: 3 2 'G' 2 'G' 2 'G'
135: 4 3 'G' 3 'G' 3 'G'
136: 5 [ 4 'G' ] [ 3 'G' ]
137: 6 [ 5 'G' ] [ 4 'G' ]
138: */
139: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
140: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift);
141: STGetShift(eps->st,&sigmar);
142: STGetOperators(eps->st,0,&A);
144: if (isSinv) {
145: /* shift-and-invert mode */
146: iparam[6] = 3;
147: if (eps->ispositive) bmat[0] = 'G';
148: else bmat[0] = 'I';
149: } else if (isShift && eps->ispositive) {
150: /* generalized shift mode with B positive definite */
151: iparam[6] = 2;
152: bmat[0] = 'G';
153: } else {
154: /* regular mode */
155: if (eps->ishermitian && eps->isgeneralized)
156: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
157: iparam[6] = 1;
158: bmat[0] = 'I';
159: }
161: #if !defined(PETSC_USE_COMPLEX)
162: if (eps->ishermitian) {
163: switch (eps->which) {
164: case EPS_TARGET_MAGNITUDE:
165: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
166: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
167: case EPS_TARGET_REAL:
168: case EPS_LARGEST_REAL: which = "LA"; break;
169: case EPS_SMALLEST_REAL: which = "SA"; break;
170: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
171: }
172: } else {
173: #endif
174: switch (eps->which) {
175: case EPS_TARGET_MAGNITUDE:
176: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
177: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
178: case EPS_TARGET_REAL:
179: case EPS_LARGEST_REAL: which = "LR"; break;
180: case EPS_SMALLEST_REAL: which = "SR"; break;
181: case EPS_TARGET_IMAGINARY:
182: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
183: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
184: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
185: }
186: #if !defined(PETSC_USE_COMPLEX)
187: }
188: #endif
190: do {
192: #if !defined(PETSC_USE_COMPLEX)
193: if (eps->ishermitian) {
194: PetscStackCall("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
195: } else {
196: PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
197: }
198: #else
199: 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));
200: #endif
202: if (ido == -1 || ido == 1 || ido == 2) {
203: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
204: /* special case for shift-and-invert with B semi-positive definite*/
205: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
206: } else {
207: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
208: }
209: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
211: if (ido == -1) {
212: /* Y = OP * X for for the initialization phase to
213: force the starting vector into the range of OP */
214: STApply(eps->st,x,y);
215: } else if (ido == 2) {
216: /* Y = B * X */
217: IPApplyMatrix(eps->ip,x,y);
218: } else { /* ido == 1 */
219: if (iparam[6] == 3 && bmat[0] == 'G') {
220: /* Y = OP * X for shift-and-invert with B semi-positive definite */
221: STMatSolve(eps->st,1,x,y);
222: } else if (iparam[6] == 2) {
223: /* X=A*X Y=B^-1*X for shift with B positive definite */
224: MatMult(A,x,y);
225: if (sigmar != 0.0) {
226: IPApplyMatrix(eps->ip,x,w);
227: VecAXPY(y,sigmar,w);
228: }
229: VecCopy(y,x);
230: STMatSolve(eps->st,1,x,y);
231: } else {
232: /* Y = OP * X */
233: STApply(eps->st,x,y);
234: }
235: IPOrthogonalize(eps->ip,0,NULL,eps->nds,NULL,eps->defl,y,NULL,NULL,NULL);
236: }
238: VecResetArray(x);
239: VecResetArray(y);
240: } else if (ido != 99) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in ARPACK reverse comunication interface (ido=%d)",ido);
242: } while (ido != 99);
244: eps->nconv = iparam[4];
245: eps->its = iparam[2];
247: 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");
248: else if (info!=0 && info!=1) SETERRQ1(PetscObjectComm((PetscObject)eps),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: 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));
257: } else {
258: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
259: 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));
260: }
261: #else
262: EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
263: 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));
264: #endif
265: if (info!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info);
266: }
268: VecRestoreArray(eps->V[0],&pV);
269: VecRestoreArray(eps->work[1],&resid);
270: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
271: else eps->reason = EPS_DIVERGED_ITS;
273: if (eps->ishermitian) {
274: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
275: } else {
276: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
277: }
279: VecDestroy(&x);
280: VecDestroy(&y);
281: return(0);
282: }
286: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
287: {
289: PetscBool isSinv;
292: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
293: if (!isSinv) {
294: EPSBackTransform_Default(eps);
295: }
296: return(0);
297: }
301: PetscErrorCode EPSReset_ARPACK(EPS eps)
302: {
304: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
307: PetscFree(ar->workev);
308: PetscFree(ar->workl);
309: PetscFree(ar->select);
310: PetscFree(ar->workd);
311: #if defined(PETSC_USE_COMPLEX)
312: PetscFree(ar->rwork);
313: #endif
314: EPSReset_Default(eps);
315: return(0);
316: }
320: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
321: {
325: PetscFree(eps->data);
326: return(0);
327: }
331: PETSC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
332: {
336: PetscNewLog(eps,EPS_ARPACK,&eps->data);
337: eps->ops->setup = EPSSetUp_ARPACK;
338: eps->ops->destroy = EPSDestroy_ARPACK;
339: eps->ops->reset = EPSReset_ARPACK;
340: eps->ops->backtransform = EPSBackTransform_ARPACK;
341: eps->ops->computevectors = EPSComputeVectors_Default;
342: return(0);
343: }