Actual source code: extchemfield.c

  1: static const char help[] = "Integrate chemistry using TChem.\n";

  3: #include <petscts.h>
  4: #include <petscdmda.h>

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

 20:     This is an extension of extchem.c to solve the reaction equations independently in each cell of a one dimensional field

 22:     Obtain the three files into this directory

 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 .

 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

 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

 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

 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

 41:         Save the images in a .gif (movie) file
 42:         -draw_save -draw_save_single_file

 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

 47:         Turn off diffusion
 48:         -diffusion no

 50:         Turn off reactions
 51:         -reactions no

 53:     The solution for component i = 0 is the temperature.

 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

 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]))

 61:     FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species.

 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: };

 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*);


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

101:   PetscInitialize(&argc,&argv,(char*)0,help);
102:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");
103:   PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL);
104:   PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL);
105:   user.pressure = 1.01325e5;    /* Pascal */
106:   PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL);
107:   user.Tini   = 1550;
108:   PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL);
109:   user.diffus = 100;
110:   PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL);
111:   PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL);
112:   user.diffusion = PETSC_TRUE;
113:   PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL);
114:   user.reactions = PETSC_TRUE;
115:   PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL);
116:   PetscOptionsEnd();

118:   TC_initChem(chemfile, thermofile, 0, 1.0);
119:   user.Nspec = TC_getNspec();
120:   user.Nreac = TC_getNreac();

122:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm);
123:   DMSetFromOptions(user.dm);
124:   DMSetUp(user.dm);
125:   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:   DMDASetUniformCoordinates(user.dm,0.0,1.0,0.0,1.0,0.0,1.0);

129:   /* set the names of each field in the DMDA based on the species name */
130:   PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names);
131:   PetscStrcpy(names,"Temp");
132:   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
133:   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:   DMDASetFieldNames(user.dm,(const char * const *)snames);
137:   PetscFree(snames);
138:   PetscFree(names);

140:   DMCreateMatrix(user.dm,&J);
141:   DMCreateGlobalVector(user.dm,&X);

143:   PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows);

145:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146:      Create timestepping solver context
147:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148:   TSCreate(PETSC_COMM_WORLD,&ts);
149:   TSSetDM(ts,user.dm);
150:   TSSetType(ts,TSARKIMEX);
151:   TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
152:   TSARKIMEXSetType(ts,TSARKIMEX4);
153:   TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
154:   TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user);

156:   ftime = 1.0;
157:   TSSetMaxTime(ts,ftime);
158:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161:      Set initial conditions
162:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163:   FormInitialSolution(ts,X,&user);
164:   TSSetSolution(ts,X);
165:   dt   = 1e-10;                 /* Initial time step */
166:   TSSetTimeStep(ts,dt);
167:   TSGetAdapt(ts,&adapt);
168:   TSAdaptSetStepLimits(adapt,1e-12,1e-4); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
169:   TSSetMaxSNESFailures(ts,-1);            /* Retry step an unlimited number of times */

171:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172:      Pass information to graphical monitoring routine
173:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174:   if (showsolutions) {
175:     DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL);
176:     for (i=xs;i<xs+xm;i++) {
177:       MonitorCell(ts,&user,i);
178:     }
179:   }

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:      Set runtime options
183:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184:   TSSetFromOptions(ts);

186:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187:      Set final conditions for sensitivities
188:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189:   DMCreateGlobalVector(user.dm,&lambda);
190:   TSSetCostGradients(ts,1,&lambda,NULL);
191:   VecSetValue(lambda,0,1.0,INSERT_VALUES);
192:   VecAssemblyBegin(lambda);
193:   VecAssemblyEnd(lambda);

195:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196:      Solve ODE
197:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198:   TSSolve(ts,X);
199:   TSGetSolveTime(ts,&ftime);
200:   TSGetStepNumber(ts,&steps);
201:   TSGetConvergedReason(ts,&reason);
202:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);

204:   {
205:     Vec                max;
206:     const char * const *names;
207:     PetscInt           i;
208:     const PetscReal    *bmax;

210:     TSMonitorEnvelopeGetBounds(ts,&max,NULL);
211:     if (max) {
212:       TSMonitorLGGetVariableNames(ts,&names);
213:       if (names) {
214:         VecGetArrayRead(max,&bmax);
215:         PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n");
216:         for (i=1; i<user.Nspec; i++) {
217:           if (bmax[i] > .01) PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i]);
218:         }
219:         VecRestoreArrayRead(max,&bmax);
220:       }
221:     }
222:   }

224:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225:      Free work space.
226:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227:   TC_reset();
228:   DMDestroy(&user.dm);
229:   MatDestroy(&J);
230:   VecDestroy(&X);
231:   VecDestroy(&lambda);
232:   TSDestroy(&ts);
233:   PetscFree3(user.tchemwork,user.Jdense,user.rows);
234:   PetscFinalize();
235:   return 0;
236: }

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;

