Actual source code: arpack.c
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: This file implements a wrapper to the ARPACK package
12: */
14: #include <slepc/private/epsimpl.h>
15: #include "arpack.h"
17: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
18: {
20: PetscInt ncv;
21: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
24: EPSCheckDefinite(eps);
25: if (eps->ncv!=PETSC_DEFAULT) {
26: if (eps->ncv<eps->nev+2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The value of ncv must be at least nev+2");
27: } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n); /* set default value of ncv */
28: if (eps->mpd!=PETSC_DEFAULT) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
29: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(300,(PetscInt)(2*eps->n/eps->ncv));
30: if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
31: if (eps->which==EPS_ALL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
32: if (eps->which==EPS_WHICH_USER) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support user-defined ordering of eigenvalues");
33: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);
34: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
36: ncv = eps->ncv;
37: #if defined(PETSC_USE_COMPLEX)
38: PetscFree(ar->rwork);
39: PetscMalloc1(ncv,&ar->rwork);
40: PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscReal));
41: ar->lworkl = 3*ncv*ncv+5*ncv;
42: PetscFree(ar->workev);
43: PetscMalloc1(3*ncv,&ar->workev);
44: PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
45: #else
46: if (eps->ishermitian) {
47: ar->lworkl = ncv*(ncv+8);
48: } else {
49: ar->lworkl = 3*ncv*ncv+6*ncv;
50: PetscFree(ar->workev);
51: PetscMalloc1(3*ncv,&ar->workev);
52: PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
53: }
54: #endif
55: PetscFree(ar->workl);
56: PetscMalloc1(ar->lworkl,&ar->workl);
57: PetscLogObjectMemory((PetscObject)eps,ar->lworkl*sizeof(PetscScalar));
58: PetscFree(ar->select);
59: PetscMalloc1(ncv,&ar->select);
60: PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscInt));
61: PetscFree(ar->workd);
62: PetscMalloc1(3*eps->nloc,&ar->workd);
63: PetscLogObjectMemory((PetscObject)eps,3*eps->nloc*sizeof(PetscScalar));
65: EPSAllocateSolution(eps,0);
66: EPS_SetInnerProduct(eps);
67: EPSSetWorkVecs(eps,2);
68: return(0);
69: }
71: PetscErrorCode EPSSolve_ARPACK(EPS eps)
72: {
74: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
75: char bmat[1],howmny[] = "A";
76: const char *which;
77: PetscInt n,iparam[11],ipntr[14],ido,info,nev,ncv,rvec;
78: #if !defined(PETSC_HAVE_MPIUNI) && !defined(PETSC_HAVE_MSMPI)
79: MPI_Fint fcomm;
80: #endif
81: PetscScalar sigmar,*pV,*resid;
82: Vec x,y,w = eps->work[0];
83: Mat A;
84: PetscBool isSinv,isShift;
85: #if !defined(PETSC_USE_COMPLEX)
86: PetscScalar sigmai = 0.0;
87: #endif
90: nev = eps->nev;
91: ncv = eps->ncv;
92: #if !defined(PETSC_HAVE_MPIUNI) && !defined(PETSC_HAVE_MSMPI)
93: fcomm = MPI_Comm_c2f(PetscObjectComm((PetscObject)eps));
94: #endif
95: n = eps->nloc;
96: EPSGetStartVector(eps,0,NULL);
97: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
98: BVCopyVec(eps->V,0,eps->work[1]);
99: BVGetArray(eps->V,&pV);
100: VecGetArray(eps->work[1],&resid);
102: ido = 0; /* first call to reverse communication interface */
103: info = 1; /* indicates an initial vector is provided */
104: iparam[0] = 1; /* use exact shifts */
105: iparam[2] = eps->max_it; /* max Arnoldi iterations */
106: iparam[3] = 1; /* blocksize */
107: iparam[4] = 0; /* number of converged Ritz values */
109: /*
110: Computational modes ([]=not supported):
111: symmetric non-symmetric complex
112: 1 1 'I' 1 'I' 1 'I'
113: 2 3 'I' 3 'I' 3 'I'
114: 3 2 'G' 2 'G' 2 'G'
115: 4 3 'G' 3 'G' 3 'G'
116: 5 [ 4 'G' ] [ 3 'G' ]
117: 6 [ 5 'G' ] [ 4 'G' ]
118: */
119: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
120: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift);
121: STGetShift(eps->st,&sigmar);
122: STGetMatrix(eps->st,0,&A);
123: MatCreateVecsEmpty(A,&x,&y);
125: if (isSinv) {
126: /* shift-and-invert mode */
127: iparam[6] = 3;
128: if (eps->ispositive) bmat[0] = 'G';
129: else bmat[0] = 'I';
130: } else if (isShift && eps->ispositive) {
131: /* generalized shift mode with B positive definite */
132: iparam[6] = 2;
133: bmat[0] = 'G';
134: } else {
135: /* regular mode */
136: if (eps->ishermitian && eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
137: iparam[6] = 1;
138: bmat[0] = 'I';
139: }
141: #if !defined(PETSC_USE_COMPLEX)
142: if (eps->ishermitian) {
143: switch (eps->which) {
144: case EPS_TARGET_MAGNITUDE:
145: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
146: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
147: case EPS_TARGET_REAL:
148: case EPS_LARGEST_REAL: which = "LA"; break;
149: case EPS_SMALLEST_REAL: which = "SA"; break;
150: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
151: }
152: } else {
153: #endif
154: switch (eps->which) {
155: case EPS_TARGET_MAGNITUDE:
156: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
157: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
158: case EPS_TARGET_REAL:
159: case EPS_LARGEST_REAL: which = "LR"; break;
160: case EPS_SMALLEST_REAL: which = "SR"; break;
161: case EPS_TARGET_IMAGINARY:
162: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
163: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
164: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
165: }
166: #if !defined(PETSC_USE_COMPLEX)
167: }
168: #endif
170: do {
172: #if !defined(PETSC_USE_COMPLEX)
173: if (eps->ishermitian) {
174: PetscStackCall("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
175: } else {
176: PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
177: }
178: #else
179: 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));
180: #endif
182: if (ido == -1 || ido == 1 || ido == 2) {
183: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
184: /* special case for shift-and-invert with B semi-positive definite*/
185: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
186: } else {
187: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
188: }
189: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
191: if (ido == -1) {
192: /* Y = OP * X for for the initialization phase to
193: force the starting vector into the range of OP */
194: STApply(eps->st,x,y);
195: } else if (ido == 2) {
196: /* Y = B * X */
197: BVApplyMatrix(eps->V,x,y);
198: } else { /* ido == 1 */
199: if (iparam[6] == 3 && bmat[0] == 'G') {
200: /* Y = OP * X for shift-and-invert with B semi-positive definite */
201: STMatSolve(eps->st,x,y);
202: } else if (iparam[6] == 2) {
203: /* X=A*X Y=B^-1*X for shift with B positive definite */
204: MatMult(A,x,y);
205: if (sigmar != 0.0) {
206: BVApplyMatrix(eps->V,x,w);
207: VecAXPY(y,sigmar,w);
208: }
209: VecCopy(y,x);
210: STMatSolve(eps->st,x,y);
211: } else {
212: /* Y = OP * X */
213: STApply(eps->st,x,y);
214: }
215: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
216: }
218: VecResetArray(x);
219: VecResetArray(y);
220: } else if (ido != 99) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in ARPACK reverse comunication interface (ido=%d)",ido);
222: } while (ido != 99);
224: eps->nconv = iparam[4];
225: eps->its = iparam[2];
227: 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");
228: else if (info!=0 && info!=1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",(int)info);
230: rvec = PETSC_TRUE;
232: if (eps->nconv > 0) {
233: #if !defined(PETSC_USE_COMPLEX)
234: if (eps->ishermitian) {
235: 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));
236: } else {
237: 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));
238: }
239: #else
240: 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));
241: #endif
242: if (info!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",(int)info);
243: }
245: BVRestoreArray(eps->V,&pV);
246: VecRestoreArray(eps->work[1],&resid);
247: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
248: else eps->reason = EPS_DIVERGED_ITS;
250: VecDestroy(&x);
251: VecDestroy(&y);
252: return(0);
253: }
255: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
256: {
258: PetscBool isSinv;
261: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
262: if (!isSinv) {
263: EPSBackTransform_Default(eps);
264: }
265: return(0);
266: }
268: PetscErrorCode EPSReset_ARPACK(EPS eps)
269: {
271: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
274: PetscFree(ar->workev);
275: PetscFree(ar->workl);
276: PetscFree(ar->select);
277: PetscFree(ar->workd);
278: #if defined(PETSC_USE_COMPLEX)
279: PetscFree(ar->rwork);
280: #endif
281: return(0);
282: }
284: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
285: {
289: PetscFree(eps->data);
290: return(0);
291: }
293: SLEPC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
294: {
295: EPS_ARPACK *ctx;
299: PetscNewLog(eps,&ctx);
300: eps->data = (void*)ctx;
302: eps->ops->solve = EPSSolve_ARPACK;
303: eps->ops->setup = EPSSetUp_ARPACK;
304: eps->ops->setupsort = EPSSetUpSort_Basic;
305: eps->ops->destroy = EPSDestroy_ARPACK;
306: eps->ops->reset = EPSReset_ARPACK;
307: eps->ops->backtransform = EPSBackTransform_ARPACK;
308: return(0);
309: }