Actual source code: ex1.c
petsc-3.8.4 2018-03-24
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_use_markers 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;
77: PetscErrorCode IFunctionView(AppCtx *ctx,PetscViewer v)
78: {
82: PetscViewerBinaryWrite (v,&ctx->k,1,PETSC_SCALAR,PETSC_FALSE );
83: return (0);
84: }
86: PetscErrorCode IFunctionLoad(AppCtx **ctx,PetscViewer v)
87: {
91: PetscNew (ctx);
92: PetscViewerBinaryRead (v,&(*ctx)->k,1,NULL,PETSC_SCALAR);
93: return (0);
94: }
96: /*
97: Defines the ODE passed to the ODE solver
98: */
99: PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
100: {
101: PetscErrorCode ierr;
102: PetscScalar *f;
103: const PetscScalar *u,*udot;
106: /* The next three lines allow us to access the entries of the vectors directly */
107: VecGetArrayRead (U,&u);
108: VecGetArrayRead (Udot,&udot);
109: VecGetArray (F,&f);
110: f[0] = udot[0] + ctx->k*u[0]*u[1];
111: f[1] = udot[1] + ctx->k*u[0]*u[1];
112: f[2] = udot[2] - ctx->k*u[0]*u[1];
113: VecRestoreArrayRead (U,&u);
114: VecRestoreArrayRead (Udot,&udot);
115: VecRestoreArray (F,&f);
116: return (0);
117: }
119: /*
120: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian () for the meaning of a and the Jacobian.
121: */
122: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
123: {
124: PetscErrorCode ierr;
125: PetscInt rowcol[] = {0,1,2};
126: PetscScalar J[3][3];
127: const PetscScalar *u,*udot;
130: VecGetArrayRead (U,&u);
131: VecGetArrayRead (Udot,&udot);
132: J[0][0] = a + ctx->k*u[1]; J[0][1] = ctx->k*u[0]; J[0][2] = 0.0;
133: J[1][0] = ctx->k*u[1]; J[1][1] = a + ctx->k*u[0]; J[1][2] = 0.0;
134: J[2][0] = -ctx->k*u[1]; J[2][1] = -ctx->k*u[0]; J[2][2] = a;
135: MatSetValues (B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES );
136: VecRestoreArrayRead (U,&u);
137: VecRestoreArrayRead (Udot,&udot);
139: MatAssemblyBegin (A,MAT_FINAL_ASSEMBLY );
140: MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY );
141: if (A != B) {
142: MatAssemblyBegin (B,MAT_FINAL_ASSEMBLY );
143: MatAssemblyEnd (B,MAT_FINAL_ASSEMBLY );
144: }
145: return (0);
146: }
148: /*
149: Defines the exact (analytic) solution to the ODE
150: */
151: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
152: {
153: PetscErrorCode ierr;
154: const PetscScalar *uinit;
155: PetscScalar *u,d0,q;
158: VecGetArrayRead (ctx->initialsolution,&uinit);
159: VecGetArray (U,&u);
160: d0 = uinit[0] - uinit[1];
161: if (d0 == 0.0) q = ctx->k*t;
162: else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0;
163: u[0] = uinit[0]/(1.0 + uinit[1]*q);
164: u[1] = u[0] - d0;
165: u[2] = uinit[1] + uinit[2] - u[1];
166: VecRestoreArray (U,&u);
167: VecRestoreArrayRead (ctx->initialsolution,&uinit);
168: return (0);
169: }
171: int main(int argc,char **argv)
172: {
173: TS ts; /* ODE integrator */
174: Vec U; /* solution will be stored here */
175: Mat A; /* Jacobian matrix */
177: PetscMPIInt size;
178: PetscInt n = 3;
179: AppCtx ctx;
180: PetscScalar *u;
181: const char * const names[] = {"U1" ,"U2" ,"U3" ,NULL};
183: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184: Initialize program
185: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186: PetscInitialize (&argc,&argv,(char*)0,help);if (ierr) return ierr;
187: MPI_Comm_size (PETSC_COMM_WORLD ,&size);
188: if (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
190: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191: Create necessary matrix and vectors
192: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193: MatCreate (PETSC_COMM_WORLD ,&A);
194: MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
195: MatSetFromOptions (A);
196: MatSetUp (A);
198: MatCreateVecs (A,&U,NULL);
200: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201: Set runtime options
202: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203: ctx.k = .9;
204: PetscOptionsGetScalar (NULL,NULL,"-k" ,&ctx.k,NULL);
205: VecDuplicate (U,&ctx.initialsolution);
206: VecGetArray (ctx.initialsolution,&u);
207: u[0] = 1;
208: u[1] = .7;
209: u[2] = 0;
210: VecRestoreArray (ctx.initialsolution,&u);
211: PetscOptionsGetVec(NULL,NULL,"-initial" ,ctx.initialsolution,NULL);
213: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: Create timestepping solver context
215: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216: TSCreate (PETSC_COMM_WORLD ,&ts);
217: TSSetProblemType (ts,TS_NONLINEAR );
218: TSSetType (ts,TSROSW );
219: TSSetIFunction (ts,NULL,(TSIFunction) IFunction,&ctx);
220: TSSetIJacobian (ts,A,A,(TSIJacobian)IJacobian,&ctx);
221: TSSetSolutionFunction (ts,(TSSolutionFunction)Solution,&ctx);
223: {
224: DM dm;
225: void *ptr;
226: TSGetDM (ts,&dm);
227: PetscDLSym (NULL,"IFunctionView" ,&ptr);
228: PetscDLSym (NULL,"IFunctionLoad" ,&ptr);
229: DMTSSetIFunctionSerialize (dm,(PetscErrorCode (*)(void*,PetscViewer ))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer ))IFunctionLoad);
230: DMTSSetIJacobianSerialize (dm,(PetscErrorCode (*)(void*,PetscViewer ))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer ))IFunctionLoad);
231: }
233: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234: Set initial conditions
235: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236: Solution(ts,0,U,&ctx);
237: TSSetSolution (ts,U);
239: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
240: Set solver options
241: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
242: TSSetTimeStep (ts,.001);
243: TSSetMaxSteps (ts,1000);
244: TSSetMaxTime (ts,20.0);
245: TSSetExactFinalTime (ts,TS_EXACTFINALTIME_STEPOVER );
246: TSSetFromOptions (ts);
247: TSMonitorLGSetVariableNames (ts,names);
249: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250: Solve nonlinear system
251: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
252: TSSolve (ts,U);
254: TSView (ts,PETSC_VIEWER_BINARY_WORLD );
256: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257: Free work space. All PETSc objects should be destroyed when they are no longer needed.
258: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
259: VecDestroy (&ctx.initialsolution);
260: MatDestroy (&A);
261: VecDestroy (&U);
262: TSDestroy (&ts);
264: PetscFinalize ();
265: return ierr;
266: }