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