1c4762a1bSJed Brown static const char help[] = "Integrate chemistry using TChem.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscts.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown #if defined(PETSC_HAVE_TCHEM) 6c4762a1bSJed Brown #if defined(MAX) 7c4762a1bSJed Brown #undef MAX 8c4762a1bSJed Brown #endif 9c4762a1bSJed Brown #if defined(MIN) 10c4762a1bSJed Brown #undef MIN 11c4762a1bSJed Brown #endif 12c4762a1bSJed Brown #include <TC_params.h> 13c4762a1bSJed Brown #include <TC_interface.h> 14c4762a1bSJed Brown #else 15c4762a1bSJed Brown #error TChem is required for this example. Reconfigure PETSc using --download-tchem. 16c4762a1bSJed Brown #endif 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown See extchem.example.1 for how to run an example 19c4762a1bSJed Brown 20c4762a1bSJed Brown See also h2_10sp.inp for another example 21c4762a1bSJed Brown 22a5b23f4aSJose E. Roman Determine sensitivity of final temperature on each variables initial conditions 23c4762a1bSJed Brown -ts_dt 1.e-5 -ts_type cn -ts_adjoint_solve -ts_adjoint_view_solution draw 24c4762a1bSJed Brown 25c4762a1bSJed Brown The solution for component i = 0 is the temperature. 26c4762a1bSJed Brown 27c4762a1bSJed Brown The solution, i > 0, is the mass fraction, massf[i], of species i, i.e. mass of species i/ total mass of all species 28c4762a1bSJed Brown 29c4762a1bSJed Brown The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species 30c4762a1bSJed Brown Define M[i] = mass per mole of species i then 31c4762a1bSJed Brown molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j])) 32c4762a1bSJed Brown 33c4762a1bSJed Brown FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species. 34c4762a1bSJed Brown 35c4762a1bSJed Brown These are other data sets for other possible runs 36c4762a1bSJed Brown https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/n_heptane_v3.1_therm.dat 37c4762a1bSJed Brown https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/nc7_ver3.1_mech.txt 38c4762a1bSJed Brown 39c4762a1bSJed Brown */ 40c4762a1bSJed Brown typedef struct _User *User; 41c4762a1bSJed Brown struct _User { 42c4762a1bSJed Brown PetscReal pressure; 43c4762a1bSJed Brown int Nspec; 44c4762a1bSJed Brown int Nreac; 45c4762a1bSJed Brown PetscReal Tini; 46c4762a1bSJed Brown double *tchemwork; 47c4762a1bSJed Brown double *Jdense; /* Dense array workspace where Tchem computes the Jacobian */ 48c4762a1bSJed Brown PetscInt *rows; 49c4762a1bSJed Brown char **snames; 50c4762a1bSJed Brown }; 51c4762a1bSJed Brown 52c4762a1bSJed Brown static PetscErrorCode PrintSpecies(User, Vec); 53c4762a1bSJed Brown static PetscErrorCode MassFractionToMoleFraction(User, Vec, Vec *); 54c4762a1bSJed Brown static PetscErrorCode MoleFractionToMassFraction(User, Vec, Vec *); 55c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS, PetscReal, Vec, Vec, void *); 56c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 57c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS, Vec, void *); 58c4762a1bSJed Brown static PetscErrorCode ComputeMassConservation(Vec, PetscReal *, void *); 59c4762a1bSJed Brown static PetscErrorCode MonitorMassConservation(TS, PetscInt, PetscReal, Vec, void *); 60c4762a1bSJed Brown static PetscErrorCode MonitorTempature(TS, PetscInt, PetscReal, Vec, void *); 61c4762a1bSJed Brown 629371c9d4SSatish Balay #define PetscCallTC(ierr) \ 639371c9d4SSatish Balay do { PetscCheck(!ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in TChem library, return code %d", ierr); } while (0) 64c4762a1bSJed Brown 659371c9d4SSatish Balay int main(int argc, char **argv) { 66c4762a1bSJed Brown TS ts; /* time integrator */ 67c4762a1bSJed Brown TSAdapt adapt; 68c4762a1bSJed Brown Vec X, lambda; /* solution vector */ 69c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 70c4762a1bSJed Brown PetscInt steps; 71c4762a1bSJed Brown PetscReal ftime, dt; 72c4762a1bSJed Brown char chemfile[PETSC_MAX_PATH_LEN], thermofile[PETSC_MAX_PATH_LEN], lchemfile[PETSC_MAX_PATH_LEN], lthermofile[PETSC_MAX_PATH_LEN], lperiodic[PETSC_MAX_PATH_LEN]; 73c4762a1bSJed Brown const char *periodic = "file://${PETSC_DIR}/${PETSC_ARCH}/share/periodictable.dat"; 74c4762a1bSJed Brown struct _User user; /* user-defined work context */ 75c4762a1bSJed Brown TSConvergedReason reason; 76c4762a1bSJed Brown char **snames, *names; 77c4762a1bSJed Brown PetscInt i; 78c4762a1bSJed Brown TSTrajectory tj; 79c4762a1bSJed Brown PetscBool flg = PETSC_FALSE, tflg = PETSC_FALSE, found; 80c4762a1bSJed Brown 81327415f7SBarry Smith PetscFunctionBeginUser; 829566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 83d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Chemistry solver options", ""); 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-chem", "CHEMKIN input file", "", chemfile, chemfile, sizeof(chemfile), NULL)); 859566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD, chemfile, lchemfile, PETSC_MAX_PATH_LEN, &found)); 863c633725SBarry Smith PetscCheck(found, PETSC_COMM_WORLD, PETSC_ERR_FILE_OPEN, "Cannot download %s and no local version %s", chemfile, lchemfile); 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-thermo", "NASA thermo input file", "", thermofile, thermofile, sizeof(thermofile), NULL)); 889566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD, thermofile, lthermofile, PETSC_MAX_PATH_LEN, &found)); 893c633725SBarry Smith PetscCheck(found, PETSC_COMM_WORLD, PETSC_ERR_FILE_OPEN, "Cannot download %s and no local version %s", thermofile, lthermofile); 90c4762a1bSJed Brown user.pressure = 1.01325e5; /* Pascal */ 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pressure", "Pressure of reaction [Pa]", "", user.pressure, &user.pressure, NULL)); 92c4762a1bSJed Brown user.Tini = 1000; /* Kelvin */ 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Tini", "Initial temperature [K]", "", user.Tini, &user.Tini, NULL)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitor_mass", "Monitor the total mass at each timestep", "", flg, &flg, NULL)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-monitor_temp", "Monitor the temperature each timestep", "", tflg, &tflg, NULL)); 96d0609cedSBarry Smith PetscOptionsEnd(); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* tchem requires periodic table in current directory */ 999566063dSJacob Faibussowitsch PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD, periodic, lperiodic, PETSC_MAX_PATH_LEN, &found)); 1003c633725SBarry Smith PetscCheck(found, PETSC_COMM_WORLD, PETSC_ERR_FILE_OPEN, "Cannot located required periodic table %s or local version %s", periodic, lperiodic); 101c4762a1bSJed Brown 1029566063dSJacob Faibussowitsch PetscCallTC(TC_initChem(lchemfile, lthermofile, 0, 1.0)); 103c4762a1bSJed Brown TC_setThermoPres(user.pressure); 104c4762a1bSJed Brown user.Nspec = TC_getNspec(); 105c4762a1bSJed Brown user.Nreac = TC_getNreac(); 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown Get names of all species in easy to use array 108c4762a1bSJed Brown */ 1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((user.Nspec + 1) * LENGTHOFSPECNAME, &names)); 1109566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(names, "Temp")); 1112f613bf5SBarry Smith TC_getSnames(user.Nspec, names + LENGTHOFSPECNAME); 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((user.Nspec + 2), &snames)); 113c4762a1bSJed Brown for (i = 0; i < user.Nspec + 1; i++) snames[i] = names + i * LENGTHOFSPECNAME; 114c4762a1bSJed Brown snames[user.Nspec + 1] = NULL; 1159566063dSJacob Faibussowitsch PetscCall(PetscStrArrayallocpy((const char *const *)snames, &user.snames)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(snames)); 1179566063dSJacob Faibussowitsch PetscCall(PetscFree(names)); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(user.Nspec + 1, &user.tchemwork, PetscSqr(user.Nspec + 1), &user.Jdense, user.Nspec + 1, &user.rows)); 1209566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, user.Nspec + 1, &X)); 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch /* PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,PETSC_DECIDE,NULL,&J)); */ 1239566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, user.Nspec + 1, user.Nspec + 1, NULL, &J)); 1249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Create timestepping solver context 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1309566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX)); 1319566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE)); 1329566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts, TSARKIMEX4)); 1339566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user)); 1349566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, J, J, FormRHSJacobian, &user)); 135c4762a1bSJed Brown 1361baa6e33SBarry Smith if (flg) PetscCall(TSMonitorSet(ts, MonitorMassConservation, NULL, NULL)); 137*48a46eb9SPierre Jolivet if (tflg) PetscCall(TSMonitorSet(ts, MonitorTempature, &user, NULL)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown ftime = 1.0; 1409566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1419566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Set initial conditions 145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1469566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts, X, &user)); 1479566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X)); 148c4762a1bSJed Brown dt = 1e-10; /* Initial time step */ 1499566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 1509566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 1519566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt, 1e-12, 1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */ 1529566063dSJacob Faibussowitsch PetscCall(TSSetMaxSNESFailures(ts, -1)); /* Retry step an unlimited number of times */ 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155c4762a1bSJed Brown Set runtime options 156c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1579566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown Set final conditions for sensitivities 161c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &lambda)); 1639566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, &lambda, NULL)); 1649566063dSJacob Faibussowitsch PetscCall(VecSetValue(lambda, 0, 1.0, INSERT_VALUES)); 1659566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(lambda)); 1669566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(lambda)); 167c4762a1bSJed Brown 1689566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj)); 169c4762a1bSJed Brown if (tj) { 1709566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetVariableNames(tj, (const char *const *)user.snames)); 1719566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetTransform(tj, (PetscErrorCode(*)(void *, Vec, Vec *))MassFractionToMoleFraction, NULL, &user)); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175c4762a1bSJed Brown Pass information to graphical monitoring routine 176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1779566063dSJacob Faibussowitsch PetscCall(TSMonitorLGSetVariableNames(ts, (const char *const *)user.snames)); 1789566063dSJacob Faibussowitsch PetscCall(TSMonitorLGSetTransform(ts, (PetscErrorCode(*)(void *, Vec, Vec *))MassFractionToMoleFraction, NULL, &user)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Solve ODE 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1839566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 1849566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 1859566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 1869566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 18763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* { 190c4762a1bSJed Brown Vec max; 191c4762a1bSJed Brown PetscInt i; 192c4762a1bSJed Brown const PetscReal *bmax; 193c4762a1bSJed Brown 1949566063dSJacob Faibussowitsch PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL)); 195c4762a1bSJed Brown if (max) { 1969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(max,&bmax)); 1979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n")); 198c4762a1bSJed Brown for (i=1; i<user.Nspec; i++) { 1999566063dSJacob Faibussowitsch if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",user.snames[i],(double)bmax[i])); 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(max,&bmax)); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown Vec y; 2069566063dSJacob Faibussowitsch PetscCall(MassFractionToMoleFraction(&user,X,&y)); 2079566063dSJacob Faibussowitsch PetscCall(PrintSpecies(&user,y)); 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); */ 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211c4762a1bSJed Brown Free work space. 212c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 213c4762a1bSJed Brown TC_reset(); 2149566063dSJacob Faibussowitsch PetscCall(PetscStrArrayDestroy(&user.snames)); 2159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 2169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 2179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda)); 2189566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2199566063dSJacob Faibussowitsch PetscCall(PetscFree3(user.tchemwork, user.Jdense, user.rows)); 2209566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 221b122ec5aSJacob Faibussowitsch return 0; 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 2249371c9d4SSatish Balay static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) { 225c4762a1bSJed Brown User user = (User)ptr; 226c4762a1bSJed Brown PetscScalar *f; 227c4762a1bSJed Brown const PetscScalar *x; 228c4762a1bSJed Brown 229c4762a1bSJed Brown PetscFunctionBeginUser; 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2319566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 232c4762a1bSJed Brown 2339566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(user->tchemwork, x, user->Nspec + 1)); 234c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize */ 2359566063dSJacob Faibussowitsch PetscCallTC(TC_getSrc(user->tchemwork, user->Nspec + 1, f)); 236c4762a1bSJed Brown f[0] /= user->Tini; /* Non-dimensionalize */ 237c4762a1bSJed Brown 2389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 2399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 240c4762a1bSJed Brown PetscFunctionReturn(0); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 2439371c9d4SSatish Balay static PetscErrorCode FormRHSJacobian(TS ts, PetscReal t, Vec X, Mat Amat, Mat Pmat, void *ptr) { 244c4762a1bSJed Brown User user = (User)ptr; 245c4762a1bSJed Brown const PetscScalar *x; 246c4762a1bSJed Brown PetscInt M = user->Nspec + 1, i; 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionBeginUser; 2499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 2509566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(user->tchemwork, x, user->Nspec + 1)); 2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 252c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */ 2539566063dSJacob Faibussowitsch PetscCall(TC_getJacTYN(user->tchemwork, user->Nspec, user->Jdense, 1)); 254c4762a1bSJed Brown 255c4762a1bSJed Brown for (i = 0; i < M; i++) user->Jdense[i + 0 * M] /= user->Tini; /* Non-dimensionalize first column */ 256c4762a1bSJed Brown for (i = 0; i < M; i++) user->Jdense[0 + i * M] /= user->Tini; /* Non-dimensionalize first row */ 257c4762a1bSJed Brown for (i = 0; i < M; i++) user->rows[i] = i; 2589566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pmat, MAT_ROW_ORIENTED, PETSC_FALSE)); 2599566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pmat, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 2609566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Pmat)); 2619566063dSJacob Faibussowitsch PetscCall(MatSetValues(Pmat, M, user->rows, M, user->rows, user->Jdense, INSERT_VALUES)); 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pmat, MAT_FINAL_ASSEMBLY)); 2649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pmat, MAT_FINAL_ASSEMBLY)); 265c4762a1bSJed Brown if (Amat != Pmat) { 2669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 2679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 2729371c9d4SSatish Balay PetscErrorCode FormInitialSolution(TS ts, Vec X, void *ctx) { 273c4762a1bSJed Brown PetscScalar *x; 274c4762a1bSJed Brown PetscInt i; 275c4762a1bSJed Brown Vec y; 276c4762a1bSJed Brown const PetscInt maxspecies = 10; 277c4762a1bSJed Brown PetscInt smax = maxspecies, mmax = maxspecies; 278c4762a1bSJed Brown char *names[maxspecies]; 279c4762a1bSJed Brown PetscReal molefracs[maxspecies], sum; 280c4762a1bSJed Brown PetscBool flg; 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscFunctionBeginUser; 2839566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 2849566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 285c4762a1bSJed Brown x[0] = 1.0; /* Non-dimensionalized by user->Tini */ 286c4762a1bSJed Brown 2879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-initial_species", names, &smax, &flg)); 2883c633725SBarry Smith PetscCheck(smax >= 2, PETSC_COMM_SELF, PETSC_ERR_USER, "Must provide at least two initial species"); 2899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-initial_mole", molefracs, &mmax, &flg)); 29063a3b9bcSJacob Faibussowitsch PetscCheck(smax == mmax, PETSC_COMM_SELF, PETSC_ERR_USER, "Must provide same number of initial species %" PetscInt_FMT " as initial moles %" PetscInt_FMT, smax, mmax); 291c4762a1bSJed Brown sum = 0; 292c4762a1bSJed Brown for (i = 0; i < smax; i++) sum += molefracs[i]; 293c4762a1bSJed Brown for (i = 0; i < smax; i++) molefracs[i] = molefracs[i] / sum; 294c4762a1bSJed Brown for (i = 0; i < smax; i++) { 295c4762a1bSJed Brown int ispec = TC_getSpos(names[i], strlen(names[i])); 2963c633725SBarry Smith PetscCheck(ispec >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "Could not find species %s", names[i]); 29763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Species %" PetscInt_FMT ": %s %g\n", i, names[i], (double)molefracs[i])); 298c4762a1bSJed Brown x[1 + ispec] = molefracs[i]; 299c4762a1bSJed Brown } 300*48a46eb9SPierre Jolivet for (i = 0; i < smax; i++) PetscCall(PetscFree(names[i])); 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 3029566063dSJacob Faibussowitsch /* PetscCall(PrintSpecies((User)ctx,X)); */ 3039566063dSJacob Faibussowitsch PetscCall(MoleFractionToMassFraction((User)ctx, X, &y)); 3049566063dSJacob Faibussowitsch PetscCall(VecCopy(y, X)); 3059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 306c4762a1bSJed Brown PetscFunctionReturn(0); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* 310c4762a1bSJed Brown Converts the input vector which is in mass fractions (used by tchem) to mole fractions 311c4762a1bSJed Brown */ 3129371c9d4SSatish Balay PetscErrorCode MassFractionToMoleFraction(User user, Vec massf, Vec *molef) { 313c4762a1bSJed Brown PetscScalar *mof; 314c4762a1bSJed Brown const PetscScalar *maf; 315c4762a1bSJed Brown 3167510d9b0SBarry Smith PetscFunctionBeginUser; 3179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(massf, molef)); 3189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(massf, &maf)); 3199566063dSJacob Faibussowitsch PetscCall(VecGetArray(*molef, &mof)); 320c4762a1bSJed Brown mof[0] = maf[0]; /* copy over temperature */ 321c4762a1bSJed Brown TC_getMs2Ml((double *)maf + 1, user->Nspec, mof + 1); 3229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*molef, &mof)); 3239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(massf, &maf)); 324c4762a1bSJed Brown PetscFunctionReturn(0); 325c4762a1bSJed Brown } 326c4762a1bSJed Brown 327c4762a1bSJed Brown /* 328c4762a1bSJed Brown Converts the input vector which is in mole fractions to mass fractions (used by tchem) 329c4762a1bSJed Brown */ 3309371c9d4SSatish Balay PetscErrorCode MoleFractionToMassFraction(User user, Vec molef, Vec *massf) { 331c4762a1bSJed Brown const PetscScalar *mof; 332c4762a1bSJed Brown PetscScalar *maf; 333c4762a1bSJed Brown 3347510d9b0SBarry Smith PetscFunctionBeginUser; 3359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(molef, massf)); 3369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(molef, &mof)); 3379566063dSJacob Faibussowitsch PetscCall(VecGetArray(*massf, &maf)); 338c4762a1bSJed Brown maf[0] = mof[0]; /* copy over temperature */ 339c4762a1bSJed Brown TC_getMl2Ms((double *)mof + 1, user->Nspec, maf + 1); 3409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(molef, &mof)); 3419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*massf, &maf)); 342c4762a1bSJed Brown PetscFunctionReturn(0); 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 3459371c9d4SSatish Balay PetscErrorCode ComputeMassConservation(Vec x, PetscReal *mass, void *ctx) { 3467510d9b0SBarry Smith PetscFunctionBeginUser; 3479566063dSJacob Faibussowitsch PetscCall(VecSum(x, mass)); 348c4762a1bSJed Brown PetscFunctionReturn(0); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 3519371c9d4SSatish Balay PetscErrorCode MonitorMassConservation(TS ts, PetscInt step, PetscReal time, Vec x, void *ctx) { 352c4762a1bSJed Brown const PetscScalar *T; 353c4762a1bSJed Brown PetscReal mass; 354c4762a1bSJed Brown 3557510d9b0SBarry Smith PetscFunctionBeginUser; 3569566063dSJacob Faibussowitsch PetscCall(ComputeMassConservation(x, &mass, ctx)); 3579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &T)); 358c4762a1bSJed Brown mass -= PetscAbsScalar(T[0]); 3599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &T)); 36063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep %" PetscInt_FMT " time %g percent mass lost or gained %g\n", step, (double)time, (double)(100. * (1.0 - mass)))); 361c4762a1bSJed Brown PetscFunctionReturn(0); 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 3649371c9d4SSatish Balay PetscErrorCode MonitorTempature(TS ts, PetscInt step, PetscReal time, Vec x, void *ctx) { 365c4762a1bSJed Brown User user = (User)ctx; 366c4762a1bSJed Brown const PetscScalar *T; 367c4762a1bSJed Brown 3687510d9b0SBarry Smith PetscFunctionBeginUser; 3699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &T)); 37063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep %" PetscInt_FMT " time %g temperature %g\n", step, (double)time, (double)(T[0] * user->Tini))); 3719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &T)); 372c4762a1bSJed Brown PetscFunctionReturn(0); 373c4762a1bSJed Brown } 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* 376c4762a1bSJed Brown Prints out each species with its name 377c4762a1bSJed Brown */ 3789371c9d4SSatish Balay PETSC_UNUSED PetscErrorCode PrintSpecies(User user, Vec molef) { 379c4762a1bSJed Brown const PetscScalar *mof; 380c4762a1bSJed Brown PetscInt i, *idx, n = user->Nspec + 1; 381c4762a1bSJed Brown 3827510d9b0SBarry Smith PetscFunctionBeginUser; 3839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &idx)); 384c4762a1bSJed Brown for (i = 0; i < n; i++) idx[i] = i; 3859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(molef, &mof)); 3869566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithPermutation(n, mof, idx)); 387*48a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%6s %g\n", user->snames[idx[n - i - 1]], (double)PetscRealPart(mof[idx[n - i - 1]]))); 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 3899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(molef, &mof)); 390c4762a1bSJed Brown PetscFunctionReturn(0); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown 393c4762a1bSJed Brown /*TEST 394c4762a1bSJed Brown build: 395c4762a1bSJed Brown requires: tchem 396c4762a1bSJed Brown 397c4762a1bSJed Brown test: 398c4762a1bSJed Brown args: -chem http://combustion.berkeley.edu/gri_mech/version30/files30/grimech30.dat -thermo http://combustion.berkeley.edu/gri_mech/version30/files30/thermo30.dat -initial_species CH4,O2,N2,AR -initial_mole 0.0948178320887,0.189635664177,0.706766236705,0.00878026702874 -Tini 1500 -Tini 1500 -ts_arkimex_fully_implicit -ts_max_snes_failures -1 -ts_adapt_dt_max 1e-4 -ts_arkimex_type 4 -ts_max_time .005 399c4762a1bSJed Brown requires: !single 400c4762a1bSJed Brown filter: grep -v iterations 401c4762a1bSJed Brown 402c4762a1bSJed Brown TEST*/ 403