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