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-2010, 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/external/arpack/arpackp.h
25: #include private/stimpl.h
27: PetscErrorCode EPSSolve_ARPACK(EPS);
31: PetscErrorCode EPSSetUp_ARPACK(EPS eps)
32: {
34: PetscInt ncv;
35: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
38: if (eps->ncv) {
39: if (eps->ncv<eps->nev+2) SETERRQ(1,"The value of ncv must be at least nev+2");
40: } else /* set default value of ncv */
41: eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),eps->n);
42: if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
43: if (!eps->max_it) eps->max_it = PetscMax(300,(PetscInt)(2*eps->n/eps->ncv));
44: if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
46: ncv = eps->ncv;
47: #if defined(PETSC_USE_COMPLEX)
48: PetscFree(ar->rwork);
49: PetscMalloc(ncv*sizeof(PetscReal),&ar->rwork);
50: ar->lworkl = PetscBLASIntCast(3*ncv*ncv+5*ncv);
51: PetscFree(ar->workev);
52: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
53: #else
54: if( eps->ishermitian ) {
55: ar->lworkl = PetscBLASIntCast(ncv*(ncv+8));
56: } else {
57: ar->lworkl = PetscBLASIntCast(3*ncv*ncv+6*ncv);
58: PetscFree(ar->workev);
59: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
60: }
61: #endif
62: PetscFree(ar->workl);
63: PetscMalloc(ar->lworkl*sizeof(PetscScalar),&ar->workl);
64: PetscFree(ar->select);
65: PetscMalloc(ncv*sizeof(PetscTruth),&ar->select);
66: PetscFree(ar->workd);
67: PetscMalloc(3*eps->nloc*sizeof(PetscScalar),&ar->workd);
69: if (eps->extraction) {
70: PetscInfo(eps,"Warning: extraction type ignored\n");
71: }
73: if (eps->balance!=EPS_BALANCE_NONE)
74: SETERRQ(PETSC_ERR_SUP,"Balancing not supported in the Arpack interface");
76: EPSAllocateSolution(eps);
77: EPSDefaultGetWork(eps,2);
79: /* dispatch solve method */
80: if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
81: eps->ops->solve = EPSSolve_ARPACK;
82: return(0);
83: }
87: PetscErrorCode EPSSolve_ARPACK(EPS eps)
88: {
90: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
91: char bmat[1], howmny[] = "A";
92: const char *which;
93: PetscBLASInt n, iparam[11], ipntr[14], ido, info,
94: nev, ncv;
95: PetscScalar sigmar, *pV, *resid;
96: Vec x, y, w = eps->work[0];
97: Mat A;
98: PetscTruth isSinv, isShift, rvec;
99: PetscBLASInt fcomm;
100: #if !defined(PETSC_USE_COMPLEX)
101: PetscScalar sigmai = 0.0;
102: #endif
104:
105: nev = PetscBLASIntCast(eps->nev);
106: ncv = PetscBLASIntCast(eps->ncv);
107: fcomm = PetscBLASIntCast(MPI_Comm_c2f(((PetscObject)eps)->comm));
108: n = PetscBLASIntCast(eps->nloc);
109: VecCreateMPIWithArray(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,PETSC_NULL,&x);
110: VecCreateMPIWithArray(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,PETSC_NULL,&y);
111: VecGetArray(eps->V[0],&pV);
112: EPSGetStartVector(eps,0,eps->work[1],PETSC_NULL);
113: VecGetArray(eps->work[1],&resid);
114:
115: ido = 0; /* first call to reverse communication interface */
116: info = 1; /* indicates a initial vector is provided */
117: iparam[0] = 1; /* use exact shifts */
118: iparam[2] = PetscBLASIntCast(eps->max_it); /* maximum number of Arnoldi update iterations */
119: iparam[3] = 1; /* blocksize */
120: iparam[4] = 0; /* number of converged Ritz values */
121:
122: /*
123: Computational modes ([]=not supported):
124: symmetric non-symmetric complex
125: 1 1 'I' 1 'I' 1 'I'
126: 2 3 'I' 3 'I' 3 'I'
127: 3 2 'G' 2 'G' 2 'G'
128: 4 3 'G' 3 'G' 3 'G'
129: 5 [ 4 'G' ] [ 3 'G' ]
130: 6 [ 5 'G' ] [ 4 'G' ]
131: */
132: PetscTypeCompare((PetscObject)eps->OP,STSINVERT,&isSinv);
133: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
134: STGetShift(eps->OP,&sigmar);
135: STGetOperators(eps->OP,&A,PETSC_NULL);
137: if (isSinv) {
138: /* shift-and-invert mode */
139: iparam[6] = 3;
140: if (eps->ispositive) bmat[0] = 'G';
141: else bmat[0] = 'I';
142: } else if (isShift && eps->ispositive) {
143: /* generalized shift mode with B positive definite */
144: iparam[6] = 2;
145: bmat[0] = 'G';
146: } else {
147: /* regular mode */
148: if (eps->ishermitian && eps->isgeneralized)
149: SETERRQ(PETSC_ERR_SUP,"Spectral transformation not supported by ARPACK hermitian solver");
150: iparam[6] = 1;
151: bmat[0] = 'I';
152: }
153:
154: #if !defined(PETSC_USE_COMPLEX)
155: if (eps->ishermitian) {
156: switch(eps->which) {
157: case EPS_TARGET_MAGNITUDE:
158: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
159: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
160: case EPS_TARGET_REAL:
161: case EPS_LARGEST_REAL: which = "LA"; break;
162: case EPS_SMALLEST_REAL: which = "SA"; break;
163: default: SETERRQ(1,"Wrong value of eps->which");
164: }
165: } else {
166: #endif
167: switch(eps->which) {
168: case EPS_TARGET_MAGNITUDE:
169: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
170: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
171: case EPS_TARGET_REAL:
172: case EPS_LARGEST_REAL: which = "LR"; break;
173: case EPS_SMALLEST_REAL: which = "SR"; break;
174: case EPS_TARGET_IMAGINARY:
175: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
176: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
177: default: SETERRQ(1,"Wrong value of eps->which");
178: }
179: #if !defined(PETSC_USE_COMPLEX)
180: }
181: #endif
183: do {
185: #if !defined(PETSC_USE_COMPLEX)
186: if (eps->ishermitian) {
187: ARsaupd_( &fcomm, &ido, bmat, &n, which, &nev, &eps->tol,
188: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
189: ar->workl, &ar->lworkl, &info, 1, 2 );
190: }
191: else {
192: ARnaupd_( &fcomm, &ido, bmat, &n, which, &nev, &eps->tol,
193: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
194: ar->workl, &ar->lworkl, &info, 1, 2 );
195: }
196: #else
197: ARnaupd_( &fcomm, &ido, bmat, &n, which, &nev, &eps->tol,
198: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
199: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 2 );
200: #endif
201:
202: if (ido == -1 || ido == 1 || ido == 2) {
203: if (ido == 1 && iparam[6] == 3 && bmat[0] == 'G') {
204: /* special case for shift-and-invert with B semi-positive definite*/
205: VecPlaceArray(x,&ar->workd[ipntr[2]-1]);
206: } else {
207: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
208: }
209: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
210:
211: if (ido == -1) {
212: /* Y = OP * X for for the initialization phase to
213: force the starting vector into the range of OP */
214: STApply(eps->OP,x,y);
215: } else if (ido == 2) {
216: /* Y = B * X */
217: IPApplyMatrix(eps->ip,x,y);
218: } else { /* ido == 1 */
219: if (iparam[6] == 3 && bmat[0] == 'G') {
220: /* Y = OP * X for shift-and-invert with B semi-positive definite */
221: STAssociatedKSPSolve(eps->OP,x,y);
222: } else if (iparam[6] == 2) {
223: /* X=A*X Y=B^-1*X for shift with B positive definite */
224: MatMult(A,x,y);
225: if (sigmar != 0.0) {
226: IPApplyMatrix(eps->ip,x,w);
227: VecAXPY(y,sigmar,w);
228: }
229: VecCopy(y,x);
230: STAssociatedKSPSolve(eps->OP,x,y);
231: } else {
232: /* Y = OP * X */
233: STApply(eps->OP,x,y);
234: }
235: IPOrthogonalize(eps->ip,0,PETSC_NULL,eps->nds,PETSC_NULL,eps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
236: }
237:
238: VecResetArray(x);
239: VecResetArray(y);
240: } else if (ido != 99) {
241: SETERRQ1(1,"Internal error in ARPACK reverse comunication interface (ido=%i)\n",ido);
242: }
243:
244: } while (ido != 99);
246: eps->nconv = iparam[4];
247: eps->its = iparam[2];
248:
249: if (info==3) { SETERRQ(1,"No shift could be applied in xxAUPD.\n"
250: "Try increasing the size of NCV relative to NEV."); }
251: else if (info!=0 && info!=1) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);}
253: rvec = PETSC_TRUE;
255: if (eps->nconv > 0) {
256: #if !defined(PETSC_USE_COMPLEX)
257: if (eps->ishermitian) {
258: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
259: ARseupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
260: pV, &n, &sigmar,
261: bmat, &n, which, &nev, &eps->tol,
262: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
263: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
264: }
265: else {
266: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
267: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr, eps->eigi,
268: pV, &n, &sigmar, &sigmai, ar->workev,
269: bmat, &n, which, &nev, &eps->tol,
270: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
271: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
272: }
273: #else
274: EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
275: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
276: pV, &n, &sigmar, ar->workev,
277: bmat, &n, which, &nev, &eps->tol,
278: resid, &ncv, pV, &n, iparam, ipntr, ar->workd,
279: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 1, 2 );
280: #endif
281: if (info!=0) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info); }
282: }
284: VecRestoreArray( eps->V[0], &pV );
285: VecRestoreArray( eps->work[1], &resid );
286: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
287: else eps->reason = EPS_DIVERGED_ITS;
289: if (eps->ishermitian) {
290: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
291: } else {
292: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
293: }
295: VecDestroy(x);
296: VecDestroy(y);
298: return(0);
299: }
303: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
304: {
306: PetscTruth isSinv;
309: PetscTypeCompare((PetscObject)eps->OP,STSINVERT,&isSinv);
310: if (!isSinv) {
311: EPSBackTransform_Default(eps);
312: }
313: return(0);
314: }
318: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
319: {
321: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
325: PetscFree(ar->workev);
326: PetscFree(ar->workl);
327: PetscFree(ar->select);
328: PetscFree(ar->workd);
329: #if defined(PETSC_USE_COMPLEX)
330: PetscFree(ar->rwork);
331: #endif
332: PetscFree(eps->data);
333: EPSDefaultFreeWork(eps);
334: EPSFreeSolution(eps);
335: return(0);
336: }
341: PetscErrorCode EPSCreate_ARPACK(EPS eps)
342: {
344: EPS_ARPACK *arpack;
347: PetscNew(EPS_ARPACK,&arpack);
348: PetscLogObjectMemory(eps,sizeof(EPS_ARPACK));
349: eps->data = (void *) arpack;
350: eps->ops->setup = EPSSetUp_ARPACK;
351: eps->ops->destroy = EPSDestroy_ARPACK;
352: eps->ops->backtransform = EPSBackTransform_ARPACK;
353: eps->ops->computevectors = EPSComputeVectors_Default;
354: return(0);
355: }