Actual source code: arpack.c
slepc-3.6.3 2016-03-29
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: EPS_SetInnerProduct(eps);
82: EPSSetWorkVecs(eps,2);
84: PetscObjectTypeCompare((PetscObject)eps->V,BVVECS,&flg);
85: if (flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires a BV with contiguous storage");
86: RGIsTrivial(eps->rg,&istrivial);
87: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
89: /* dispatch solve method */
90: eps->ops->solve = EPSSolve_ARPACK;
91: return(0);
92: }
96: PetscErrorCode EPSSolve_ARPACK(EPS eps)
97: {
99: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
100: char bmat[1],howmny[] = "A";
101: const char *which;
102: PetscBLASInt n,iparam[11],ipntr[14],ido,info,nev,ncv;
103: #if !defined(PETSC_HAVE_MPIUNI)
104: PetscBLASInt fcomm;
105: #endif
106: PetscScalar sigmar,*pV,*resid;
107: Vec v0,x,y,w = eps->work[0];
108: Mat A;
109: PetscBool isSinv,isShift,rvec;
110: #if !defined(PETSC_USE_COMPLEX)
111: PetscScalar sigmai = 0.0;
112: #endif
115: PetscBLASIntCast(eps->nev,&nev);
116: PetscBLASIntCast(eps->ncv,&ncv);
117: #if !defined(PETSC_HAVE_MPIUNI)
118: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&fcomm);
119: #endif
120: PetscBLASIntCast(eps->nloc,&n);
121: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&x);
122: VecCreateMPIWithArray(PetscObjectComm((PetscObject)eps),1,eps->nloc,PETSC_DECIDE,NULL,&y);
123: EPSGetStartVector(eps,0,NULL);
124: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
125: BVGetColumn(eps->V,0,&v0);
126: VecCopy(v0,eps->work[1]);
127: VecGetArray(v0,&pV);
128: VecGetArray(eps->work[1],&resid);
130: ido = 0; /* first call to reverse communication interface */
131: info = 1; /* indicates a initial vector is provided */
132: iparam[0] = 1; /* use exact shifts */
133: PetscBLASIntCast(eps->max_it,&iparam[2]); /* max Arnoldi iterations */
134: iparam[3] = 1; /* blocksize */
135: iparam[4] = 0; /* number of converged Ritz values */
137: /*
138: Computational modes ([]=not supported):
139: symmetric non-symmetric complex
140: 1 1 'I' 1 'I' 1 'I'
141: 2 3 'I' 3 'I' 3 'I'
142: 3 2 'G' 2 'G' 2 'G'
143: 4 3 'G' 3 'G' 3 'G'
144: 5 [ 4 'G' ] [ 3 'G' ]
145: 6 [ 5 'G' ] [ 4 'G' ]
146: */
147: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
148: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift);
149: STGetShift(eps->st,&sigmar);
150: STGetOperators(eps->st,0,&A);
152: if (isSinv) {
153: /* shift-and-invert mode */
154: iparam[6] = 3;
155: if (eps->ispositive) bmat[0] = 'G';
156: else bmat[0] = 'I';
157: } else if (isShift && eps->ispositive) {
158: /* generalized shift mode with B positive definite */
159: iparam[6] = 2;
160: bmat[0] = 'G';
161: } else {
162: /* regular mode */
163: if (eps->ishermitian && eps->isgeneralized)
164: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
165: iparam[6] = 1;
166: bmat[0] = 'I';
167: }
169: #if !defined(PETSC_USE_COMPLEX)
170: if (eps->ishermitian) {
171: switch (eps->which) {
172: case EPS_TARGET_MAGNITUDE:
173: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
174: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
175: case EPS_TARGET_REAL:
176: case EPS_LARGEST_REAL: which = "LA"; break;
177: case EPS_SMALLEST_REAL: which = "SA"; break;
178: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
179: }
180: } else {
181: #endif
182: switch (eps->which) {
183: case EPS_TARGET_MAGNITUDE:
184: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
185: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
186: case EPS_TARGET_REAL:
187: case EPS_LARGEST_REAL: which = "LR"; break;
188: case EPS_SMALLEST_REAL: which = "SR"; break;
189: case EPS_TARGET_IMAGINARY:
190: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
191: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
192: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
193: }
194: #if !defined(PETSC_USE_COMPLEX)
195: }
196: #endif
198: do {
200: #if !defined(PETSC_USE_COMPLEX)
201: if (eps->ishermitian) {
202: PetscStackCall("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
203: } else {
204: PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
205: }
206: #else
207: 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));
208: #endif
210: if (ido == -1 || ido == 1 || ido == 2) {
211: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
212: /* special case for shift-and-invert with B semi-positive definite*/
213: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
214: } else {
215: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
216: }
217: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
219: if (ido == -1) {
220: /* Y = OP * X for for the initialization phase to
221: force the starting vector into the range of OP */
222: STApply(eps->st,x,y);
223: } else if (ido == 2) {
224: /* Y = B * X */
225: BVApplyMatrix(eps->V,x,y);
226: } else { /* ido == 1 */
227: if (iparam[6] == 3 && bmat[0] == 'G') {
228: /* Y = OP * X for shift-and-invert with B semi-positive definite */
229: STMatSolve(eps->st,x,y);
230: } else if (iparam[6] == 2) {
231: /* X=A*X Y=B^-1*X for shift with B positive definite */
232: MatMult(A,x,y);
233: if (sigmar != 0.0) {
234: BVApplyMatrix(eps->V,x,w);
235: VecAXPY(y,sigmar,w);
236: }
237: VecCopy(y,x);
238: STMatSolve(eps->st,x,y);
239: } else {
240: /* Y = OP * X */
241: STApply(eps->st,x,y);
242: }
243: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
244: }
246: VecResetArray(x);
247: VecResetArray(y);
248: } else if (ido != 99) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in ARPACK reverse comunication interface (ido=%d)",ido);
250: } while (ido != 99);
252: eps->nconv = iparam[4];
253: eps->its = iparam[2];
255: 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");
256: else if (info!=0 && info!=1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);
258: rvec = PETSC_TRUE;
260: if (eps->nconv > 0) {
261: #if !defined(PETSC_USE_COMPLEX)
262: if (eps->ishermitian) {
263: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
264: 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));
265: } else {
266: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
267: 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));
268: }
269: #else
270: EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
271: 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));
272: #endif
273: if (info!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info);
274: }
276: VecRestoreArray(v0,&pV);
277: BVRestoreColumn(eps->V,0,&v0);
278: VecRestoreArray(eps->work[1],&resid);
279: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
280: else eps->reason = EPS_DIVERGED_ITS;
282: if (eps->ishermitian) {
283: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
284: } else {
285: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
286: }
288: VecDestroy(&x);
289: VecDestroy(&y);
290: return(0);
291: }
295: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
296: {
298: PetscBool isSinv;
301: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
302: if (!isSinv) {
303: EPSBackTransform_Default(eps);
304: }
305: return(0);
306: }
310: PetscErrorCode EPSReset_ARPACK(EPS eps)
311: {
313: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
316: PetscFree(ar->workev);
317: PetscFree(ar->workl);
318: PetscFree(ar->select);
319: PetscFree(ar->workd);
320: #if defined(PETSC_USE_COMPLEX)
321: PetscFree(ar->rwork);
322: #endif
323: return(0);
324: }
328: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
329: {
333: PetscFree(eps->data);
334: return(0);
335: }
339: PETSC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
340: {
341: EPS_ARPACK *ctx;
345: PetscNewLog(eps,&ctx);
346: eps->data = (void*)ctx;
348: eps->ops->setup = EPSSetUp_ARPACK;
349: eps->ops->destroy = EPSDestroy_ARPACK;
350: eps->ops->reset = EPSReset_ARPACK;
351: eps->ops->backtransform = EPSBackTransform_ARPACK;
352: return(0);
353: }