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