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