1c4762a1bSJed Brown static const char help[] = "Integrate chemistry using TChem.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscts.h> 4c4762a1bSJed Brown #include <petscdmda.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #if defined(PETSC_HAVE_TCHEM) 7c4762a1bSJed Brown #if defined(MAX) 8c4762a1bSJed Brown #undef MAX 9c4762a1bSJed Brown #endif 10c4762a1bSJed Brown #if defined(MIN) 11c4762a1bSJed Brown #undef MIN 12c4762a1bSJed Brown #endif 13c4762a1bSJed Brown # include <TC_params.h> 14c4762a1bSJed Brown # include <TC_interface.h> 15c4762a1bSJed Brown #else 16c4762a1bSJed Brown # error TChem is required for this example. Reconfigure PETSc using --download-tchem. 17c4762a1bSJed Brown #endif 18c4762a1bSJed Brown /* 19c4762a1bSJed Brown 20c4762a1bSJed Brown This is an extension of extchem.c to solve the reaction equations independently in each cell of a one dimensional field 21c4762a1bSJed Brown 22c4762a1bSJed Brown Obtain the three files into this directory 23c4762a1bSJed Brown 24c4762a1bSJed Brown curl http://combustion.berkeley.edu/gri_mech/version30/files30/grimech30.dat > chem.inp 25c4762a1bSJed Brown curl http://combustion.berkeley.edu/gri_mech/version30/files30/thermo30.dat > therm.dat 26c4762a1bSJed Brown cp $PETSC_DIR/$PETSC_ARCH/externalpackages/tchem/data/periodictable.dat . 27c4762a1bSJed Brown 28c4762a1bSJed Brown Run with 29c4762a1bSJed Brown ./extchemfield -ts_arkimex_fully_implicit -ts_max_snes_failures -1 -ts_adapt_monitor -ts_adapt_dt_max 1e-4 -ts_arkimex_type 4 -ts_max_time .005 30c4762a1bSJed Brown 31c4762a1bSJed Brown Options for visualizing the solution: 32c4762a1bSJed Brown Watch certain variables in each cell evolve with time 33c4762a1bSJed Brown -draw_solution 1 -ts_monitor_lg_solution_variables Temp,H2,O2,H2O,CH4,CO,CO2,C2H2,N2 -lg_use_markers false -draw_pause -2 34c4762a1bSJed Brown 35c4762a1bSJed Brown Watch certain variables in all cells evolve with time 36c4762a1bSJed Brown -da_refine 4 -ts_monitor_draw_solution -draw_fields_by_name Temp,H2 -draw_vec_mark_points -draw_pause -2 37c4762a1bSJed Brown 38c4762a1bSJed Brown Keep the initial temperature distribution as one monitors the current temperature distribution 39c4762a1bSJed Brown -ts_monitor_draw_solution_initial -draw_bounds .9,1.7 -draw_fields_by_name Temp 40c4762a1bSJed Brown 41c4762a1bSJed Brown Save the images in a .gif (movie) file 42c4762a1bSJed Brown -draw_save -draw_save_single_file 43c4762a1bSJed Brown 44a5b23f4aSJose E. Roman Compute the sensitivies of the solution of the first temperature on the initial conditions 45c4762a1bSJed Brown -ts_adjoint_solve -ts_dt 1.e-5 -ts_type cn -ts_adjoint_view_solution draw 46c4762a1bSJed Brown 47c4762a1bSJed Brown Turn off diffusion 48c4762a1bSJed Brown -diffusion no 49c4762a1bSJed Brown 50c4762a1bSJed Brown Turn off reactions 51c4762a1bSJed Brown -reactions no 52c4762a1bSJed Brown 53c4762a1bSJed Brown The solution for component i = 0 is the temperature. 54c4762a1bSJed Brown 55c4762a1bSJed Brown The solution, i > 0, is the mass fraction, massf[i], of species i, i.e. mass of species i/ total mass of all species 56c4762a1bSJed Brown 57c4762a1bSJed Brown The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species 58c4762a1bSJed Brown Define M[i] = mass per mole of species i then 59c4762a1bSJed Brown molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j])) 60c4762a1bSJed Brown 61c4762a1bSJed Brown FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species. 62c4762a1bSJed Brown 63c4762a1bSJed Brown */ 64c4762a1bSJed Brown typedef struct _User *User; 65c4762a1bSJed Brown struct _User { 66c4762a1bSJed Brown PetscReal pressure; 67c4762a1bSJed Brown int Nspec; 68c4762a1bSJed Brown int Nreac; 69c4762a1bSJed Brown PetscReal Tini,dx; 70c4762a1bSJed Brown PetscReal diffus; 71c4762a1bSJed Brown DM dm; 72c4762a1bSJed Brown PetscBool diffusion,reactions; 73c4762a1bSJed Brown double *tchemwork; 74c4762a1bSJed Brown double *Jdense; /* Dense array workspace where Tchem computes the Jacobian */ 75c4762a1bSJed Brown PetscInt *rows; 76c4762a1bSJed Brown }; 77c4762a1bSJed Brown 78c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS,User,PetscInt); 79c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 80c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 81c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch #define PetscCallTC(ierr) do {PetscCheck(!ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in TChem library, return code %d",ierr);} while (0) 84c4762a1bSJed Brown 85c4762a1bSJed Brown int main(int argc,char **argv) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown TS ts; /* time integrator */ 88c4762a1bSJed Brown TSAdapt adapt; 89c4762a1bSJed Brown Vec X; /* solution vector */ 90c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 91c4762a1bSJed Brown PetscInt steps,ncells,xs,xm,i; 92c4762a1bSJed Brown PetscReal ftime,dt; 93c4762a1bSJed Brown char chemfile[PETSC_MAX_PATH_LEN] = "chem.inp",thermofile[PETSC_MAX_PATH_LEN] = "therm.dat"; 94c4762a1bSJed Brown struct _User user; 95c4762a1bSJed Brown TSConvergedReason reason; 96c4762a1bSJed Brown PetscBool showsolutions = PETSC_FALSE; 97c4762a1bSJed Brown char **snames,*names; 98c4762a1bSJed Brown Vec lambda; /* used with TSAdjoint for sensitivities */ 99c4762a1bSJed Brown 1009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 101*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options",""); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL)); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL)); 104c4762a1bSJed Brown user.pressure = 1.01325e5; /* Pascal */ 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL)); 106c4762a1bSJed Brown user.Tini = 1550; 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL)); 108c4762a1bSJed Brown user.diffus = 100; 1099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL)); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL)); 111c4762a1bSJed Brown user.diffusion = PETSC_TRUE; 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL)); 113c4762a1bSJed Brown user.reactions = PETSC_TRUE; 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL)); 115*d0609cedSBarry Smith PetscOptionsEnd(); 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCallTC(TC_initChem(chemfile, thermofile, 0, 1.0)); 118c4762a1bSJed Brown user.Nspec = TC_getNspec(); 119c4762a1bSJed Brown user.Nreac = TC_getNreac(); 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm)); 1229566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 1239566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 1249566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user.dm,NULL,&ncells,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 125c4762a1bSJed Brown user.dx = 1.0/ncells; /* Set the coordinates of the cell centers; note final ghost cell is at x coordinate 1.0 */ 1269566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(user.dm,0.0,1.0,0.0,1.0,0.0,1.0)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* set the names of each field in the DMDA based on the species name */ 1299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names)); 1309566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(names,"Temp")); 1312f613bf5SBarry Smith TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME); 1329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((user.Nspec+2),&snames)); 133c4762a1bSJed Brown for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME; 134c4762a1bSJed Brown snames[user.Nspec+1] = NULL; 1359566063dSJacob Faibussowitsch PetscCall(DMDASetFieldNames(user.dm,(const char * const *)snames)); 1369566063dSJacob Faibussowitsch PetscCall(PetscFree(snames)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFree(names)); 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm,&J)); 1409566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm,&X)); 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145c4762a1bSJed Brown Create timestepping solver context 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1479566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1489566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,user.dm)); 1499566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSARKIMEX)); 1509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE)); 1519566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts,TSARKIMEX4)); 1529566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user)); 1539566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown ftime = 1.0; 1569566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ftime)); 1579566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown Set initial conditions 161c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1629566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts,X,&user)); 1639566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,X)); 164c4762a1bSJed Brown dt = 1e-10; /* Initial time step */ 1659566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 1669566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 1679566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt,1e-12,1e-4)); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */ 1689566063dSJacob Faibussowitsch PetscCall(TSSetMaxSNESFailures(ts,-1)); /* Retry step an unlimited number of times */ 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171c4762a1bSJed Brown Pass information to graphical monitoring routine 172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 173c4762a1bSJed Brown if (showsolutions) { 1749566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL)); 175c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 1769566063dSJacob Faibussowitsch PetscCall(MonitorCell(ts,&user,i)); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Set runtime options 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1839566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186c4762a1bSJed Brown Set final conditions for sensitivities 187c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1889566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm,&lambda)); 1899566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts,1,&lambda,NULL)); 1909566063dSJacob Faibussowitsch PetscCall(VecSetValue(lambda,0,1.0,INSERT_VALUES)); 1919566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(lambda)); 1929566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(lambda)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 195c4762a1bSJed Brown Solve ODE 196c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1979566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 1989566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 1999566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 2009566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts,&reason)); 2019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown { 204c4762a1bSJed Brown Vec max; 205c4762a1bSJed Brown const char * const *names; 206c4762a1bSJed Brown PetscInt i; 207c4762a1bSJed Brown const PetscReal *bmax; 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(TSMonitorEnvelopeGetBounds(ts,&max,NULL)); 210c4762a1bSJed Brown if (max) { 2119566063dSJacob Faibussowitsch PetscCall(TSMonitorLGGetVariableNames(ts,&names)); 212c4762a1bSJed Brown if (names) { 2139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(max,&bmax)); 2149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n")); 215c4762a1bSJed Brown for (i=1; i<user.Nspec; i++) { 2169566063dSJacob Faibussowitsch if (bmax[i] > .01) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i])); 217c4762a1bSJed Brown } 2189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(max,&bmax)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221c4762a1bSJed Brown } 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 224c4762a1bSJed Brown Free work space. 225c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 226c4762a1bSJed Brown TC_reset(); 2279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 2289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 2299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda)); 2319566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2329566063dSJacob Faibussowitsch PetscCall(PetscFree3(user.tchemwork,user.Jdense,user.rows)); 2339566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 234b122ec5aSJacob Faibussowitsch return 0; 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* 238c4762a1bSJed Brown Applies the second order centered difference diffusion operator on a one dimensional periodic domain 239c4762a1bSJed Brown */ 240c4762a1bSJed Brown static PetscErrorCode FormDiffusionFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 241c4762a1bSJed Brown { 242c4762a1bSJed Brown User user = (User)ptr; 243c4762a1bSJed Brown PetscScalar **f; 244c4762a1bSJed Brown const PetscScalar **x; 245c4762a1bSJed Brown DM dm; 246c4762a1bSJed Brown PetscInt i,xs,xm,j,dof; 247c4762a1bSJed Brown Vec Xlocal; 248c4762a1bSJed Brown PetscReal idx; 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscFunctionBeginUser; 2519566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 2529566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 2539566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xlocal)); 2549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal)); 2559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal)); 2569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(dm,Xlocal,&x)); 2579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(dm,F,&f)); 2589566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL)); 259c4762a1bSJed Brown 260c4762a1bSJed Brown idx = 1.0*user->diffus/user->dx; 261c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 262c4762a1bSJed Brown for (j=0; j<dof; j++) { 263c4762a1bSJed Brown f[i][j] += idx*(x[i+1][j] - 2.0*x[i][j] + x[i-1][j]); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown } 2669566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x)); 2679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(dm,F,&f)); 2689566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xlocal)); 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown Produces the second order centered difference diffusion operator on a one dimensional periodic domain 274c4762a1bSJed Brown */ 275c4762a1bSJed Brown static PetscErrorCode FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 276c4762a1bSJed Brown { 277c4762a1bSJed Brown User user = (User)ptr; 278c4762a1bSJed Brown DM dm; 279c4762a1bSJed Brown PetscInt i,xs,xm,j,dof; 280c4762a1bSJed Brown PetscReal idx,values[3]; 281c4762a1bSJed Brown MatStencil row,col[3]; 282c4762a1bSJed Brown 283c4762a1bSJed Brown PetscFunctionBeginUser; 2849566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 2859566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 2869566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL)); 287c4762a1bSJed Brown 288c4762a1bSJed Brown idx = 1.0*user->diffus/user->dx; 289c4762a1bSJed Brown values[0] = idx; 290c4762a1bSJed Brown values[1] = -2.0*idx; 291c4762a1bSJed Brown values[2] = idx; 292c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 293c4762a1bSJed Brown for (j=0; j<dof; j++) { 294c4762a1bSJed Brown row.i = i; row.c = j; 295c4762a1bSJed Brown col[0].i = i-1; col[0].c = j; 296c4762a1bSJed Brown col[1].i = i; col[1].c = j; 297c4762a1bSJed Brown col[2].i = i+1; col[2].c = j; 2989566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES)); 299c4762a1bSJed Brown } 300c4762a1bSJed Brown } 3019566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY)); 3029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY)); 303c4762a1bSJed Brown PetscFunctionReturn(0); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 307c4762a1bSJed Brown { 308c4762a1bSJed Brown User user = (User)ptr; 309c4762a1bSJed Brown PetscScalar **f; 310c4762a1bSJed Brown const PetscScalar **x; 311c4762a1bSJed Brown DM dm; 312c4762a1bSJed Brown PetscInt i,xs,xm; 313c4762a1bSJed Brown 314c4762a1bSJed Brown PetscFunctionBeginUser; 315c4762a1bSJed Brown if (user->reactions) { 3169566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 3179566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(dm,X,&x)); 3189566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(dm,F,&f)); 3199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL)); 320c4762a1bSJed Brown 321c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 3229566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1)); 323c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize */ 3249566063dSJacob Faibussowitsch PetscCallTC(TC_getSrc(user->tchemwork,user->Nspec+1,f[i])); 325c4762a1bSJed Brown f[i][0] /= user->Tini; /* Non-dimensionalize */ 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 3289566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(dm,X,&x)); 3299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(dm,F,&f)); 330c4762a1bSJed Brown } else { 3319566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown if (user->diffusion) { 3349566063dSJacob Faibussowitsch PetscCall(FormDiffusionFunction(ts,t,X,F,ptr)); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown PetscFunctionReturn(0); 337c4762a1bSJed Brown } 338c4762a1bSJed Brown 339c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 340c4762a1bSJed Brown { 341c4762a1bSJed Brown User user = (User)ptr; 342c4762a1bSJed Brown const PetscScalar **x; 343c4762a1bSJed Brown PetscInt M = user->Nspec+1,i,j,xs,xm; 344c4762a1bSJed Brown DM dm; 345c4762a1bSJed Brown 346c4762a1bSJed Brown PetscFunctionBeginUser; 347c4762a1bSJed Brown if (user->reactions) { 3489566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 3499566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Pmat)); 3509566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE)); 3519566063dSJacob Faibussowitsch PetscCall(MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 3529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(dm,X,&x)); 3539566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL)); 354c4762a1bSJed Brown 355c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 3569566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(user->tchemwork,x[i],user->Nspec+1)); 357c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */ 3589566063dSJacob Faibussowitsch PetscCall(TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1)); 359c4762a1bSJed Brown 360c4762a1bSJed Brown for (j=0; j<M; j++) user->Jdense[j + 0*M] /= user->Tini; /* Non-dimensionalize first column */ 361c4762a1bSJed Brown for (j=0; j<M; j++) user->Jdense[0 + j*M] /= user->Tini; /* Non-dimensionalize first row */ 362c4762a1bSJed Brown for (j=0; j<M; j++) user->rows[j] = i*M+j; 3639566063dSJacob Faibussowitsch PetscCall(MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES)); 364c4762a1bSJed Brown } 3659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(dm,X,&x)); 3669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY)); 3679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY)); 368c4762a1bSJed Brown } else { 3699566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Pmat)); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown if (user->diffusion) { 3729566063dSJacob Faibussowitsch PetscCall(FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr)); 373c4762a1bSJed Brown } 374c4762a1bSJed Brown if (Amat != Pmat) { 3759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY)); 3769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY)); 377c4762a1bSJed Brown } 378c4762a1bSJed Brown PetscFunctionReturn(0); 379c4762a1bSJed Brown } 380c4762a1bSJed Brown 381c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 382c4762a1bSJed Brown { 383c4762a1bSJed Brown PetscScalar **x,*xc; 384c4762a1bSJed Brown struct {const char *name; PetscReal massfrac;} initial[] = { 385c4762a1bSJed Brown {"CH4", 0.0948178320887}, 386c4762a1bSJed Brown {"O2", 0.189635664177}, 387c4762a1bSJed Brown {"N2", 0.706766236705}, 388c4762a1bSJed Brown {"AR", 0.00878026702874} 389c4762a1bSJed Brown }; 390c4762a1bSJed Brown PetscInt i,j,xs,xm; 391c4762a1bSJed Brown DM dm; 392c4762a1bSJed Brown 393c4762a1bSJed Brown PetscFunctionBeginUser; 3949566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 3959566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 3969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL)); 397c4762a1bSJed Brown 3989566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateArray(dm,&xc)); 3999566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(dm,X,&x)); 400c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 401c4762a1bSJed Brown x[i][0] = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[i]); /* Non-dimensionalized by user->Tini */ 402dd39110bSPierre Jolivet for (j=0; j<PETSC_STATIC_ARRAY_LENGTH(initial); j++) { 403c4762a1bSJed Brown int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name)); 4043c633725SBarry Smith PetscCheck(ispec >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",initial[j].name); 4059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",j,initial[j].name,initial[j].massfrac)); 406c4762a1bSJed Brown x[i][1+ispec] = initial[j].massfrac; 407c4762a1bSJed Brown } 408c4762a1bSJed Brown } 4099566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(dm,X,&x)); 4109566063dSJacob Faibussowitsch PetscCall(DMDARestoreCoordinateArray(dm,&xc)); 411c4762a1bSJed Brown PetscFunctionReturn(0); 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* 415c4762a1bSJed Brown Routines for displaying the solutions 416c4762a1bSJed Brown */ 417c4762a1bSJed Brown typedef struct { 418c4762a1bSJed Brown PetscInt cell; 419c4762a1bSJed Brown User user; 420c4762a1bSJed Brown } UserLGCtx; 421c4762a1bSJed Brown 422c4762a1bSJed Brown static PetscErrorCode FormMoleFraction(UserLGCtx *ctx,Vec massf,Vec *molef) 423c4762a1bSJed Brown { 424c4762a1bSJed Brown User user = ctx->user; 425c4762a1bSJed Brown PetscReal *M,tM=0; 426c4762a1bSJed Brown PetscInt i,n = user->Nspec+1; 427c4762a1bSJed Brown PetscScalar *mof; 428c4762a1bSJed Brown const PetscScalar **maf; 429c4762a1bSJed Brown 430c4762a1bSJed Brown PetscFunctionBegin; 4319566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,molef)); 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->Nspec,&M)); 433c4762a1bSJed Brown TC_getSmass(user->Nspec, M); 4349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOFRead(user->dm,massf,&maf)); 4359566063dSJacob Faibussowitsch PetscCall(VecGetArray(*molef,&mof)); 436c4762a1bSJed Brown mof[0] = maf[ctx->cell][0]; /* copy over temperature */ 437c4762a1bSJed Brown for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1]; 438c4762a1bSJed Brown for (i=1; i<n; i++) { 439c4762a1bSJed Brown mof[i] = maf[ctx->cell][i]/(M[i-1]*tM); 440c4762a1bSJed Brown } 4419566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf)); 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*molef,&mof)); 4439566063dSJacob Faibussowitsch PetscCall(PetscFree(M)); 444c4762a1bSJed Brown PetscFunctionReturn(0); 445c4762a1bSJed Brown } 446c4762a1bSJed Brown 447c4762a1bSJed Brown static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx) 448c4762a1bSJed Brown { 449c4762a1bSJed Brown PetscFunctionBegin; 4509566063dSJacob Faibussowitsch PetscCall(PetscFree(uctx)); 451c4762a1bSJed Brown PetscFunctionReturn(0); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454c4762a1bSJed Brown /* 455c4762a1bSJed Brown Use TSMonitorLG to monitor the reactions in a particular cell 456c4762a1bSJed Brown */ 457c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell) 458c4762a1bSJed Brown { 459c4762a1bSJed Brown TSMonitorLGCtx ctx; 460c4762a1bSJed Brown char **snames; 461c4762a1bSJed Brown UserLGCtx *uctx; 462c4762a1bSJed Brown char label[128]; 463c4762a1bSJed Brown PetscReal temp,*xc; 464c4762a1bSJed Brown PetscMPIInt rank; 465c4762a1bSJed Brown 466c4762a1bSJed Brown PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateArray(user->dm,&xc)); 468c4762a1bSJed Brown temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]); /* Non-dimensionalized by user->Tini */ 4699566063dSJacob Faibussowitsch PetscCall(DMDARestoreCoordinateArray(user->dm,&xc)); 4709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 4719566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank)); 4729566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx)); 4739566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldNames(user->dm,(const char * const **)&snames)); 4749566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames)); 4759566063dSJacob Faibussowitsch PetscCall(PetscNew(&uctx)); 476c4762a1bSJed Brown uctx->cell = cell; 477c4762a1bSJed Brown uctx->user = user; 4789566063dSJacob Faibussowitsch PetscCall(TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx)); 4799566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy)); 480c4762a1bSJed Brown PetscFunctionReturn(0); 481c4762a1bSJed Brown } 482