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: } else /* set default value of ncv */
21: eps->ncv = PetscMin(PetscMax(20,2*eps->nev+1),N);
22: if (!eps->max_it) eps->max_it = PetscMax(300,(int)(2*N/eps->ncv));
24: ncv = eps->ncv;
25: #if defined(PETSC_USE_COMPLEX)
26: PetscFree(ar->rwork);
27: PetscMalloc(ncv*sizeof(PetscReal),&ar->rwork);
28: ar->lworkl = 3*ncv*ncv+5*ncv;
29: PetscFree(ar->workev);
30: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
31: #else
32: if( eps->ishermitian ) {
33: ar->lworkl = ncv*(ncv+8);
34: } else {
35: ar->lworkl = 3*ncv*ncv+6*ncv;
36: PetscFree(ar->workev);
37: PetscMalloc(3*ncv*sizeof(PetscScalar),&ar->workev);
38: }
39: #endif
40: PetscFree(ar->workl);
41: PetscMalloc(ar->lworkl*sizeof(PetscScalar),&ar->workl);
42: PetscFree(ar->select);
43: PetscMalloc(ncv*sizeof(PetscTruth),&ar->select);
44: VecGetLocalSize(eps->vec_initial,&n);
45: PetscFree(ar->workd);
46: PetscMalloc(3*n*sizeof(PetscScalar),&ar->workd);
48: EPSDefaultGetWork(eps,1);
49: EPSAllocateSolutionContiguous(eps);
51: return(0);
52: }
56: PetscErrorCode EPSSolve_ARPACK(EPS eps)
57: {
59: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
60: char bmat[1], howmny[] = "A";
61: const char *which;
62: PetscInt nn;
63: int n, iparam[11], ipntr[14], ido, info;
64: PetscScalar sigmar = 0.0, sigmai, *pV, *resid;
65: Vec x, y, w;
66: Mat A,B;
67: PetscTruth isSinv,isShift,rvec;
68: MPI_Fint fcomm;
69:
72: fcomm = MPI_Comm_c2f(eps->comm);
73: VecGetLocalSize(eps->vec_initial,&nn);
74: n = nn;
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] = 'I';
98: iparam[6] = 1;
99: if (eps->ishermitian && eps->isgeneralized) {
100: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
101: PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);
102: if (isSinv) {
103: bmat[0] = 'G';
104: iparam[6] = 3;
105: STGetShift(eps->OP,&sigmar);
106: sigmai = 0.0;
107: } else if (isShift) {
108: bmat[0] = 'G';
109: iparam[6] = 2;
110: }
111: }
112:
113: #if !defined(PETSC_USE_COMPLEX)
114: if (eps->ishermitian) {
115: switch(eps->which) {
116: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
117: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
118: case EPS_LARGEST_REAL: which = "LA"; break;
119: case EPS_SMALLEST_REAL: which = "SA"; break;
120: default: SETERRQ(1,"Wrong value of eps->which");
121: }
122: } else {
123: #endif
124: switch(eps->which) {
125: case EPS_LARGEST_MAGNITUDE: which = "LM"; break;
126: case EPS_SMALLEST_MAGNITUDE: which = "SM"; break;
127: case EPS_LARGEST_REAL: which = "LR"; break;
128: case EPS_SMALLEST_REAL: which = "SR"; break;
129: case EPS_LARGEST_IMAGINARY: which = "LI"; break;
130: case EPS_SMALLEST_IMAGINARY: which = "SI"; break;
131: default: SETERRQ(1,"Wrong value of eps->which");
132: }
133: #if !defined(PETSC_USE_COMPLEX)
134: }
135: #endif
137: do {
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: }
145: else {
146: ARnaupd_( &fcomm, &ido, bmat, &n, which, &eps->nev, &eps->tol,
147: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
148: ar->workl, &ar->lworkl, &info, 1, 2 );
149: }
150: #else
151: ARnaupd_( &fcomm, &ido, bmat, &n, which, &eps->nev, &eps->tol,
152: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
153: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 2 );
154: #endif
155:
156: if (ido >= -1 && ido <= 2) {
157: VecPlaceArray(x,&ar->workd[ipntr[0]-1]);
158: VecPlaceArray(y,&ar->workd[ipntr[1]-1]);
159: if (ido == 1 || ido == -1) { /* Y=OP*X */
160: STApply(eps->OP,x,y);
161: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,y,PETSC_NULL,PETSC_NULL,PETSC_NULL);
162: if (ido == 1 && iparam[6] == 2) { /* X=A*X */
163: w = eps->work[0];
164: STGetOperators(eps->OP,&A,PETSC_NULL);
165: MatMult(A,x,w);
166: VecCopy(w,x);
167: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,x,PETSC_NULL,PETSC_NULL,PETSC_NULL);
168: }
169: } else if (ido == 2) { /* Y=B*X */
170: STGetOperators(eps->OP,PETSC_NULL,&B);
171: MatMult(B,x,y);
172: }
173: VecResetArray(x);
174: VecResetArray(y);
175: } else if (ido != 99) {
176: SETERRQ1(1,"Internal error in ARPACK reverse comunication interface (ido=%i)\n",ido);
177: }
178:
179: } while (ido != 99);
181: eps->nconv = iparam[4];
182: eps->its = iparam[2];
183:
184: if (info==3) { SETERRQ(1,"No shift could be applied in xxAUPD.\n"
185: "Try increasing the size of NCV relative to NEV."); }
186: else if (info!=0 && info!=1) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxAUPD (%d)",info);}
188: rvec = PETSC_TRUE;
190: if (eps->nconv > 0) {
191: #if !defined(PETSC_USE_COMPLEX)
192: if (eps->ishermitian) {
193: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,&ar->workl[ipntr[6]-1],eps->ncv);
194: ARseupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
195: pV, &n, &sigmar,
196: bmat, &n, which, &eps->nev, &eps->tol,
197: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
198: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
199: }
200: else {
201: EPSMonitor(eps,iparam[2],iparam[4],&ar->workl[ipntr[5]-1],&ar->workl[ipntr[6]-1],&ar->workl[ipntr[7]-1],eps->ncv);
202: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr, eps->eigi,
203: pV, &n, &sigmar, &sigmai, ar->workev,
204: bmat, &n, which, &eps->nev, &eps->tol,
205: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
206: ar->workl, &ar->lworkl, &info, 1, 1, 2 );
207: }
208: #else
209: EPSMonitor(eps,eps->its,iparam[4],&ar->workl[ipntr[5]-1],eps->eigi,(PetscReal*)&ar->workl[ipntr[7]-1],eps->ncv);
210: ARneupd_ ( &fcomm, &rvec, howmny, ar->select, eps->eigr,
211: pV, &n, &sigmar, ar->workev,
212: bmat, &n, which, &eps->nev, &eps->tol,
213: resid, &eps->ncv, pV, &n, iparam, ipntr, ar->workd,
214: ar->workl, &ar->lworkl, ar->rwork, &info, 1, 1, 2 );
215: #endif
216: if (info!=0) { SETERRQ1(PETSC_ERR_LIB,"Error reported by ARPACK subroutine xxEUPD (%d)",info); }
217: }
219: VecRestoreArray( eps->V[0], &pV );
220: VecRestoreArray( eps->vec_initial, &resid );
221: if( eps->nconv >= eps->nev ) eps->reason = EPS_CONVERGED_TOL;
222: else eps->reason = EPS_DIVERGED_ITS;
224: if (eps->ishermitian) {
225: PetscMemcpy(eps->errest,&ar->workl[ipntr[8]-1],eps->nconv);
226: } else {
227: PetscMemcpy(eps->errest,&ar->workl[ipntr[10]-1],eps->nconv);
228: }
230: VecDestroy(x);
231: VecDestroy(y);
233: return(0);
234: }
238: PetscErrorCode EPSBackTransform_ARPACK(EPS eps)
239: {
241: PetscTruth isShift,isSinv;
244: if (eps->ishermitian && eps->isgeneralized) {
245: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
246: PetscTypeCompare((PetscObject)eps->OP,STSINV,&isSinv);
247: if (isSinv || isShift) return(0);
248: }
249: EPSBackTransform_Default(eps);
250: return(0);
251: }
255: PetscErrorCode EPSDestroy_ARPACK(EPS eps)
256: {
258: EPS_ARPACK *ar = (EPS_ARPACK *)eps->data;
262: PetscFree(ar->workev);
263: PetscFree(ar->workl);
264: PetscFree(ar->select);
265: PetscFree(ar->workd);
266: #if defined(PETSC_USE_COMPLEX)
267: PetscFree(ar->rwork);
268: #endif
269: PetscFree(eps->data);
270: EPSDefaultFreeWork(eps);
271: EPSFreeSolutionContiguous(eps);
272: return(0);
273: }
278: PetscErrorCode EPSCreate_ARPACK(EPS eps)
279: {
281: EPS_ARPACK *arpack;
284: PetscNew(EPS_ARPACK,&arpack);
285: PetscLogObjectMemory(eps,sizeof(EPS_ARPACK));
286: eps->data = (void *) arpack;
287: eps->ops->solve = EPSSolve_ARPACK;
288: eps->ops->setup = EPSSetUp_ARPACK;
289: eps->ops->destroy = EPSDestroy_ARPACK;
290: eps->ops->backtransform = EPSBackTransform_ARPACK;
291: eps->ops->computevectors = EPSComputeVectors_Default;
292: return(0);
293: }