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