xref: /petsc/src/ts/tutorials/extchem.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown static const char help[] = "Integrate chemistry using TChem.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #if defined(PETSC_HAVE_TCHEM)
6c4762a1bSJed Brown #if defined(MAX)
7c4762a1bSJed Brown #undef MAX
8c4762a1bSJed Brown #endif
9c4762a1bSJed Brown #if defined(MIN)
10c4762a1bSJed Brown #undef MIN
11c4762a1bSJed Brown #endif
12c4762a1bSJed Brown #  include <TC_params.h>
13c4762a1bSJed Brown #  include <TC_interface.h>
14c4762a1bSJed Brown #else
15c4762a1bSJed Brown #  error TChem is required for this example.  Reconfigure PETSc using --download-tchem.
16c4762a1bSJed Brown #endif
17c4762a1bSJed Brown /*
18c4762a1bSJed Brown     See extchem.example.1 for how to run an example
19c4762a1bSJed Brown 
20c4762a1bSJed Brown     See also h2_10sp.inp for another example
21c4762a1bSJed Brown 
22a5b23f4aSJose E. Roman     Determine sensitivity of final temperature on each variables initial conditions
23c4762a1bSJed Brown     -ts_dt 1.e-5 -ts_type cn -ts_adjoint_solve -ts_adjoint_view_solution draw
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     The solution for component i = 0 is the temperature.
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     The solution, i > 0, is the mass fraction, massf[i], of species i, i.e. mass of species i/ total mass of all species
28c4762a1bSJed Brown 
29c4762a1bSJed Brown     The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species
30c4762a1bSJed Brown         Define M[i] = mass per mole of species i then
31c4762a1bSJed Brown         molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j]))
32c4762a1bSJed Brown 
33c4762a1bSJed Brown     FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species.
34c4762a1bSJed Brown 
35c4762a1bSJed Brown     These are other data sets for other possible runs
36c4762a1bSJed Brown        https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/n_heptane_v3.1_therm.dat
37c4762a1bSJed Brown        https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/nc7_ver3.1_mech.txt
38c4762a1bSJed Brown 
39c4762a1bSJed Brown */
40c4762a1bSJed Brown typedef struct _User *User;
41c4762a1bSJed Brown struct _User {
42c4762a1bSJed Brown   PetscReal pressure;
43c4762a1bSJed Brown   int       Nspec;
44c4762a1bSJed Brown   int       Nreac;
45c4762a1bSJed Brown   PetscReal Tini;
46c4762a1bSJed Brown   double    *tchemwork;
47c4762a1bSJed Brown   double    *Jdense;        /* Dense array workspace where Tchem computes the Jacobian */
48c4762a1bSJed Brown   PetscInt  *rows;
49c4762a1bSJed Brown   char      **snames;
50c4762a1bSJed Brown };
51c4762a1bSJed Brown 
52c4762a1bSJed Brown static PetscErrorCode PrintSpecies(User,Vec);
53c4762a1bSJed Brown static PetscErrorCode MassFractionToMoleFraction(User,Vec,Vec*);
54c4762a1bSJed Brown static PetscErrorCode MoleFractionToMassFraction(User,Vec,Vec*);
55c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
56c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
57c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
58c4762a1bSJed Brown static PetscErrorCode ComputeMassConservation(Vec,PetscReal*,void*);
59c4762a1bSJed Brown static PetscErrorCode MonitorMassConservation(TS,PetscInt,PetscReal,Vec,void*);
60c4762a1bSJed Brown static PetscErrorCode MonitorTempature(TS,PetscInt,PetscReal,Vec,void*);
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch #define PetscCallTC(ierr) do {PetscCheck(!ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0)
63c4762a1bSJed Brown 
64c4762a1bSJed Brown int main(int argc,char **argv)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   TS                ts;         /* time integrator */
67c4762a1bSJed Brown   TSAdapt           adapt;
68c4762a1bSJed Brown   Vec               X,lambda;          /* solution vector */
69c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
70c4762a1bSJed Brown   PetscInt          steps;
71c4762a1bSJed Brown   PetscReal         ftime,dt;
72c4762a1bSJed Brown   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];
73c4762a1bSJed Brown   const char        *periodic = "file://${PETSC_DIR}/${PETSC_ARCH}/share/periodictable.dat";
74c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
75c4762a1bSJed Brown   TSConvergedReason reason;
76c4762a1bSJed Brown   char              **snames,*names;
77c4762a1bSJed Brown   PetscInt          i;
78c4762a1bSJed Brown   TSTrajectory      tj;
79c4762a1bSJed Brown   PetscBool         flg = PETSC_FALSE,tflg = PETSC_FALSE,found;
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
82*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,chemfile,lchemfile,PETSC_MAX_PATH_LEN,&found));
853c633725SBarry Smith   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",chemfile,lchemfile);
869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL));
879566063dSJacob Faibussowitsch   PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,thermofile,lthermofile,PETSC_MAX_PATH_LEN,&found));
883c633725SBarry Smith   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",thermofile,lthermofile);
89c4762a1bSJed Brown   user.pressure = 1.01325e5;    /* Pascal */
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL));
91c4762a1bSJed Brown   user.Tini = 1000;             /* Kelvin */
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL));
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitor_mass","Monitor the total mass at each timestep","",flg,&flg,NULL));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitor_temp","Monitor the temperature each timestep","",tflg,&tflg,NULL));
95*d0609cedSBarry Smith   PetscOptionsEnd();
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* tchem requires periodic table in current directory */
989566063dSJacob Faibussowitsch   PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,periodic,lperiodic,PETSC_MAX_PATH_LEN,&found));
993c633725SBarry Smith   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot located required periodic table %s or local version %s",periodic,lperiodic);
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch   PetscCallTC(TC_initChem(lchemfile, lthermofile, 0, 1.0));
102c4762a1bSJed Brown   TC_setThermoPres(user.pressure);
103c4762a1bSJed Brown   user.Nspec = TC_getNspec();
104c4762a1bSJed Brown   user.Nreac = TC_getNreac();
105c4762a1bSJed Brown   /*
106c4762a1bSJed Brown       Get names of all species in easy to use array
107c4762a1bSJed Brown   */
1089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names));
1099566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(names,"Temp"));
1102f613bf5SBarry Smith   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
1119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((user.Nspec+2),&snames));
112c4762a1bSJed Brown   for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME;
113c4762a1bSJed Brown   snames[user.Nspec+1] = NULL;
1149566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayallocpy((const char *const *)snames,&user.snames));
1159566063dSJacob Faibussowitsch   PetscCall(PetscFree(snames));
1169566063dSJacob Faibussowitsch   PetscCall(PetscFree(names));
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows));
1199566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.Nspec+1,&X));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   /* PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,PETSC_DECIDE,NULL,&J)); */
1229566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,NULL,&J));
1239566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Create timestepping solver context
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1289566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1299566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSARKIMEX));
1309566063dSJacob Faibussowitsch   PetscCall(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
1319566063dSJacob Faibussowitsch   PetscCall(TSARKIMEXSetType(ts,TSARKIMEX4));
1329566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
1339566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   if (flg) {
1369566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,MonitorMassConservation,NULL,NULL));
137c4762a1bSJed Brown   }
138c4762a1bSJed Brown   if (tflg) {
1399566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,MonitorTempature,&user,NULL));
140c4762a1bSJed Brown   }
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   ftime = 1.0;
1439566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
1449566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown      Set initial conditions
148c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1499566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts,X,&user));
1509566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,X));
151c4762a1bSJed Brown   dt   = 1e-10;                 /* Initial time step */
1529566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
1539566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
1549566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
1559566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts,-1));            /* Retry step an unlimited number of times */
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown      Set runtime options
159c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1609566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown      Set final conditions for sensitivities
164c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X,&lambda));
1669566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,&lambda,NULL));
1679566063dSJacob Faibussowitsch   PetscCall(VecSetValue(lambda,0,1.0,INSERT_VALUES));
1689566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(lambda));
1699566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(lambda));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts,&tj));
172c4762a1bSJed Brown   if (tj) {
1739566063dSJacob Faibussowitsch     PetscCall(TSTrajectorySetVariableNames(tj,(const char * const *)user.snames));
1749566063dSJacob Faibussowitsch     PetscCall(TSTrajectorySetTransform(tj,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user));
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown      Pass information to graphical monitoring routine
179c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1809566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGSetVariableNames(ts,(const char * const *)user.snames));
1819566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGSetTransform(ts,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184c4762a1bSJed Brown      Solve ODE
185c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1869566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,X));
1879566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
1889566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
1899566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts,&reason));
1909566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* {
193c4762a1bSJed Brown     Vec                max;
194c4762a1bSJed Brown     PetscInt           i;
195c4762a1bSJed Brown     const PetscReal    *bmax;
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch     PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL));
198c4762a1bSJed Brown     if (max) {
1999566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(max,&bmax));
2009566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n"));
201c4762a1bSJed Brown       for (i=1; i<user.Nspec; i++) {
2029566063dSJacob Faibussowitsch         if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",user.snames[i],(double)bmax[i]));
203c4762a1bSJed Brown       }
2049566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(max,&bmax));
205c4762a1bSJed Brown     }
206c4762a1bSJed Brown   }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   Vec y;
2099566063dSJacob Faibussowitsch   PetscCall(MassFractionToMoleFraction(&user,X,&y));
2109566063dSJacob Faibussowitsch   PetscCall(PrintSpecies(&user,y));
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y)); */
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown      Free work space.
215c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216c4762a1bSJed Brown   TC_reset();
2179566063dSJacob Faibussowitsch   PetscCall(PetscStrArrayDestroy(&user.snames));
2189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
2199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
2209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda));
2219566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2229566063dSJacob Faibussowitsch   PetscCall(PetscFree3(user.tchemwork,user.Jdense,user.rows));
2239566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
224b122ec5aSJacob Faibussowitsch   return 0;
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
227c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
228c4762a1bSJed Brown {
229c4762a1bSJed Brown   User              user = (User)ptr;
230c4762a1bSJed Brown   PetscScalar       *f;
231c4762a1bSJed Brown   const PetscScalar *x;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   PetscFunctionBeginUser;
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
236c4762a1bSJed Brown 
2379566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(user->tchemwork,x,user->Nspec+1));
238c4762a1bSJed Brown   user->tchemwork[0] *= user->Tini; /* Dimensionalize */
2399566063dSJacob Faibussowitsch   PetscCallTC(TC_getSrc(user->tchemwork,user->Nspec+1,f));
240c4762a1bSJed Brown   f[0] /= user->Tini;           /* Non-dimensionalize */
241c4762a1bSJed Brown 
2429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
2439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
244c4762a1bSJed Brown   PetscFunctionReturn(0);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
248c4762a1bSJed Brown {
249c4762a1bSJed Brown   User              user = (User)ptr;
250c4762a1bSJed Brown   const PetscScalar *x;
251c4762a1bSJed Brown   PetscInt          M = user->Nspec+1,i;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   PetscFunctionBeginUser;
2549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2559566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(user->tchemwork,x,user->Nspec+1));
2569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
257c4762a1bSJed Brown   user->tchemwork[0] *= user->Tini;  /* Dimensionalize temperature (first row) because that is what Tchem wants */
2589566063dSJacob Faibussowitsch   PetscCall(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   for (i=0; i<M; i++) user->Jdense[i + 0*M] /= user->Tini; /* Non-dimensionalize first column */
261c4762a1bSJed Brown   for (i=0; i<M; i++) user->Jdense[0 + i*M] /= user->Tini; /* Non-dimensionalize first row */
262c4762a1bSJed Brown   for (i=0; i<M; i++) user->rows[i] = i;
2639566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE));
2649566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
2659566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Pmat));
2669566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES));
267c4762a1bSJed Brown 
2689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
2699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
270c4762a1bSJed Brown   if (Amat != Pmat) {
2719566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY));
2729566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY));
273c4762a1bSJed Brown   }
274c4762a1bSJed Brown   PetscFunctionReturn(0);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
277c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
278c4762a1bSJed Brown {
279c4762a1bSJed Brown   PetscScalar    *x;
280c4762a1bSJed Brown   PetscInt       i;
281c4762a1bSJed Brown   Vec            y;
282c4762a1bSJed Brown   const PetscInt maxspecies = 10;
283c4762a1bSJed Brown   PetscInt       smax = maxspecies,mmax = maxspecies;
284c4762a1bSJed Brown   char           *names[maxspecies];
285c4762a1bSJed Brown   PetscReal      molefracs[maxspecies],sum;
286c4762a1bSJed Brown   PetscBool      flg;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionBeginUser;
2899566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
2909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
291c4762a1bSJed Brown   x[0] = 1.0;  /* Non-dimensionalized by user->Tini */
292c4762a1bSJed Brown 
2939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-initial_species",names,&smax,&flg));
2943c633725SBarry Smith   PetscCheck(smax >= 2,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide at least two initial species");
2959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-initial_mole",molefracs,&mmax,&flg));
2963c633725SBarry Smith   PetscCheck(smax == mmax,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide same number of initial species %D as initial moles %D",smax,mmax);
297c4762a1bSJed Brown   sum = 0;
298c4762a1bSJed Brown   for (i=0; i<smax; i++) sum += molefracs[i];
299c4762a1bSJed Brown   for (i=0; i<smax; i++) molefracs[i] = molefracs[i]/sum;
300c4762a1bSJed Brown   for (i=0; i<smax; i++) {
301c4762a1bSJed Brown     int ispec = TC_getSpos(names[i], strlen(names[i]));
3023c633725SBarry Smith     PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",names[i]);
3039566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",i,names[i],molefracs[i]));
304c4762a1bSJed Brown     x[1+ispec] = molefracs[i];
305c4762a1bSJed Brown   }
306c4762a1bSJed Brown   for (i=0; i<smax; i++) {
3079566063dSJacob Faibussowitsch     PetscCall(PetscFree(names[i]));
308c4762a1bSJed Brown   }
3099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
3109566063dSJacob Faibussowitsch   /* PetscCall(PrintSpecies((User)ctx,X)); */
3119566063dSJacob Faibussowitsch   PetscCall(MoleFractionToMassFraction((User)ctx,X,&y));
3129566063dSJacob Faibussowitsch   PetscCall(VecCopy(y,X));
3139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
314c4762a1bSJed Brown   PetscFunctionReturn(0);
315c4762a1bSJed Brown }
316c4762a1bSJed Brown 
317c4762a1bSJed Brown /*
318c4762a1bSJed Brown    Converts the input vector which is in mass fractions (used by tchem) to mole fractions
319c4762a1bSJed Brown */
320c4762a1bSJed Brown PetscErrorCode MassFractionToMoleFraction(User user,Vec massf,Vec *molef)
321c4762a1bSJed Brown {
322c4762a1bSJed Brown   PetscScalar       *mof;
323c4762a1bSJed Brown   const PetscScalar *maf;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   PetscFunctionBegin;
3269566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(massf,molef));
3279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(massf,&maf));
3289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*molef,&mof));
329c4762a1bSJed Brown   mof[0] = maf[0]; /* copy over temperature */
330c4762a1bSJed Brown   TC_getMs2Ml((double*)maf+1,user->Nspec,mof+1);
3319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*molef,&mof));
3329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(massf,&maf));
333c4762a1bSJed Brown   PetscFunctionReturn(0);
334c4762a1bSJed Brown }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown /*
337c4762a1bSJed Brown    Converts the input vector which is in mole fractions to mass fractions (used by tchem)
338c4762a1bSJed Brown */
339c4762a1bSJed Brown PetscErrorCode MoleFractionToMassFraction(User user,Vec molef,Vec *massf)
340c4762a1bSJed Brown {
341c4762a1bSJed Brown   const PetscScalar *mof;
342c4762a1bSJed Brown   PetscScalar       *maf;
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   PetscFunctionBegin;
3459566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(molef,massf));
3469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(molef,&mof));
3479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*massf,&maf));
348c4762a1bSJed Brown   maf[0] = mof[0]; /* copy over temperature */
349c4762a1bSJed Brown   TC_getMl2Ms((double*)mof+1,user->Nspec,maf+1);
3509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(molef,&mof));
3519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*massf,&maf));
352c4762a1bSJed Brown   PetscFunctionReturn(0);
353c4762a1bSJed Brown }
354c4762a1bSJed Brown 
355c4762a1bSJed Brown PetscErrorCode ComputeMassConservation(Vec x,PetscReal *mass,void* ctx)
356c4762a1bSJed Brown {
357c4762a1bSJed Brown   PetscFunctionBegin;
3589566063dSJacob Faibussowitsch   PetscCall(VecSum(x,mass));
359c4762a1bSJed Brown   PetscFunctionReturn(0);
360c4762a1bSJed Brown }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown PetscErrorCode MonitorMassConservation(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
363c4762a1bSJed Brown {
364c4762a1bSJed Brown   const PetscScalar  *T;
365c4762a1bSJed Brown   PetscReal          mass;
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCall(ComputeMassConservation(x,&mass,ctx));
3699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&T));
370c4762a1bSJed Brown   mass -= PetscAbsScalar(T[0]);
3719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&T));
3729566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g percent mass lost or gained %g\n",step,(double)time,(double)100.*(1.0 - mass)));
373c4762a1bSJed Brown   PetscFunctionReturn(0);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown 
376c4762a1bSJed Brown PetscErrorCode MonitorTempature(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
377c4762a1bSJed Brown {
378c4762a1bSJed Brown   User               user = (User) ctx;
379c4762a1bSJed Brown   const PetscScalar  *T;
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   PetscFunctionBegin;
3829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&T));
3839566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g temperature %g\n",step,(double)time,(double)T[0]*user->Tini));
3849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&T));
385c4762a1bSJed Brown   PetscFunctionReturn(0);
386c4762a1bSJed Brown }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown /*
389c4762a1bSJed Brown    Prints out each species with its name
390c4762a1bSJed Brown */
391c4762a1bSJed Brown PETSC_UNUSED PetscErrorCode PrintSpecies(User user,Vec molef)
392c4762a1bSJed Brown {
393c4762a1bSJed Brown   const PetscScalar *mof;
394c4762a1bSJed Brown   PetscInt          i,*idx,n = user->Nspec+1;
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&idx));
398c4762a1bSJed Brown   for (i=0; i<n;i++) idx[i] = i;
3999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(molef,&mof));
4009566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithPermutation(n,mof,idx));
401c4762a1bSJed Brown   for (i=0; i<n; i++) {
4029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"%6s %g\n",user->snames[idx[n-i-1]],mof[idx[n-i-1]]));
403c4762a1bSJed Brown   }
4049566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
4059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(molef,&mof));
406c4762a1bSJed Brown   PetscFunctionReturn(0);
407c4762a1bSJed Brown }
408c4762a1bSJed Brown 
409c4762a1bSJed Brown /*TEST
410c4762a1bSJed Brown     build:
411c4762a1bSJed Brown       requires: tchem
412c4762a1bSJed Brown 
413c4762a1bSJed Brown     test:
414c4762a1bSJed Brown       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
415c4762a1bSJed Brown       requires: !single
416c4762a1bSJed Brown       filter: grep -v iterations
417c4762a1bSJed Brown 
418c4762a1bSJed Brown TEST*/
419