Actual source code: ex1.c
petsc-3.5.4 2015-05-23
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,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: return (0);
152: }
156: /*
157: Defines the exact (analytic) solution to the ODE
158: */
159: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
160: {
161: PetscErrorCode ierr;
162: const PetscScalar *uinit;
163: PetscScalar *u,d0,q;
166: VecGetArrayRead (ctx->initialsolution,&uinit);
167: VecGetArray (U,&u);
168: d0 = uinit[0] - uinit[1];
169: if (d0 == 0.0) q = ctx->k*t;
170: else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0;
171: u[0] = uinit[0]/(1.0 + uinit[1]*q);
172: u[1] = u[0] - d0;
173: u[2] = uinit[1] + uinit[2] - u[1];
174: VecRestoreArray (U,&u);
175: VecRestoreArrayRead (ctx->initialsolution,&uinit);
176: return (0);
177: }
181: int main(int argc,char **argv)
182: {
183: TS ts; /* ODE integrator */
184: Vec U; /* solution will be stored here */
185: Mat A; /* Jacobian matrix */
187: PetscMPIInt size;
188: PetscInt n = 3;
189: AppCtx ctx;
190: PetscScalar *u;
192: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: Initialize program
194: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195: PetscInitialize (&argc,&argv,(char*)0,help);
196: MPI_Comm_size (PETSC_COMM_WORLD ,&size);
197: if (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
199: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200: Create necessary matrix and vectors
201: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202: MatCreate (PETSC_COMM_WORLD ,&A);
203: MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
204: MatSetFromOptions (A);
205: MatSetUp (A);
207: MatGetVecs (A,&U,NULL);
209: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210: Set runtime options
211: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212: PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Reaction options" ,"" );
213: {
214: ctx.k = .9;
215: PetscOptionsScalar ("-k" ,"Reaction coefficient" ,"" ,ctx.k,&ctx.k,NULL);
216: VecDuplicate (U,&ctx.initialsolution);
217: VecGetArray (ctx.initialsolution,&u);
218: u[0] = 1;
219: u[1] = .7;
220: u[2] = 0;
221: VecRestoreArray (ctx.initialsolution,&u);
222: PetscOptionsVec("-initial" ,"Initial values" ,"" ,ctx.initialsolution,NULL);
223: }
224: PetscOptionsEnd ();
226: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227: Create timestepping solver context
228: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229: TSCreate (PETSC_COMM_WORLD ,&ts);
230: TSSetProblemType (ts,TS_NONLINEAR);
231: TSSetType (ts,TSROSW );
232: TSSetIFunction (ts,NULL,(TSIFunction) IFunction,&ctx);
233: TSSetIJacobian (ts,A,A,(TSIJacobian)IJacobian,&ctx);
234: TSSetSolutionFunction (ts,(TSSolutionFunction)Solution,&ctx);
236: {
237: DM dm;
238: void *ptr;
239: TSGetDM (ts,&dm);
240: PetscDLSym (NULL,"IFunctionView" ,&ptr);
241: PetscDLSym (NULL,"IFunctionLoad" ,&ptr);
242: DMTSSetIFunctionSerialize (dm,(PetscErrorCode (*)(void*,PetscViewer ))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer ))IFunctionLoad);
243: DMTSSetIJacobianSerialize (dm,(PetscErrorCode (*)(void*,PetscViewer ))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer ))IFunctionLoad);
244: }
246: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247: Set initial conditions
248: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
249: Solution(ts,0,U,&ctx);
250: TSSetSolution (ts,U);
252: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253: Set solver options
254: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255: TSSetDuration (ts,1000,20.0);
256: TSSetInitialTimeStep (ts,0.0,.001);
257: TSSetFromOptions (ts);
259: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260: Solve nonlinear system
261: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262: TSSolve (ts,U);
264: TSView (ts,PETSC_VIEWER_BINARY_WORLD );
266: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267: Free work space. All PETSc objects should be destroyed when they are no longer needed.
268: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269: VecDestroy (&ctx.initialsolution);
270: MatDestroy (&A);
271: VecDestroy (&U);
272: TSDestroy (&ts);
274: PetscFinalize ();
275: return (0);
276: }