xref: /petsc/src/ts/tutorials/extchem.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
623c633725SBarry Smith #define CHKERRTC(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   PetscErrorCode    ierr;
72c4762a1bSJed Brown   PetscReal         ftime,dt;
73c4762a1bSJed 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];
74c4762a1bSJed Brown   const char        *periodic = "file://${PETSC_DIR}/${PETSC_ARCH}/share/periodictable.dat";
75c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
76c4762a1bSJed Brown   TSConvergedReason reason;
77c4762a1bSJed Brown   char              **snames,*names;
78c4762a1bSJed Brown   PetscInt          i;
79c4762a1bSJed Brown   TSTrajectory      tj;
80c4762a1bSJed Brown   PetscBool         flg = PETSC_FALSE,tflg = PETSC_FALSE,found;
81c4762a1bSJed Brown 
82*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
83c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");CHKERRQ(ierr);
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFileRetrieve(PETSC_COMM_WORLD,chemfile,lchemfile,PETSC_MAX_PATH_LEN,&found));
863c633725SBarry Smith   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",chemfile,lchemfile);
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFileRetrieve(PETSC_COMM_WORLD,thermofile,lthermofile,PETSC_MAX_PATH_LEN,&found));
893c633725SBarry Smith   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",thermofile,lthermofile);
90c4762a1bSJed Brown   user.pressure = 1.01325e5;    /* Pascal */
915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL));
92c4762a1bSJed Brown   user.Tini = 1000;             /* Kelvin */
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-monitor_mass","Monitor the total mass at each timestep","",flg,&flg,NULL));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-monitor_temp","Monitor the temperature each timestep","",tflg,&tflg,NULL));
96c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* tchem requires periodic table in current directory */
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFileRetrieve(PETSC_COMM_WORLD,periodic,lperiodic,PETSC_MAX_PATH_LEN,&found));
1003c633725SBarry Smith   PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot located required periodic table %s or local version %s",periodic,lperiodic);
101c4762a1bSJed Brown 
1025f80ce2aSJacob Faibussowitsch   CHKERRTC(TC_initChem(lchemfile, lthermofile, 0, 1.0));
103c4762a1bSJed Brown   TC_setThermoPres(user.pressure);
104c4762a1bSJed Brown   user.Nspec = TC_getNspec();
105c4762a1bSJed Brown   user.Nreac = TC_getNreac();
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown       Get names of all species in easy to use array
108c4762a1bSJed Brown   */
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcpy(names,"Temp"));
1112f613bf5SBarry Smith   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1((user.Nspec+2),&snames));
113c4762a1bSJed Brown   for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME;
114c4762a1bSJed Brown   snames[user.Nspec+1] = NULL;
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrArrayallocpy((const char *const *)snames,&user.snames));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(snames));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(names));
118c4762a1bSJed Brown 
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,user.Nspec+1,&X));
121c4762a1bSJed Brown 
1225f80ce2aSJacob Faibussowitsch   /* CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,PETSC_DECIDE,NULL,&J)); */
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,NULL,&J));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Create timestepping solver context
128c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSARKIMEX));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSARKIMEXSetType(ts,TSARKIMEX4));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   if (flg) {
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,MonitorMassConservation,NULL,NULL));
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown   if (tflg) {
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,MonitorTempature,&user,NULL));
141c4762a1bSJed Brown   }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   ftime = 1.0;
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,ftime));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown      Set initial conditions
149c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialSolution(ts,X,&user));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,X));
152c4762a1bSJed Brown   dt   = 1e-10;                 /* Initial time step */
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetAdapt(ts,&adapt));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSNESFailures(ts,-1));            /* Retry step an unlimited number of times */
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159c4762a1bSJed Brown      Set runtime options
160c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown      Set final conditions for sensitivities
165c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&lambda));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,&lambda,NULL));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValue(lambda,0,1.0,INSERT_VALUES));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(lambda));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(lambda));
171c4762a1bSJed Brown 
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTrajectory(ts,&tj));
173c4762a1bSJed Brown   if (tj) {
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectorySetVariableNames(tj,(const char * const *)user.snames));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(TSTrajectorySetTransform(tj,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user));
176c4762a1bSJed Brown   }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179c4762a1bSJed Brown      Pass information to graphical monitoring routine
180c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitorLGSetVariableNames(ts,(const char * const *)user.snames));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitorLGSetTransform(ts,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185c4762a1bSJed Brown      Solve ODE
186c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,X));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetConvergedReason(ts,&reason));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /* {
194c4762a1bSJed Brown     Vec                max;
195c4762a1bSJed Brown     PetscInt           i;
196c4762a1bSJed Brown     const PetscReal    *bmax;
197c4762a1bSJed Brown 
1985f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorEnvelopeGetBounds(ts,&max,NULL));
199c4762a1bSJed Brown     if (max) {
2005f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArrayRead(max,&bmax));
2015f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n"));
202c4762a1bSJed Brown       for (i=1; i<user.Nspec; i++) {
2035f80ce2aSJacob Faibussowitsch         if (bmax[i] > .01) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",user.snames[i],(double)bmax[i]));
204c4762a1bSJed Brown       }
2055f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArrayRead(max,&bmax));
206c4762a1bSJed Brown     }
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   Vec y;
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(MassFractionToMoleFraction(&user,X,&y));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(PrintSpecies(&user,y));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y)); */
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215c4762a1bSJed Brown      Free work space.
216c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217c4762a1bSJed Brown   TC_reset();
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrArrayDestroy(&user.snames));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&lambda));
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(user.tchemwork,user.Jdense,user.rows));
224*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
225*b122ec5aSJacob Faibussowitsch   return 0;
226c4762a1bSJed Brown }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
229c4762a1bSJed Brown {
230c4762a1bSJed Brown   User              user = (User)ptr;
231c4762a1bSJed Brown   PetscErrorCode    ierr;
232c4762a1bSJed Brown   PetscScalar       *f;
233c4762a1bSJed Brown   const PetscScalar *x;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   PetscFunctionBeginUser;
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
238c4762a1bSJed Brown 
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(user->tchemwork,x,user->Nspec+1));
240c4762a1bSJed Brown   user->tchemwork[0] *= user->Tini; /* Dimensionalize */
2415f80ce2aSJacob Faibussowitsch   CHKERRTC(TC_getSrc(user->tchemwork,user->Nspec+1,f));
242c4762a1bSJed Brown   f[0] /= user->Tini;           /* Non-dimensionalize */
243c4762a1bSJed Brown 
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
246c4762a1bSJed Brown   PetscFunctionReturn(0);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
250c4762a1bSJed Brown {
251c4762a1bSJed Brown   User              user = (User)ptr;
252c4762a1bSJed Brown   PetscErrorCode    ierr;
253c4762a1bSJed Brown   const PetscScalar *x;
254c4762a1bSJed Brown   PetscInt          M = user->Nspec+1,i;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   PetscFunctionBeginUser;
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(user->tchemwork,x,user->Nspec+1));
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
260c4762a1bSJed Brown   user->tchemwork[0] *= user->Tini;  /* Dimensionalize temperature (first row) because that is what Tchem wants */
2615f80ce2aSJacob Faibussowitsch   CHKERRQ(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   for (i=0; i<M; i++) user->Jdense[i + 0*M] /= user->Tini; /* Non-dimensionalize first column */
264c4762a1bSJed Brown   for (i=0; i<M; i++) user->Jdense[0 + i*M] /= user->Tini; /* Non-dimensionalize first row */
265c4762a1bSJed Brown   for (i=0; i<M; i++) user->rows[i] = i;
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(Pmat));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES));
270c4762a1bSJed Brown 
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
273c4762a1bSJed Brown   if (Amat != Pmat) {
2745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY));
2755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY));
276c4762a1bSJed Brown   }
277c4762a1bSJed Brown   PetscFunctionReturn(0);
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
281c4762a1bSJed Brown {
282c4762a1bSJed Brown   PetscScalar    *x;
283c4762a1bSJed Brown   PetscErrorCode ierr;
284c4762a1bSJed Brown   PetscInt       i;
285c4762a1bSJed Brown   Vec            y;
286c4762a1bSJed Brown   const PetscInt maxspecies = 10;
287c4762a1bSJed Brown   PetscInt       smax = maxspecies,mmax = maxspecies;
288c4762a1bSJed Brown   char           *names[maxspecies];
289c4762a1bSJed Brown   PetscReal      molefracs[maxspecies],sum;
290c4762a1bSJed Brown   PetscBool      flg;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   PetscFunctionBeginUser;
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(X));
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(X,&x));
295c4762a1bSJed Brown   x[0] = 1.0;  /* Non-dimensionalized by user->Tini */
296c4762a1bSJed Brown 
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetStringArray(NULL,NULL,"-initial_species",names,&smax,&flg));
2983c633725SBarry Smith   PetscCheck(smax >= 2,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide at least two initial species");
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetRealArray(NULL,NULL,"-initial_mole",molefracs,&mmax,&flg));
3003c633725SBarry Smith   PetscCheck(smax == mmax,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide same number of initial species %D as initial moles %D",smax,mmax);
301c4762a1bSJed Brown   sum = 0;
302c4762a1bSJed Brown   for (i=0; i<smax; i++) sum += molefracs[i];
303c4762a1bSJed Brown   for (i=0; i<smax; i++) molefracs[i] = molefracs[i]/sum;
304c4762a1bSJed Brown   for (i=0; i<smax; i++) {
305c4762a1bSJed Brown     int ispec = TC_getSpos(names[i], strlen(names[i]));
3063c633725SBarry Smith     PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",names[i]);
3075f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",i,names[i],molefracs[i]));
308c4762a1bSJed Brown     x[1+ispec] = molefracs[i];
309c4762a1bSJed Brown   }
310c4762a1bSJed Brown   for (i=0; i<smax; i++) {
3115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(names[i]));
312c4762a1bSJed Brown   }
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(X,&x));
3145f80ce2aSJacob Faibussowitsch   /* CHKERRQ(PrintSpecies((User)ctx,X)); */
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(MoleFractionToMassFraction((User)ctx,X,&y));
3165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(y,X));
3175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
318c4762a1bSJed Brown   PetscFunctionReturn(0);
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321c4762a1bSJed Brown /*
322c4762a1bSJed Brown    Converts the input vector which is in mass fractions (used by tchem) to mole fractions
323c4762a1bSJed Brown */
324c4762a1bSJed Brown PetscErrorCode MassFractionToMoleFraction(User user,Vec massf,Vec *molef)
325c4762a1bSJed Brown {
326c4762a1bSJed Brown   PetscErrorCode    ierr;
327c4762a1bSJed Brown   PetscScalar       *mof;
328c4762a1bSJed Brown   const PetscScalar *maf;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBegin;
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(massf,molef));
3325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(massf,&maf));
3335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(*molef,&mof));
334c4762a1bSJed Brown   mof[0] = maf[0]; /* copy over temperature */
335c4762a1bSJed Brown   TC_getMs2Ml((double*)maf+1,user->Nspec,mof+1);
3365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(*molef,&mof));
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(massf,&maf));
338c4762a1bSJed Brown   PetscFunctionReturn(0);
339c4762a1bSJed Brown }
340c4762a1bSJed Brown 
341c4762a1bSJed Brown /*
342c4762a1bSJed Brown    Converts the input vector which is in mole fractions to mass fractions (used by tchem)
343c4762a1bSJed Brown */
344c4762a1bSJed Brown PetscErrorCode MoleFractionToMassFraction(User user,Vec molef,Vec *massf)
345c4762a1bSJed Brown {
346c4762a1bSJed Brown   PetscErrorCode    ierr;
347c4762a1bSJed Brown   const PetscScalar *mof;
348c4762a1bSJed Brown   PetscScalar       *maf;
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   PetscFunctionBegin;
3515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(molef,massf));
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(molef,&mof));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(*massf,&maf));
354c4762a1bSJed Brown   maf[0] = mof[0]; /* copy over temperature */
355c4762a1bSJed Brown   TC_getMl2Ms((double*)mof+1,user->Nspec,maf+1);
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(molef,&mof));
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(*massf,&maf));
358c4762a1bSJed Brown   PetscFunctionReturn(0);
359c4762a1bSJed Brown }
360c4762a1bSJed Brown 
361c4762a1bSJed Brown PetscErrorCode ComputeMassConservation(Vec x,PetscReal *mass,void* ctx)
362c4762a1bSJed Brown {
363c4762a1bSJed Brown   PetscErrorCode ierr;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   PetscFunctionBegin;
3665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSum(x,mass));
367c4762a1bSJed Brown   PetscFunctionReturn(0);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown 
370c4762a1bSJed Brown PetscErrorCode MonitorMassConservation(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
371c4762a1bSJed Brown {
372c4762a1bSJed Brown   const PetscScalar  *T;
373c4762a1bSJed Brown   PetscReal          mass;
374c4762a1bSJed Brown   PetscErrorCode     ierr;
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   PetscFunctionBegin;
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(ComputeMassConservation(x,&mass,ctx));
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&T));
379c4762a1bSJed Brown   mass -= PetscAbsScalar(T[0]);
3805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&T));
3815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g percent mass lost or gained %g\n",step,(double)time,(double)100.*(1.0 - mass)));
382c4762a1bSJed Brown   PetscFunctionReturn(0);
383c4762a1bSJed Brown }
384c4762a1bSJed Brown 
385c4762a1bSJed Brown PetscErrorCode MonitorTempature(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx)
386c4762a1bSJed Brown {
387c4762a1bSJed Brown   User               user = (User) ctx;
388c4762a1bSJed Brown   const PetscScalar  *T;
389c4762a1bSJed Brown   PetscErrorCode     ierr;
390c4762a1bSJed Brown 
391c4762a1bSJed Brown   PetscFunctionBegin;
3925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&T));
3935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Timestep %D time %g temperature %g\n",step,(double)time,(double)T[0]*user->Tini));
3945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&T));
395c4762a1bSJed Brown   PetscFunctionReturn(0);
396c4762a1bSJed Brown }
397c4762a1bSJed Brown 
398c4762a1bSJed Brown /*
399c4762a1bSJed Brown    Prints out each species with its name
400c4762a1bSJed Brown */
401c4762a1bSJed Brown PETSC_UNUSED PetscErrorCode PrintSpecies(User user,Vec molef)
402c4762a1bSJed Brown {
403c4762a1bSJed Brown   PetscErrorCode    ierr;
404c4762a1bSJed Brown   const PetscScalar *mof;
405c4762a1bSJed Brown   PetscInt          i,*idx,n = user->Nspec+1;
406c4762a1bSJed Brown 
407c4762a1bSJed Brown   PetscFunctionBegin;
4085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&idx));
409c4762a1bSJed Brown   for (i=0; i<n;i++) idx[i] = i;
4105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(molef,&mof));
4115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortRealWithPermutation(n,mof,idx));
412c4762a1bSJed Brown   for (i=0; i<n; i++) {
4135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%6s %g\n",user->snames[idx[n-i-1]],mof[idx[n-i-1]]));
414c4762a1bSJed Brown   }
4155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
4165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(molef,&mof));
417c4762a1bSJed Brown   PetscFunctionReturn(0);
418c4762a1bSJed Brown }
419c4762a1bSJed Brown 
420c4762a1bSJed Brown /*TEST
421c4762a1bSJed Brown     build:
422c4762a1bSJed Brown       requires: tchem
423c4762a1bSJed Brown 
424c4762a1bSJed Brown     test:
425c4762a1bSJed 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
426c4762a1bSJed Brown       requires: !single
427c4762a1bSJed Brown       filter: grep -v iterations
428c4762a1bSJed Brown 
429c4762a1bSJed Brown TEST*/
430