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