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