Actual source code: arpack.c
slepc-3.19.2 2023-09-05
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: {
19: PetscInt ncv;
20: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
22: PetscFunctionBegin;
23: EPSCheckDefinite(eps);
24: if (eps->ncv!=PETSC_DEFAULT) {
25: PetscCheck(eps->ncv>=eps->nev+2,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The value of ncv must be at least nev+2");
26: } else eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n); /* set default value of ncv */
27: if (eps->mpd!=PETSC_DEFAULT) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
28: if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(300,(PetscInt)(2*eps->n/eps->ncv));
29: if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
30: PetscCheck(eps->which!=EPS_ALL,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
31: PetscCheck(eps->which!=EPS_WHICH_USER,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support user-defined ordering of eigenvalues");
32: EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE | EPS_FEATURE_STOPPING | EPS_FEATURE_TWOSIDED);
33: EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
35: ncv = eps->ncv;
36: #if defined(PETSC_USE_COMPLEX)
37: PetscCall(PetscFree(ar->rwork));
38: PetscCall(PetscMalloc1(ncv,&ar->rwork));
39: ar->lworkl = 3*ncv*ncv+5*ncv;
40: PetscCall(PetscFree(ar->workev));
41: PetscCall(PetscMalloc1(3*ncv,&ar->workev));
42: #else
43: if (eps->ishermitian) {
44: ar->lworkl = ncv*(ncv+8);
45: } else {
46: ar->lworkl = 3*ncv*ncv+6*ncv;
47: PetscCall(PetscFree(ar->workev));
48: PetscCall(PetscMalloc1(3*ncv,&ar->workev));
49: }
50: #endif
51: PetscCall(PetscFree(ar->workl));
52: PetscCall(PetscMalloc1(ar->lworkl,&ar->workl));
53: PetscCall(PetscFree(ar->select));
54: PetscCall(PetscMalloc1(ncv,&ar->select));
55: PetscCall(PetscFree(ar->workd));
56: PetscCall(PetscMalloc1(3*eps->nloc,&ar->workd));
58: PetscCall(EPSAllocateSolution(eps,0));
59: PetscCall(EPS_SetInnerProduct(eps));
60: PetscCall(EPSSetWorkVecs(eps,2));
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: PetscErrorCode EPSSolve_ARPACK(EPS eps)
65: {
66: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
67: char bmat[1],howmny[] = "A";
68: const char *which;
69: PetscInt n,iparam[11],ipntr[14],ido,info,nev,ncv,rvec;
70: #if !defined(PETSC_HAVE_MPIUNI) && !defined(PETSC_HAVE_MSMPI)
71: MPI_Fint fcomm;
72: #endif
73: PetscScalar sigmar,*pV,*resid;
74: Vec x,y,w = eps->work[0];
75: Mat A;
76: PetscBool isSinv,isShift;
77: #if !defined(PETSC_USE_COMPLEX)
78: PetscScalar sigmai = 0.0;
79: #endif
81: PetscFunctionBegin;
82: nev = eps->nev;
83: ncv = eps->ncv;
84: #if !defined(PETSC_HAVE_MPIUNI) && !defined(PETSC_HAVE_MSMPI)
85: fcomm = MPI_Comm_c2f(PetscObjectComm((PetscObject)eps));
86: #endif
87: n = eps->nloc;
88: PetscCall(EPSGetStartVector(eps,0,NULL));
89: PetscCall(BVSetActiveColumns(eps->V,0,0)); /* just for deflation space */
90: PetscCall(BVCopyVec(eps->V,0,eps->work[1]));
91: PetscCall(BVGetArray(eps->V,&pV));
92: PetscCall(VecGetArray(eps->work[1],&resid));
94: ido = 0; /* first call to reverse communication interface */
95: info = 1; /* indicates an initial vector is provided */
96: iparam[0] = 1; /* use exact shifts */
97: iparam[2] = eps->max_it; /* max Arnoldi iterations */
98: iparam[3] = 1; /* blocksize */
99: iparam[4] = 0; /* number of converged Ritz values */
101: /*
102: Computational modes ([]=not supported):
103: symmetric non-symmetric complex
104: 1 1 'I' 1 'I' 1 'I'
105: 2 3 'I' 3 'I' 3 'I'
106: 3 2 'G' 2 'G' 2 'G'
107: 4 3 'G' 3 'G' 3 'G'
108: 5 [ 4 'G' ] [ 3 'G' ]
109: 6 [ 5 'G' ] [ 4 'G' ]
110: */
111: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv));
112: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift));
113: PetscCall(STGetShift(eps->st,&sigmar));
114: PetscCall(STGetMatrix(eps->st,0,&A));
115: PetscCall(MatCreateVecsEmpty(A,&x,&y));
117: if (isSinv) {
118: /* shift-and-invert mode */
119: iparam[6] = 3;
120: if (eps->ispositive) bmat[0] = 'G';
121: else bmat[0] = 'I';
122: } else if (isShift && eps->ispositive) {
123: /* generalized shift mode with B positive definite */
124: iparam[6] = 2;
125: bmat[0] = 'G';
126: } else {
127: /* regular mode */
128: PetscCheck(!eps->ishermitian || !eps->isgeneralized,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
129: iparam[6] = 1;
130: bmat[0] = 'I';
131: }
133: #if !defined(PETSC_USE_COMPLEX)
134: if (eps->ishermitian) {
135: switch (eps->which) {
136: case EPS_TARGET_MAGNITUDE:
137: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
138: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
139: case EPS_TARGET_REAL:
140: case EPS_LARGEST_REAL: which = "LA"; break;
141: case EPS_SMALLEST_REAL: which = "SA"; break;
142: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
143: }
144: } else {
145: #endif
146: switch (eps->which) {
147: case EPS_TARGET_MAGNITUDE:
148: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
149: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
150: case EPS_TARGET_REAL:
151: case EPS_LARGEST_REAL: which = "LR"; break;
152: case EPS_SMALLEST_REAL: which = "SR"; break;
153: case EPS_TARGET_IMAGINARY:
154: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
155: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
156: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
157: }
158: #if !defined(PETSC_USE_COMPLEX)
159: }
160: #endif
162: do {
164: #if !defined(PETSC_USE_COMPLEX)
165: if (eps->ishermitian) {
166: PetscStackCallExternalVoid("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
167: } else {
168: PetscStackCallExternalVoid("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
169: }
170: #else
171: PetscStackCallExternalVoid("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));
172: #endif
174: if (ido == -1 || ido == 1 || ido == 2) {
175: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') PetscCall(VecPlaceArray(x,&ar->workd[ipntr[2]-1])); /* special case for shift-and-invert with B semi-positive definite*/
176: else PetscCall(VecPlaceArray(x,&ar->workd[ipntr[0]-1]));
177: PetscCall(VecPlaceArray(y,&ar->workd[ipntr[1]-1]));
179: if (ido == -1) {
180: /* Y = OP * X for for the initialization phase to
181: force the starting vector into the range of OP */
182: PetscCall(STApply(eps->st,x,y));
183: } else if (ido == 2) {
184: /* Y = B * X */
185: PetscCall(BVApplyMatrix(eps->V,x,y));
186: } else { /* ido == 1 */
187: if (iparam[6] == 3 && bmat[0] == 'G') {
188: /* Y = OP * X for shift-and-invert with B semi-positive definite */
189: PetscCall(STMatSolve(eps->st,x,y));
190: } else if (iparam[6] == 2) {
191: /* X=A*X Y=B^-1*X for shift with B positive definite */
192: PetscCall(MatMult(A,x,y));
193: if (sigmar != 0.0) {
194: PetscCall(BVApplyMatrix(eps->V,x,w));
195: PetscCall(VecAXPY(y,sigmar,w));
196: }
197: PetscCall(VecCopy(y,x));
198: PetscCall(STMatSolve(eps->st,x,y));
199: } else {
200: /* Y = OP * X */
201: PetscCall(STApply(eps->st,x,y));
202: }
203: PetscCall(BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL));
204: }
206: PetscCall(VecResetArray(x));
207: PetscCall(VecResetArray(y));
208: } else PetscCheck(ido==99,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in ARPACK reverse communication interface (ido=%" PetscInt_FMT ")",ido);
210: } while (ido != 99);
212: eps->nconv = iparam[4];
213: eps->its = iparam[2];
215: PetscCheck(info!=3,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"No shift could be applied in xxAUPD.\nTry increasing the size of NCV relative to NEV");
216: PetscCheck(info==0 || info==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%" PetscInt_FMT ")",info);
218: rvec = PETSC_TRUE;
220: if (eps->nconv > 0) {
221: #if !defined(PETSC_USE_COMPLEX)
222: if (eps->ishermitian) {
223: PetscStackCallExternalVoid("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));
224: } else {
225: PetscStackCallExternalVoid("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));
226: }
227: #else
228: PetscStackCallExternalVoid("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));
229: #endif
230: PetscCheck(info==0,PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%" PetscInt_FMT ")",info);
231: }
233: PetscCall(BVRestoreArray(eps->V,&pV));
234: PetscCall(VecRestoreArray(eps->work[1],&resid));
235: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
236: else eps->reason = EPS_DIVERGED_ITS;
238: PetscCall(VecDestroy(&x));
239: PetscCall(VecDestroy(&y));
240: PetscFunctionReturn(PETSC_SUCCESS);
241: }
243: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
244: {
245: PetscBool isSinv;
247: PetscFunctionBegin;
248: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv));
249: if (!isSinv) PetscCall(EPSBackTransform_Default(eps));
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: PetscErrorCode EPSReset_ARPACK(EPS eps)
254: {
255: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
257: PetscFunctionBegin;
258: PetscCall(PetscFree(ar->workev));
259: PetscCall(PetscFree(ar->workl));
260: PetscCall(PetscFree(ar->select));
261: PetscCall(PetscFree(ar->workd));
262: #if defined(PETSC_USE_COMPLEX)
263: PetscCall(PetscFree(ar->rwork));
264: #endif
265: PetscFunctionReturn(PETSC_SUCCESS);
266: }
268: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
269: {
270: PetscFunctionBegin;
271: PetscCall(PetscFree(eps->data));
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: SLEPC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
276: {
277: EPS_ARPACK *ctx;
279: PetscFunctionBegin;
280: PetscCall(PetscNew(&ctx));
281: eps->data = (void*)ctx;
283: eps->ops->solve = EPSSolve_ARPACK;
284: eps->ops->setup = EPSSetUp_ARPACK;
285: eps->ops->setupsort = EPSSetUpSort_Basic;
286: eps->ops->destroy = EPSDestroy_ARPACK;
287: eps->ops->reset = EPSReset_ARPACK;
288: eps->ops->backtransform = EPSBackTransform_ARPACK;
289: PetscFunctionReturn(PETSC_SUCCESS);
290: }