xref: /petsc/src/ts/tutorials/extchemfield.c (revision 6aad120caa16b1027d343a5f30f73d01448e4dc0)
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 PetscCallTC(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   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
102   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");PetscCall(ierr);
103   PetscCall(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL));
104   PetscCall(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL));
105   user.pressure = 1.01325e5;    /* Pascal */
106   PetscCall(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL));
107   user.Tini   = 1550;
108   PetscCall(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL));
109   user.diffus = 100;
110   PetscCall(PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL));
111   PetscCall(PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL));
112   user.diffusion = PETSC_TRUE;
113   PetscCall(PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL));
114   user.reactions = PETSC_TRUE;
115   PetscCall(PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL));
116   ierr = PetscOptionsEnd();PetscCall(ierr);
117 
118   PetscCallTC(TC_initChem(chemfile, thermofile, 0, 1.0));
119   user.Nspec = TC_getNspec();
120   user.Nreac = TC_getNreac();
121 
122   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm));
123   PetscCall(DMSetFromOptions(user.dm));
124   PetscCall(DMSetUp(user.dm));
125   PetscCall(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   PetscCall(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   PetscCall(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names));
131   PetscCall(PetscStrcpy(names,"Temp"));
132   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
133   PetscCall(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   PetscCall(DMDASetFieldNames(user.dm,(const char * const *)snames));
137   PetscCall(PetscFree(snames));
138   PetscCall(PetscFree(names));
139 
140   PetscCall(DMCreateMatrix(user.dm,&J));
141   PetscCall(DMCreateGlobalVector(user.dm,&X));
142 
143   PetscCall(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   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
149   PetscCall(TSSetDM(ts,user.dm));
150   PetscCall(TSSetType(ts,TSARKIMEX));
151   PetscCall(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
152   PetscCall(TSARKIMEXSetType(ts,TSARKIMEX4));
153   PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
154   PetscCall(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user));
155 
156   ftime = 1.0;
157   PetscCall(TSSetMaxTime(ts,ftime));
158   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
159 
160   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161      Set initial conditions
162    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163   PetscCall(FormInitialSolution(ts,X,&user));
164   PetscCall(TSSetSolution(ts,X));
165   dt   = 1e-10;                 /* Initial time step */
166   PetscCall(TSSetTimeStep(ts,dt));
167   PetscCall(TSGetAdapt(ts,&adapt));
168   PetscCall(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
169   PetscCall(TSSetMaxSNESFailures(ts,-1));            /* Retry step an unlimited number of times */
170 
171   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172      Pass information to graphical monitoring routine
173    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174   if (showsolutions) {
175     PetscCall(DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL));
176     for (i=xs;i<xs+xm;i++) {
177       PetscCall(MonitorCell(ts,&user,i));
178     }
179   }
180 
181   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182      Set runtime options
183    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184   PetscCall(TSSetFromOptions(ts));
185 
186   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187      Set final conditions for sensitivities
188    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189   PetscCall(DMCreateGlobalVector(user.dm,&lambda));
190   PetscCall(TSSetCostGradients(ts,1,&lambda,NULL));
191   PetscCall(VecSetValue(lambda,0,1.0,INSERT_VALUES));
192   PetscCall(VecAssemblyBegin(lambda));
193   PetscCall(VecAssemblyEnd(lambda));
194 
195   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196      Solve ODE
197      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198   PetscCall(TSSolve(ts,X));
199   PetscCall(TSGetSolveTime(ts,&ftime));
200   PetscCall(TSGetStepNumber(ts,&steps));
201   PetscCall(TSGetConvergedReason(ts,&reason));
202   PetscCall(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     PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL));
211     if (max) {
212       PetscCall(TSMonitorLGGetVariableNames(ts,&names));
213       if (names) {
214         PetscCall(VecGetArrayRead(max,&bmax));
215         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n"));
216         for (i=1; i<user.Nspec; i++) {
217           if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i]));
218         }
219         PetscCall(VecRestoreArrayRead(max,&bmax));
220       }
221     }
222   }
223 
224   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225      Free work space.
226    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227   TC_reset();
228   PetscCall(DMDestroy(&user.dm));
229   PetscCall(MatDestroy(&J));
230   PetscCall(VecDestroy(&X));
231   PetscCall(VecDestroy(&lambda));
232   PetscCall(TSDestroy(&ts));
233   PetscCall(PetscFree3(user.tchemwork,user.Jdense,user.rows));
234   PetscCall(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   PetscScalar       **f;
245   const PetscScalar **x;
246   DM                dm;
247   PetscInt          i,xs,xm,j,dof;
248   Vec               Xlocal;
249   PetscReal         idx;
250 
251   PetscFunctionBeginUser;
252   PetscCall(TSGetDM(ts,&dm));
253   PetscCall(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
254   PetscCall(DMGetLocalVector(dm,&Xlocal));
255   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal));
256   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal));
257   PetscCall(DMDAVecGetArrayDOFRead(dm,Xlocal,&x));
258   PetscCall(DMDAVecGetArrayDOF(dm,F,&f));
259   PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
260 
261   idx = 1.0*user->diffus/user->dx;
262   for (i=xs; i<xs+xm; i++) {
263     for (j=0; j<dof; j++) {
264       f[i][j] += idx*(x[i+1][j] - 2.0*x[i][j] + x[i-1][j]);
265     }
266   }
267   PetscCall(DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x));
268   PetscCall(DMDAVecRestoreArrayDOF(dm,F,&f));
269   PetscCall(DMRestoreLocalVector(dm,&Xlocal));
270   PetscFunctionReturn(0);
271 }
272 
273 /*
274    Produces the second order centered difference diffusion operator on a one dimensional periodic domain
275 */
276 static PetscErrorCode FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
277 {
278   User              user = (User)ptr;
279   DM                dm;
280   PetscInt          i,xs,xm,j,dof;
281   PetscReal         idx,values[3];
282   MatStencil        row,col[3];
283 
284   PetscFunctionBeginUser;
285   PetscCall(TSGetDM(ts,&dm));
286   PetscCall(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
287   PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
288 
289   idx = 1.0*user->diffus/user->dx;
290   values[0] = idx;
291   values[1] = -2.0*idx;
292   values[2] = idx;
293   for (i=xs; i<xs+xm; i++) {
294     for (j=0; j<dof; j++) {
295       row.i = i;      row.c = j;
296       col[0].i = i-1; col[0].c = j;
297       col[1].i = i;   col[1].c = j;
298       col[2].i = i+1; col[2].c = j;
299       PetscCall(MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES));
300     }
301   }
302   PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
303   PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
304   PetscFunctionReturn(0);
305 }
306 
307 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
308 {
309   User              user = (User)ptr;
310   PetscScalar       **f;
311   const PetscScalar **x;
312   DM                dm;
313   PetscInt          i,xs,xm;
314 
315   PetscFunctionBeginUser;
316   if (user->reactions) {
317     PetscCall(TSGetDM(ts,&dm));
318     PetscCall(DMDAVecGetArrayDOFRead(dm,X,&x));
319     PetscCall(DMDAVecGetArrayDOF(dm,F,&f));
320     PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
321 
322     for (i=xs; i<xs+xm; i++) {
323       PetscCall(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1));
324       user->tchemwork[0] *= user->Tini; /* Dimensionalize */
325       PetscCallTC(TC_getSrc(user->tchemwork,user->Nspec+1,f[i]));
326       f[i][0] /= user->Tini;           /* Non-dimensionalize */
327     }
328 
329     PetscCall(DMDAVecRestoreArrayDOFRead(dm,X,&x));
330     PetscCall(DMDAVecRestoreArrayDOF(dm,F,&f));
331   } else {
332     PetscCall(VecZeroEntries(F));
333   }
334   if (user->diffusion) {
335     PetscCall(FormDiffusionFunction(ts,t,X,F,ptr));
336   }
337   PetscFunctionReturn(0);
338 }
339 
340 static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
341 {
342   User              user = (User)ptr;
343   const PetscScalar **x;
344   PetscInt          M = user->Nspec+1,i,j,xs,xm;
345   DM                dm;
346 
347   PetscFunctionBeginUser;
348   if (user->reactions) {
349     PetscCall(TSGetDM(ts,&dm));
350     PetscCall(MatZeroEntries(Pmat));
351     PetscCall(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE));
352     PetscCall(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
353     PetscCall(DMDAVecGetArrayDOFRead(dm,X,&x));
354     PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
355 
356     for (i=xs; i<xs+xm; i++) {
357       PetscCall(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1));
358       user->tchemwork[0] *= user->Tini;  /* Dimensionalize temperature (first row) because that is what Tchem wants */
359       PetscCall(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1));
360 
361       for (j=0; j<M; j++) user->Jdense[j + 0*M] /= user->Tini; /* Non-dimensionalize first column */
362       for (j=0; j<M; j++) user->Jdense[0 + j*M] /= user->Tini; /* Non-dimensionalize first row */
363       for (j=0; j<M; j++) user->rows[j] = i*M+j;
364       PetscCall(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES));
365     }
366     PetscCall(DMDAVecRestoreArrayDOFRead(dm,X,&x));
367     PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
368     PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
369   } else {
370     PetscCall(MatZeroEntries(Pmat));
371   }
372   if (user->diffusion) {
373     PetscCall(FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr));
374   }
375   if (Amat != Pmat) {
376     PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY));
377     PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY));
378   }
379   PetscFunctionReturn(0);
380 }
381 
382 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
383 {
384   PetscScalar    **x,*xc;
385   struct {const char *name; PetscReal massfrac;} initial[] = {
386     {"CH4", 0.0948178320887},
387     {"O2", 0.189635664177},
388     {"N2", 0.706766236705},
389     {"AR", 0.00878026702874}
390   };
391   PetscInt       i,j,xs,xm;
392   DM             dm;
393 
394   PetscFunctionBeginUser;
395   PetscCall(VecZeroEntries(X));
396   PetscCall(TSGetDM(ts,&dm));
397   PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
398 
399   PetscCall(DMDAGetCoordinateArray(dm,&xc));
400   PetscCall(DMDAVecGetArrayDOF(dm,X,&x));
401   for (i=xs; i<xs+xm; i++) {
402     x[i][0] = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[i]);  /* Non-dimensionalized by user->Tini */
403     for (j=0; j<sizeof(initial)/sizeof(initial[0]); j++) {
404       int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name));
405       PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",initial[j].name);
406       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",j,initial[j].name,initial[j].massfrac));
407       x[i][1+ispec] = initial[j].massfrac;
408     }
409   }
410   PetscCall(DMDAVecRestoreArrayDOF(dm,X,&x));
411   PetscCall(DMDARestoreCoordinateArray(dm,&xc));
412   PetscFunctionReturn(0);
413 }
414 
415 /*
416     Routines for displaying the solutions
417 */
418 typedef struct {
419   PetscInt cell;
420   User     user;
421 } UserLGCtx;
422 
423 static PetscErrorCode FormMoleFraction(UserLGCtx *ctx,Vec massf,Vec *molef)
424 {
425   User              user = ctx->user;
426   PetscReal         *M,tM=0;
427   PetscInt          i,n = user->Nspec+1;
428   PetscScalar       *mof;
429   const PetscScalar **maf;
430 
431   PetscFunctionBegin;
432   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,molef));
433   PetscCall(PetscMalloc1(user->Nspec,&M));
434   TC_getSmass(user->Nspec, M);
435   PetscCall(DMDAVecGetArrayDOFRead(user->dm,massf,&maf));
436   PetscCall(VecGetArray(*molef,&mof));
437   mof[0] = maf[ctx->cell][0]; /* copy over temperature */
438   for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1];
439   for (i=1; i<n; i++) {
440     mof[i] = maf[ctx->cell][i]/(M[i-1]*tM);
441   }
442   PetscCall(DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf));
443   PetscCall(VecRestoreArray(*molef,&mof));
444   PetscCall(PetscFree(M));
445   PetscFunctionReturn(0);
446 }
447 
448 static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx)
449 {
450   PetscFunctionBegin;
451   PetscCall(PetscFree(uctx));
452   PetscFunctionReturn(0);
453 }
454 
455 /*
456    Use TSMonitorLG to monitor the reactions in a particular cell
457 */
458 static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell)
459 {
460   TSMonitorLGCtx ctx;
461   char           **snames;
462   UserLGCtx      *uctx;
463   char           label[128];
464   PetscReal      temp,*xc;
465   PetscMPIInt    rank;
466 
467   PetscFunctionBegin;
468   PetscCall(DMDAGetCoordinateArray(user->dm,&xc));
469   temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]);  /* Non-dimensionalized by user->Tini */
470   PetscCall(DMDARestoreCoordinateArray(user->dm,&xc));
471   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
472   PetscCall(PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank));
473   PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx));
474   PetscCall(DMDAGetFieldNames(user->dm,(const char * const **)&snames));
475   PetscCall(TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames));
476   PetscCall(PetscNew(&uctx));
477   uctx->cell = cell;
478   uctx->user = user;
479   PetscCall(TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx));
480   PetscCall(TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy));
481   PetscFunctionReturn(0);
482 }
483