xref: /petsc/src/ts/tutorials/extchemfield.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254) !
1 static const char help[] = "Integrate chemistry using TChem.\n";
2 
3 #include <petscts.h>
4 #include <petscdmda.h>
5 
6 #if defined(PETSC_HAVE_TCHEM)
7 #if defined(MAX)
8 #undef MAX
9 #endif
10 #if defined(MIN)
11 #undef MIN
12 #endif
13 #  include <TC_params.h>
14 #  include <TC_interface.h>
15 #else
16 #  error TChem is required for this example.  Reconfigure PETSc using --download-tchem.
17 #endif
18 /*
19 
20     This is an extension of extchem.c to solve the reaction equations independently in each cell of a one dimensional field
21 
22     Obtain the three files into this directory
23 
24        curl http://combustion.berkeley.edu/gri_mech/version30/files30/grimech30.dat > chem.inp
25        curl http://combustion.berkeley.edu/gri_mech/version30/files30/thermo30.dat > therm.dat
26        cp $PETSC_DIR/$PETSC_ARCH/externalpackages/tchem/data/periodictable.dat .
27 
28     Run with
29      ./extchemfield  -ts_arkimex_fully_implicit -ts_max_snes_failures -1 -ts_adapt_monitor -ts_adapt_dt_max 1e-4 -ts_arkimex_type 4 -ts_max_time .005
30 
31      Options for visualizing the solution:
32         Watch certain variables in each cell evolve with time
33         -draw_solution 1 -ts_monitor_lg_solution_variables Temp,H2,O2,H2O,CH4,CO,CO2,C2H2,N2 -lg_use_markers false  -draw_pause -2
34 
35         Watch certain variables in all cells evolve with time
36         -da_refine 4 -ts_monitor_draw_solution -draw_fields_by_name Temp,H2 -draw_vec_mark_points  -draw_pause -2
37 
38         Keep the initial temperature distribution as one monitors the current temperature distribution
39         -ts_monitor_draw_solution_initial -draw_bounds .9,1.7 -draw_fields_by_name Temp
40 
41         Save the images in a .gif (movie) file
42         -draw_save -draw_save_single_file
43 
44         Compute the sensitivies of the solution of the first temperature on the initial conditions
45         -ts_adjoint_solve  -ts_dt 1.e-5 -ts_type cn -ts_adjoint_view_solution draw
46 
47         Turn off diffusion
48         -diffusion no
49 
50         Turn off reactions
51         -reactions no
52 
53     The solution for component i = 0 is the temperature.
54 
55     The solution, i > 0, is the mass fraction, massf[i], of species i, i.e. mass of species i/ total mass of all species
56 
57     The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species
58         Define M[i] = mass per mole of species i then
59         molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j]))
60 
61     FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species.
62 
63 */
64 typedef struct _User *User;
65 struct _User {
66   PetscReal pressure;
67   int       Nspec;
68   int       Nreac;
69   PetscReal Tini,dx;
70   PetscReal diffus;
71   DM        dm;
72   PetscBool diffusion,reactions;
73   double    *tchemwork;
74   double    *Jdense;        /* Dense array workspace where Tchem computes the Jacobian */
75   PetscInt  *rows;
76 };
77 
78 static PetscErrorCode MonitorCell(TS,User,PetscInt);
79 static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
80 static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
81 static PetscErrorCode FormInitialSolution(TS,Vec,void*);
82 
83 #define CHKERRTC(ierr) do {PetscCheck(!ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0)
84 
85 int main(int argc,char **argv)
86 {
87   TS                ts;         /* time integrator */
88   TSAdapt           adapt;
89   Vec               X;          /* solution vector */
90   Mat               J;          /* Jacobian matrix */
91   PetscInt          steps,ncells,xs,xm,i;
92   PetscErrorCode    ierr;
93   PetscReal         ftime,dt;
94   char              chemfile[PETSC_MAX_PATH_LEN] = "chem.inp",thermofile[PETSC_MAX_PATH_LEN] = "therm.dat";
95   struct _User      user;
96   TSConvergedReason reason;
97   PetscBool         showsolutions = PETSC_FALSE;
98   char              **snames,*names;
99   Vec               lambda;     /* used with TSAdjoint for sensitivities */
100 
101   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
102   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");CHKERRQ(ierr);
103   CHKERRQ(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL));
104   CHKERRQ(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL));
105   user.pressure = 1.01325e5;    /* Pascal */
106   CHKERRQ(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL));
107   user.Tini   = 1550;
108   CHKERRQ(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL));
109   user.diffus = 100;
110   CHKERRQ(PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL));
111   CHKERRQ(PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL));
112   user.diffusion = PETSC_TRUE;
113   CHKERRQ(PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL));
114   user.reactions = PETSC_TRUE;
115   CHKERRQ(PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL));
116   ierr = PetscOptionsEnd();CHKERRQ(ierr);
117 
118   CHKERRTC(TC_initChem(chemfile, thermofile, 0, 1.0));
119   user.Nspec = TC_getNspec();
120   user.Nreac = TC_getNreac();
121 
122   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm));
123   CHKERRQ(DMSetFromOptions(user.dm));
124   CHKERRQ(DMSetUp(user.dm));
125   CHKERRQ(DMDAGetInfo(user.dm,NULL,&ncells,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
126   user.dx = 1.0/ncells;  /* Set the coordinates of the cell centers; note final ghost cell is at x coordinate 1.0 */
127   CHKERRQ(DMDASetUniformCoordinates(user.dm,0.0,1.0,0.0,1.0,0.0,1.0));
128 
129   /* set the names of each field in the DMDA based on the species name */
130   CHKERRQ(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names));
131   CHKERRQ(PetscStrcpy(names,"Temp"));
132   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
133   CHKERRQ(PetscMalloc1((user.Nspec+2),&snames));
134   for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME;
135   snames[user.Nspec+1] = NULL;
136   CHKERRQ(DMDASetFieldNames(user.dm,(const char * const *)snames));
137   CHKERRQ(PetscFree(snames));
138   CHKERRQ(PetscFree(names));
139 
140   CHKERRQ(DMCreateMatrix(user.dm,&J));
141   CHKERRQ(DMCreateGlobalVector(user.dm,&X));
142 
143   CHKERRQ(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows));
144 
145   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146      Create timestepping solver context
147      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
149   CHKERRQ(TSSetDM(ts,user.dm));
150   CHKERRQ(TSSetType(ts,TSARKIMEX));
151   CHKERRQ(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
152   CHKERRQ(TSARKIMEXSetType(ts,TSARKIMEX4));
153   CHKERRQ(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
154   CHKERRQ(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user));
155 
156   ftime = 1.0;
157   CHKERRQ(TSSetMaxTime(ts,ftime));
158   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
159 
160   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161      Set initial conditions
162    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163   CHKERRQ(FormInitialSolution(ts,X,&user));
164   CHKERRQ(TSSetSolution(ts,X));
165   dt   = 1e-10;                 /* Initial time step */
166   CHKERRQ(TSSetTimeStep(ts,dt));
167   CHKERRQ(TSGetAdapt(ts,&adapt));
168   CHKERRQ(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
169   CHKERRQ(TSSetMaxSNESFailures(ts,-1));            /* Retry step an unlimited number of times */
170 
171   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172      Pass information to graphical monitoring routine
173    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174   if (showsolutions) {
175     CHKERRQ(DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL));
176     for (i=xs;i<xs+xm;i++) {
177       CHKERRQ(MonitorCell(ts,&user,i));
178     }
179   }
180 
181   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182      Set runtime options
183    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184   CHKERRQ(TSSetFromOptions(ts));
185 
186   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187      Set final conditions for sensitivities
188    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189   CHKERRQ(DMCreateGlobalVector(user.dm,&lambda));
190   CHKERRQ(TSSetCostGradients(ts,1,&lambda,NULL));
191   CHKERRQ(VecSetValue(lambda,0,1.0,INSERT_VALUES));
192   CHKERRQ(VecAssemblyBegin(lambda));
193   CHKERRQ(VecAssemblyEnd(lambda));
194 
195   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196      Solve ODE
197      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198   CHKERRQ(TSSolve(ts,X));
199   CHKERRQ(TSGetSolveTime(ts,&ftime));
200   CHKERRQ(TSGetStepNumber(ts,&steps));
201   CHKERRQ(TSGetConvergedReason(ts,&reason));
202   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
203 
204   {
205     Vec                max;
206     const char * const *names;
207     PetscInt           i;
208     const PetscReal    *bmax;
209 
210     CHKERRQ(TSMonitorEnvelopeGetBounds(ts,&max,NULL));
211     if (max) {
212       CHKERRQ(TSMonitorLGGetVariableNames(ts,&names));
213       if (names) {
214         CHKERRQ(VecGetArrayRead(max,&bmax));
215         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n"));
216         for (i=1; i<user.Nspec; i++) {
217           if (bmax[i] > .01) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i]));
218         }
219         CHKERRQ(VecRestoreArrayRead(max,&bmax));
220       }
221     }
222   }
223 
224   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225      Free work space.
226    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227   TC_reset();
228   CHKERRQ(DMDestroy(&user.dm));
229   CHKERRQ(MatDestroy(&J));
230   CHKERRQ(VecDestroy(&X));
231   CHKERRQ(VecDestroy(&lambda));
232   CHKERRQ(TSDestroy(&ts));
233   CHKERRQ(PetscFree3(user.tchemwork,user.Jdense,user.rows));
234   CHKERRQ(PetscFinalize());
235   return 0;
236 }
237 
238 /*
239    Applies the second order centered difference diffusion operator on a one dimensional periodic domain
240 */
241 static PetscErrorCode FormDiffusionFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
242 {
243   User              user = (User)ptr;
244   PetscErrorCode    ierr;
245   PetscScalar       **f;
246   const PetscScalar **x;
247   DM                dm;
248   PetscInt          i,xs,xm,j,dof;
249   Vec               Xlocal;
250   PetscReal         idx;
251 
252   PetscFunctionBeginUser;
253   CHKERRQ(TSGetDM(ts,&dm));
254   CHKERRQ(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
255   CHKERRQ(DMGetLocalVector(dm,&Xlocal));
256   CHKERRQ(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal));
257   CHKERRQ(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal));
258   CHKERRQ(DMDAVecGetArrayDOFRead(dm,Xlocal,&x));
259   CHKERRQ(DMDAVecGetArrayDOF(dm,F,&f));
260   CHKERRQ(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
261 
262   idx = 1.0*user->diffus/user->dx;
263   for (i=xs; i<xs+xm; i++) {
264     for (j=0; j<dof; j++) {
265       f[i][j] += idx*(x[i+1][j] - 2.0*x[i][j] + x[i-1][j]);
266     }
267   }
268   CHKERRQ(DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x));
269   CHKERRQ(DMDAVecRestoreArrayDOF(dm,F,&f));
270   CHKERRQ(DMRestoreLocalVector(dm,&Xlocal));
271   PetscFunctionReturn(0);
272 }
273 
274 /*
275    Produces the second order centered difference diffusion operator on a one dimensional periodic domain
276 */
277 static PetscErrorCode FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
278 {
279   User              user = (User)ptr;
280   PetscErrorCode    ierr;
281   DM                dm;
282   PetscInt          i,xs,xm,j,dof;
283   PetscReal         idx,values[3];
284   MatStencil        row,col[3];
285 
286   PetscFunctionBeginUser;
287   CHKERRQ(TSGetDM(ts,&dm));
288   CHKERRQ(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
289   CHKERRQ(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
290 
291   idx = 1.0*user->diffus/user->dx;
292   values[0] = idx;
293   values[1] = -2.0*idx;
294   values[2] = idx;
295   for (i=xs; i<xs+xm; i++) {
296     for (j=0; j<dof; j++) {
297       row.i = i;      row.c = j;
298       col[0].i = i-1; col[0].c = j;
299       col[1].i = i;   col[1].c = j;
300       col[2].i = i+1; col[2].c = j;
301       CHKERRQ(MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES));
302     }
303   }
304   CHKERRQ(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
305   CHKERRQ(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
306   PetscFunctionReturn(0);
307 }
308 
309 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
310 {
311   User              user = (User)ptr;
312   PetscErrorCode    ierr;
313   PetscScalar       **f;
314   const PetscScalar **x;
315   DM                dm;
316   PetscInt          i,xs,xm;
317 
318   PetscFunctionBeginUser;
319   if (user->reactions) {
320     CHKERRQ(TSGetDM(ts,&dm));
321     CHKERRQ(DMDAVecGetArrayDOFRead(dm,X,&x));
322     CHKERRQ(DMDAVecGetArrayDOF(dm,F,&f));
323     CHKERRQ(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
324 
325     for (i=xs; i<xs+xm; i++) {
326       CHKERRQ(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1));
327       user->tchemwork[0] *= user->Tini; /* Dimensionalize */
328       CHKERRTC(TC_getSrc(user->tchemwork,user->Nspec+1,f[i]));
329       f[i][0] /= user->Tini;           /* Non-dimensionalize */
330     }
331 
332     CHKERRQ(DMDAVecRestoreArrayDOFRead(dm,X,&x));
333     CHKERRQ(DMDAVecRestoreArrayDOF(dm,F,&f));
334   } else {
335     CHKERRQ(VecZeroEntries(F));
336   }
337   if (user->diffusion) {
338     CHKERRQ(FormDiffusionFunction(ts,t,X,F,ptr));
339   }
340   PetscFunctionReturn(0);
341 }
342 
343 static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
344 {
345   User              user = (User)ptr;
346   PetscErrorCode    ierr;
347   const PetscScalar **x;
348   PetscInt          M = user->Nspec+1,i,j,xs,xm;
349   DM                dm;
350 
351   PetscFunctionBeginUser;
352   if (user->reactions) {
353     CHKERRQ(TSGetDM(ts,&dm));
354     CHKERRQ(MatZeroEntries(Pmat));
355     CHKERRQ(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE));
356     CHKERRQ(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
357     CHKERRQ(DMDAVecGetArrayDOFRead(dm,X,&x));
358     CHKERRQ(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
359 
360     for (i=xs; i<xs+xm; i++) {
361       CHKERRQ(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1));
362       user->tchemwork[0] *= user->Tini;  /* Dimensionalize temperature (first row) because that is what Tchem wants */
363       CHKERRQ(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1));
364 
365       for (j=0; j<M; j++) user->Jdense[j + 0*M] /= user->Tini; /* Non-dimensionalize first column */
366       for (j=0; j<M; j++) user->Jdense[0 + j*M] /= user->Tini; /* Non-dimensionalize first row */
367       for (j=0; j<M; j++) user->rows[j] = i*M+j;
368       CHKERRQ(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES));
369     }
370     CHKERRQ(DMDAVecRestoreArrayDOFRead(dm,X,&x));
371     CHKERRQ(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
372     CHKERRQ(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
373   } else {
374     CHKERRQ(MatZeroEntries(Pmat));
375   }
376   if (user->diffusion) {
377     CHKERRQ(FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr));
378   }
379   if (Amat != Pmat) {
380     CHKERRQ(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY));
381     CHKERRQ(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY));
382   }
383   PetscFunctionReturn(0);
384 }
385 
386 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
387 {
388   PetscScalar    **x,*xc;
389   PetscErrorCode ierr;
390   struct {const char *name; PetscReal massfrac;} initial[] = {
391     {"CH4", 0.0948178320887},
392     {"O2", 0.189635664177},
393     {"N2", 0.706766236705},
394     {"AR", 0.00878026702874}
395   };
396   PetscInt       i,j,xs,xm;
397   DM             dm;
398 
399   PetscFunctionBeginUser;
400   CHKERRQ(VecZeroEntries(X));
401   CHKERRQ(TSGetDM(ts,&dm));
402   CHKERRQ(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
403 
404   CHKERRQ(DMDAGetCoordinateArray(dm,&xc));
405   CHKERRQ(DMDAVecGetArrayDOF(dm,X,&x));
406   for (i=xs; i<xs+xm; i++) {
407     x[i][0] = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[i]);  /* Non-dimensionalized by user->Tini */
408     for (j=0; j<sizeof(initial)/sizeof(initial[0]); j++) {
409       int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name));
410       PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",initial[j].name);
411       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",j,initial[j].name,initial[j].massfrac));
412       x[i][1+ispec] = initial[j].massfrac;
413     }
414   }
415   CHKERRQ(DMDAVecRestoreArrayDOF(dm,X,&x));
416   CHKERRQ(DMDARestoreCoordinateArray(dm,&xc));
417   PetscFunctionReturn(0);
418 }
419 
420 /*
421     Routines for displaying the solutions
422 */
423 typedef struct {
424   PetscInt cell;
425   User     user;
426 } UserLGCtx;
427 
428 static PetscErrorCode FormMoleFraction(UserLGCtx *ctx,Vec massf,Vec *molef)
429 {
430   User              user = ctx->user;
431   PetscErrorCode    ierr;
432   PetscReal         *M,tM=0;
433   PetscInt          i,n = user->Nspec+1;
434   PetscScalar       *mof;
435   const PetscScalar **maf;
436 
437   PetscFunctionBegin;
438   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,molef));
439   CHKERRQ(PetscMalloc1(user->Nspec,&M));
440   TC_getSmass(user->Nspec, M);
441   CHKERRQ(DMDAVecGetArrayDOFRead(user->dm,massf,&maf));
442   CHKERRQ(VecGetArray(*molef,&mof));
443   mof[0] = maf[ctx->cell][0]; /* copy over temperature */
444   for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1];
445   for (i=1; i<n; i++) {
446     mof[i] = maf[ctx->cell][i]/(M[i-1]*tM);
447   }
448   CHKERRQ(DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf));
449   CHKERRQ(VecRestoreArray(*molef,&mof));
450   CHKERRQ(PetscFree(M));
451   PetscFunctionReturn(0);
452 }
453 
454 static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx)
455 {
456   PetscErrorCode ierr;
457 
458   PetscFunctionBegin;
459   CHKERRQ(PetscFree(uctx));
460   PetscFunctionReturn(0);
461 }
462 
463 /*
464    Use TSMonitorLG to monitor the reactions in a particular cell
465 */
466 static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell)
467 {
468   PetscErrorCode ierr;
469   TSMonitorLGCtx ctx;
470   char           **snames;
471   UserLGCtx      *uctx;
472   char           label[128];
473   PetscReal      temp,*xc;
474   PetscMPIInt    rank;
475 
476   PetscFunctionBegin;
477   CHKERRQ(DMDAGetCoordinateArray(user->dm,&xc));
478   temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]);  /* Non-dimensionalized by user->Tini */
479   CHKERRQ(DMDARestoreCoordinateArray(user->dm,&xc));
480   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
481   CHKERRQ(PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank));
482   CHKERRQ(TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx));
483   CHKERRQ(DMDAGetFieldNames(user->dm,(const char * const **)&snames));
484   CHKERRQ(TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames));
485   CHKERRQ(PetscNew(&uctx));
486   uctx->cell = cell;
487   uctx->user = user;
488   CHKERRQ(TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx));
489   CHKERRQ(TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy));
490   PetscFunctionReturn(0);
491 }
492