1 static const char help[] = "Integrate chemistry using TChem.\n"; 2 3 #include <petscts.h> 4 5 #if defined(PETSC_HAVE_TCHEM) 6 #if defined(MAX) 7 #undef MAX 8 #endif 9 #if defined(MIN) 10 #undef MIN 11 #endif 12 # include <TC_params.h> 13 # include <TC_interface.h> 14 #else 15 # error TChem is required for this example. Reconfigure PETSc using --download-tchem. 16 #endif 17 /* 18 See extchem.example.1 for how to run an example 19 20 See also h2_10sp.inp for another example 21 22 Determine sensitivity of final temperature on each variables initial conditions 23 -ts_dt 1.e-5 -ts_type cn -ts_adjoint_solve -ts_adjoint_view_solution draw 24 25 The solution for component i = 0 is the temperature. 26 27 The solution, i > 0, is the mass fraction, massf[i], of species i, i.e. mass of species i/ total mass of all species 28 29 The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species 30 Define M[i] = mass per mole of species i then 31 molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j])) 32 33 FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species. 34 35 These are other data sets for other possible runs 36 https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/n_heptane_v3.1_therm.dat 37 https://www-pls.llnl.gov/data/docs/science_and_technology/chemistry/combustion/nc7_ver3.1_mech.txt 38 39 */ 40 typedef struct _User *User; 41 struct _User { 42 PetscReal pressure; 43 int Nspec; 44 int Nreac; 45 PetscReal Tini; 46 double *tchemwork; 47 double *Jdense; /* Dense array workspace where Tchem computes the Jacobian */ 48 PetscInt *rows; 49 char **snames; 50 }; 51 52 static PetscErrorCode PrintSpecies(User,Vec); 53 static PetscErrorCode MassFractionToMoleFraction(User,Vec,Vec*); 54 static PetscErrorCode MoleFractionToMassFraction(User,Vec,Vec*); 55 static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 56 static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 57 static PetscErrorCode FormInitialSolution(TS,Vec,void*); 58 static PetscErrorCode ComputeMassConservation(Vec,PetscReal*,void*); 59 static PetscErrorCode MonitorMassConservation(TS,PetscInt,PetscReal,Vec,void*); 60 static PetscErrorCode MonitorTempature(TS,PetscInt,PetscReal,Vec,void*); 61 62 #define PetscCallTC(ierr) do {PetscCheck(!ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0) 63 64 int main(int argc,char **argv) 65 { 66 TS ts; /* time integrator */ 67 TSAdapt adapt; 68 Vec X,lambda; /* solution vector */ 69 Mat J; /* Jacobian matrix */ 70 PetscInt steps; 71 PetscReal ftime,dt; 72 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]; 73 const char *periodic = "file://${PETSC_DIR}/${PETSC_ARCH}/share/periodictable.dat"; 74 struct _User user; /* user-defined work context */ 75 TSConvergedReason reason; 76 char **snames,*names; 77 PetscInt i; 78 TSTrajectory tj; 79 PetscBool flg = PETSC_FALSE,tflg = PETSC_FALSE,found; 80 81 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 82 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options",""); 83 PetscCall(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL)); 84 PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,chemfile,lchemfile,PETSC_MAX_PATH_LEN,&found)); 85 PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",chemfile,lchemfile); 86 PetscCall(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL)); 87 PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,thermofile,lthermofile,PETSC_MAX_PATH_LEN,&found)); 88 PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot download %s and no local version %s",thermofile,lthermofile); 89 user.pressure = 1.01325e5; /* Pascal */ 90 PetscCall(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL)); 91 user.Tini = 1000; /* Kelvin */ 92 PetscCall(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL)); 93 PetscCall(PetscOptionsBool("-monitor_mass","Monitor the total mass at each timestep","",flg,&flg,NULL)); 94 PetscCall(PetscOptionsBool("-monitor_temp","Monitor the temperature each timestep","",tflg,&tflg,NULL)); 95 PetscOptionsEnd(); 96 97 /* tchem requires periodic table in current directory */ 98 PetscCall(PetscFileRetrieve(PETSC_COMM_WORLD,periodic,lperiodic,PETSC_MAX_PATH_LEN,&found)); 99 PetscCheck(found,PETSC_COMM_WORLD,PETSC_ERR_FILE_OPEN,"Cannot located required periodic table %s or local version %s",periodic,lperiodic); 100 101 PetscCallTC(TC_initChem(lchemfile, lthermofile, 0, 1.0)); 102 TC_setThermoPres(user.pressure); 103 user.Nspec = TC_getNspec(); 104 user.Nreac = TC_getNreac(); 105 /* 106 Get names of all species in easy to use array 107 */ 108 PetscCall(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names)); 109 PetscCall(PetscStrcpy(names,"Temp")); 110 TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME); 111 PetscCall(PetscMalloc1((user.Nspec+2),&snames)); 112 for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME; 113 snames[user.Nspec+1] = NULL; 114 PetscCall(PetscStrArrayallocpy((const char *const *)snames,&user.snames)); 115 PetscCall(PetscFree(snames)); 116 PetscCall(PetscFree(names)); 117 118 PetscCall(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows)); 119 PetscCall(VecCreateSeq(PETSC_COMM_SELF,user.Nspec+1,&X)); 120 121 /* PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,PETSC_DECIDE,NULL,&J)); */ 122 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,user.Nspec+1,user.Nspec+1,NULL,&J)); 123 PetscCall(MatSetFromOptions(J)); 124 125 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126 Create timestepping solver context 127 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 128 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 129 PetscCall(TSSetType(ts,TSARKIMEX)); 130 PetscCall(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE)); 131 PetscCall(TSARKIMEXSetType(ts,TSARKIMEX4)); 132 PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 133 PetscCall(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user)); 134 135 if (flg) PetscCall(TSMonitorSet(ts,MonitorMassConservation,NULL,NULL)); 136 if (tflg) { 137 PetscCall(TSMonitorSet(ts,MonitorTempature,&user,NULL)); 138 } 139 140 ftime = 1.0; 141 PetscCall(TSSetMaxTime(ts,ftime)); 142 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 143 144 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145 Set initial conditions 146 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147 PetscCall(FormInitialSolution(ts,X,&user)); 148 PetscCall(TSSetSolution(ts,X)); 149 dt = 1e-10; /* Initial time step */ 150 PetscCall(TSSetTimeStep(ts,dt)); 151 PetscCall(TSGetAdapt(ts,&adapt)); 152 PetscCall(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */ 153 PetscCall(TSSetMaxSNESFailures(ts,-1)); /* Retry step an unlimited number of times */ 154 155 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156 Set runtime options 157 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158 PetscCall(TSSetFromOptions(ts)); 159 160 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161 Set final conditions for sensitivities 162 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 163 PetscCall(VecDuplicate(X,&lambda)); 164 PetscCall(TSSetCostGradients(ts,1,&lambda,NULL)); 165 PetscCall(VecSetValue(lambda,0,1.0,INSERT_VALUES)); 166 PetscCall(VecAssemblyBegin(lambda)); 167 PetscCall(VecAssemblyEnd(lambda)); 168 169 PetscCall(TSGetTrajectory(ts,&tj)); 170 if (tj) { 171 PetscCall(TSTrajectorySetVariableNames(tj,(const char * const *)user.snames)); 172 PetscCall(TSTrajectorySetTransform(tj,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user)); 173 } 174 175 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176 Pass information to graphical monitoring routine 177 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 178 PetscCall(TSMonitorLGSetVariableNames(ts,(const char * const *)user.snames)); 179 PetscCall(TSMonitorLGSetTransform(ts,(PetscErrorCode (*)(void*,Vec,Vec*))MassFractionToMoleFraction,NULL,&user)); 180 181 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182 Solve ODE 183 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 184 PetscCall(TSSolve(ts,X)); 185 PetscCall(TSGetSolveTime(ts,&ftime)); 186 PetscCall(TSGetStepNumber(ts,&steps)); 187 PetscCall(TSGetConvergedReason(ts,&reason)); 188 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %" PetscInt_FMT " steps\n",TSConvergedReasons[reason],(double)ftime,steps)); 189 190 /* { 191 Vec max; 192 PetscInt i; 193 const PetscReal *bmax; 194 195 PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL)); 196 if (max) { 197 PetscCall(VecGetArrayRead(max,&bmax)); 198 PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n")); 199 for (i=1; i<user.Nspec; i++) { 200 if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",user.snames[i],(double)bmax[i])); 201 } 202 PetscCall(VecRestoreArrayRead(max,&bmax)); 203 } 204 } 205 206 Vec y; 207 PetscCall(MassFractionToMoleFraction(&user,X,&y)); 208 PetscCall(PrintSpecies(&user,y)); 209 PetscCall(VecDestroy(&y)); */ 210 211 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212 Free work space. 213 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 214 TC_reset(); 215 PetscCall(PetscStrArrayDestroy(&user.snames)); 216 PetscCall(MatDestroy(&J)); 217 PetscCall(VecDestroy(&X)); 218 PetscCall(VecDestroy(&lambda)); 219 PetscCall(TSDestroy(&ts)); 220 PetscCall(PetscFree3(user.tchemwork,user.Jdense,user.rows)); 221 PetscCall(PetscFinalize()); 222 return 0; 223 } 224 225 static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 226 { 227 User user = (User)ptr; 228 PetscScalar *f; 229 const PetscScalar *x; 230 231 PetscFunctionBeginUser; 232 PetscCall(VecGetArrayRead(X,&x)); 233 PetscCall(VecGetArray(F,&f)); 234 235 PetscCall(PetscArraycpy(user->tchemwork,x,user->Nspec+1)); 236 user->tchemwork[0] *= user->Tini; /* Dimensionalize */ 237 PetscCallTC(TC_getSrc(user->tchemwork,user->Nspec+1,f)); 238 f[0] /= user->Tini; /* Non-dimensionalize */ 239 240 PetscCall(VecRestoreArrayRead(X,&x)); 241 PetscCall(VecRestoreArray(F,&f)); 242 PetscFunctionReturn(0); 243 } 244 245 static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 246 { 247 User user = (User)ptr; 248 const PetscScalar *x; 249 PetscInt M = user->Nspec+1,i; 250 251 PetscFunctionBeginUser; 252 PetscCall(VecGetArrayRead(X,&x)); 253 PetscCall(PetscArraycpy(user->tchemwork,x,user->Nspec+1)); 254 PetscCall(VecRestoreArrayRead(X,&x)); 255 user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */ 256 PetscCall(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1)); 257 258 for (i=0; i<M; i++) user->Jdense[i + 0*M] /= user->Tini; /* Non-dimensionalize first column */ 259 for (i=0; i<M; i++) user->Jdense[0 + i*M] /= user->Tini; /* Non-dimensionalize first row */ 260 for (i=0; i<M; i++) user->rows[i] = i; 261 PetscCall(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE)); 262 PetscCall(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 263 PetscCall(MatZeroEntries(Pmat)); 264 PetscCall(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES)); 265 266 PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY)); 267 PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY)); 268 if (Amat != Pmat) { 269 PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY)); 270 PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY)); 271 } 272 PetscFunctionReturn(0); 273 } 274 275 PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 276 { 277 PetscScalar *x; 278 PetscInt i; 279 Vec y; 280 const PetscInt maxspecies = 10; 281 PetscInt smax = maxspecies,mmax = maxspecies; 282 char *names[maxspecies]; 283 PetscReal molefracs[maxspecies],sum; 284 PetscBool flg; 285 286 PetscFunctionBeginUser; 287 PetscCall(VecZeroEntries(X)); 288 PetscCall(VecGetArray(X,&x)); 289 x[0] = 1.0; /* Non-dimensionalized by user->Tini */ 290 291 PetscCall(PetscOptionsGetStringArray(NULL,NULL,"-initial_species",names,&smax,&flg)); 292 PetscCheck(smax >= 2,PETSC_COMM_SELF,PETSC_ERR_USER,"Must provide at least two initial species"); 293 PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-initial_mole",molefracs,&mmax,&flg)); 294 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); 295 sum = 0; 296 for (i=0; i<smax; i++) sum += molefracs[i]; 297 for (i=0; i<smax; i++) molefracs[i] = molefracs[i]/sum; 298 for (i=0; i<smax; i++) { 299 int ispec = TC_getSpos(names[i], strlen(names[i])); 300 PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",names[i]); 301 PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species %" PetscInt_FMT ": %s %g\n",i,names[i],(double)molefracs[i])); 302 x[1+ispec] = molefracs[i]; 303 } 304 for (i=0; i<smax; i++) { 305 PetscCall(PetscFree(names[i])); 306 } 307 PetscCall(VecRestoreArray(X,&x)); 308 /* PetscCall(PrintSpecies((User)ctx,X)); */ 309 PetscCall(MoleFractionToMassFraction((User)ctx,X,&y)); 310 PetscCall(VecCopy(y,X)); 311 PetscCall(VecDestroy(&y)); 312 PetscFunctionReturn(0); 313 } 314 315 /* 316 Converts the input vector which is in mass fractions (used by tchem) to mole fractions 317 */ 318 PetscErrorCode MassFractionToMoleFraction(User user,Vec massf,Vec *molef) 319 { 320 PetscScalar *mof; 321 const PetscScalar *maf; 322 323 PetscFunctionBegin; 324 PetscCall(VecDuplicate(massf,molef)); 325 PetscCall(VecGetArrayRead(massf,&maf)); 326 PetscCall(VecGetArray(*molef,&mof)); 327 mof[0] = maf[0]; /* copy over temperature */ 328 TC_getMs2Ml((double*)maf+1,user->Nspec,mof+1); 329 PetscCall(VecRestoreArray(*molef,&mof)); 330 PetscCall(VecRestoreArrayRead(massf,&maf)); 331 PetscFunctionReturn(0); 332 } 333 334 /* 335 Converts the input vector which is in mole fractions to mass fractions (used by tchem) 336 */ 337 PetscErrorCode MoleFractionToMassFraction(User user,Vec molef,Vec *massf) 338 { 339 const PetscScalar *mof; 340 PetscScalar *maf; 341 342 PetscFunctionBegin; 343 PetscCall(VecDuplicate(molef,massf)); 344 PetscCall(VecGetArrayRead(molef,&mof)); 345 PetscCall(VecGetArray(*massf,&maf)); 346 maf[0] = mof[0]; /* copy over temperature */ 347 TC_getMl2Ms((double*)mof+1,user->Nspec,maf+1); 348 PetscCall(VecRestoreArrayRead(molef,&mof)); 349 PetscCall(VecRestoreArray(*massf,&maf)); 350 PetscFunctionReturn(0); 351 } 352 353 PetscErrorCode ComputeMassConservation(Vec x,PetscReal *mass,void* ctx) 354 { 355 PetscFunctionBegin; 356 PetscCall(VecSum(x,mass)); 357 PetscFunctionReturn(0); 358 } 359 360 PetscErrorCode MonitorMassConservation(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx) 361 { 362 const PetscScalar *T; 363 PetscReal mass; 364 365 PetscFunctionBegin; 366 PetscCall(ComputeMassConservation(x,&mass,ctx)); 367 PetscCall(VecGetArrayRead(x,&T)); 368 mass -= PetscAbsScalar(T[0]); 369 PetscCall(VecRestoreArrayRead(x,&T)); 370 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)))); 371 PetscFunctionReturn(0); 372 } 373 374 PetscErrorCode MonitorTempature(TS ts,PetscInt step,PetscReal time,Vec x,void* ctx) 375 { 376 User user = (User) ctx; 377 const PetscScalar *T; 378 379 PetscFunctionBegin; 380 PetscCall(VecGetArrayRead(x,&T)); 381 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Timestep %" PetscInt_FMT " time %g temperature %g\n",step,(double)time,(double)(T[0]*user->Tini))); 382 PetscCall(VecRestoreArrayRead(x,&T)); 383 PetscFunctionReturn(0); 384 } 385 386 /* 387 Prints out each species with its name 388 */ 389 PETSC_UNUSED PetscErrorCode PrintSpecies(User user,Vec molef) 390 { 391 const PetscScalar *mof; 392 PetscInt i,*idx,n = user->Nspec+1; 393 394 PetscFunctionBegin; 395 PetscCall(PetscMalloc1(n,&idx)); 396 for (i=0; i<n;i++) idx[i] = i; 397 PetscCall(VecGetArrayRead(molef,&mof)); 398 PetscCall(PetscSortRealWithPermutation(n,mof,idx)); 399 for (i=0; i<n; i++) { 400 PetscCall(PetscPrintf(PETSC_COMM_SELF,"%6s %g\n",user->snames[idx[n-i-1]],(double)PetscRealPart(mof[idx[n-i-1]]))); 401 } 402 PetscCall(PetscFree(idx)); 403 PetscCall(VecRestoreArrayRead(molef,&mof)); 404 PetscFunctionReturn(0); 405 } 406 407 /*TEST 408 build: 409 requires: tchem 410 411 test: 412 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 413 requires: !single 414 filter: grep -v iterations 415 416 TEST*/ 417