Actual source code: arpack.c
1: /*
2: This file implements a wrapper to the ARPACK package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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 src/eps/impls/arpack/arpackp.h
28: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
29: {
31: PetscInt N, n;
32: PetscInt ncv;
33: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
36: VecGetSize(eps->vec_initial,&N);
37: if (eps->ncv) {
38: if (eps->ncv<eps->nev+2) SETERRQ(1,"The value of ncv must be at least nev+2");
39: } else /* set default value of ncv */
40: eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),N);
41: if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
42: if (!eps->max_it) eps->max_it = PetscMax(300,(PetscInt)(2*N/eps->ncv));
44: ncv = eps->ncv;
45: #if defined(PETSC_USE_COMPLEX)
46: PetscFree(ar->rwork);
47: PetscMalloc(ncv*sizeof(PetscReal),&ar->rwork);
48: ar->lworkl = PetscBLASIntCast(3*ncv*ncv+5*ncv);
49: PetscFree(ar->workev);
50: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
51: #else
52: if( eps->ishermitian ) {
53: ar->lworkl = PetscBLASIntCast(ncv*(ncv+8));
54: } else {
55: ar->lworkl = PetscBLASIntCast(3*ncv*ncv+6*ncv);
56: PetscFree(ar->workev);
57: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
58: }
59: #endif
60: PetscFree(ar->workl);
61: PetscMalloc(ar->lworkl*sizeof(PetscScalar),&ar->workl);
62: PetscFree(ar->select);
63: PetscMalloc(ncv*sizeof(PetscTruth),&ar->select);
64: VecGetLocalSize(eps->vec_initial,&n);
65: PetscFree(ar->workd);
66: PetscMalloc(3*n*sizeof(PetscScalar),&ar->workd);
68: if (eps->extraction) {
69: PetscInfo(eps,"Warning: extraction type ignored\n");
70: }
72: EPSDefaultGetWork(eps,2);
73: EPSAllocateSolution(eps);
75: return(0);
76: }
80: PetscErrorCode EPSSolve_ARPACK(EPS eps)
81: {
83: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
84: char bmat[1], howmny[] = "A";
85: const char *which;
86: PetscInt nn;
87: PetscBLASInt n, iparam[11], ipntr[14], ido, info,
88: nev, ncv;
89: PetscScalar sigmar, *pV, *resid;
90: Vec x, y, w = eps->work[0];
91: Mat A;
92: PetscTruth isSinv, isShift, rvec;
93: PetscBLASInt fcomm;
94: #if !defined(PETSC_USE_COMPLEX)
95: PetscScalar sigmai = 0.0;
96: #endif
98:
99: nev = PetscBLASIntCast(eps->nev);
100: ncv = PetscBLASIntCast(eps->ncv);
101: fcomm = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm));
102: VecGetLocalSize(eps->vec_initial,&nn);
103: n = PetscBLASIntCast(nn);
104: VecCreateMPIWithArray(((PetscObject)eps)->comm,n,PETSC_DECIDE,PETSC_NULL,&x);
105: VecCreateMPIWithArray(((PetscObject)eps)->comm,n,PETSC_DECIDE,PETSC_NULL,&y);
106: VecGetArray(eps->V[0],&pV);
107: VecCopy(eps->vec_initial,eps->work[1]);
108: VecGetArray(eps->work[1],&resid);
109:
110: ido = 0; /* first call to reverse communication interface */
111: info = 1; /* indicates a initial vector is provided */
112: iparam[0] = 1; /* use exact shifts */
113: iparam[2] = PetscBLASIntCast(eps->max_it); /* maximum number of Arnoldi update iterations */
114: iparam[3] = 1; /* blocksize */
115: iparam[4] = 0; /* number of converged Ritz values */
116:
117: /*
118: Computational modes ([]=not supported):
119: symmetric non-symmetric complex
120: 1 1 'I' 1 'I' 1 'I'
121: 2 3 'I' 3 'I' 3 'I'
122: 3 2 'G' 2 'G' 2 'G'
123: 4 3 'G' 3 'G' 3 'G'
124: 5 [ 4 'G' ] [ 3 'G' ]
125: 6 [ 5 'G' ] [ 4 'G' ]
126: */
127: PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);
128: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
129: STGetShift(eps->OP,&sigmar);
130: STGetOperators(eps->OP,&A,PETSC_NULL);
132: if (isSinv) {
133: /* shift-and-invert mode */
134: iparam[6] = 3;
135: if (eps->ispositive) bmat[0] = 'G';
136: else bmat[0] = 'I';
137: } else if (isShift && eps->ispositive) {
138: /* generalized shift mode with B positive definite */
139: iparam[6] = 2;
140: bmat[0] = 'G';
141: } else {
142: /* regular mode */
143: if (eps->ishermitian && eps->isgeneralized)
144: SETERRQ(PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
145: iparam[6] = 1;
146: bmat[0] = 'I';
147: }
148:
149: #if !defined(PETSC_USE_COMPLEX)
150: if (eps->ishermitian) {
151: switch(eps->which) {
152: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
153: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
154: case EPS_LARGEST_REAL: which = "LA"; break;
155: case EPS_SMALLEST_REAL: which = "SA"; break;
156: default: SETERRQ(1,"Wrong value of eps->which");
157: }
158: } else {
159: #endif
160: switch(eps->which) {
161: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
162: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
163: case EPS_LARGEST_REAL: which = "LR"; break;
164: case EPS_SMALLEST_REAL: which = "SR"; break;
165: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
166: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
167: default: SETERRQ(1,"Wrong value of eps->which");
168: }
169: #if !defined(PETSC_USE_COMPLEX)
170: }
171: #endif
173: do {
175: #if !defined(PETSC_USE_COMPLEX)
176: if (eps->ishermitian) {
177: ARsaupd_( &fcomm, &ido, bmat, &n, which, &nev, &eps->tol,
178: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
179: ar->workl, &ar->lworkl, &info, 1, 2 );
180: }
181: else {
182: ARnaupd_( &fcomm, &ido, bmat, &n, which, &nev, &eps->tol,
183: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
184: ar->workl, &ar->lworkl, &info, 1, 2 );
185: }
186: #else
187: ARnaupd_( &fcomm, &ido, bmat, &n, which, &nev, &eps->tol,
188: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
189: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 2 );
190: #endif
191:
192: if (ido == -1 || ido == 1 || ido == 2) {
193: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
194: /* special case for shift-and-invert with B semi-positive definite*/
195: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
196: } else {
197: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
198: }
199: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
200:
201: if (ido == -1) {
202: /* Y = OP * X for for the initialization phase to
203: force the starting vector into the range of OP */
204: STApply(eps->OP,x,y);
205: } else if (ido == 2) {
206: /* Y = B * X */
207: IPApplyMatrix(eps->ip,x,y);
208: } else { /* ido == 1 */
209: if (iparam[6] == 3 && bmat[0] == 'G') {
210: /* Y = OP * X for shift-and-invert with B semi-positive definite */
211: STAssociatedKSPSolve(eps->OP,x,y);
212: } else if (iparam[6] == 2) {
213: /* X=A*X Y=B^-1*X for shift with B positive definite */
214: MatMult(A,x,y);
215: if (sigmar != 0.0) {
216: IPApplyMatrix(eps->ip,x,w);
217: VecAXPY(y,sigmar,w);
218: }
219: VecCopy(y,x);
220: STAssociatedKSPSolve(eps->OP,x,y);
221: } else {
222: /* Y = OP * X */
223: STApply(eps->OP,x,y);
224: }
225: IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL,w,PETSC_NULL);
226: }
227:
228: VecResetArray(x);
229: VecResetArray(y);
230: } else if (ido != 99) {
231: SETERRQ1(1,"Internal error in ARPACK reverse comunication interface (ido=%i)\n",ido);
232: }
233:
234: } while (ido != 99);
236: eps->nconv = iparam[4];
237: eps->its = iparam[2];
238:
239: if (info==3) { SETERRQ(1,"No shift could be applied in xxAUPD.\n"
240: "Try increasing the size of NCV relative to NEV."); }
241: else if (info!=0 && info!=1) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);}
243: rvec = PETSC_TRUE;
245: if (eps->nconv > 0) {
246: #if !defined(PETSC_USE_COMPLEX)
247: if (eps->ishermitian) {
248: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
249: ARseupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
250: pV, &n, &sigmar,
251: bmat, &n, which, &nev, &eps->tol,
252: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
253: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
254: }
255: else {
256: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
257: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr, eps->eigi,
258: pV, &n, &sigmar, &sigmai, ar->workev,
259: bmat, &n, which, &nev, &eps->tol,
260: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
261: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
262: }
263: #else
264: EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
265: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
266: pV, &n, &sigmar, ar->workev,
267: bmat, &n, which, &nev, &eps->tol,
268: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
269: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 1, 2 );
270: #endif
271: if (info!=0) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info); }
272: }
274: VecRestoreArray( eps->V[0], &pV );
275: VecRestoreArray( eps->work[1], &resid );
276: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
277: else eps->reason = EPS_DIVERGED_ITS;
279: if (eps->ishermitian) {
280: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
281: } else {
282: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
283: }
285: VecDestroy(x);
286: VecDestroy(y);
288: return(0);
289: }
293: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
294: {
296: PetscTruth isSinv;
299: PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);
300: if (!isSinv) {
301: EPSBackTransform_Default(eps);
302: }
303: return(0);
304: }
308: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
309: {
311: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
315: PetscFree(ar->workev);
316: PetscFree(ar->workl);
317: PetscFree(ar->select);
318: PetscFree(ar->workd);
319: #if defined(PETSC_USE_COMPLEX)
320: PetscFree(ar->rwork);
321: #endif
322: PetscFree(eps->data);
323: EPSDefaultFreeWork(eps);
324: EPSFreeSolution(eps);
325: return(0);
326: }
331: PetscErrorCode EPSCreate_ARPACK(EPS eps)
332: {
334: EPS_ARPACK *arpack;
337: PetscNew(EPS_ARPACK,&arpack);
338: PetscLogObjectMemory(eps,sizeof(EPS_ARPACK));
339: eps->data = (void *) arpack;
340: eps->ops->solve = EPSSolve_ARPACK;
341: eps->ops->setup = EPSSetUp_ARPACK;
342: eps->ops->destroy = EPSDestroy_ARPACK;
343: eps->ops->backtransform = EPSBackTransform_ARPACK;
344: eps->ops->computevectors = EPSComputeVectors_Default;
345: return(0);
346: }