Actual source code: ex1.c


  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: */

 71: #include <petscts.h>

 73: typedef struct {
 74:   PetscScalar k;
 75:   Vec         initialsolution;
 76: } AppCtx;

 78: PetscErrorCode IFunctionView(AppCtx *ctx,PetscViewer v)
 79: {
 80:   PetscViewerBinaryWrite(v,&ctx->k,1,PETSC_SCALAR);
 81:   return 0;
 82: }

 84: PetscErrorCode IFunctionLoad(AppCtx **ctx,PetscViewer v)
 85: {
 86:   PetscNew(ctx);
 87:   PetscViewerBinaryRead(v,&(*ctx)->k,1,NULL,PETSC_SCALAR);
 88:   return 0;
 89: }

 91: /*
 92:      Defines the ODE passed to the ODE solver
 93: */
 94: PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 95: {
 96:   PetscScalar       *f;
 97:   const PetscScalar *u,*udot;

 99:   /*  The next three lines allow us to access the entries of the vectors directly */
100:   VecGetArrayRead(U,&u);
101:   VecGetArrayRead(Udot,&udot);
102:   VecGetArrayWrite(F,&f);
103:   f[0] = udot[0] + ctx->k*u[0]*u[1];
104:   f[1] = udot[1] + ctx->k*u[0]*u[1];
105:   f[2] = udot[2] - ctx->k*u[0]*u[1];
106:   VecRestoreArrayRead(U,&u);
107:   VecRestoreArrayRead(Udot,&udot);
108:   VecRestoreArrayWrite(F,&f);
109:   return 0;
110: }

112: /*
113:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
114: */
115: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
116: {
117:   PetscInt          rowcol[] = {0,1,2};
118:   PetscScalar       J[3][3];
119:   const PetscScalar *u,*udot;

121:   VecGetArrayRead(U,&u);
122:   VecGetArrayRead(Udot,&udot);
123:   J[0][0] = a + ctx->k*u[1];   J[0][1] = ctx->k*u[0];       J[0][2] = 0.0;
124:   J[1][0] = ctx->k*u[1];       J[1][1] = a + ctx->k*u[0];   J[1][2] = 0.0;
125:   J[2][0] = -ctx->k*u[1];      J[2][1] = -ctx->k*u[0];      J[2][2] = a;
126:   MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES);
127:   VecRestoreArrayRead(U,&u);
128:   VecRestoreArrayRead(Udot,&udot);

130:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
131:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
132:   if (A != B) {
133:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
134:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
135:   }
136:   return 0;
137: }

139: /*
140:      Defines the exact (analytic) solution to the ODE
141: */
142: static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
143: {
144:   const PetscScalar *uinit;
145:   PetscScalar       *u,d0,q;

147:   VecGetArrayRead(ctx->initialsolution,&uinit);
148:   VecGetArrayWrite(U,&u);
149:   d0   = uinit[0] - uinit[1];
150:   if (d0 == 0.0) q = ctx->k*t;
151:   else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0;
152:   u[0] = uinit[0]/(1.0 + uinit[1]*q);
153:   u[1] = u[0] - d0;
154:   u[2] = uinit[1] + uinit[2] - u[1];
155:   VecRestoreArrayWrite(U,&u);
156:   VecRestoreArrayRead(ctx->initialsolution,&uinit);
157:   return 0;
158: }

160: int main(int argc,char **argv)
161: {
162:   TS             ts;            /* ODE integrator */
163:   Vec            U;             /* solution will be stored here */
164:   Mat            A;             /* Jacobian matrix */
165:   PetscMPIInt    size;
166:   PetscInt       n = 3;
167:   AppCtx         ctx;
168:   PetscScalar    *u;
169:   const char     * const names[] = {"U1","U2","U3",NULL};

171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Initialize program
173:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174:   PetscInitialize(&argc,&argv,(char*)0,help);
175:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

178:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179:     Create necessary matrix and vectors
180:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181:   MatCreate(PETSC_COMM_WORLD,&A);
182:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
183:   MatSetFromOptions(A);
184:   MatSetUp(A);

186:   MatCreateVecs(A,&U,NULL);

188:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189:     Set runtime options
190:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191:   ctx.k = .9;
192:   PetscOptionsGetScalar(NULL,NULL,"-k",&ctx.k,NULL);
193:   VecDuplicate(U,&ctx.initialsolution);
194:   VecGetArrayWrite(ctx.initialsolution,&u);
195:   u[0]  = 1;
196:   u[1]  = .7;
197:   u[2]  = 0;
198:   VecRestoreArrayWrite(ctx.initialsolution,&u);
199:   PetscOptionsGetVec(NULL,NULL,"-initial",ctx.initialsolution,NULL);

201:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202:      Create timestepping solver context
203:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204:   TSCreate(PETSC_COMM_WORLD,&ts);
205:   TSSetProblemType(ts,TS_NONLINEAR);
206:   TSSetType(ts,TSROSW);
207:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
208:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
209:   TSSetSolutionFunction(ts,(TSSolutionFunction)Solution,&ctx);

211:   {
212:     DM   dm;
213:     void *ptr;
214:     TSGetDM(ts,&dm);
215:     PetscDLSym(NULL,"IFunctionView",&ptr);
216:     PetscDLSym(NULL,"IFunctionLoad",&ptr);
217:     DMTSSetIFunctionSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad);
218:     DMTSSetIJacobianSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad);
219:   }

221:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222:      Set initial conditions
223:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224:   Solution(ts,0,U,&ctx);
225:   TSSetSolution(ts,U);

227:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228:      Set solver options
229:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230:   TSSetTimeStep(ts,.001);
231:   TSSetMaxSteps(ts,1000);
232:   TSSetMaxTime(ts,20.0);
233:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
234:   TSSetFromOptions(ts);
235:   TSMonitorLGSetVariableNames(ts,names);

237:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238:      Solve nonlinear system
239:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240:   TSSolve(ts,U);

242:   TSView(ts,PETSC_VIEWER_BINARY_WORLD);

244:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
246:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
247:   VecDestroy(&ctx.initialsolution);
248:   MatDestroy(&A);
249:   VecDestroy(&U);
250:   TSDestroy(&ts);

252:   PetscFinalize();
253:   return 0;
254: }

256: /*TEST

258:    test:
259:      args: -ts_view
260:      requires: dlsym defined(PETSC_HAVE_DYNAMIC_LIBRARIES)

262:    test:
263:      suffix: 2
264:      args: -ts_monitor_lg_error -ts_monitor_lg_solution  -ts_view
265:      requires: x dlsym defined(PETSC_HAVE_DYNAMIC_LIBRARIES)
266:      output_file: output/ex1_1.out

268: TEST*/