Actual source code: arpack.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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 <../src/eps/impls/external/arpack/arpackp.h>
17: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
18: {
20: PetscInt ncv;
21: PetscBool istrivial;
22: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
25: if (eps->ncv) {
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) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
29: if (!eps->max_it) eps->max_it = PetscMax(300,(PetscInt)(2*eps->n/eps->ncv));
30: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
32: ncv = eps->ncv;
33: #if defined(PETSC_USE_COMPLEX)
34: PetscFree(ar->rwork);
35: PetscMalloc1(ncv,&ar->rwork);
36: PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscReal));
37: PetscBLASIntCast(3*ncv*ncv+5*ncv,&ar->lworkl);
38: PetscFree(ar->workev);
39: PetscMalloc1(3*ncv,&ar->workev);
40: PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
41: #else
42: if (eps->ishermitian) {
43: PetscBLASIntCast(ncv*(ncv+8),&ar->lworkl);
44: } else {
45: PetscBLASIntCast(3*ncv*ncv+6*ncv,&ar->lworkl);
46: PetscFree(ar->workev);
47: PetscMalloc1(3*ncv,&ar->workev);
48: PetscLogObjectMemory((PetscObject)eps,3*ncv*sizeof(PetscScalar));
49: }
50: #endif
51: PetscFree(ar->workl);
52: PetscMalloc1(ar->lworkl,&ar->workl);
53: PetscLogObjectMemory((PetscObject)eps,ar->lworkl*sizeof(PetscScalar));
54: PetscFree(ar->select);
55: PetscMalloc1(ncv,&ar->select);
56: PetscLogObjectMemory((PetscObject)eps,ncv*sizeof(PetscBool));
57: PetscFree(ar->workd);
58: PetscMalloc1(3*eps->nloc,&ar->workd);
59: PetscLogObjectMemory((PetscObject)eps,3*eps->nloc*sizeof(PetscScalar));
61: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
63: if (eps->balance!=EPS_BALANCE_NONE) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Balancing not supported in the Arpack interface");
64: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
65: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
67: EPSAllocateSolution(eps,0);
68: EPS_SetInnerProduct(eps);
69: EPSSetWorkVecs(eps,2);
71: RGIsTrivial(eps->rg,&istrivial);
72: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
73: return(0);
74: }
76: PetscErrorCode EPSSolve_ARPACK(EPS eps)
77: {
79: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
80: char bmat[1],howmny[] = "A";
81: const char *which;
82: PetscBLASInt n,iparam[11],ipntr[14],ido,info,nev,ncv;
83: #if !defined(PETSC_HAVE_MPIUNI)
84: PetscBLASInt fcomm;
85: #endif
86: PetscScalar sigmar,*pV,*resid;
87: Vec x,y,w = eps->work[0];
88: Mat A;
89: PetscBool isSinv,isShift,rvec;
90: #if !defined(PETSC_USE_COMPLEX)
91: PetscScalar sigmai = 0.0;
92: #endif
95: PetscBLASIntCast(eps->nev,&nev);
96: PetscBLASIntCast(eps->ncv,&ncv);
97: #if !defined(PETSC_HAVE_MPIUNI)
98: PetscBLASIntCast(MPI_Comm_c2f(PetscObjectComm((PetscObject)eps)),&fcomm);
99: #endif
100: PetscBLASIntCast(eps->nloc,&n);
101: EPSGetStartVector(eps,0,NULL);
102: BVSetActiveColumns(eps->V,0,0); /* just for deflation space */
103: BVCopyVec(eps->V,0,eps->work[1]);
104: BVGetArray(eps->V,&pV);
105: VecGetArray(eps->work[1],&resid);
107: ido = 0; /* first call to reverse communication interface */
108: info = 1; /* indicates an initial vector is provided */
109: iparam[0] = 1; /* use exact shifts */
110: PetscBLASIntCast(eps->max_it,&iparam[2]); /* max Arnoldi iterations */
111: iparam[3] = 1; /* blocksize */
112: iparam[4] = 0; /* number of converged Ritz values */
114: /*
115: Computational modes ([]=not supported):
116: symmetric non-symmetric complex
117: 1 1 'I' 1 'I' 1 'I'
118: 2 3 'I' 3 'I' 3 'I'
119: 3 2 'G' 2 'G' 2 'G'
120: 4 3 'G' 3 'G' 3 'G'
121: 5 [ 4 'G' ] [ 3 'G' ]
122: 6 [ 5 'G' ] [ 4 'G' ]
123: */
124: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
125: PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&isShift);
126: STGetShift(eps->st,&sigmar);
127: STGetMatrix(eps->st,0,&A);
128: MatCreateVecsEmpty(A,&x,&y);
130: if (isSinv) {
131: /* shift-and-invert mode */
132: iparam[6] = 3;
133: if (eps->ispositive) bmat[0] = 'G';
134: else bmat[0] = 'I';
135: } else if (isShift && eps->ispositive) {
136: /* generalized shift mode with B positive definite */
137: iparam[6] = 2;
138: bmat[0] = 'G';
139: } else {
140: /* regular mode */
141: if (eps->ishermitian && eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
142: iparam[6] = 1;
143: bmat[0] = 'I';
144: }
146: #if !defined(PETSC_USE_COMPLEX)
147: if (eps->ishermitian) {
148: switch (eps->which) {
149: case EPS_TARGET_MAGNITUDE:
150: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
151: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
152: case EPS_TARGET_REAL:
153: case EPS_LARGEST_REAL: which = "LA"; break;
154: case EPS_SMALLEST_REAL: which = "SA"; break;
155: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
156: }
157: } else {
158: #endif
159: switch (eps->which) {
160: case EPS_TARGET_MAGNITUDE:
161: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
162: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
163: case EPS_TARGET_REAL:
164: case EPS_LARGEST_REAL: which = "LR"; break;
165: case EPS_SMALLEST_REAL: which = "SR"; break;
166: case EPS_TARGET_IMAGINARY:
167: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
168: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
169: default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Wrong value of eps->which");
170: }
171: #if !defined(PETSC_USE_COMPLEX)
172: }
173: #endif
175: do {
177: #if !defined(PETSC_USE_COMPLEX)
178: if (eps->ishermitian) {
179: PetscStackCall("ARPACKsaupd",ARPACKsaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
180: } else {
181: PetscStackCall("ARPACKnaupd",ARPACKnaupd_(&fcomm,&ido,bmat,&n,which,&nev,&eps->tol,resid,&ncv,pV,&n,iparam,ipntr,ar->workd,ar->workl,&ar->lworkl,&info));
182: }
183: #else
184: 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));
185: #endif
187: if (ido == -1 || ido == 1 || ido == 2) {
188: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
189: /* special case for shift-and-invert with B semi-positive definite*/
190: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
191: } else {
192: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
193: }
194: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
196: if (ido == -1) {
197: /* Y = OP * X for for the initialization phase to
198: force the starting vector into the range of OP */
199: STApply(eps->st,x,y);
200: } else if (ido == 2) {
201: /* Y = B * X */
202: BVApplyMatrix(eps->V,x,y);
203: } else { /* ido == 1 */
204: if (iparam[6] == 3 && bmat[0] == 'G') {
205: /* Y = OP * X for shift-and-invert with B semi-positive definite */
206: STMatSolve(eps->st,x,y);
207: } else if (iparam[6] == 2) {
208: /* X=A*X Y=B^-1*X for shift with B positive definite */
209: MatMult(A,x,y);
210: if (sigmar != 0.0) {
211: BVApplyMatrix(eps->V,x,w);
212: VecAXPY(y,sigmar,w);
213: }
214: VecCopy(y,x);
215: STMatSolve(eps->st,x,y);
216: } else {
217: /* Y = OP * X */
218: STApply(eps->st,x,y);
219: }
220: BVOrthogonalizeVec(eps->V,y,NULL,NULL,NULL);
221: }
223: VecResetArray(x);
224: VecResetArray(y);
225: } else if (ido != 99) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Internal error in ARPACK reverse comunication interface (ido=%d)",ido);
227: } while (ido != 99);
229: eps->nconv = iparam[4];
230: eps->its = iparam[2];
232: 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");
233: else if (info!=0 && info!=1) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",(int)info);
235: rvec = PETSC_TRUE;
237: if (eps->nconv > 0) {
238: #if !defined(PETSC_USE_COMPLEX)
239: if (eps->ishermitian) {
240: 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));
241: } else {
242: 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));
243: }
244: #else
245: 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));
246: #endif
247: if (info!=0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",(int)info);
248: }
250: BVRestoreArray(eps->V,&pV);
251: VecRestoreArray(eps->work[1],&resid);
252: if (eps->nconv >= eps->nev) eps->reason = EPS_CONVERGED_TOL;
253: else eps->reason = EPS_DIVERGED_ITS;
255: VecDestroy(&x);
256: VecDestroy(&y);
257: return(0);
258: }
260: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
261: {
263: PetscBool isSinv;
266: PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&isSinv);
267: if (!isSinv) {
268: EPSBackTransform_Default(eps);
269: }
270: return(0);
271: }
273: PetscErrorCode EPSReset_ARPACK(EPS eps)
274: {
276: EPS_ARPACK *ar = (EPS_ARPACK*)eps->data;
279: PetscFree(ar->workev);
280: PetscFree(ar->workl);
281: PetscFree(ar->select);
282: PetscFree(ar->workd);
283: #if defined(PETSC_USE_COMPLEX)
284: PetscFree(ar->rwork);
285: #endif
286: return(0);
287: }
289: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
290: {
294: PetscFree(eps->data);
295: return(0);
296: }
298: SLEPC_EXTERN PetscErrorCode EPSCreate_ARPACK(EPS eps)
299: {
300: EPS_ARPACK *ctx;
304: PetscNewLog(eps,&ctx);
305: eps->data = (void*)ctx;
307: eps->ops->solve = EPSSolve_ARPACK;
308: eps->ops->setup = EPSSetUp_ARPACK;
309: eps->ops->destroy = EPSDestroy_ARPACK;
310: eps->ops->reset = EPSReset_ARPACK;
311: eps->ops->backtransform = EPSBackTransform_ARPACK;
312: return(0);
313: }