xref: /petsc/src/ts/tutorials/extchemfield.c (revision dd39110ba674af879e763a5939b18b2b35b4c87f)
1c4762a1bSJed Brown static const char help[] = "Integrate chemistry using TChem.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #if defined(PETSC_HAVE_TCHEM)
7c4762a1bSJed Brown #if defined(MAX)
8c4762a1bSJed Brown #undef MAX
9c4762a1bSJed Brown #endif
10c4762a1bSJed Brown #if defined(MIN)
11c4762a1bSJed Brown #undef MIN
12c4762a1bSJed Brown #endif
13c4762a1bSJed Brown #  include <TC_params.h>
14c4762a1bSJed Brown #  include <TC_interface.h>
15c4762a1bSJed Brown #else
16c4762a1bSJed Brown #  error TChem is required for this example.  Reconfigure PETSc using --download-tchem.
17c4762a1bSJed Brown #endif
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown 
20c4762a1bSJed Brown     This is an extension of extchem.c to solve the reaction equations independently in each cell of a one dimensional field
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     Obtain the three files into this directory
23c4762a1bSJed Brown 
24c4762a1bSJed Brown        curl http://combustion.berkeley.edu/gri_mech/version30/files30/grimech30.dat > chem.inp
25c4762a1bSJed Brown        curl http://combustion.berkeley.edu/gri_mech/version30/files30/thermo30.dat > therm.dat
26c4762a1bSJed Brown        cp $PETSC_DIR/$PETSC_ARCH/externalpackages/tchem/data/periodictable.dat .
27c4762a1bSJed Brown 
28c4762a1bSJed Brown     Run with
29c4762a1bSJed Brown      ./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
30c4762a1bSJed Brown 
31c4762a1bSJed Brown      Options for visualizing the solution:
32c4762a1bSJed Brown         Watch certain variables in each cell evolve with time
33c4762a1bSJed Brown         -draw_solution 1 -ts_monitor_lg_solution_variables Temp,H2,O2,H2O,CH4,CO,CO2,C2H2,N2 -lg_use_markers false  -draw_pause -2
34c4762a1bSJed Brown 
35c4762a1bSJed Brown         Watch certain variables in all cells evolve with time
36c4762a1bSJed Brown         -da_refine 4 -ts_monitor_draw_solution -draw_fields_by_name Temp,H2 -draw_vec_mark_points  -draw_pause -2
37c4762a1bSJed Brown 
38c4762a1bSJed Brown         Keep the initial temperature distribution as one monitors the current temperature distribution
39c4762a1bSJed Brown         -ts_monitor_draw_solution_initial -draw_bounds .9,1.7 -draw_fields_by_name Temp
40c4762a1bSJed Brown 
41c4762a1bSJed Brown         Save the images in a .gif (movie) file
42c4762a1bSJed Brown         -draw_save -draw_save_single_file
43c4762a1bSJed Brown 
44a5b23f4aSJose E. Roman         Compute the sensitivies of the solution of the first temperature on the initial conditions
45c4762a1bSJed Brown         -ts_adjoint_solve  -ts_dt 1.e-5 -ts_type cn -ts_adjoint_view_solution draw
46c4762a1bSJed Brown 
47c4762a1bSJed Brown         Turn off diffusion
48c4762a1bSJed Brown         -diffusion no
49c4762a1bSJed Brown 
50c4762a1bSJed Brown         Turn off reactions
51c4762a1bSJed Brown         -reactions no
52c4762a1bSJed Brown 
53c4762a1bSJed Brown     The solution for component i = 0 is the temperature.
54c4762a1bSJed Brown 
55c4762a1bSJed 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
56c4762a1bSJed Brown 
57c4762a1bSJed Brown     The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species
58c4762a1bSJed Brown         Define M[i] = mass per mole of species i then
59c4762a1bSJed Brown         molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j]))
60c4762a1bSJed Brown 
61c4762a1bSJed Brown     FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species.
62c4762a1bSJed Brown 
63c4762a1bSJed Brown */
64c4762a1bSJed Brown typedef struct _User *User;
65c4762a1bSJed Brown struct _User {
66c4762a1bSJed Brown   PetscReal pressure;
67c4762a1bSJed Brown   int       Nspec;
68c4762a1bSJed Brown   int       Nreac;
69c4762a1bSJed Brown   PetscReal Tini,dx;
70c4762a1bSJed Brown   PetscReal diffus;
71c4762a1bSJed Brown   DM        dm;
72c4762a1bSJed Brown   PetscBool diffusion,reactions;
73c4762a1bSJed Brown   double    *tchemwork;
74c4762a1bSJed Brown   double    *Jdense;        /* Dense array workspace where Tchem computes the Jacobian */
75c4762a1bSJed Brown   PetscInt  *rows;
76c4762a1bSJed Brown };
77c4762a1bSJed Brown 
78c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS,User,PetscInt);
79c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
80c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
81c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch #define PetscCallTC(ierr) do {PetscCheck(!ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0)
84c4762a1bSJed Brown 
85c4762a1bSJed Brown int main(int argc,char **argv)
86c4762a1bSJed Brown {
87c4762a1bSJed Brown   TS                ts;         /* time integrator */
88c4762a1bSJed Brown   TSAdapt           adapt;
89c4762a1bSJed Brown   Vec               X;          /* solution vector */
90c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
91c4762a1bSJed Brown   PetscInt          steps,ncells,xs,xm,i;
92c4762a1bSJed Brown   PetscErrorCode    ierr;
93c4762a1bSJed Brown   PetscReal         ftime,dt;
94c4762a1bSJed Brown   char              chemfile[PETSC_MAX_PATH_LEN] = "chem.inp",thermofile[PETSC_MAX_PATH_LEN] = "therm.dat";
95c4762a1bSJed Brown   struct _User      user;
96c4762a1bSJed Brown   TSConvergedReason reason;
97c4762a1bSJed Brown   PetscBool         showsolutions = PETSC_FALSE;
98c4762a1bSJed Brown   char              **snames,*names;
99c4762a1bSJed Brown   Vec               lambda;     /* used with TSAdjoint for sensitivities */
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
1029566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");PetscCall(ierr);
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL));
105c4762a1bSJed Brown   user.pressure = 1.01325e5;    /* Pascal */
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL));
107c4762a1bSJed Brown   user.Tini   = 1550;
1089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL));
109c4762a1bSJed Brown   user.diffus = 100;
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL));
112c4762a1bSJed Brown   user.diffusion = PETSC_TRUE;
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL));
114c4762a1bSJed Brown   user.reactions = PETSC_TRUE;
1159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL));
1169566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCallTC(TC_initChem(chemfile, thermofile, 0, 1.0));
119c4762a1bSJed Brown   user.Nspec = TC_getNspec();
120c4762a1bSJed Brown   user.Nreac = TC_getNreac();
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm));
1239566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
1249566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
1259566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user.dm,NULL,&ncells,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
126c4762a1bSJed Brown   user.dx = 1.0/ncells;  /* Set the coordinates of the cell centers; note final ghost cell is at x coordinate 1.0 */
1279566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(user.dm,0.0,1.0,0.0,1.0,0.0,1.0));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* set the names of each field in the DMDA based on the species name */
1309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names));
1319566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(names,"Temp"));
1322f613bf5SBarry Smith   TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);
1339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((user.Nspec+2),&snames));
134c4762a1bSJed Brown   for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME;
135c4762a1bSJed Brown   snames[user.Nspec+1] = NULL;
1369566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldNames(user.dm,(const char * const *)snames));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFree(snames));
1389566063dSJacob Faibussowitsch   PetscCall(PetscFree(names));
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.dm,&J));
1419566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm,&X));
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146c4762a1bSJed Brown      Create timestepping solver context
147c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1489566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1499566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,user.dm));
1509566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSARKIMEX));
1519566063dSJacob Faibussowitsch   PetscCall(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
1529566063dSJacob Faibussowitsch   PetscCall(TSARKIMEXSetType(ts,TSARKIMEX4));
1539566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
1549566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   ftime = 1.0;
1579566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
1589566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161c4762a1bSJed Brown      Set initial conditions
162c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1639566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts,X,&user));
1649566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,X));
165c4762a1bSJed Brown   dt   = 1e-10;                 /* Initial time step */
1669566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
1679566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
1689566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
1699566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts,-1));            /* Retry step an unlimited number of times */
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172c4762a1bSJed Brown      Pass information to graphical monitoring routine
173c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174c4762a1bSJed Brown   if (showsolutions) {
1759566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL));
176c4762a1bSJed Brown     for (i=xs;i<xs+xm;i++) {
1779566063dSJacob Faibussowitsch       PetscCall(MonitorCell(ts,&user,i));
178c4762a1bSJed Brown     }
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182c4762a1bSJed Brown      Set runtime options
183c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1849566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown      Set final conditions for sensitivities
188c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1899566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm,&lambda));
1909566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,&lambda,NULL));
1919566063dSJacob Faibussowitsch   PetscCall(VecSetValue(lambda,0,1.0,INSERT_VALUES));
1929566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(lambda));
1939566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(lambda));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196c4762a1bSJed Brown      Solve ODE
197c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1989566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,X));
1999566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
2009566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
2019566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts,&reason));
2029566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   {
205c4762a1bSJed Brown     Vec                max;
206c4762a1bSJed Brown     const char * const *names;
207c4762a1bSJed Brown     PetscInt           i;
208c4762a1bSJed Brown     const PetscReal    *bmax;
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch     PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL));
211c4762a1bSJed Brown     if (max) {
2129566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGGetVariableNames(ts,&names));
213c4762a1bSJed Brown       if (names) {
2149566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(max,&bmax));
2159566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n"));
216c4762a1bSJed Brown         for (i=1; i<user.Nspec; i++) {
2179566063dSJacob Faibussowitsch           if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i]));
218c4762a1bSJed Brown         }
2199566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(max,&bmax));
220c4762a1bSJed Brown       }
221c4762a1bSJed Brown     }
222c4762a1bSJed Brown   }
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225c4762a1bSJed Brown      Free work space.
226c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227c4762a1bSJed Brown   TC_reset();
2289566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
2299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
2309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
2319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda));
2329566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2339566063dSJacob Faibussowitsch   PetscCall(PetscFree3(user.tchemwork,user.Jdense,user.rows));
2349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
235b122ec5aSJacob Faibussowitsch   return 0;
236c4762a1bSJed Brown }
237c4762a1bSJed Brown 
238c4762a1bSJed Brown /*
239c4762a1bSJed Brown    Applies the second order centered difference diffusion operator on a one dimensional periodic domain
240c4762a1bSJed Brown */
241c4762a1bSJed Brown static PetscErrorCode FormDiffusionFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
242c4762a1bSJed Brown {
243c4762a1bSJed Brown   User              user = (User)ptr;
244c4762a1bSJed Brown   PetscScalar       **f;
245c4762a1bSJed Brown   const PetscScalar **x;
246c4762a1bSJed Brown   DM                dm;
247c4762a1bSJed Brown   PetscInt          i,xs,xm,j,dof;
248c4762a1bSJed Brown   Vec               Xlocal;
249c4762a1bSJed Brown   PetscReal         idx;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   PetscFunctionBeginUser;
2529566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
2539566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
2549566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm,&Xlocal));
2559566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal));
2569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal));
2579566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOFRead(dm,Xlocal,&x));
2589566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(dm,F,&f));
2599566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   idx = 1.0*user->diffus/user->dx;
262c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
263c4762a1bSJed Brown     for (j=0; j<dof; j++) {
264c4762a1bSJed Brown       f[i][j] += idx*(x[i+1][j] - 2.0*x[i][j] + x[i-1][j]);
265c4762a1bSJed Brown     }
266c4762a1bSJed Brown   }
2679566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x));
2689566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(dm,F,&f));
2699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm,&Xlocal));
270c4762a1bSJed Brown   PetscFunctionReturn(0);
271c4762a1bSJed Brown }
272c4762a1bSJed Brown 
273c4762a1bSJed Brown /*
274c4762a1bSJed Brown    Produces the second order centered difference diffusion operator on a one dimensional periodic domain
275c4762a1bSJed Brown */
276c4762a1bSJed Brown static PetscErrorCode FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
277c4762a1bSJed Brown {
278c4762a1bSJed Brown   User              user = (User)ptr;
279c4762a1bSJed Brown   DM                dm;
280c4762a1bSJed Brown   PetscInt          i,xs,xm,j,dof;
281c4762a1bSJed Brown   PetscReal         idx,values[3];
282c4762a1bSJed Brown   MatStencil        row,col[3];
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   PetscFunctionBeginUser;
2859566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
2869566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
2879566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   idx = 1.0*user->diffus/user->dx;
290c4762a1bSJed Brown   values[0] = idx;
291c4762a1bSJed Brown   values[1] = -2.0*idx;
292c4762a1bSJed Brown   values[2] = idx;
293c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
294c4762a1bSJed Brown     for (j=0; j<dof; j++) {
295c4762a1bSJed Brown       row.i = i;      row.c = j;
296c4762a1bSJed Brown       col[0].i = i-1; col[0].c = j;
297c4762a1bSJed Brown       col[1].i = i;   col[1].c = j;
298c4762a1bSJed Brown       col[2].i = i+1; col[2].c = j;
2999566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES));
300c4762a1bSJed Brown     }
301c4762a1bSJed Brown   }
3029566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
3039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
304c4762a1bSJed Brown   PetscFunctionReturn(0);
305c4762a1bSJed Brown }
306c4762a1bSJed Brown 
307c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
308c4762a1bSJed Brown {
309c4762a1bSJed Brown   User              user = (User)ptr;
310c4762a1bSJed Brown   PetscScalar       **f;
311c4762a1bSJed Brown   const PetscScalar **x;
312c4762a1bSJed Brown   DM                dm;
313c4762a1bSJed Brown   PetscInt          i,xs,xm;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   PetscFunctionBeginUser;
316c4762a1bSJed Brown   if (user->reactions) {
3179566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts,&dm));
3189566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayDOFRead(dm,X,&x));
3199566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayDOF(dm,F,&f));
3209566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
321c4762a1bSJed Brown 
322c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
3239566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1));
324c4762a1bSJed Brown       user->tchemwork[0] *= user->Tini; /* Dimensionalize */
3259566063dSJacob Faibussowitsch       PetscCallTC(TC_getSrc(user->tchemwork,user->Nspec+1,f[i]));
326c4762a1bSJed Brown       f[i][0] /= user->Tini;           /* Non-dimensionalize */
327c4762a1bSJed Brown     }
328c4762a1bSJed Brown 
3299566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayDOFRead(dm,X,&x));
3309566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayDOF(dm,F,&f));
331c4762a1bSJed Brown   } else {
3329566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown   if (user->diffusion) {
3359566063dSJacob Faibussowitsch     PetscCall(FormDiffusionFunction(ts,t,X,F,ptr));
336c4762a1bSJed Brown   }
337c4762a1bSJed Brown   PetscFunctionReturn(0);
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr)
341c4762a1bSJed Brown {
342c4762a1bSJed Brown   User              user = (User)ptr;
343c4762a1bSJed Brown   const PetscScalar **x;
344c4762a1bSJed Brown   PetscInt          M = user->Nspec+1,i,j,xs,xm;
345c4762a1bSJed Brown   DM                dm;
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   PetscFunctionBeginUser;
348c4762a1bSJed Brown   if (user->reactions) {
3499566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts,&dm));
3509566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(Pmat));
3519566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE));
3529566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
3539566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayDOFRead(dm,X,&x));
3549566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
3579566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1));
358c4762a1bSJed Brown       user->tchemwork[0] *= user->Tini;  /* Dimensionalize temperature (first row) because that is what Tchem wants */
3599566063dSJacob Faibussowitsch       PetscCall(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown       for (j=0; j<M; j++) user->Jdense[j + 0*M] /= user->Tini; /* Non-dimensionalize first column */
362c4762a1bSJed Brown       for (j=0; j<M; j++) user->Jdense[0 + j*M] /= user->Tini; /* Non-dimensionalize first row */
363c4762a1bSJed Brown       for (j=0; j<M; j++) user->rows[j] = i*M+j;
3649566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES));
365c4762a1bSJed Brown     }
3669566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayDOFRead(dm,X,&x));
3679566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY));
3689566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY));
369c4762a1bSJed Brown   } else {
3709566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(Pmat));
371c4762a1bSJed Brown   }
372c4762a1bSJed Brown   if (user->diffusion) {
3739566063dSJacob Faibussowitsch     PetscCall(FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr));
374c4762a1bSJed Brown   }
375c4762a1bSJed Brown   if (Amat != Pmat) {
3769566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY));
3779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY));
378c4762a1bSJed Brown   }
379c4762a1bSJed Brown   PetscFunctionReturn(0);
380c4762a1bSJed Brown }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
383c4762a1bSJed Brown {
384c4762a1bSJed Brown   PetscScalar    **x,*xc;
385c4762a1bSJed Brown   struct {const char *name; PetscReal massfrac;} initial[] = {
386c4762a1bSJed Brown     {"CH4", 0.0948178320887},
387c4762a1bSJed Brown     {"O2", 0.189635664177},
388c4762a1bSJed Brown     {"N2", 0.706766236705},
389c4762a1bSJed Brown     {"AR", 0.00878026702874}
390c4762a1bSJed Brown   };
391c4762a1bSJed Brown   PetscInt       i,j,xs,xm;
392c4762a1bSJed Brown   DM             dm;
393c4762a1bSJed Brown 
394c4762a1bSJed Brown   PetscFunctionBeginUser;
3959566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
3969566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&dm));
3979566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL));
398c4762a1bSJed Brown 
3999566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateArray(dm,&xc));
4009566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(dm,X,&x));
401c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
402c4762a1bSJed Brown     x[i][0] = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[i]);  /* Non-dimensionalized by user->Tini */
403*dd39110bSPierre Jolivet     for (j=0; j<PETSC_STATIC_ARRAY_LENGTH(initial); j++) {
404c4762a1bSJed Brown       int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name));
4053c633725SBarry Smith       PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",initial[j].name);
4069566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",j,initial[j].name,initial[j].massfrac));
407c4762a1bSJed Brown       x[i][1+ispec] = initial[j].massfrac;
408c4762a1bSJed Brown     }
409c4762a1bSJed Brown   }
4109566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(dm,X,&x));
4119566063dSJacob Faibussowitsch   PetscCall(DMDARestoreCoordinateArray(dm,&xc));
412c4762a1bSJed Brown   PetscFunctionReturn(0);
413c4762a1bSJed Brown }
414c4762a1bSJed Brown 
415c4762a1bSJed Brown /*
416c4762a1bSJed Brown     Routines for displaying the solutions
417c4762a1bSJed Brown */
418c4762a1bSJed Brown typedef struct {
419c4762a1bSJed Brown   PetscInt cell;
420c4762a1bSJed Brown   User     user;
421c4762a1bSJed Brown } UserLGCtx;
422c4762a1bSJed Brown 
423c4762a1bSJed Brown static PetscErrorCode FormMoleFraction(UserLGCtx *ctx,Vec massf,Vec *molef)
424c4762a1bSJed Brown {
425c4762a1bSJed Brown   User              user = ctx->user;
426c4762a1bSJed Brown   PetscReal         *M,tM=0;
427c4762a1bSJed Brown   PetscInt          i,n = user->Nspec+1;
428c4762a1bSJed Brown   PetscScalar       *mof;
429c4762a1bSJed Brown   const PetscScalar **maf;
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   PetscFunctionBegin;
4329566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,molef));
4339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->Nspec,&M));
434c4762a1bSJed Brown   TC_getSmass(user->Nspec, M);
4359566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOFRead(user->dm,massf,&maf));
4369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*molef,&mof));
437c4762a1bSJed Brown   mof[0] = maf[ctx->cell][0]; /* copy over temperature */
438c4762a1bSJed Brown   for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1];
439c4762a1bSJed Brown   for (i=1; i<n; i++) {
440c4762a1bSJed Brown     mof[i] = maf[ctx->cell][i]/(M[i-1]*tM);
441c4762a1bSJed Brown   }
4429566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf));
4439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*molef,&mof));
4449566063dSJacob Faibussowitsch   PetscCall(PetscFree(M));
445c4762a1bSJed Brown   PetscFunctionReturn(0);
446c4762a1bSJed Brown }
447c4762a1bSJed Brown 
448c4762a1bSJed Brown static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx)
449c4762a1bSJed Brown {
450c4762a1bSJed Brown   PetscFunctionBegin;
4519566063dSJacob Faibussowitsch   PetscCall(PetscFree(uctx));
452c4762a1bSJed Brown   PetscFunctionReturn(0);
453c4762a1bSJed Brown }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown /*
456c4762a1bSJed Brown    Use TSMonitorLG to monitor the reactions in a particular cell
457c4762a1bSJed Brown */
458c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell)
459c4762a1bSJed Brown {
460c4762a1bSJed Brown   TSMonitorLGCtx ctx;
461c4762a1bSJed Brown   char           **snames;
462c4762a1bSJed Brown   UserLGCtx      *uctx;
463c4762a1bSJed Brown   char           label[128];
464c4762a1bSJed Brown   PetscReal      temp,*xc;
465c4762a1bSJed Brown   PetscMPIInt    rank;
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateArray(user->dm,&xc));
469c4762a1bSJed Brown   temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]);  /* Non-dimensionalized by user->Tini */
4709566063dSJacob Faibussowitsch   PetscCall(DMDARestoreCoordinateArray(user->dm,&xc));
4719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
4729566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank));
4739566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx));
4749566063dSJacob Faibussowitsch   PetscCall(DMDAGetFieldNames(user->dm,(const char * const **)&snames));
4759566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames));
4769566063dSJacob Faibussowitsch   PetscCall(PetscNew(&uctx));
477c4762a1bSJed Brown   uctx->cell = cell;
478c4762a1bSJed Brown   uctx->user = user;
4799566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx));
4809566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy));
481c4762a1bSJed Brown   PetscFunctionReturn(0);
482c4762a1bSJed Brown }
483