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 839371c9d4SSatish Balay #define PetscCallTC(ierr) \ 849371c9d4SSatish Balay do { PetscCheck(!ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in TChem library, return code %d", ierr); } while (0) 85c4762a1bSJed Brown 869371c9d4SSatish Balay int main(int argc, char **argv) { 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 PetscReal ftime, dt; 93c4762a1bSJed Brown char chemfile[PETSC_MAX_PATH_LEN] = "chem.inp", thermofile[PETSC_MAX_PATH_LEN] = "therm.dat"; 94c4762a1bSJed Brown struct _User user; 95c4762a1bSJed Brown TSConvergedReason reason; 96c4762a1bSJed Brown PetscBool showsolutions = PETSC_FALSE; 97c4762a1bSJed Brown char **snames, *names; 98c4762a1bSJed Brown Vec lambda; /* used with TSAdjoint for sensitivities */ 99c4762a1bSJed Brown 100327415f7SBarry Smith PetscFunctionBeginUser; 1019566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 102d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Chemistry solver options", ""); 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)); 116d0609cedSBarry Smith PetscOptionsEnd(); 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)); 176*48a46eb9SPierre Jolivet for (i = xs; i < xs + xm; i++) PetscCall(MonitorCell(ts, &user, i)); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180c4762a1bSJed Brown Set runtime options 181c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1829566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185c4762a1bSJed Brown Set final conditions for sensitivities 186c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1879566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &lambda)); 1889566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, &lambda, NULL)); 1899566063dSJacob Faibussowitsch PetscCall(VecSetValue(lambda, 0, 1.0, INSERT_VALUES)); 1909566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(lambda)); 1919566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(lambda)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 194c4762a1bSJed Brown Solve ODE 195c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1969566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 1979566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 1989566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 1999566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 20063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown { 203c4762a1bSJed Brown Vec max; 204c4762a1bSJed Brown const char *const *names; 205c4762a1bSJed Brown PetscInt i; 206c4762a1bSJed Brown const PetscReal *bmax; 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(TSMonitorEnvelopeGetBounds(ts, &max, NULL)); 209c4762a1bSJed Brown if (max) { 2109566063dSJacob Faibussowitsch PetscCall(TSMonitorLGGetVariableNames(ts, &names)); 211c4762a1bSJed Brown if (names) { 2129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(max, &bmax)); 2139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Species - maximum mass fraction\n")); 214c4762a1bSJed Brown for (i = 1; i < user.Nspec; i++) { 21563a3b9bcSJacob Faibussowitsch if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%s %g\n", names[i], (double)bmax[i])); 216c4762a1bSJed Brown } 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(max, &bmax)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 223c4762a1bSJed Brown Free work space. 224c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 225c4762a1bSJed Brown TC_reset(); 2269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 2299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda)); 2309566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2319566063dSJacob Faibussowitsch PetscCall(PetscFree3(user.tchemwork, user.Jdense, user.rows)); 2329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 233b122ec5aSJacob Faibussowitsch return 0; 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* 237c4762a1bSJed Brown Applies the second order centered difference diffusion operator on a one dimensional periodic domain 238c4762a1bSJed Brown */ 2399371c9d4SSatish Balay static PetscErrorCode FormDiffusionFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) { 240c4762a1bSJed Brown User user = (User)ptr; 241c4762a1bSJed Brown PetscScalar **f; 242c4762a1bSJed Brown const PetscScalar **x; 243c4762a1bSJed Brown DM dm; 244c4762a1bSJed Brown PetscInt i, xs, xm, j, dof; 245c4762a1bSJed Brown Vec Xlocal; 246c4762a1bSJed Brown PetscReal idx; 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionBeginUser; 2499566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 2509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 2519566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xlocal)); 2529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xlocal)); 2539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xlocal)); 2549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(dm, Xlocal, &x)); 2559566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(dm, F, &f)); 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown idx = 1.0 * user->diffus / user->dx; 259c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 2609371c9d4SSatish Balay for (j = 0; j < dof; j++) { f[i][j] += idx * (x[i + 1][j] - 2.0 * x[i][j] + x[i - 1][j]); } 261c4762a1bSJed Brown } 2629566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(dm, Xlocal, &x)); 2639566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(dm, F, &f)); 2649566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xlocal)); 265c4762a1bSJed Brown PetscFunctionReturn(0); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* 269c4762a1bSJed Brown Produces the second order centered difference diffusion operator on a one dimensional periodic domain 270c4762a1bSJed Brown */ 2719371c9d4SSatish Balay static PetscErrorCode FormDiffusionJacobian(TS ts, PetscReal t, Vec X, Mat Amat, Mat Pmat, void *ptr) { 272c4762a1bSJed Brown User user = (User)ptr; 273c4762a1bSJed Brown DM dm; 274c4762a1bSJed Brown PetscInt i, xs, xm, j, dof; 275c4762a1bSJed Brown PetscReal idx, values[3]; 276c4762a1bSJed Brown MatStencil row, col[3]; 277c4762a1bSJed Brown 278c4762a1bSJed Brown PetscFunctionBeginUser; 2799566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 2809566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 2819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL)); 282c4762a1bSJed Brown 283c4762a1bSJed Brown idx = 1.0 * user->diffus / user->dx; 284c4762a1bSJed Brown values[0] = idx; 285c4762a1bSJed Brown values[1] = -2.0 * idx; 286c4762a1bSJed Brown values[2] = idx; 287c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 288c4762a1bSJed Brown for (j = 0; j < dof; j++) { 2899371c9d4SSatish Balay row.i = i; 2909371c9d4SSatish Balay row.c = j; 2919371c9d4SSatish Balay col[0].i = i - 1; 2929371c9d4SSatish Balay col[0].c = j; 2939371c9d4SSatish Balay col[1].i = i; 2949371c9d4SSatish Balay col[1].c = j; 2959371c9d4SSatish Balay col[2].i = i + 1; 2969371c9d4SSatish Balay col[2].c = j; 2979566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Pmat, 1, &row, 3, col, values, ADD_VALUES)); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown } 3009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pmat, MAT_FINAL_ASSEMBLY)); 3019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pmat, MAT_FINAL_ASSEMBLY)); 302c4762a1bSJed Brown PetscFunctionReturn(0); 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 3059371c9d4SSatish Balay static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) { 306c4762a1bSJed Brown User user = (User)ptr; 307c4762a1bSJed Brown PetscScalar **f; 308c4762a1bSJed Brown const PetscScalar **x; 309c4762a1bSJed Brown DM dm; 310c4762a1bSJed Brown PetscInt i, xs, xm; 311c4762a1bSJed Brown 312c4762a1bSJed Brown PetscFunctionBeginUser; 313c4762a1bSJed Brown if (user->reactions) { 3149566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(dm, X, &x)); 3169566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(dm, F, &f)); 3179566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL)); 318c4762a1bSJed Brown 319c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 3209566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(user->tchemwork, x[i], user->Nspec + 1)); 321c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize */ 3229566063dSJacob Faibussowitsch PetscCallTC(TC_getSrc(user->tchemwork, user->Nspec + 1, f[i])); 323c4762a1bSJed Brown f[i][0] /= user->Tini; /* Non-dimensionalize */ 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 3269566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(dm, X, &x)); 3279566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(dm, F, &f)); 328c4762a1bSJed Brown } else { 3299566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 330c4762a1bSJed Brown } 3311baa6e33SBarry Smith if (user->diffusion) PetscCall(FormDiffusionFunction(ts, t, X, F, ptr)); 332c4762a1bSJed Brown PetscFunctionReturn(0); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 3359371c9d4SSatish Balay static PetscErrorCode FormRHSJacobian(TS ts, PetscReal t, Vec X, Mat Amat, Mat Pmat, void *ptr) { 336c4762a1bSJed Brown User user = (User)ptr; 337c4762a1bSJed Brown const PetscScalar **x; 338c4762a1bSJed Brown PetscInt M = user->Nspec + 1, i, j, xs, xm; 339c4762a1bSJed Brown DM dm; 340c4762a1bSJed Brown 341c4762a1bSJed Brown PetscFunctionBeginUser; 342c4762a1bSJed Brown if (user->reactions) { 3439566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3449566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Pmat)); 3459566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pmat, MAT_ROW_ORIENTED, PETSC_FALSE)); 3469566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pmat, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 3479566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(dm, X, &x)); 3489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL)); 349c4762a1bSJed Brown 350c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 3519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(user->tchemwork, x[i], user->Nspec + 1)); 352c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */ 3539566063dSJacob Faibussowitsch PetscCall(TC_getJacTYN(user->tchemwork, user->Nspec, user->Jdense, 1)); 354c4762a1bSJed Brown 355c4762a1bSJed Brown for (j = 0; j < M; j++) user->Jdense[j + 0 * M] /= user->Tini; /* Non-dimensionalize first column */ 356c4762a1bSJed Brown for (j = 0; j < M; j++) user->Jdense[0 + j * M] /= user->Tini; /* Non-dimensionalize first row */ 357c4762a1bSJed Brown for (j = 0; j < M; j++) user->rows[j] = i * M + j; 3589566063dSJacob Faibussowitsch PetscCall(MatSetValues(Pmat, M, user->rows, M, user->rows, user->Jdense, INSERT_VALUES)); 359c4762a1bSJed Brown } 3609566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(dm, X, &x)); 3619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pmat, MAT_FINAL_ASSEMBLY)); 3629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pmat, MAT_FINAL_ASSEMBLY)); 363c4762a1bSJed Brown } else { 3649566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Pmat)); 365c4762a1bSJed Brown } 3661baa6e33SBarry Smith if (user->diffusion) PetscCall(FormDiffusionJacobian(ts, t, X, Amat, Pmat, ptr)); 367c4762a1bSJed Brown if (Amat != Pmat) { 3689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 3699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown PetscFunctionReturn(0); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 3749371c9d4SSatish Balay PetscErrorCode FormInitialSolution(TS ts, Vec X, void *ctx) { 375c4762a1bSJed Brown PetscScalar **x, *xc; 3769371c9d4SSatish Balay struct { 3779371c9d4SSatish Balay const char *name; 3789371c9d4SSatish Balay PetscReal massfrac; 3799371c9d4SSatish Balay } initial[] = { 380c4762a1bSJed Brown {"CH4", 0.0948178320887 }, 381c4762a1bSJed Brown {"O2", 0.189635664177 }, 382c4762a1bSJed Brown {"N2", 0.706766236705 }, 383c4762a1bSJed Brown {"AR", 0.00878026702874} 384c4762a1bSJed Brown }; 385c4762a1bSJed Brown PetscInt i, j, xs, xm; 386c4762a1bSJed Brown DM dm; 387c4762a1bSJed Brown 388c4762a1bSJed Brown PetscFunctionBeginUser; 3899566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 3909566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 3919566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm, &xs, NULL, NULL, &xm, NULL, NULL)); 392c4762a1bSJed Brown 3939566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateArray(dm, &xc)); 3949566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(dm, X, &x)); 395c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 396c4762a1bSJed Brown x[i][0] = 1.0 + .05 * PetscSinScalar(2. * PETSC_PI * xc[i]); /* Non-dimensionalized by user->Tini */ 397dd39110bSPierre Jolivet for (j = 0; j < PETSC_STATIC_ARRAY_LENGTH(initial); j++) { 398c4762a1bSJed Brown int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name)); 3993c633725SBarry Smith PetscCheck(ispec >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Could not find species %s", initial[j].name); 40063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Species %d: %s %g\n", j, initial[j].name, (double)initial[j].massfrac)); 401c4762a1bSJed Brown x[i][1 + ispec] = initial[j].massfrac; 402c4762a1bSJed Brown } 403c4762a1bSJed Brown } 4049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(dm, X, &x)); 4059566063dSJacob Faibussowitsch PetscCall(DMDARestoreCoordinateArray(dm, &xc)); 406c4762a1bSJed Brown PetscFunctionReturn(0); 407c4762a1bSJed Brown } 408c4762a1bSJed Brown 409c4762a1bSJed Brown /* 410c4762a1bSJed Brown Routines for displaying the solutions 411c4762a1bSJed Brown */ 412c4762a1bSJed Brown typedef struct { 413c4762a1bSJed Brown PetscInt cell; 414c4762a1bSJed Brown User user; 415c4762a1bSJed Brown } UserLGCtx; 416c4762a1bSJed Brown 4179371c9d4SSatish Balay static PetscErrorCode FormMoleFraction(UserLGCtx *ctx, Vec massf, Vec *molef) { 418c4762a1bSJed Brown User user = ctx->user; 419c4762a1bSJed Brown PetscReal *M, tM = 0; 420c4762a1bSJed Brown PetscInt i, n = user->Nspec + 1; 421c4762a1bSJed Brown PetscScalar *mof; 422c4762a1bSJed Brown const PetscScalar **maf; 423c4762a1bSJed Brown 4247510d9b0SBarry Smith PetscFunctionBeginUser; 4259566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, molef)); 4269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->Nspec, &M)); 427c4762a1bSJed Brown TC_getSmass(user->Nspec, M); 4289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(user->dm, massf, &maf)); 4299566063dSJacob Faibussowitsch PetscCall(VecGetArray(*molef, &mof)); 430c4762a1bSJed Brown mof[0] = maf[ctx->cell][0]; /* copy over temperature */ 431c4762a1bSJed Brown for (i = 1; i < n; i++) tM += maf[ctx->cell][i] / M[i - 1]; 4329371c9d4SSatish Balay for (i = 1; i < n; i++) { mof[i] = maf[ctx->cell][i] / (M[i - 1] * tM); } 4339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(user->dm, massf, &maf)); 4349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*molef, &mof)); 4359566063dSJacob Faibussowitsch PetscCall(PetscFree(M)); 436c4762a1bSJed Brown PetscFunctionReturn(0); 437c4762a1bSJed Brown } 438c4762a1bSJed Brown 4399371c9d4SSatish Balay static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx) { 4407510d9b0SBarry Smith PetscFunctionBeginUser; 4419566063dSJacob Faibussowitsch PetscCall(PetscFree(uctx)); 442c4762a1bSJed Brown PetscFunctionReturn(0); 443c4762a1bSJed Brown } 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown Use TSMonitorLG to monitor the reactions in a particular cell 447c4762a1bSJed Brown */ 4489371c9d4SSatish Balay static PetscErrorCode MonitorCell(TS ts, User user, PetscInt cell) { 449c4762a1bSJed Brown TSMonitorLGCtx ctx; 450c4762a1bSJed Brown char **snames; 451c4762a1bSJed Brown UserLGCtx *uctx; 452c4762a1bSJed Brown char label[128]; 453c4762a1bSJed Brown PetscReal temp, *xc; 454c4762a1bSJed Brown PetscMPIInt rank; 455c4762a1bSJed Brown 4567510d9b0SBarry Smith PetscFunctionBeginUser; 4579566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateArray(user->dm, &xc)); 458c4762a1bSJed Brown temp = 1.0 + .05 * PetscSinScalar(2. * PETSC_PI * xc[cell]); /* Non-dimensionalized by user->Tini */ 4599566063dSJacob Faibussowitsch PetscCall(DMDARestoreCoordinateArray(user->dm, &xc)); 4609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 4619566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(label, sizeof(label), "Initial Temperature %g Cell %d Rank %d", (double)user->Tini * temp, (int)cell, rank)); 4629566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF, NULL, label, PETSC_DECIDE, PETSC_DECIDE, 600, 400, 1, &ctx)); 4639566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldNames(user->dm, (const char *const **)&snames)); 4649566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetVariableNames(ctx, (const char *const *)snames)); 4659566063dSJacob Faibussowitsch PetscCall(PetscNew(&uctx)); 466c4762a1bSJed Brown uctx->cell = cell; 467c4762a1bSJed Brown uctx->user = user; 4689566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetTransform(ctx, (PetscErrorCode(*)(void *, Vec, Vec *))FormMoleFraction, (PetscErrorCode(*)(void *))MonitorCellDestroy, uctx)); 4699566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, TSMonitorLGSolution, ctx, (PetscErrorCode(*)(void **))TSMonitorLGCtxDestroy)); 470c4762a1bSJed Brown PetscFunctionReturn(0); 471c4762a1bSJed Brown } 472