Actual source code: arpack.c
2: /*
3: This file implements a wrapper to the ARPACK package
4: */
5: #include src/eps/impls/arpack/arpackp.h
9: static int EPSSetUp_ARPACK(EPS eps)
10: {
11: int ierr, n, ncv;
12: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
15: ncv = eps->ncv;
16: #if defined(PETSC_USE_COMPLEX)
17: PetscMalloc(ncv*sizeof(PetscReal),&ar->rwork);
18: ar->lworkl = 3*ncv*ncv+5*ncv;
19: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
20: #else
21: if( eps->ishermitian ) {
22: ar->lworkl = ncv*(ncv+8);
23: }
24: else {
25: ar->lworkl = 3*ncv*ncv+6*ncv;
26: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
27: }
28: #endif
29: PetscMalloc(ar->lworkl*sizeof(PetscScalar),&ar->workl);
30: PetscMalloc(ncv*sizeof(PetscTruth),&ar->select);
31: VecGetLocalSize(eps->vec_initial,&n);
32: PetscMalloc(3*n*sizeof(PetscScalar),&ar->workd);
34: EPSDefaultGetWork(eps,1);
36: return(0);
37: }
41: static int EPSSetDefaults_ARPACK(EPS eps)
42: {
43: int ierr, N;
46: VecGetSize(eps->vec_initial,&N);
47: if (eps->ncv) {
48: if (eps->ncv<eps->nev+2) SETERRQ(1,"The value of ncv must be at least nev+2");
49: if (eps->ncv>N) SETERRQ(1,"The value of ncv cannot be larger than N");
50: }
51: else /* set default value of ncv */
52: eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),N);
54: if (!eps->max_it) eps->max_it = PetscMax(300,(int)ceil(2*N/eps->ncv));
55: if (!eps->tol) eps->tol = 1.e-7;
56: return(0);
57: }
61: static int EPSSolve_ARPACK(EPS eps)
62: {
63: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
64: char bmat[1], *which, howmny[] = "A";
65: int i, n, iparam[11], ipntr[14], ido, info, ierr;
66: PetscScalar sigmar, sigmai, *pV, *resid;
67: Vec x, y, w;
68: Mat A;
69: PetscTruth isSinv,rvec;
70: MPI_Fint fcomm;
71:
74: fcomm = MPI_Comm_c2f(eps->comm);
75: VecGetLocalSize(eps->vec_initial,&n);
76: VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&x);
77: VecCreateMPIWithArray(eps->comm,n,PETSC_DECIDE,PETSC_NULL,&y);
78: VecGetArray(eps->V[0],&pV);
79: VecGetArray(eps->vec_initial,&resid);
80:
81: ido = 0; /* first call to reverse communication interface */
82: info = 1; /* indicates a initial vector is provided */
83: iparam[0] = 1; /* use exact shifts */
84: iparam[2] = eps->max_it; /* maximum number of Arnoldi update iterations */
85: iparam[3] = 1; /* blocksize */
86: iparam[4] = 0; /* number of converged Ritz values */
87:
88: /*
89: Computational modes ([]=not supported):
90: symmetric non-symmetric complex
91: 1 1 'I' 1 'I' 1 'I'
92: 2 3 'I' 3 'I' 3 'I'
93: 3 2 'G' 2 'G' 2 'G'
94: 4 3 'G' 3 'G' 3 'G'
95: 5 [ 4 'G' ] [ 3 'G' ]
96: 6 [ 5 'G' ] [ 4 'G' ]
97: */
98: bmat[0] = eps->isgeneralized? 'G': 'I';
99: PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);
100: if(isSinv) {
101: iparam[6] = 3;
102: STGetShift(eps->OP,&sigmar);
103: sigmai = 0.0;
104: }
105: else iparam[6] = eps->isgeneralized? 2: 1;
106:
107: #if !defined(PETSC_USE_COMPLEX)
108: if (eps->ishermitian) {
109: switch(eps->which) {
110: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
111: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
112: case EPS_LARGEST_REAL: which = "LA"; break;
113: case EPS_SMALLEST_REAL: which = "SA"; break;
114: default: SETERRQ(1,"Wrong value of eps->which");
115: }
116: } else {
117: #endif
118: switch(eps->which) {
119: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
120: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
121: case EPS_LARGEST_REAL: which = "LR"; break;
122: case EPS_SMALLEST_REAL: which = "SR"; break;
123: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
124: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
125: default: SETERRQ(1,"Wrong value of eps->which");
126: }
127: #if !defined(PETSC_USE_COMPLEX)
128: }
129: #endif
131: for(;;) {
133: #if !defined(PETSC_USE_COMPLEX)
134: if (eps->ishermitian) {
135: ARsaupd_( &fcomm, &ido, bmat, &n, which, &eps->nev, &eps->tol,
136: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
137: ar->workl, &ar->lworkl, &info, 1, 2 );
138: EPSMonitorEstimates(eps,iparam[2],iparam[4],&ar->workl[ipntr[6]-1],eps->ncv);
139: EPSMonitorValues(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],PETSC_NULL,eps->ncv);
140: }
141: else {
142: ARnaupd_( &fcomm, &ido, bmat, &n, which, &eps->nev, &eps->tol,
143: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
144: ar->workl, &ar->lworkl, &info, 1, 2 );
145: EPSMonitorEstimates(eps,iparam[2],iparam[4],&ar->workl[ipntr[7]-1],eps->ncv);
146: EPSMonitorValues(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],eps->ncv);
147: }
148: #else
149: ARnaupd_( &fcomm, &ido, bmat, &n, which, &eps->nev, &eps->tol,
150: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
151: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 2 );
152: EPSMonitorEstimates(eps,iparam[2],iparam[4],(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
153: EPSMonitorValues(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],PETSC_NULL,eps->ncv);
154: #endif
156: if( ido == -1 ) {
157: VecPlaceArray( x, &ar->workd[ipntr[0]-1] );
158: VecPlaceArray( y, &ar->workd[ipntr[1]-1] );
159: STApply( eps->OP, x, y );
160: if (iparam[6]==2) { /* regular inverse mode */
161: w = eps->work[0];
162: STGetOperators( eps->OP, &A, PETSC_NULL );
163: MatMult(A,x,w);
164: VecCopy( w, x );
165: }
166: }
167: else if( ido == 1 ) {
168: VecPlaceArray( y, &ar->workd[ipntr[1]-1] );
169: if (isSinv && eps->isgeneralized) {
170: VecPlaceArray( x, &ar->workd[ipntr[2]-1] );
171: STApplyNoB( eps->OP, x, y );
172: }
173: else {
174: VecPlaceArray( x, &ar->workd[ipntr[0]-1] );
175: STApply( eps->OP, x, y );
176: }
177: if (iparam[6]==2) { /* regular inverse mode */
178: w = eps->work[0];
179: STGetOperators( eps->OP, &A, PETSC_NULL );
180: MatMult(A,x,w);
181: VecCopy( w, x );
182: }
183: }
184: else if( ido == 2 ) {
185: VecPlaceArray( x, &ar->workd[ipntr[0]-1] );
186: VecPlaceArray( y, &ar->workd[ipntr[1]-1] );
187: STApplyB( eps->OP, x, y );
188: }
189: else break;
190: }
192: eps->nconv = iparam[4];
193: eps->its = iparam[2];
194:
195: if (info==1) { SETERRQ(0,"Exceeded maximum number of iterations in xxAUPD"); }
196: else if(info==3) { SETERRQ(1,"No shift could be applied in xxAUPD.\n"
197: "Try increasing the size of NCV relative to NEV."); }
198: else if (info!=0) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);}
200: if (eps->dropvectors) rvec = PETSC_FALSE;
201: else rvec = PETSC_TRUE;
203: #if !defined(PETSC_USE_COMPLEX)
204: if (eps->ishermitian) {
205: ARseupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
206: pV, &n, &sigmar,
207: bmat, &n, which, &eps->nev, &eps->tol,
208: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
209: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
210: for (i=0;i<eps->nconv;i++) eps->eigi[i]=0.0;
211: }
212: else {
213: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr, eps->eigi,
214: pV, &n, &sigmar, &sigmai, ar->workev,
215: bmat, &n, which, &eps->nev, &eps->tol,
216: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
217: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
218: }
219: #else
220: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
221: pV, &n, &sigmar, ar->workev,
222: bmat, &n, which, &eps->nev, &eps->tol,
223: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
224: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 1, 2 );
225: #endif
227: if (info!=0) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info); }
229: VecRestoreArray( eps->V[0], &pV );
230: VecRestoreArray( eps->vec_initial, &resid );
231: eps->reason = EPS_CONVERGED_TOL;
233: if (eps->ishermitian) {
234: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
235: } else {
236: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
237: }
239: VecDestroy(x);
240: VecDestroy(y);
242: return(0);
243: }
247: int EPSBackTransform_ARPACK(EPS eps)
248: {
249: ST st;
250: int ierr,i;
251: PetscTruth isSinv;
255: EPSGetST(eps,&st);
256: PetscTypeCompare((PetscObject)st,STSINV,&isSinv);
257: if (!isSinv)
258: for (i=0;i<eps->nconv;i++) {
259: STBackTransform(st,&eps->eigr[i],&eps->eigi[i]);
260: }
262: return(0);
263: }
267: int EPSDestroy_ARPACK(EPS eps)
268: {
269: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
270: int ierr;
274: if (ar->workev) { PetscFree(ar->workev); }
275: if (ar->workl) { PetscFree(ar->workl); }
276: if (ar->select) { PetscFree(ar->select); }
277: if (ar->workd) { PetscFree(ar->workd); }
278: #if defined(PETSC_USE_COMPLEX)
279: if (ar->rwork) { PetscFree(ar->rwork); }
280: #endif
281: EPSDefaultDestroy(eps);
282: return(0);
283: }
285: EXTERN_C_BEGIN
288: int EPSCreate_ARPACK(EPS eps)
289: {
290: EPS_ARPACK *arpack;
291: int ierr;
294: PetscNew(EPS_ARPACK,&arpack);
295: PetscMemzero(arpack,sizeof(EPS_ARPACK));
296: PetscLogObjectMemory(eps,sizeof(EPS_ARPACK));
297: eps->data = (void *) arpack;
298: eps->ops->setup = EPSSetUp_ARPACK;
299: eps->ops->setdefaults = EPSSetDefaults_ARPACK;
300: eps->ops->solve = EPSSolve_ARPACK;
301: eps->ops->destroy = EPSDestroy_ARPACK;
302: eps->ops->view = 0;
303: eps->ops->backtransform = EPSBackTransform_ARPACK;
304: return(0);
305: }
306: EXTERN_C_END