xref: /petsc/src/ts/tutorials/extchemfield.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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
2909cb0f53SBarry Smith      ./extchemfield  -ts_arkimex_fully_implicit -ts_max_snes_failures unlimited -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 
44baca6076SPierre Jolivet         Compute the sensitivities of the solution of the first temperature on the initial conditions
45188af4bfSBarry Smith         -ts_adjoint_solve  -ts_time_step 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 
839371c9d4SSatish Balay #define PetscCallTC(ierr) \
84d71ae5a4SJacob Faibussowitsch   do { \
85d71ae5a4SJacob Faibussowitsch     PetscCheck(!ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in TChem library, return code %d", ierr); \
86d71ae5a4SJacob Faibussowitsch   } while (0)
87c4762a1bSJed Brown 
main(int argc,char ** argv)88d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
89d71ae5a4SJacob Faibussowitsch {
90c4762a1bSJed Brown   TS                ts; /* time integrator */
91c4762a1bSJed Brown   TSAdapt           adapt;
92c4762a1bSJed Brown   Vec               X; /* solution vector */
93c4762a1bSJed Brown   Mat               J; /* Jacobian matrix */
94c4762a1bSJed Brown   PetscInt          steps, ncells, xs, xm, i;
95c4762a1bSJed Brown   PetscReal         ftime, dt;
96c4762a1bSJed Brown   char              chemfile[PETSC_MAX_PATH_LEN] = "chem.inp", thermofile[PETSC_MAX_PATH_LEN] = "therm.dat";
97c4762a1bSJed Brown   struct _User      user;
98c4762a1bSJed Brown   TSConvergedReason reason;
99c4762a1bSJed Brown   PetscBool         showsolutions = PETSC_FALSE;
100c4762a1bSJed Brown   char            **snames, *names;
101c4762a1bSJed Brown   Vec               lambda; /* used with TSAdjoint for sensitivities */
102c4762a1bSJed Brown 
103327415f7SBarry Smith   PetscFunctionBeginUser;
104c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
105d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Chemistry solver options", "");
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-chem", "CHEMKIN input file", "", chemfile, chemfile, sizeof(chemfile), NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-thermo", "NASA thermo input file", "", thermofile, thermofile, sizeof(thermofile), NULL));
108c4762a1bSJed Brown   user.pressure = 1.01325e5; /* Pascal */
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pressure", "Pressure of reaction [Pa]", "", user.pressure, &user.pressure, NULL));
110c4762a1bSJed Brown   user.Tini = 1550;
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-Tini", "Initial temperature [K]", "", user.Tini, &user.Tini, NULL));
112c4762a1bSJed Brown   user.diffus = 100;
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-diffus", "Diffusion constant", "", user.diffus, &user.diffus, NULL));
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-draw_solution", "Plot the solution for each cell", "", showsolutions, &showsolutions, NULL));
115c4762a1bSJed Brown   user.diffusion = PETSC_TRUE;
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-diffusion", "Have diffusion", "", user.diffusion, &user.diffusion, NULL));
117c4762a1bSJed Brown   user.reactions = PETSC_TRUE;
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-reactions", "Have reactions", "", user.reactions, &user.reactions, NULL));
119d0609cedSBarry Smith   PetscOptionsEnd();
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCallTC(TC_initChem(chemfile, thermofile, 0, 1.0));
122c4762a1bSJed Brown   user.Nspec = TC_getNspec();
123c4762a1bSJed Brown   user.Nreac = TC_getNreac();
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, user.Nspec + 1, 1, NULL, &user.dm));
1269566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
1279566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
1289566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user.dm, NULL, &ncells, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
129c4762a1bSJed Brown   user.dx = 1.0 / ncells; /* Set the coordinates of the cell centers; note final ghost cell is at x coordinate 1.0 */
1309566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(user.dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* set the names of each field in the DMDA based on the species name */
1339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((user.Nspec + 1) * LENGTHOFSPECNAME, &names));
1346e21546bSBarry Smith   PetscCall(PetscStrncpy(names, "Temp", (user.Nspec + 1) * LENGTHOFSPECNAME));
1352f613bf5SBarry Smith   TC_getSnames(user.Nspec, names + LENGTHOFSPECNAME);
13657508eceSPierre Jolivet   PetscCall(PetscMalloc1(user.Nspec + 2, &snames));
137c4762a1bSJed Brown   for (i = 0; i < user.Nspec + 1; i++) snames[i] = names + i * LENGTHOFSPECNAME;
138c4762a1bSJed Brown   snames[user.Nspec + 1] = NULL;
1399566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldNames(user.dm, (const char *const *)snames));
1409566063dSJacob Faibussowitsch   PetscCall(PetscFree(snames));
1419566063dSJacob Faibussowitsch   PetscCall(PetscFree(names));
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.dm, &J));
1449566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm, &X));
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(user.Nspec + 1, &user.tchemwork, PetscSqr(user.Nspec + 1), &user.Jdense, user.Nspec + 1, &user.rows));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown      Create timestepping solver context
150c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1519566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1529566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, user.dm));
1539566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
1549566063dSJacob Faibussowitsch   PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE));
1559566063dSJacob Faibussowitsch   PetscCall(TSARKIMEXSetType(ts, TSARKIMEX4));
1569566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user));
1579566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, J, J, FormRHSJacobian, &user));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   ftime = 1.0;
1609566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1619566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown      Set initial conditions
165c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1669566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts, X, &user));
1679566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
168c4762a1bSJed Brown   dt = 1e-10; /* Initial time step */
1699566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
1709566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
1719566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt, 1e-12, 1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */
1729566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts, -1));             /* Retry step an unlimited number of times */
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Pass information to graphical monitoring routine
176c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177c4762a1bSJed Brown   if (showsolutions) {
1789566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(user.dm, &xs, NULL, NULL, &xm, NULL, NULL));
17948a46eb9SPierre Jolivet     for (i = xs; i < xs + xm; i++) PetscCall(MonitorCell(ts, &user, i));
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183c4762a1bSJed Brown      Set runtime options
184c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1859566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown      Set final conditions for sensitivities
189c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1909566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm, &lambda));
1919566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, &lambda, NULL));
1929566063dSJacob Faibussowitsch   PetscCall(VecSetValue(lambda, 0, 1.0, INSERT_VALUES));
1939566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(lambda));
1949566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(lambda));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197c4762a1bSJed Brown      Solve ODE
198c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1999566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
2009566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2019566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
2029566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
20363a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   {
206c4762a1bSJed Brown     Vec                max;
207c4762a1bSJed Brown     const char *const *names;
208c4762a1bSJed Brown     PetscInt           i;
209c4762a1bSJed Brown     const PetscReal   *bmax;
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch     PetscCall(TSMonitorEnvelopeGetBounds(ts, &max, NULL));
212c4762a1bSJed Brown     if (max) {
2139566063dSJacob Faibussowitsch       PetscCall(TSMonitorLGGetVariableNames(ts, &names));
214c4762a1bSJed Brown       if (names) {
2159566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(max, &bmax));
2169566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Species - maximum mass fraction\n"));
217c4762a1bSJed Brown         for (i = 1; i < user.Nspec; i++) {
21863a3b9bcSJacob Faibussowitsch           if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s %g\n", names[i], (double)bmax[i]));
219c4762a1bSJed Brown         }
2209566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(max, &bmax));
221c4762a1bSJed Brown       }
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226c4762a1bSJed Brown      Free work space.
227c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228c4762a1bSJed Brown   TC_reset();
2299566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
2309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
2319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
2329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda));
2339566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2349566063dSJacob Faibussowitsch   PetscCall(PetscFree3(user.tchemwork, user.Jdense, user.rows));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
236b122ec5aSJacob Faibussowitsch   return 0;
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown /*
240c4762a1bSJed Brown    Applies the second order centered difference diffusion operator on a one dimensional periodic domain
241c4762a1bSJed Brown */
FormDiffusionFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)242d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormDiffusionFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
243d71ae5a4SJacob Faibussowitsch {
244c4762a1bSJed Brown   User                user = (User)ptr;
245c4762a1bSJed Brown   PetscScalar       **f;
246c4762a1bSJed Brown   const PetscScalar **x;
247c4762a1bSJed Brown   DM                  dm;
248c4762a1bSJed Brown   PetscInt            i, xs, xm, j, dof;
249c4762a1bSJed Brown   Vec                 Xlocal;
250c4762a1bSJed Brown   PetscReal           idx;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBeginUser;
2539566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2549566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
2559566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &Xlocal));
2569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xlocal));
2579566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xlocal));
2589566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOFRead(dm, Xlocal, &x));
2599566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(dm, F, &f));
2609566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   idx = 1.0 * user->diffus / user->dx;
263c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
264ad540459SPierre Jolivet     for (j = 0; j < dof; j++) f[i][j] += idx * (x[i + 1][j] - 2.0 * x[i][j] + x[i - 1][j]);
265c4762a1bSJed Brown   }
2669566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOFRead(dm, Xlocal, &x));
2679566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(dm, F, &f));
2689566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &Xlocal));
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270c4762a1bSJed Brown }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown /*
273c4762a1bSJed Brown    Produces the second order centered difference diffusion operator on a one dimensional periodic domain
274c4762a1bSJed Brown */
FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void * ptr)275d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormDiffusionJacobian(TS ts, PetscReal t, Vec X, Mat Amat, Mat Pmat, void *ptr)
276d71ae5a4SJacob Faibussowitsch {
277c4762a1bSJed Brown   User       user = (User)ptr;
278c4762a1bSJed Brown   DM         dm;
279c4762a1bSJed Brown   PetscInt   i, xs, xm, j, dof;
280c4762a1bSJed Brown   PetscReal  idx, values[3];
281c4762a1bSJed Brown   MatStencil row, col[3];
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   PetscFunctionBeginUser;
2849566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2859566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
2869566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL));
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   idx       = 1.0 * user->diffus / user->dx;
289c4762a1bSJed Brown   values[0] = idx;
290c4762a1bSJed Brown   values[1] = -2.0 * idx;
291c4762a1bSJed Brown   values[2] = idx;
292c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
293c4762a1bSJed Brown     for (j = 0; j < dof; j++) {
2949371c9d4SSatish Balay       row.i    = i;
2959371c9d4SSatish Balay       row.c    = j;
2969371c9d4SSatish Balay       col[0].i = i - 1;
2979371c9d4SSatish Balay       col[0].c = j;
2989371c9d4SSatish Balay       col[1].i = i;
2999371c9d4SSatish Balay       col[1].c = j;
3009371c9d4SSatish Balay       col[2].i = i + 1;
3019371c9d4SSatish Balay       col[2].c = j;
3029566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Pmat, 1, &row, 3, col, values, ADD_VALUES));
303c4762a1bSJed Brown     }
304c4762a1bSJed Brown   }
3059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pmat, MAT_FINAL_ASSEMBLY));
3069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pmat, MAT_FINAL_ASSEMBLY));
3073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
308c4762a1bSJed Brown }
309c4762a1bSJed Brown 
FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)310d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
311d71ae5a4SJacob Faibussowitsch {
312c4762a1bSJed Brown   User                user = (User)ptr;
313c4762a1bSJed Brown   PetscScalar       **f;
314c4762a1bSJed Brown   const PetscScalar **x;
315c4762a1bSJed Brown   DM                  dm;
316c4762a1bSJed Brown   PetscInt            i, xs, xm;
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   PetscFunctionBeginUser;
319c4762a1bSJed Brown   if (user->reactions) {
3209566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
3219566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayDOFRead(dm, X, &x));
3229566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayDOF(dm, F, &f));
3239566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL));
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
3269566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(user->tchemwork, x[i], user->Nspec + 1));
327c4762a1bSJed Brown       user->tchemwork[0] *= user->Tini; /* Dimensionalize */
3289566063dSJacob Faibussowitsch       PetscCallTC(TC_getSrc(user->tchemwork, user->Nspec + 1, f[i]));
329c4762a1bSJed Brown       f[i][0] /= user->Tini; /* Non-dimensionalize */
330c4762a1bSJed Brown     }
331c4762a1bSJed Brown 
3329566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayDOFRead(dm, X, &x));
3339566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayDOF(dm, F, &f));
334c4762a1bSJed Brown   } else {
3359566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
336c4762a1bSJed Brown   }
3371baa6e33SBarry Smith   if (user->diffusion) PetscCall(FormDiffusionFunction(ts, t, X, F, ptr));
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339c4762a1bSJed Brown }
340c4762a1bSJed Brown 
FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void * ptr)341d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormRHSJacobian(TS ts, PetscReal t, Vec X, Mat Amat, Mat Pmat, void *ptr)
342d71ae5a4SJacob Faibussowitsch {
343c4762a1bSJed Brown   User                user = (User)ptr;
344c4762a1bSJed Brown   const PetscScalar **x;
345c4762a1bSJed Brown   PetscInt            M = user->Nspec + 1, i, j, xs, xm;
346c4762a1bSJed Brown   DM                  dm;
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   PetscFunctionBeginUser;
349c4762a1bSJed Brown   if (user->reactions) {
3509566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
3519566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(Pmat));
3529566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Pmat, MAT_ROW_ORIENTED, PETSC_FALSE));
3539566063dSJacob Faibussowitsch     PetscCall(MatSetOption(Pmat, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
3549566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayDOFRead(dm, X, &x));
3559566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
3589566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(user->tchemwork, x[i], user->Nspec + 1));
359c4762a1bSJed Brown       user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */
3609566063dSJacob Faibussowitsch       PetscCall(TC_getJacTYN(user->tchemwork, user->Nspec, user->Jdense, 1));
361c4762a1bSJed Brown 
362c4762a1bSJed Brown       for (j = 0; j < M; j++) user->Jdense[j + 0 * M] /= user->Tini; /* Non-dimensionalize first column */
363c4762a1bSJed Brown       for (j = 0; j < M; j++) user->Jdense[0 + j * M] /= user->Tini; /* Non-dimensionalize first row */
364c4762a1bSJed Brown       for (j = 0; j < M; j++) user->rows[j] = i * M + j;
3659566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Pmat, M, user->rows, M, user->rows, user->Jdense, INSERT_VALUES));
366c4762a1bSJed Brown     }
3679566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayDOFRead(dm, X, &x));
3689566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Pmat, MAT_FINAL_ASSEMBLY));
3699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Pmat, MAT_FINAL_ASSEMBLY));
370c4762a1bSJed Brown   } else {
3719566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(Pmat));
372c4762a1bSJed Brown   }
3731baa6e33SBarry Smith   if (user->diffusion) PetscCall(FormDiffusionJacobian(ts, t, X, Amat, Pmat, ptr));
374c4762a1bSJed Brown   if (Amat != Pmat) {
3759566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
3769566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
377c4762a1bSJed Brown   }
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
379c4762a1bSJed Brown }
380c4762a1bSJed Brown 
FormInitialSolution(TS ts,Vec X,PetscCtx ctx)381*2a8381b2SBarry Smith PetscErrorCode FormInitialSolution(TS ts, Vec X, PetscCtx ctx)
382d71ae5a4SJacob Faibussowitsch {
383c4762a1bSJed Brown   PetscScalar **x, *xc;
3849371c9d4SSatish Balay   struct {
3859371c9d4SSatish Balay     const char *name;
3869371c9d4SSatish Balay     PetscReal   massfrac;
3879371c9d4SSatish Balay   } initial[] = {
388c4762a1bSJed Brown     {"CH4", 0.0948178320887 },
389c4762a1bSJed Brown     {"O2",  0.189635664177  },
390c4762a1bSJed Brown     {"N2",  0.706766236705  },
391c4762a1bSJed Brown     {"AR",  0.00878026702874}
392c4762a1bSJed Brown   };
393c4762a1bSJed Brown   PetscInt i, j, xs, xm;
394c4762a1bSJed Brown   DM       dm;
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   PetscFunctionBeginUser;
3979566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
3989566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
3999566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL));
400c4762a1bSJed Brown 
4019566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateArray(dm, &xc));
4029566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(dm, X, &x));
403c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
404c4762a1bSJed Brown     x[i][0] = 1.0 + .05 * PetscSinScalar(2. * PETSC_PI * xc[i]); /* Non-dimensionalized by user->Tini */
405dd39110bSPierre Jolivet     for (j = 0; j < PETSC_STATIC_ARRAY_LENGTH(initial); j++) {
4066e21546bSBarry Smith       int ispec = TC_getSpos(initial[j].name, (int)strlen(initial[j].name));
4073c633725SBarry Smith       PetscCheck(ispec >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Could not find species %s", initial[j].name);
40863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Species %d: %s %g\n", j, initial[j].name, (double)initial[j].massfrac));
409c4762a1bSJed Brown       x[i][1 + ispec] = initial[j].massfrac;
410c4762a1bSJed Brown     }
411c4762a1bSJed Brown   }
4129566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(dm, X, &x));
4139566063dSJacob Faibussowitsch   PetscCall(DMDARestoreCoordinateArray(dm, &xc));
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415c4762a1bSJed Brown }
416c4762a1bSJed Brown 
417c4762a1bSJed Brown /*
418c4762a1bSJed Brown     Routines for displaying the solutions
419c4762a1bSJed Brown */
420c4762a1bSJed Brown typedef struct {
421c4762a1bSJed Brown   PetscInt cell;
422c4762a1bSJed Brown   User     user;
423c4762a1bSJed Brown } UserLGCtx;
424c4762a1bSJed Brown 
FormMoleFraction(UserLGCtx * ctx,Vec massf,Vec * molef)425d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormMoleFraction(UserLGCtx *ctx, Vec massf, Vec *molef)
426d71ae5a4SJacob Faibussowitsch {
427c4762a1bSJed Brown   User                user = ctx->user;
428c4762a1bSJed Brown   PetscReal          *M, tM = 0;
429c4762a1bSJed Brown   PetscInt            i, n  = user->Nspec + 1;
430c4762a1bSJed Brown   PetscScalar        *mof;
431c4762a1bSJed Brown   const PetscScalar **maf;
432c4762a1bSJed Brown 
4337510d9b0SBarry Smith   PetscFunctionBeginUser;
4349566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, molef));
4359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->Nspec, &M));
436c4762a1bSJed Brown   TC_getSmass(user->Nspec, M);
4379566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOFRead(user->dm, massf, &maf));
4389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(*molef, &mof));
439c4762a1bSJed Brown   mof[0] = maf[ctx->cell][0]; /* copy over temperature */
440c4762a1bSJed Brown   for (i = 1; i < n; i++) tM += maf[ctx->cell][i] / M[i - 1];
441ad540459SPierre Jolivet   for (i = 1; i < n; i++) mof[i] = maf[ctx->cell][i] / (M[i - 1] * tM);
4429566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOFRead(user->dm, massf, &maf));
4439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(*molef, &mof));
4449566063dSJacob Faibussowitsch   PetscCall(PetscFree(M));
4453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
446c4762a1bSJed Brown }
447c4762a1bSJed Brown 
MonitorCellDestroy(UserLGCtx ** uctx)4486e21546bSBarry Smith static PetscErrorCode MonitorCellDestroy(UserLGCtx **uctx)
449d71ae5a4SJacob Faibussowitsch {
4507510d9b0SBarry Smith   PetscFunctionBeginUser;
4516e21546bSBarry Smith   PetscCall(PetscFree(*uctx));
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
453c4762a1bSJed Brown }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown /*
456c4762a1bSJed Brown    Use TSMonitorLG to monitor the reactions in a particular cell
457c4762a1bSJed Brown */
MonitorCell(TS ts,User user,PetscInt cell)458d71ae5a4SJacob Faibussowitsch static PetscErrorCode MonitorCell(TS ts, User user, PetscInt cell)
459d71ae5a4SJacob Faibussowitsch {
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 
4677510d9b0SBarry Smith   PetscFunctionBeginUser;
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;
4796e21546bSBarry Smith   PetscCall(TSMonitorLGCtxSetTransform(ctx, (PetscErrorCode (*)(void *, Vec, Vec *))FormMoleFraction, (PetscCtxDestroyFn *)MonitorCellDestroy, uctx));
48049abdd8aSBarry Smith   PetscCall(TSMonitorSet(ts, TSMonitorLGSolution, ctx, (PetscCtxDestroyFn *)TSMonitorLGCtxDestroy));
4813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
482c4762a1bSJed Brown }
483