Actual source code: extchem.c
1: static const char help[] = "Integrate chemistry using TChem.\n";
3: #include <petscts.h>
5: #if defined(PETSC_HAVE_TCHEM)
6: #if defined(MAX)
7: #undef MAX
8: #endif
9: #if defined(MIN)
10: #undef MIN
11: #endif
12: # include <TC_params.h>
13: # include <TC_interface.h>
14: #else
15: # error TChem is required for this example. Reconfigure PETSc using --download-tchem.
16: #endif
17: /*
18: See extchem.example.1 for how to run an example
20: See also h2_10sp.inp for another example
22: Determine sensitivity of final temperature on each variables initial conditions
23: -ts_dt 1.e-5 -ts_type cn -ts_adjoint_solve -ts_adjoint_view_solution draw
25: The solution for component i = 0 is the temperature.
27: The solution, i > 0, is the mass fraction, massf[i], of species i, i.e. mass of species i/ total mass of all species
29: The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species
30: Define M[i] = mass per mole of species i then
31: molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j]))
33: FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species.
35: These are other data sets for other possible runs
36: https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/n_heptane_v3.1_therm.dat
37: https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/nc7_ver3.1_mech.txt
39: */
40: typedef struct _User *User;
41: struct _User {
42: PetscReal pressure;
43: int Nspec;
44: int Nreac;
45: PetscReal Tini;
46: double *tchemwork;
47: double *Jdense; /* Dense array workspace where Tchem computes the Jacobian */
48: PetscInt *rows;
49: char **snames;
50: };
52: static PetscErrorCode PrintSpecies(User,Vec);
53: static PetscErrorCode MassFractionToMoleFraction(User,Vec,Vec*);
54: static PetscErrorCode MoleFractionToMassFraction(User,Vec,Vec*);
55: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
56: static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
57: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
58: static PetscErrorCode ComputeMassConservation(Vec,PetscReal*,void*);
59: static PetscErrorCode MonitorMassConservation(TS,PetscInt,PetscReal,Vec,void*);
60: static PetscErrorCode MonitorTempature(TS,PetscInt,PetscReal,Vec,void*);
62: #define TCCHKERRQ(ierr) do {if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0)
64: int main(int argc,char **argv)
65: {
66: TS ts; /* time integrator */
67: TSAdapt adapt;
68: Vec X,lambda; /* solution vector */
69: Mat J; /* Jacobian matrix */
70: PetscInt steps;
71: PetscErrorCode ierr;
72: PetscReal ftime,dt;
73: char chemfile[PETSC_MAX_PATH_LEN],thermofile[PETSC_MAX_PATH_LEN],lchemfile[PETSC_MAX_PATH_LEN],lthermofile[PETSC_MAX_PATH_LEN],lperiodic[PETSC_MAX_PATH_LEN];
74: const char *periodic = "file://${PETSC_DIR}/${PETSC_ARCH}/share/periodictable.dat";
75: struct _User user; /* user-defined work context */
76: TSConvergedReason reason;
77: char **snames,*names;
78: PetscInt i;
79: TSTrajectory tj;
80: PetscBool flg = PETSC_FALSE,tflg = PETSC_FALSE,found;
82: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
83: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");
84: PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL);
85: PetscFileRetrieve(PETSC_COMM_WORLD,chemfile,lchemfile,PETSC_MAX_PATH_LEN,&found);
86: if (!found) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",chemfile,lchemfile);
87: PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL);
88: PetscFileRetrieve(PETSC_COMM_WORLD,thermofile,lthermofile,PETSC_MAX_PATH_LEN,&found);
89: if (!found) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",thermofile,lthermofile);
90: user.pressure = 1.01325e5; /* Pascal */
91: PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL);
92: user.Tini = 1000; /* Kelvin */
93: PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL);
94: PetscOptionsBool("-monitor_mass","Monitor the total mass at each timestep","",flg,&flg,NULL);
95: PetscOptionsBool("-monitor_temp","Monitor the temperature each timestep","",tflg,&tflg,NULL);
96: PetscOptionsEnd();
98: /* tchem requires periodic table in current directory */
99: PetscFileRetrieve(PETSC_COMM_WORLD,periodic,lperiodic,PETSC_MAX_PATH_LEN,&found);
100: if (!found) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot located required periodic table %s or local version %s",periodic,lperiodic);
102: TC_initChem(lchemfile, lthermofile, 0, 1.0);TC
103: TC_setThermoPres(user.pressure);
104: user.Nspec = TC_getNspec();
105: user.Nreac = TC_getNreac();
106: /*
107: Get names of all species in easy to use array
108: */
109: PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names);
110: PetscStrcpy(names,"Temp");
111: TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
112: PetscMalloc1((user.Nspec+2),&snames);
113: for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME;
114: snames[user.Nspec+1] = NULL;
115: PetscStrArrayallocpy((const char *const *)snames,&user.snames);
116: PetscFree(snames);
117: PetscFree(names);
119: PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows);
120: VecCreateSeq(PETSC_COMM_SELF,user.Nspec+1,&X);
122: /* MatCreateSeqAIJ(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,PETSC_DECIDE,NULL,&J); */
123: MatCreateSeqDense(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,NULL,&J);
124: MatSetFromOptions(J);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Create timestepping solver context
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: TSCreate(PETSC_COMM_WORLD,&ts);
130: TSSetType(ts,TSARKIMEX);
131: TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
132: TSARKIMEXSetType(ts,TSARKIMEX4);
133: TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
134: TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user);
136: if (flg) {
137: TSMonitorSet(ts,MonitorMassConservation,NULL,NULL);
138: }
139: if (tflg) {
140: TSMonitorSet(ts,MonitorTempature,&user,NULL);
141: }
143: ftime = 1.0;
144: TSSetMaxTime(ts,ftime);
145: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
147: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: Set initial conditions
149: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150: FormInitialSolution(ts,X,&user);
151: TSSetSolution(ts,X);
152: dt = 1e-10; /* Initial time step */
153: TSSetTimeStep(ts,dt);
154: TSGetAdapt(ts,&adapt);
155: TSAdaptSetStepLimits(adapt,1e-12,1e-4); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
156: TSSetMaxSNESFailures(ts,-1); /* Retry step an unlimited number of times */
158: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159: Set runtime options
160: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: TSSetFromOptions(ts);
163: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: Set final conditions for sensitivities
165: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166: VecDuplicate(X,&lambda);
167: TSSetCostGradients(ts,1,&lambda,NULL);
168: VecSetValue(lambda,0,1.0,INSERT_VALUES);
169: VecAssemblyBegin(lambda);
170: VecAssemblyEnd(lambda);
172: TSGetTrajectory(ts,&tj);
173: if (tj) {
174: TSTrajectorySetVariableNames(tj,(const char * const *)user.snames);
175: TSTrajectorySetTransform(tj,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user);
176: }
178: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179: Pass information to graphical monitoring routine
180: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181: TSMonitorLGSetVariableNames(ts,(const char * const *)user.snames);
182: TSMonitorLGSetTransform(ts,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user);
184: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: Solve ODE
186: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187: TSSolve(ts,X);
188: TSGetSolveTime(ts,&ftime);
189: TSGetStepNumber(ts,&steps);
190: TSGetConvergedReason(ts,&reason);
191: PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);
193: /* {
194: Vec max;
195: PetscInt i;
196: const PetscReal *bmax;
198: TSMonitorEnvelopeGetBounds(ts,&max,NULL);
199: if (max) {
200: VecGetArrayRead(max,&bmax);
201: PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n");
202: for (i=1; i<user.Nspec; i++) {
203: if (bmax[i] > .01) {PetscPrintf(PETSC_COMM_SELF,"%s %g\n",user.snames[i],(double)bmax[i]);}
204: }
205: VecRestoreArrayRead(max,&bmax);
206: }
207: }
209: Vec y;
210: MassFractionToMoleFraction(&user,X,&y);
211: PrintSpecies(&user,y);
212: VecDestroy(&y); */
214: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215: Free work space.
216: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217: TC_reset();
218: PetscStrArrayDestroy(&user.snames);
219: MatDestroy(&J);
220: VecDestroy(&X);
221: VecDestroy(&lambda);
222: TSDestroy(&ts);
223: PetscFree3(user.tchemwork,user.Jdense,user.rows);
224: PetscFinalize();
225: return ierr;
226: }
228: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
229: {
230: User user = (User)ptr;
231: PetscErrorCode ierr;
232: PetscScalar *f;
233: const PetscScalar *x;
236: VecGetArrayRead(X,&x);
237: VecGetArray(F,&f);
239: PetscArraycpy(user->tchemwork,x,user->Nspec+1);
240: user->tchemwork[0] *= user->Tini; /* Dimensionalize */
241: TC_getSrc(user->tchemwork,user->Nspec+1,f);TC
242: f[0] /= user->Tini; /* Non-dimensionalize */
244: VecRestoreArrayRead(X,&x);
245: VecRestoreArray(F,&f);
246: return(0);
247: }
249: static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
250: {
251: User user = (User)ptr;
252: PetscErrorCode ierr;
253: const PetscScalar *x;
254: PetscInt M = user->Nspec+1,i;
257: VecGetArrayRead(X,&x);
258: PetscArraycpy(user->tchemwork,x,user->Nspec+1);
259: VecRestoreArrayRead(X,&x);
260: user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */
261: TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1);
263: for (i=0; i<M; i++) user->Jdense[i + 0*M] /= user->Tini; /* Non-dimensionalize first column */
264: for (i=0; i<M; i++) user->Jdense[0 + i*M] /= user->Tini; /* Non-dimensionalize first row */
265: for (i=0; i<M; i++) user->rows[i] = i;
266: MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE);
267: MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
268: MatZeroEntries(Pmat);
269: MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES);
271: MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);
272: MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);
273: if (Amat != Pmat) {
274: MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
275: MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
276: }
277: return(0);
278: }
280: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
281: {
282: PetscScalar *x;
284: PetscInt i;
285: Vec y;
286: const PetscInt maxspecies = 10;
287: PetscInt smax = maxspecies,mmax = maxspecies;
288: char *names[maxspecies];
289: PetscReal molefracs[maxspecies],sum;
290: PetscBool flg;
293: VecZeroEntries(X);
294: VecGetArray(X,&x);
295: x[0] = 1.0; /* Non-dimensionalized by user->Tini */
297: PetscOptionsGetStringArray(NULL,NULL,"-initial_species",names,&smax,&flg);
298: if (smax < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide at least two initial species");
299: PetscOptionsGetRealArray(NULL,NULL,"-initial_mole",molefracs,&mmax,&flg);
300: if (smax != mmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide same number of initial species %D as initial moles %D",smax,mmax);
301: sum = 0;
302: for (i=0; i<smax; i++) sum += molefracs[i];
303: for (i=0; i<smax; i++) molefracs[i] = molefracs[i]/sum;
304: for (i=0; i<smax; i++) {
305: int ispec = TC_getSpos(names[i], strlen(names[i]));
306: if (ispec < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",names[i]);
307: PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",i,names[i],molefracs[i]);
308: x[1+ispec] = molefracs[i];
309: }
310: for (i=0; i<smax; i++) {
311: PetscFree(names[i]);
312: }
313: VecRestoreArray(X,&x);
314: /* PrintSpecies((User)ctx,X); */
315: MoleFractionToMassFraction((User)ctx,X,&y);
316: VecCopy(y,X);
317: VecDestroy(&y);
318: return(0);
319: }
321: /*
322: Converts the input vector which is in mass fractions (used by tchem) to mole fractions
323: */
324: PetscErrorCode MassFractionToMoleFraction(User user,Vec massf,Vec *molef)
325: {
326: PetscErrorCode ierr;
327: PetscScalar *mof;
328: const PetscScalar *maf;
331: VecDuplicate(massf,molef);
332: VecGetArrayRead(massf,&maf);
333: VecGetArray(*molef,&mof);
334: mof[0] = maf[0]; /* copy over temperature */
335: TC_getMs2Ml((double*)maf+1,user->Nspec,mof+1);
336: VecRestoreArray(*molef,&mof);
337: VecRestoreArrayRead(massf,&maf);
338: return(0);
339: }
341: /*
342: Converts the input vector which is in mole fractions to mass fractions (used by tchem)
343: */
344: PetscErrorCode MoleFractionToMassFraction(User user,Vec molef,Vec *massf)
345: {
346: PetscErrorCode ierr;
347: const PetscScalar *mof;
348: PetscScalar *maf;
351: VecDuplicate(molef,massf);
352: VecGetArrayRead(molef,&mof);
353: VecGetArray(*massf,&maf);
354: maf[0] = mof[0]; /* copy over temperature */
355: TC_getMl2Ms((double*)mof+1,user->Nspec,maf+1);
356: VecRestoreArrayRead(molef,&mof);
357: VecRestoreArray(*massf,&maf);
358: return(0);
359: }
361: PetscErrorCode ComputeMassConservation(Vec x,PetscReal *mass,void* ctx)
362: {
366: VecSum(x,mass);
367: return(0);
368: }
370: PetscErrorCode MonitorMassConservation(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
371: {
372: const PetscScalar *T;
373: PetscReal mass;
374: PetscErrorCode ierr;
377: ComputeMassConservation(x,&mass,ctx);
378: VecGetArrayRead(x,&T);
379: mass -= PetscAbsScalar(T[0]);
380: VecRestoreArrayRead(x,&T);
381: PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g percent mass lost or gained %g\n",step,(double)time,(double)100.*(1.0 - mass));
382: return(0);
383: }
385: PetscErrorCode MonitorTempature(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
386: {
387: User user = (User) ctx;
388: const PetscScalar *T;
389: PetscErrorCode ierr;
392: VecGetArrayRead(x,&T);
393: PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g temperature %g\n",step,(double)time,(double)T[0]*user->Tini);
394: VecRestoreArrayRead(x,&T);
395: return(0);
396: }
398: /*
399: Prints out each species with its name
400: */
401: PETSC_UNUSED PetscErrorCode PrintSpecies(User user,Vec molef)
402: {
403: PetscErrorCode ierr;
404: const PetscScalar *mof;
405: PetscInt i,*idx,n = user->Nspec+1;
408: PetscMalloc1(n,&idx);
409: for (i=0; i<n;i++) idx[i] = i;
410: VecGetArrayRead(molef,&mof);
411: PetscSortRealWithPermutation(n,mof,idx);
412: for (i=0; i<n; i++) {
413: PetscPrintf(PETSC_COMM_SELF,"%6s %g\n",user->snames[idx[n-i-1]],mof[idx[n-i-1]]);
414: }
415: PetscFree(idx);
416: VecRestoreArrayRead(molef,&mof);
417: return(0);
418: }
420: /*TEST
421: build:
422: requires: tchem
424: test:
425: args: -chem http://combustion.berkeley.edu/gri_mech/version30/files30/grimech30.dat -thermo http://combustion.berkeley.edu/gri_mech/version30/files30/thermo30.dat -initial_species CH4,O2,N2,AR -initial_mole 0.0948178320887,0.189635664177,0.706766236705,0.00878026702874 -Tini 1500 -Tini 1500 -ts_arkimex_fully_implicit -ts_max_snes_failures -1 -ts_adapt_dt_max 1e-4 -ts_arkimex_type 4 -ts_max_time .005
426: requires: !single
427: filter: grep -v iterations
429: TEST*/