Actual source code: ex9.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Solves a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
25: " -L <L>, where <L> = bifurcation parameter.\n"
26: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
27: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
29: #include slepceps.h
31: /*
32: This example computes the eigenvalues with largest real part of the
33: following matrix
35: A = [ tau1*T+(beta-1)*I alpha^2*I
36: -beta*I tau2*T-alpha^2*I ],
38: where
40: T = tridiag{1,-2,1}
41: h = 1/(n+1)
42: tau1 = delta1/(h*L)^2
43: tau2 = delta2/(h*L)^2
44: */
47: /*
48: Matrix operations
49: */
50: PetscErrorCode MatBrussel_Mult(Mat,Vec,Vec);
51: PetscErrorCode MatBrussel_Shift(PetscScalar*,Mat);
52: PetscErrorCode MatBrussel_GetDiagonal(Mat,Vec);
54: typedef struct {
55: Mat T;
56: Vec x1, x2, y1, y2;
57: PetscScalar alpha, beta, tau1, tau2, sigma;
58: } CTX_BRUSSEL;
62: int main( int argc, char **argv )
63: {
64: Mat A; /* eigenvalue problem matrix */
65: EPS eps; /* eigenproblem solver context */
66: const EPSType type;
67: PetscReal error, tol, re, im;
68: PetscScalar delta1, delta2, L, h, kr, ki, value[3];
69: PetscInt N=30, n, i, col[3], Istart, Iend, nev, maxit, its, nconv;
70: PetscTruth FirstBlock=PETSC_FALSE, LastBlock=PETSC_FALSE;
72: CTX_BRUSSEL *ctx;
74: SlepcInitialize(&argc,&argv,(char*)0,help);
76: PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
77: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%d\n\n",N);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Generate the matrix
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: /*
84: Create shell matrix context and set default parameters
85: */
86: PetscNew(CTX_BRUSSEL,&ctx);
87: ctx->alpha = 2.0;
88: ctx->beta = 5.45;
89: delta1 = 0.008;
90: delta2 = 0.004;
91: L = 0.51302;
93: /*
94: Look the command line for user-provided parameters
95: */
96: PetscOptionsGetScalar(PETSC_NULL,"-L",&L,PETSC_NULL);
97: PetscOptionsGetScalar(PETSC_NULL,"-alpha",&ctx->alpha,PETSC_NULL);
98: PetscOptionsGetScalar(PETSC_NULL,"-beta",&ctx->beta,PETSC_NULL);
99: PetscOptionsGetScalar(PETSC_NULL,"-delta1",&delta1,PETSC_NULL);
100: PetscOptionsGetScalar(PETSC_NULL,"-delta2",&delta2,PETSC_NULL);
102: /*
103: Create matrix T
104: */
105: MatCreate(PETSC_COMM_WORLD,&ctx->T);
106: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
107: MatSetFromOptions(ctx->T);
108:
109: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
110: if (Istart==0) FirstBlock=PETSC_TRUE;
111: if (Iend==N) LastBlock=PETSC_TRUE;
112: value[0]=1.0; value[1]=-2.0; value[2]=1.0;
113: for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
114: col[0]=i-1; col[1]=i; col[2]=i+1;
115: MatSetValues(ctx->T,1,&i,3,col,value,INSERT_VALUES);
116: }
117: if (LastBlock) {
118: i=N-1; col[0]=N-2; col[1]=N-1;
119: MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
120: }
121: if (FirstBlock) {
122: i=0; col[0]=0; col[1]=1; value[0]=-2.0; value[1]=1.0;
123: MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
124: }
126: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
127: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
128: MatGetLocalSize(ctx->T,&n,PETSC_NULL);
130: /*
131: Fill the remaining information in the shell matrix context
132: and create auxiliary vectors
133: */
134: h = 1.0 / (PetscReal)(N+1);
135: ctx->tau1 = delta1 / ((h*L)*(h*L));
136: ctx->tau2 = delta2 / ((h*L)*(h*L));
137: ctx->sigma = 0.0;
138: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x1);
139: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x2);
140: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y1);
141: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y2);
143: /*
144: Create the shell matrix
145: */
146: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
147: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatBrussel_Mult);
148: MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatBrussel_Shift);
149: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatBrussel_GetDiagonal);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Create the eigensolver and set various options
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: /*
156: Create eigensolver context
157: */
158: EPSCreate(PETSC_COMM_WORLD,&eps);
160: /*
161: Set operators. In this case, it is a standard eigenvalue problem
162: */
163: EPSSetOperators(eps,A,PETSC_NULL);
164: EPSSetProblemType(eps,EPS_NHEP);
166: /*
167: Ask for the rightmost eigenvalues
168: */
169: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
171: /*
172: Set other solver options at runtime
173: */
174: EPSSetFromOptions(eps);
176: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177: Solve the eigensystem
178: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: EPSSolve(eps);
181: EPSGetIterationNumber(eps, &its);
182: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
183:
184: /*
185: Optional: Get some information from the solver and display it
186: */
187: EPSGetType(eps,&type);
188: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
189: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
190: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
191: EPSGetTolerances(eps,&tol,&maxit);
192: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Display solution and clean up
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198: /*
199: Get number of converged eigenpairs
200: */
201: EPSGetConverged(eps,&nconv);
202: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
204: if (nconv>0) {
205: /*
206: Display eigenvalues and relative errors
207: */
208: PetscPrintf(PETSC_COMM_WORLD,
209: " k ||Ax-kx||/||kx||\n"
210: " --------------------- ------------------\n" );
211: for( i=0; i<nconv; i++ ) {
212: /*
213: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
214: ki (imaginary part)
215: */
216: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
218: /*
219: Compute the relative error associated to each eigenpair
220: */
221: EPSComputeRelativeError(eps,i,&error);
223: #if defined(PETSC_USE_COMPLEX)
224: re = PetscRealPart(kr);
225: im = PetscImaginaryPart(kr);
226: #else
227: re = kr;
228: im = ki;
229: #endif
230: if( im != 0.0 ) {
231: PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
232: } else {
233: PetscPrintf(PETSC_COMM_WORLD," % 6f ",re);
234: }
235: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
236: }
237: PetscPrintf(PETSC_COMM_WORLD,"\n" );
238: }
239:
240: /*
241: Free work space
242: */
243: EPSDestroy(eps);
244: MatDestroy(A);
245: MatDestroy(ctx->T);
246: VecDestroy(ctx->x1);
247: VecDestroy(ctx->x2);
248: VecDestroy(ctx->y1);
249: VecDestroy(ctx->y2);
250: PetscFree(ctx);
251: SlepcFinalize();
252: return 0;
253: }
257: PetscErrorCode MatBrussel_Mult(Mat A,Vec x,Vec y)
258: {
260: PetscInt n;
261: PetscScalar *px, *py;
262: CTX_BRUSSEL *ctx;
265: MatShellGetContext(A,(void**)&ctx);
266: MatGetLocalSize(ctx->T,&n,PETSC_NULL);
267: VecGetArray(x,&px);
268: VecGetArray(y,&py);
269: VecPlaceArray(ctx->x1,px);
270: VecPlaceArray(ctx->x2,px+n);
271: VecPlaceArray(ctx->y1,py);
272: VecPlaceArray(ctx->y2,py+n);
274: MatMult(ctx->T,ctx->x1,ctx->y1);
275: VecScale(ctx->y1,ctx->tau1);
276: VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
277: VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);
279: MatMult(ctx->T,ctx->x2,ctx->y2);
280: VecScale(ctx->y2,ctx->tau2);
281: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
282: VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);
284: VecRestoreArray(x,&px);
285: VecRestoreArray(y,&py);
286: VecResetArray(ctx->x1);
287: VecResetArray(ctx->x2);
288: VecResetArray(ctx->y1);
289: VecResetArray(ctx->y2);
290: return(0);
291: }
295: PetscErrorCode MatBrussel_Shift(PetscScalar* a,Mat Y)
296: {
297: CTX_BRUSSEL *ctx;
301: MatShellGetContext( Y, (void**)&ctx );
302: ctx->sigma += *a;
303: return(0);
304: }
308: PetscErrorCode MatBrussel_GetDiagonal(Mat A,Vec diag)
309: {
310: Vec d1, d2;
312: PetscInt n;
313: PetscScalar *pd;
314: MPI_Comm comm;
315: CTX_BRUSSEL *ctx;
318: MatShellGetContext(A,(void**)&ctx);
319: PetscObjectGetComm((PetscObject)A,&comm);
320: MatGetLocalSize(ctx->T,&n,PETSC_NULL);
321: VecGetArray(diag,&pd);
322: VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd,&d1);
323: VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd+n,&d2);
325: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
326: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
328: VecDestroy(d1);
329: VecDestroy(d2);
330: VecRestoreArray(diag,&pd);
331: return(0);
332: }