xref: /petsc/src/ts/tutorials/extchem.c (revision 6aad120caa16b1027d343a5f30f73d01448e4dc0)
1 static const char help[] = "Integrate chemistry using TChem.\n";
2 
3 #include <petscts.h>
4 
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
19 
20     See also h2_10sp.inp for another example
21 
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
24 
25     The solution for component i = 0 is the temperature.
26 
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
28 
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]))
32 
33     FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species.
34 
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
38 
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 };
51 
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*);
61 
62 #define PetscCallTC(ierr) do {PetscCheck(!ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0)
63 
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;
81 
82   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
83   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");PetscCall(ierr);
84   PetscCall(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL));
85   PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,chemfile,lchemfile,PETSC_MAX_PATH_LEN,&found));
86   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",chemfile,lchemfile);
87   PetscCall(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL));
88   PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,thermofile,lthermofile,PETSC_MAX_PATH_LEN,&found));
89   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",thermofile,lthermofile);
90   user.pressure = 1.01325e5;    /* Pascal */
91   PetscCall(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL));
92   user.Tini = 1000;             /* Kelvin */
93   PetscCall(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL));
94   PetscCall(PetscOptionsBool("-monitor_mass","Monitor the total mass at each timestep","",flg,&flg,NULL));
95   PetscCall(PetscOptionsBool("-monitor_temp","Monitor the temperature each timestep","",tflg,&tflg,NULL));
96   ierr = PetscOptionsEnd();PetscCall(ierr);
97 
98   /* tchem requires periodic table in current directory */
99   PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,periodic,lperiodic,PETSC_MAX_PATH_LEN,&found));
100   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot located required periodic table %s or local version %s",periodic,lperiodic);
101 
102   PetscCallTC(TC_initChem(lchemfile, lthermofile, 0, 1.0));
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   PetscCall(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names));
110   PetscCall(PetscStrcpy(names,"Temp"));
111   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
112   PetscCall(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   PetscCall(PetscStrArrayallocpy((const char *const *)snames,&user.snames));
116   PetscCall(PetscFree(snames));
117   PetscCall(PetscFree(names));
118 
119   PetscCall(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows));
120   PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.Nspec+1,&X));
121 
122   /* PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,PETSC_DECIDE,NULL,&J)); */
123   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,NULL,&J));
124   PetscCall(MatSetFromOptions(J));
125 
126   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127      Create timestepping solver context
128      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
130   PetscCall(TSSetType(ts,TSARKIMEX));
131   PetscCall(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
132   PetscCall(TSARKIMEXSetType(ts,TSARKIMEX4));
133   PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
134   PetscCall(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user));
135 
136   if (flg) {
137     PetscCall(TSMonitorSet(ts,MonitorMassConservation,NULL,NULL));
138   }
139   if (tflg) {
140     PetscCall(TSMonitorSet(ts,MonitorTempature,&user,NULL));
141   }
142 
143   ftime = 1.0;
144   PetscCall(TSSetMaxTime(ts,ftime));
145   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
146 
147   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148      Set initial conditions
149    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150   PetscCall(FormInitialSolution(ts,X,&user));
151   PetscCall(TSSetSolution(ts,X));
152   dt   = 1e-10;                 /* Initial time step */
153   PetscCall(TSSetTimeStep(ts,dt));
154   PetscCall(TSGetAdapt(ts,&adapt));
155   PetscCall(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
156   PetscCall(TSSetMaxSNESFailures(ts,-1));            /* Retry step an unlimited number of times */
157 
158   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159      Set runtime options
160    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161   PetscCall(TSSetFromOptions(ts));
162 
163   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164      Set final conditions for sensitivities
165    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166   PetscCall(VecDuplicate(X,&lambda));
167   PetscCall(TSSetCostGradients(ts,1,&lambda,NULL));
168   PetscCall(VecSetValue(lambda,0,1.0,INSERT_VALUES));
169   PetscCall(VecAssemblyBegin(lambda));
170   PetscCall(VecAssemblyEnd(lambda));
171 
172   PetscCall(TSGetTrajectory(ts,&tj));
173   if (tj) {
174     PetscCall(TSTrajectorySetVariableNames(tj,(const char * const *)user.snames));
175     PetscCall(TSTrajectorySetTransform(tj,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user));
176   }
177 
178   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179      Pass information to graphical monitoring routine
180    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181   PetscCall(TSMonitorLGSetVariableNames(ts,(const char * const *)user.snames));
182   PetscCall(TSMonitorLGSetTransform(ts,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user));
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185      Solve ODE
186      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   PetscCall(TSSolve(ts,X));
188   PetscCall(TSGetSolveTime(ts,&ftime));
189   PetscCall(TSGetStepNumber(ts,&steps));
190   PetscCall(TSGetConvergedReason(ts,&reason));
191   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
192 
193   /* {
194     Vec                max;
195     PetscInt           i;
196     const PetscReal    *bmax;
197 
198     PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL));
199     if (max) {
200       PetscCall(VecGetArrayRead(max,&bmax));
201       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n"));
202       for (i=1; i<user.Nspec; i++) {
203         if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",user.snames[i],(double)bmax[i]));
204       }
205       PetscCall(VecRestoreArrayRead(max,&bmax));
206     }
207   }
208 
209   Vec y;
210   PetscCall(MassFractionToMoleFraction(&user,X,&y));
211   PetscCall(PrintSpecies(&user,y));
212   PetscCall(VecDestroy(&y)); */
213 
214   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215      Free work space.
216    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217   TC_reset();
218   PetscCall(PetscStrArrayDestroy(&user.snames));
219   PetscCall(MatDestroy(&J));
220   PetscCall(VecDestroy(&X));
221   PetscCall(VecDestroy(&lambda));
222   PetscCall(TSDestroy(&ts));
223   PetscCall(PetscFree3(user.tchemwork,user.Jdense,user.rows));
224   PetscCall(PetscFinalize());
225   return 0;
226 }
227 
228 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
229 {
230   User              user = (User)ptr;
231   PetscScalar       *f;
232   const PetscScalar *x;
233 
234   PetscFunctionBeginUser;
235   PetscCall(VecGetArrayRead(X,&x));
236   PetscCall(VecGetArray(F,&f));
237 
238   PetscCall(PetscArraycpy(user->tchemwork,x,user->Nspec+1));
239   user->tchemwork[0] *= user->Tini; /* Dimensionalize */
240   PetscCallTC(TC_getSrc(user->tchemwork,user->Nspec+1,f));
241   f[0] /= user->Tini;           /* Non-dimensionalize */
242 
243   PetscCall(VecRestoreArrayRead(X,&x));
244   PetscCall(VecRestoreArray(F,&f));
245   PetscFunctionReturn(0);
246 }
247 
248 static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
249 {
250   User              user = (User)ptr;
251   const PetscScalar *x;
252   PetscInt          M = user->Nspec+1,i;
253 
254   PetscFunctionBeginUser;
255   PetscCall(VecGetArrayRead(X,&x));
256   PetscCall(PetscArraycpy(user->tchemwork,x,user->Nspec+1));
257   PetscCall(VecRestoreArrayRead(X,&x));
258   user->tchemwork[0] *= user->Tini;  /* Dimensionalize temperature (first row) because that is what Tchem wants */
259   PetscCall(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1));
260 
261   for (i=0; i<M; i++) user->Jdense[i + 0*M] /= user->Tini; /* Non-dimensionalize first column */
262   for (i=0; i<M; i++) user->Jdense[0 + i*M] /= user->Tini; /* Non-dimensionalize first row */
263   for (i=0; i<M; i++) user->rows[i] = i;
264   PetscCall(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE));
265   PetscCall(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
266   PetscCall(MatZeroEntries(Pmat));
267   PetscCall(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES));
268 
269   PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
270   PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
271   if (Amat != Pmat) {
272     PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY));
273     PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY));
274   }
275   PetscFunctionReturn(0);
276 }
277 
278 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
279 {
280   PetscScalar    *x;
281   PetscInt       i;
282   Vec            y;
283   const PetscInt maxspecies = 10;
284   PetscInt       smax = maxspecies,mmax = maxspecies;
285   char           *names[maxspecies];
286   PetscReal      molefracs[maxspecies],sum;
287   PetscBool      flg;
288 
289   PetscFunctionBeginUser;
290   PetscCall(VecZeroEntries(X));
291   PetscCall(VecGetArray(X,&x));
292   x[0] = 1.0;  /* Non-dimensionalized by user->Tini */
293 
294   PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-initial_species",names,&smax,&flg));
295   PetscCheck(smax >= 2,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide at least two initial species");
296   PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-initial_mole",molefracs,&mmax,&flg));
297   PetscCheck(smax == mmax,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide same number of initial species %D as initial moles %D",smax,mmax);
298   sum = 0;
299   for (i=0; i<smax; i++) sum += molefracs[i];
300   for (i=0; i<smax; i++) molefracs[i] = molefracs[i]/sum;
301   for (i=0; i<smax; i++) {
302     int ispec = TC_getSpos(names[i], strlen(names[i]));
303     PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",names[i]);
304     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",i,names[i],molefracs[i]));
305     x[1+ispec] = molefracs[i];
306   }
307   for (i=0; i<smax; i++) {
308     PetscCall(PetscFree(names[i]));
309   }
310   PetscCall(VecRestoreArray(X,&x));
311   /* PetscCall(PrintSpecies((User)ctx,X)); */
312   PetscCall(MoleFractionToMassFraction((User)ctx,X,&y));
313   PetscCall(VecCopy(y,X));
314   PetscCall(VecDestroy(&y));
315   PetscFunctionReturn(0);
316 }
317 
318 /*
319    Converts the input vector which is in mass fractions (used by tchem) to mole fractions
320 */
321 PetscErrorCode MassFractionToMoleFraction(User user,Vec massf,Vec *molef)
322 {
323   PetscScalar       *mof;
324   const PetscScalar *maf;
325 
326   PetscFunctionBegin;
327   PetscCall(VecDuplicate(massf,molef));
328   PetscCall(VecGetArrayRead(massf,&maf));
329   PetscCall(VecGetArray(*molef,&mof));
330   mof[0] = maf[0]; /* copy over temperature */
331   TC_getMs2Ml((double*)maf+1,user->Nspec,mof+1);
332   PetscCall(VecRestoreArray(*molef,&mof));
333   PetscCall(VecRestoreArrayRead(massf,&maf));
334   PetscFunctionReturn(0);
335 }
336 
337 /*
338    Converts the input vector which is in mole fractions to mass fractions (used by tchem)
339 */
340 PetscErrorCode MoleFractionToMassFraction(User user,Vec molef,Vec *massf)
341 {
342   const PetscScalar *mof;
343   PetscScalar       *maf;
344 
345   PetscFunctionBegin;
346   PetscCall(VecDuplicate(molef,massf));
347   PetscCall(VecGetArrayRead(molef,&mof));
348   PetscCall(VecGetArray(*massf,&maf));
349   maf[0] = mof[0]; /* copy over temperature */
350   TC_getMl2Ms((double*)mof+1,user->Nspec,maf+1);
351   PetscCall(VecRestoreArrayRead(molef,&mof));
352   PetscCall(VecRestoreArray(*massf,&maf));
353   PetscFunctionReturn(0);
354 }
355 
356 PetscErrorCode ComputeMassConservation(Vec x,PetscReal *mass,void* ctx)
357 {
358   PetscFunctionBegin;
359   PetscCall(VecSum(x,mass));
360   PetscFunctionReturn(0);
361 }
362 
363 PetscErrorCode MonitorMassConservation(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
364 {
365   const PetscScalar  *T;
366   PetscReal          mass;
367 
368   PetscFunctionBegin;
369   PetscCall(ComputeMassConservation(x,&mass,ctx));
370   PetscCall(VecGetArrayRead(x,&T));
371   mass -= PetscAbsScalar(T[0]);
372   PetscCall(VecRestoreArrayRead(x,&T));
373   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g percent mass lost or gained %g\n",step,(double)time,(double)100.*(1.0 - mass)));
374   PetscFunctionReturn(0);
375 }
376 
377 PetscErrorCode MonitorTempature(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
378 {
379   User               user = (User) ctx;
380   const PetscScalar  *T;
381 
382   PetscFunctionBegin;
383   PetscCall(VecGetArrayRead(x,&T));
384   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g temperature %g\n",step,(double)time,(double)T[0]*user->Tini));
385   PetscCall(VecRestoreArrayRead(x,&T));
386   PetscFunctionReturn(0);
387 }
388 
389 /*
390    Prints out each species with its name
391 */
392 PETSC_UNUSED PetscErrorCode PrintSpecies(User user,Vec molef)
393 {
394   const PetscScalar *mof;
395   PetscInt          i,*idx,n = user->Nspec+1;
396 
397   PetscFunctionBegin;
398   PetscCall(PetscMalloc1(n,&idx));
399   for (i=0; i<n;i++) idx[i] = i;
400   PetscCall(VecGetArrayRead(molef,&mof));
401   PetscCall(PetscSortRealWithPermutation(n,mof,idx));
402   for (i=0; i<n; i++) {
403     PetscCall(PetscPrintf(PETSC_COMM_SELF,"%6s %g\n",user->snames[idx[n-i-1]],mof[idx[n-i-1]]));
404   }
405   PetscCall(PetscFree(idx));
406   PetscCall(VecRestoreArrayRead(molef,&mof));
407   PetscFunctionReturn(0);
408 }
409 
410 /*TEST
411     build:
412       requires: tchem
413 
414     test:
415       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
416       requires: !single
417       filter: grep -v iterations
418 
419 TEST*/
420