Actual source code: ex9.c
2: static char help[] = "Solves a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
3: "The command line options are:\n"
4: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
5: " -L <L>, where <L> = bifurcation parameter.\n"
6: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
7: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
9: #include slepceps.h
11: /*
12: This example computes the eigenvalues with largest real part of the
13: following matrix
15: A = [ tau1*T+(beta-1)*I alpha^2*I
16: -beta*I tau2*T-alpha^2*I ],
18: where
20: T = tridiag{1,-2,1}
21: h = 1/(n+1)
22: tau1 = delta1/(h*L)^2
23: tau2 = delta2/(h*L)^2
24: */
27: /*
28: Matrix operations
29: */
30: PetscErrorCode MatBrussel_Mult(Mat,Vec,Vec);
31: PetscErrorCode MatBrussel_Shift(PetscScalar*,Mat);
32: PetscErrorCode MatBrussel_GetDiagonal(Mat,Vec);
34: typedef struct {
35: Mat T;
36: Vec x1, x2, y1, y2;
37: PetscScalar alpha, beta, tau1, tau2, sigma;
38: } CTX_BRUSSEL;
42: int main( int argc, char **argv )
43: {
44: Mat A; /* eigenvalue problem matrix */
45: EPS eps; /* eigenproblem solver context */
46: EPSType type;
47: PetscReal error, tol, re, im;
48: PetscScalar delta1, delta2, L, h, kr, ki, value[3];
49: PetscInt N=30, n, i, col[3], Istart, Iend;
50: int nev, maxit, its, nconv;
51: PetscTruth FirstBlock=PETSC_FALSE, LastBlock=PETSC_FALSE;
53: CTX_BRUSSEL *ctx;
55: SlepcInitialize(&argc,&argv,(char*)0,help);
57: PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);
58: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%d\n\n",N);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Generate the matrix
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: /*
65: Create shell matrix context and set default parameters
66: */
67: PetscNew(CTX_BRUSSEL,&ctx);
68: ctx->alpha = 2.0;
69: ctx->beta = 5.45;
70: delta1 = 0.008;
71: delta2 = 0.004;
72: L = 0.51302;
74: /*
75: Look the command line for user-provided parameters
76: */
77: PetscOptionsGetScalar(PETSC_NULL,"-L",&L,PETSC_NULL);
78: PetscOptionsGetScalar(PETSC_NULL,"-alpha",&ctx->alpha,PETSC_NULL);
79: PetscOptionsGetScalar(PETSC_NULL,"-beta",&ctx->beta,PETSC_NULL);
80: PetscOptionsGetScalar(PETSC_NULL,"-delta1",&delta1,PETSC_NULL);
81: PetscOptionsGetScalar(PETSC_NULL,"-delta2",&delta2,PETSC_NULL);
83: /*
84: Create matrix T
85: */
86: MatCreate(PETSC_COMM_WORLD,&ctx->T);
87: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
88: MatSetFromOptions(ctx->T);
89:
90: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
91: if (Istart==0) FirstBlock=PETSC_TRUE;
92: if (Iend==N) LastBlock=PETSC_TRUE;
93: value[0]=1.0; value[1]=-2.0; value[2]=1.0;
94: for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
95: col[0]=i-1; col[1]=i; col[2]=i+1;
96: MatSetValues(ctx->T,1,&i,3,col,value,INSERT_VALUES);
97: }
98: if (LastBlock) {
99: i=N-1; col[0]=N-2; col[1]=N-1;
100: MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
101: }
102: if (FirstBlock) {
103: i=0; col[0]=0; col[1]=1; value[0]=-2.0; value[1]=1.0;
104: MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);
105: }
107: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
108: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
109: MatGetLocalSize(ctx->T,&n,PETSC_NULL);
111: /*
112: Fill the remaining information in the shell matrix context
113: and create auxiliary vectors
114: */
115: h = 1.0 / (double)(N+1);
116: ctx->tau1 = delta1 / ((h*L)*(h*L));
117: ctx->tau2 = delta2 / ((h*L)*(h*L));
118: ctx->sigma = 0.0;
119: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x1);
120: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x2);
121: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y1);
122: VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y2);
124: /*
125: Create the shell matrix
126: */
127: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
128: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatBrussel_Mult);
129: MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatBrussel_Shift);
130: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatBrussel_GetDiagonal);
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Create the eigensolver and set various options
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: /*
137: Create eigensolver context
138: */
139: EPSCreate(PETSC_COMM_WORLD,&eps);
141: /*
142: Set operators. In this case, it is a standard eigenvalue problem
143: */
144: EPSSetOperators(eps,A,PETSC_NULL);
145: EPSSetProblemType(eps,EPS_NHEP);
147: /*
148: Ask for the rightmost eigenvalues
149: */
150: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
152: /*
153: Set other solver options at runtime
154: */
155: EPSSetFromOptions(eps);
157: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: Solve the eigensystem
159: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: EPSSolve(eps);
162: EPSGetIterationNumber(eps, &its);
163: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
164:
165: /*
166: Optional: Get some information from the solver and display it
167: */
168: EPSGetType(eps,&type);
169: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
170: EPSGetDimensions(eps,&nev,PETSC_NULL);
171: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
172: EPSGetTolerances(eps,&tol,&maxit);
173: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
175: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176: Display solution and clean up
177: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
179: /*
180: Get number of converged eigenpairs
181: */
182: EPSGetConverged(eps,&nconv);
183: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
185: if (nconv>0) {
186: /*
187: Display eigenvalues and relative errors
188: */
189: PetscPrintf(PETSC_COMM_WORLD,
190: " k ||Ax-kx||/||kx||\n"
191: " --------------------- ------------------\n" );
192: for( i=0; i<nconv; i++ ) {
193: /*
194: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
195: ki (imaginary part)
196: */
197: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
199: /*
200: Compute the relative error associated to each eigenpair
201: */
202: EPSComputeRelativeError(eps,i,&error);
204: #if defined(PETSC_USE_COMPLEX)
205: re = PetscRealPart(kr);
206: im = PetscImaginaryPart(kr);
207: #else
208: re = kr;
209: im = ki;
210: #endif
211: if( im != 0.0 ) {
212: PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);
213: } else {
214: PetscPrintf(PETSC_COMM_WORLD," % 6f ",re);
215: }
216: PetscPrintf(PETSC_COMM_WORLD," % 12g\n",error);
217: }
218: PetscPrintf(PETSC_COMM_WORLD,"\n" );
219: }
220:
221: /*
222: Free work space
223: */
224: EPSDestroy(eps);
225: MatDestroy(A);
226: MatDestroy(ctx->T);
227: VecDestroy(ctx->x1);
228: VecDestroy(ctx->x2);
229: VecDestroy(ctx->y1);
230: VecDestroy(ctx->y2);
231: PetscFree(ctx);
232: SlepcFinalize();
233: return 0;
234: }
238: PetscErrorCode MatBrussel_Mult(Mat A,Vec x,Vec y)
239: {
241: PetscInt n;
242: PetscScalar *px, *py;
243: CTX_BRUSSEL *ctx;
246: MatShellGetContext(A,(void**)&ctx);
247: MatGetLocalSize(ctx->T,&n,PETSC_NULL);
248: VecGetArray(x,&px);
249: VecGetArray(y,&py);
250: VecPlaceArray(ctx->x1,px);
251: VecPlaceArray(ctx->x2,px+n);
252: VecPlaceArray(ctx->y1,py);
253: VecPlaceArray(ctx->y2,py+n);
255: MatMult(ctx->T,ctx->x1,ctx->y1);
256: VecScale(ctx->y1,ctx->tau1);
257: VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
258: VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);
260: MatMult(ctx->T,ctx->x2,ctx->y2);
261: VecScale(ctx->y2,ctx->tau2);
262: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
263: VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);
265: VecRestoreArray(x,&px);
266: VecRestoreArray(y,&py);
267: VecResetArray(ctx->x1);
268: VecResetArray(ctx->x2);
269: VecResetArray(ctx->y1);
270: VecResetArray(ctx->y2);
271: return(0);
272: }
276: PetscErrorCode MatBrussel_Shift( PetscScalar* a, Mat Y )
277: {
278: CTX_BRUSSEL *ctx;
282: MatShellGetContext( Y, (void**)&ctx );
283: ctx->sigma += *a;
284: return(0);
285: }
289: int MatBrussel_GetDiagonal(Mat A,Vec diag)
290: {
291: Vec d1, d2;
293: PetscInt n;
294: PetscScalar *pd;
295: MPI_Comm comm;
296: CTX_BRUSSEL *ctx;
299: MatShellGetContext(A,(void**)&ctx);
300: PetscObjectGetComm((PetscObject)A,&comm);
301: MatGetLocalSize(ctx->T,&n,PETSC_NULL);
302: VecGetArray(diag,&pd);
303: VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd,&d1);
304: VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd+n,&d2);
306: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
307: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
309: VecDestroy(d1);
310: VecDestroy(d2);
311: VecRestoreArray(diag,&pd);
312: return(0);
313: }