Actual source code: ex1.c
petsc-3.4.2 2013-07-02
2: static char help[] = "Nonlinear Reaction Problem from Chemistry.\n" ;
This directory contains examples based on the PDES/ODES given in the book
Numerical Solution of Time-Dependent Advection-Diffusion-Reaction Equations by
W. Hundsdorf and J.G. Verwer
Page 3, Section 1.1 Nonlinear Reaction Problems from Chemistry
\begin{eqnarray}
{U_1}_t - k U_1 U_2 & = & 0 \\
{U_2}_t - k U_1 U_2 & = & 0 \\
{U_3}_t - k U_1 U_2 & = & 0
\end{eqnarray}
Helpful runtime monitoring options:
-ts_view - prints information about the solver being used
-ts_monitor - prints the progess of the solver
-ts_adapt_monitor - prints the progress of the time-step adaptor
-ts_monitor_lg_timestep - plots the size of each timestep (at each time-step)
-ts_monitor_lg_solution - plots each component of the solution as a function of time (at each timestep)
-ts_monitor_lg_error - plots each component of the error in the solution as a function of time (at each timestep)
-draw_pause -2 - hold the plots a the end of the solution process, enter a mouse press in each window to end the process
-ts_monitor_lg_timestep -1 - plots the size of each timestep (at the end of the solution process)
-ts_monitor_lg_solution -1 - plots each component of the solution as a function of time (at the end of the solution process)
-ts_monitor_lg_error -1 - plots each component of the error in the solution as a function of time (at the end of the solution process)
-lg_indicate_data_points false - do NOT show the data points on the plots
-draw_save - save the timestep and solution plot as a .Gif image file
35: /*
36: Project: Generate a nicely formated HTML page using
37: 1) the HTML version of the source code and text in this file, use make html to generate the file ex1.c.html
38: 2) the images (Draw_XXX_0_0.Gif Draw_YYY_0_0.Gif Draw_ZZZ_1_0.Gif) and
39: 3) the text output (output.txt) generated by running the following commands.
40: 4) <iframe src="generated_topics.html" scrolling="no" frameborder="0" width=600 height=300></iframe>
42: rm -rf *.Gif
43: ./ex1 -ts_monitor_lg_error -1 -ts_monitor_lg_solution -1 -draw_pause -2 -ts_max_steps 100 -ts_monitor_lg_timestep -1 -draw_save -draw_virtual -ts_monitor -ts_adapt_monitor -ts_view > output.txt
45: For example something like
46: <!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN">
47: <html>
48: <head>
49: <meta http-equiv="content-type" content="text/html;charset=utf-8">
50: <title>PETSc Example -- Nonlinear Reaction Problem from Chemistry</title>
51: </head>
52: <body>
53: <iframe src="ex1.c.html" scrolling="yes" frameborder="1" width=2000 height=400></iframe>
54: <img alt="" src="Draw_0x84000000_0_0.Gif"/><img alt="" src="Draw_0x84000001_0_0.Gif"/><img alt="" src="Draw_0x84000001_1_0.Gif"/>
55: <iframe src="output.txt" scrolling="yes" frameborder="1" width=2000 height=1000></iframe>
56: </body>
57: </html>
59: */
61: /*
62: Include "petscts.h" so that we can use TS solvers. Note that this
63: file automatically includes:
64: petscsys.h - base PETSc routines petscvec.h - vectors
65: petscmat.h - matrices
66: petscis.h - index sets petscksp.h - Krylov subspace methods
67: petscviewer.h - viewers petscpc.h - preconditioners
68: petscksp.h - linear solvers
69: */
70: #include <petscts.h>
72: typedef struct {
73: PetscScalar k;
74: Vec initialsolution;
75: } AppCtx;
79: PetscErrorCode IFunctionView(AppCtx *ctx,PetscViewer v)
80: {
84: PetscViewerBinaryWrite (v,&ctx->k,1,PETSC_SCALAR,PETSC_FALSE );
85: return (0);
86: }
90: PetscErrorCode IFunctionLoad(AppCtx **ctx,PetscViewer v)
91: {
95: PetscMalloc (sizeof (AppCtx),ctx);
96: PetscViewerBinaryRead (v,&(*ctx)->k,1,PETSC_SCALAR);
97: return (0);
98: }
102: /*
103: Defines the ODE passed to the ODE solver
104: */
105: PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
106: {
108: PetscScalar *u,*udot,*f;
111: /* The next three lines allow us to access the entries of the vectors directly */
112: VecGetArray (U,&u);
113: VecGetArray (Udot,&udot);
114: VecGetArray (F,&f);
115: f[0] = udot[0] + ctx->k*u[0]*u[1];
116: f[1] = udot[1] + ctx->k*u[0]*u[1];
117: f[2] = udot[2] - ctx->k*u[0]*u[1];
118: VecRestoreArray (U,&u);
119: VecRestoreArray (Udot,&udot);
120: VecRestoreArray (F,&f);
121: return (0);
122: }
126: /*
127: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian () for the meaning of a and the Jacobian.
128: */
129: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat *A,Mat *B,MatStructure *flag,AppCtx *ctx)
130: {
132: PetscInt rowcol[] = {0,1,2};
133: PetscScalar *u,*udot,J[3][3];
136: VecGetArray (U,&u);
137: VecGetArray (Udot,&udot);
138: J[0][0] = a + ctx->k*u[1]; J[0][1] = ctx->k*u[0]; J[0][2] = 0.0;
139: J[1][0] = ctx->k*u[1]; J[1][1] = a + ctx->k*u[0]; J[1][2] = 0.0;
140: J[2][0] = -ctx->k*u[1]; J[2][1] = -ctx->k*u[0]; J[2][2] = a;
141: MatSetValues (*B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES );
142: VecRestoreArray (U,&u);
143: VecRestoreArray (Udot,&udot);
145: MatAssemblyBegin (*A,MAT_FINAL_ASSEMBLY);
146: MatAssemblyEnd (*A,MAT_FINAL_ASSEMBLY);
147: if (*A != *B) {
148: MatAssemblyBegin (*B,MAT_FINAL_ASSEMBLY);
149: MatAssemblyEnd (*B,MAT_FINAL_ASSEMBLY);
150: }
151: *flag = SAME_NONZERO_PATTERN;
152: return (0);
153: }
157: /*
158: Defines the exact (analytic) solution to the ODE
159: */
160: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
161: {
162: PetscErrorCode ierr;
163: const PetscScalar *uinit;
164: PetscScalar *u,d0,q;
167: VecGetArrayRead (ctx->initialsolution,&uinit);
168: VecGetArray (U,&u);
169: d0 = uinit[0] - uinit[1];
170: if (d0 == 0.0) q = ctx->k*t;
171: else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0;
172: u[0] = uinit[0]/(1.0 + uinit[1]*q);
173: u[1] = u[0] - d0;
174: u[2] = uinit[1] + uinit[2] - u[1];
175: VecRestoreArray (U,&u);
176: VecRestoreArrayRead (ctx->initialsolution,&uinit);
177: return (0);
178: }
182: int main(int argc,char **argv)
183: {
184: TS ts; /* ODE integrator */
185: Vec U; /* solution will be stored here */
186: Mat A; /* Jacobian matrix */
188: PetscMPIInt size;
189: PetscInt n = 3;
190: AppCtx ctx;
191: PetscScalar *u;
193: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194: Initialize program
195: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196: PetscInitialize (&argc,&argv,(char*)0,help);
197: MPI_Comm_size (PETSC_COMM_WORLD ,&size);
198: if (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
200: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201: Create necessary matrix and vectors
202: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203: MatCreate (PETSC_COMM_WORLD ,&A);
204: MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
205: MatSetFromOptions (A);
206: MatSetUp (A);
208: MatGetVecs (A,&U,NULL);
210: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: Set runtime options
212: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213: PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Reaction options" ,"" );
214: {
215: ctx.k = .9;
216: PetscOptionsScalar ("-k" ,"Reaction coefficient" ,"" ,ctx.k,&ctx.k,NULL);
217: VecDuplicate (U,&ctx.initialsolution);
218: VecGetArray (ctx.initialsolution,&u);
219: u[0] = 1;
220: u[1] = .7;
221: u[2] = 0;
222: VecRestoreArray (ctx.initialsolution,&u);
223: PetscOptionsVec("-initial" ,"Initial values" ,"" ,ctx.initialsolution,NULL);
224: }
225: PetscOptionsEnd ();
227: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228: Create timestepping solver context
229: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230: TSCreate (PETSC_COMM_WORLD ,&ts);
231: TSSetProblemType (ts,TS_NONLINEAR);
232: TSSetType (ts,TSROSW );
233: TSSetIFunction (ts,NULL,(TSIFunction) IFunction,&ctx);
234: TSSetIJacobian (ts,A,A,(TSIJacobian)IJacobian,&ctx);
235: TSSetSolutionFunction (ts,(TSSolutionFunction)Solution,&ctx);
237: {
238: DM dm;
239: void *ptr;
240: TSGetDM (ts,&dm);
241: PetscDLSym (NULL,"IFunctionView" ,&ptr);
242: PetscDLSym (NULL,"IFunctionLoad" ,&ptr);
243: DMTSSetIFunctionSerialize (dm,(PetscErrorCode (*)(void*,PetscViewer ))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer ))IFunctionLoad);
244: DMTSSetIJacobianSerialize (dm,(PetscErrorCode (*)(void*,PetscViewer ))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer ))IFunctionLoad);
245: }
247: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248: Set initial conditions
249: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250: Solution(ts,0,U,&ctx);
251: TSSetSolution (ts,U);
253: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254: Set solver options
255: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
256: TSSetDuration (ts,1000,20.0);
257: TSSetInitialTimeStep (ts,0.0,.001);
258: TSSetFromOptions (ts);
260: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261: Solve nonlinear system
262: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
263: TSSolve (ts,U);
265: TSView (ts,PETSC_VIEWER_BINARY_WORLD );
267: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268: Free work space. All PETSc objects should be destroyed when they are no longer needed.
269: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
270: VecDestroy (&ctx.initialsolution);
271: MatDestroy (&A);
272: VecDestroy (&U);
273: TSDestroy (&ts);
275: PetscFinalize ();
276: return (0);
277: }