252:   TSGetDM(ts,&dm);
253:   DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
254:   DMGetLocalVector(dm,&Xlocal);
255:   DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal);
256:   DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal);
257:   DMDAVecGetArrayDOFRead(dm,Xlocal,&x);
258:   DMDAVecGetArrayDOF(dm,F,&f);
259:   DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);

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:   DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x);
268:   DMDAVecRestoreArrayDOF(dm,F,&f);
269:   DMRestoreLocalVector(dm,&Xlocal);
270:   return 0;
271: }

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];

285:   TSGetDM(ts,&dm);
286:   DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
287:   DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);

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:       MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES);
300:     }
301:   }
302:   MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);
303:   MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);
304:   return 0;
305: }

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;

316:   if (user->reactions) {
317:     TSGetDM(ts,&dm);
318:     DMDAVecGetArrayDOFRead(dm,X,&x);
319:     DMDAVecGetArrayDOF(dm,F,&f);
320:     DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);

322:     for (i=xs; i<xs+xm; i++) {
323:       PetscArraycpy(user->tchemwork,x[i],user->Nspec+1);
324:       user->tchemwork[0] *= user->Tini; /* Dimensionalize */
325:       TC_getSrc(user->tchemwork,user->Nspec+1,f[i]);
326:       f[i][0] /= user->Tini;           /* Non-dimensionalize */
327:     }

329:     DMDAVecRestoreArrayDOFRead(dm,X,&x);
330:     DMDAVecRestoreArrayDOF(dm,F,&f);
331:   } else {
332:     VecZeroEntries(F);
333:   }
334:   if (user->diffusion) {
335:     FormDiffusionFunction(ts,t,X,F,ptr);
336:   }
337:   return 0;
338: }

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;

348:   if (user->reactions) {
349:     TSGetDM(ts,&dm);
350:     MatZeroEntries(Pmat);
351:     MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE);
352:     MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
353:     DMDAVecGetArrayDOFRead(dm,X,&x);
354:     DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);

356:     for (i=xs; i<xs+xm; i++) {
357:       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:       TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1);

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:       MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES);
365:     }
366:     DMDAVecRestoreArrayDOFRead(dm,X,&x);
367:     MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);
368:     MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);
369:   } else {
370:     MatZeroEntries(Pmat);
371:   }
372:   if (user->diffusion) {
373:     FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr);
374:   }
375:   if (Amat != Pmat) {
376:     MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
377:     MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
378:   }
379:   return 0;
380: }

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;

395:   VecZeroEntries(X);
396:   TSGetDM(ts,&dm);
397:   DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);

399:   DMDAGetCoordinateArray(dm,&xc);
400:   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));
406:       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:   DMDAVecRestoreArrayDOF(dm,X,&x);
411:   DMDARestoreCoordinateArray(dm,&xc);
412:   return 0;
413: }

415: /*
416:     Routines for displaying the solutions
417: */
418: typedef struct {
419:   PetscInt cell;
420:   User     user;
421: } UserLGCtx;

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;

431:   VecCreateSeq(PETSC_COMM_SELF,n,molef);
432:   PetscMalloc1(user->Nspec,&M);
433:   TC_getSmass(user->Nspec, M);
434:   DMDAVecGetArrayDOFRead(user->dm,massf,&maf);
435:   VecGetArray(*molef,&mof);
436:   mof[0] = maf[ctx->cell][0]; /* copy over temperature */
437:   for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1];
438:   for (i=1; i<n; i++) {
439:     mof[i] = maf[ctx->cell][i]/(M[i-1]*tM);
440:   }
441:   DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf);
442:   VecRestoreArray(*molef,&mof);
443:   PetscFree(M);
444:   return 0;
445: }

447: static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx)
448: {
449:   PetscFree(uctx);
450:   return 0;
451: }

453: /*
454:    Use TSMonitorLG to monitor the reactions in a particular cell
455: */
456: static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell)
457: {
458:   TSMonitorLGCtx ctx;
459:   char           **snames;
460:   UserLGCtx      *uctx;
461:   char           label[128];
462:   PetscReal      temp,*xc;
463:   PetscMPIInt    rank;

465:   DMDAGetCoordinateArray(user->dm,&xc);
466:   temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]);  /* Non-dimensionalized by user->Tini */
467:   DMDARestoreCoordinateArray(user->dm,&xc);
468:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
469:   PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank);
470:   TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx);
471:   DMDAGetFieldNames(user->dm,(const char * const **)&snames);
472:   TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames);
473:   PetscNew(&uctx);
474:   uctx->cell = cell;
475:   uctx->user = user;
476:   TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx);
477:   TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);
478:   return 0;
479: }