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 832e58f0efSBarry Smith #define CHKERRTC(ierr) do {if (ierr) SETERRQ1(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 PetscErrorCode ierr; 93c4762a1bSJed Brown PetscReal ftime,dt; 94c4762a1bSJed Brown char chemfile[PETSC_MAX_PATH_LEN] = "chem.inp",thermofile[PETSC_MAX_PATH_LEN] = "therm.dat"; 95c4762a1bSJed Brown struct _User user; 96c4762a1bSJed Brown TSConvergedReason reason; 97c4762a1bSJed Brown PetscBool showsolutions = PETSC_FALSE; 98c4762a1bSJed Brown char **snames,*names; 99c4762a1bSJed Brown Vec lambda; /* used with TSAdjoint for sensitivities */ 100c4762a1bSJed Brown 101c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 102c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL);CHKERRQ(ierr); 105c4762a1bSJed Brown user.pressure = 1.01325e5; /* Pascal */ 106c4762a1bSJed Brown ierr = PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL);CHKERRQ(ierr); 107c4762a1bSJed Brown user.Tini = 1550; 108c4762a1bSJed Brown ierr = PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL);CHKERRQ(ierr); 109c4762a1bSJed Brown user.diffus = 100; 110c4762a1bSJed Brown ierr = PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL);CHKERRQ(ierr); 112c4762a1bSJed Brown user.diffusion = PETSC_TRUE; 113c4762a1bSJed Brown ierr = PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL);CHKERRQ(ierr); 114c4762a1bSJed Brown user.reactions = PETSC_TRUE; 115c4762a1bSJed Brown ierr = PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 117c4762a1bSJed Brown 1182e58f0efSBarry Smith ierr = TC_initChem(chemfile, thermofile, 0, 1.0);CHKERRTC(ierr); 119c4762a1bSJed Brown user.Nspec = TC_getNspec(); 120c4762a1bSJed Brown user.Nreac = TC_getNreac(); 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = DMSetUp(user.dm);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = DMDAGetInfo(user.dm,NULL,&ncells,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 126c4762a1bSJed Brown user.dx = 1.0/ncells; /* Set the coordinates of the cell centers; note final ghost cell is at x coordinate 1.0 */ 127c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(user.dm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* set the names of each field in the DMDA based on the species name */ 130c4762a1bSJed Brown ierr = PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = PetscStrcpy(names,"Temp");CHKERRQ(ierr); 132*2f613bf5SBarry Smith TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME); 133c4762a1bSJed Brown ierr = PetscMalloc1((user.Nspec+2),&snames);CHKERRQ(ierr); 134c4762a1bSJed Brown for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME; 135c4762a1bSJed Brown snames[user.Nspec+1] = NULL; 136c4762a1bSJed Brown ierr = DMDASetFieldNames(user.dm,(const char * const *)snames);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = PetscFree(snames);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = PetscFree(names);CHKERRQ(ierr); 139c4762a1bSJed Brown 140c4762a1bSJed Brown ierr = DMCreateMatrix(user.dm,&J);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.dm,&X);CHKERRQ(ierr); 142c4762a1bSJed Brown 143c4762a1bSJed Brown ierr = PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows);CHKERRQ(ierr); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Create timestepping solver context 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = TSSetDM(ts,user.dm);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr); 152c4762a1bSJed Brown ierr = TSARKIMEXSetType(ts,TSARKIMEX4);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user);CHKERRQ(ierr); 155c4762a1bSJed Brown 156c4762a1bSJed Brown ftime = 1.0; 157c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161c4762a1bSJed Brown Set initial conditions 162c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 163c4762a1bSJed Brown ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 165c4762a1bSJed Brown dt = 1e-10; /* Initial time step */ 166c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,1e-12,1e-4);CHKERRQ(ierr); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */ 169c4762a1bSJed Brown ierr = TSSetMaxSNESFailures(ts,-1);CHKERRQ(ierr); /* Retry step an unlimited number of times */ 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172c4762a1bSJed Brown Pass information to graphical monitoring routine 173c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 174c4762a1bSJed Brown if (showsolutions) { 175c4762a1bSJed Brown ierr = DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 176c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 177c4762a1bSJed Brown ierr = MonitorCell(ts,&user,i);CHKERRQ(ierr); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182c4762a1bSJed Brown Set runtime options 183c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 184c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 187c4762a1bSJed Brown Set final conditions for sensitivities 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 189c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.dm,&lambda);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,&lambda,NULL);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = VecSetValue(lambda,0,1.0,INSERT_VALUES);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = VecAssemblyBegin(lambda);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = VecAssemblyEnd(lambda);CHKERRQ(ierr); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196c4762a1bSJed Brown Solve ODE 197c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 198c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr); 203c4762a1bSJed Brown 204c4762a1bSJed Brown { 205c4762a1bSJed Brown Vec max; 206c4762a1bSJed Brown const char * const *names; 207c4762a1bSJed Brown PetscInt i; 208c4762a1bSJed Brown const PetscReal *bmax; 209c4762a1bSJed Brown 210c4762a1bSJed Brown ierr = TSMonitorEnvelopeGetBounds(ts,&max,NULL);CHKERRQ(ierr); 211c4762a1bSJed Brown if (max) { 212c4762a1bSJed Brown ierr = TSMonitorLGGetVariableNames(ts,&names);CHKERRQ(ierr); 213c4762a1bSJed Brown if (names) { 214c4762a1bSJed Brown ierr = VecGetArrayRead(max,&bmax);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n");CHKERRQ(ierr); 216c4762a1bSJed Brown for (i=1; i<user.Nspec; i++) { 217c4762a1bSJed Brown if (bmax[i] > .01) {ierr = PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i]);CHKERRQ(ierr);} 218c4762a1bSJed Brown } 219c4762a1bSJed Brown ierr = VecRestoreArrayRead(max,&bmax);CHKERRQ(ierr); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown } 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Free work space. 226c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 227c4762a1bSJed Brown TC_reset(); 228c4762a1bSJed Brown ierr = DMDestroy(&user.dm);CHKERRQ(ierr); 229c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = VecDestroy(&lambda);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = PetscFree3(user.tchemwork,user.Jdense,user.rows);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = PetscFinalize(); 235c4762a1bSJed Brown return ierr; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Applies the second order centered difference diffusion operator on a one dimensional periodic domain 240c4762a1bSJed Brown */ 241c4762a1bSJed Brown static PetscErrorCode FormDiffusionFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 242c4762a1bSJed Brown { 243c4762a1bSJed Brown User user = (User)ptr; 244c4762a1bSJed Brown PetscErrorCode ierr; 245c4762a1bSJed Brown PetscScalar **f; 246c4762a1bSJed Brown const PetscScalar **x; 247c4762a1bSJed Brown DM dm; 248c4762a1bSJed Brown PetscInt i,xs,xm,j,dof; 249c4762a1bSJed Brown Vec Xlocal; 250c4762a1bSJed Brown PetscReal idx; 251c4762a1bSJed Brown 252c4762a1bSJed Brown PetscFunctionBeginUser; 253c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = DMGetLocalVector(dm,&Xlocal);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(dm,Xlocal,&x);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,F,&f);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 261c4762a1bSJed Brown 262c4762a1bSJed Brown idx = 1.0*user->diffus/user->dx; 263c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 264c4762a1bSJed Brown for (j=0; j<dof; j++) { 265c4762a1bSJed Brown f[i][j] += idx*(x[i+1][j] - 2.0*x[i][j] + x[i-1][j]); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown } 268c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x);CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,F,&f);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xlocal);CHKERRQ(ierr); 271c4762a1bSJed Brown PetscFunctionReturn(0); 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* 275c4762a1bSJed Brown Produces the second order centered difference diffusion operator on a one dimensional periodic domain 276c4762a1bSJed Brown */ 277c4762a1bSJed Brown static PetscErrorCode FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 278c4762a1bSJed Brown { 279c4762a1bSJed Brown User user = (User)ptr; 280c4762a1bSJed Brown PetscErrorCode ierr; 281c4762a1bSJed Brown DM dm; 282c4762a1bSJed Brown PetscInt i,xs,xm,j,dof; 283c4762a1bSJed Brown PetscReal idx,values[3]; 284c4762a1bSJed Brown MatStencil row,col[3]; 285c4762a1bSJed Brown 286c4762a1bSJed Brown PetscFunctionBeginUser; 287c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 290c4762a1bSJed Brown 291c4762a1bSJed Brown idx = 1.0*user->diffus/user->dx; 292c4762a1bSJed Brown values[0] = idx; 293c4762a1bSJed Brown values[1] = -2.0*idx; 294c4762a1bSJed Brown values[2] = idx; 295c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 296c4762a1bSJed Brown for (j=0; j<dof; j++) { 297c4762a1bSJed Brown row.i = i; row.c = j; 298c4762a1bSJed Brown col[0].i = i-1; col[0].c = j; 299c4762a1bSJed Brown col[1].i = i; col[1].c = j; 300c4762a1bSJed Brown col[2].i = i+1; col[2].c = j; 301c4762a1bSJed Brown ierr = MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES);CHKERRQ(ierr); 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown ierr = MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 305c4762a1bSJed Brown ierr = MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 306c4762a1bSJed Brown PetscFunctionReturn(0); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown 309c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 310c4762a1bSJed Brown { 311c4762a1bSJed Brown User user = (User)ptr; 312c4762a1bSJed Brown PetscErrorCode ierr; 313c4762a1bSJed Brown PetscScalar **f; 314c4762a1bSJed Brown const PetscScalar **x; 315c4762a1bSJed Brown DM dm; 316c4762a1bSJed Brown PetscInt i,xs,xm; 317c4762a1bSJed Brown 318c4762a1bSJed Brown PetscFunctionBeginUser; 319c4762a1bSJed Brown if (user->reactions) { 320c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 321c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,F,&f);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 324c4762a1bSJed Brown 325c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 326c4762a1bSJed Brown ierr = PetscArraycpy(user->tchemwork,x[i],user->Nspec+1);CHKERRQ(ierr); 327c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize */ 3282e58f0efSBarry Smith ierr = TC_getSrc(user->tchemwork,user->Nspec+1,f[i]);CHKERRTC(ierr); 329c4762a1bSJed Brown f[i][0] /= user->Tini; /* Non-dimensionalize */ 330c4762a1bSJed Brown } 331c4762a1bSJed Brown 332c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 333c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,F,&f);CHKERRQ(ierr); 334c4762a1bSJed Brown } else { 335c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 336c4762a1bSJed Brown } 337c4762a1bSJed Brown if (user->diffusion) { 338c4762a1bSJed Brown ierr = FormDiffusionFunction(ts,t,X,F,ptr);CHKERRQ(ierr); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown PetscFunctionReturn(0); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 344c4762a1bSJed Brown { 345c4762a1bSJed Brown User user = (User)ptr; 346c4762a1bSJed Brown PetscErrorCode ierr; 347c4762a1bSJed Brown const PetscScalar **x; 348c4762a1bSJed Brown PetscInt M = user->Nspec+1,i,j,xs,xm; 349c4762a1bSJed Brown DM dm; 350c4762a1bSJed Brown 351c4762a1bSJed Brown PetscFunctionBeginUser; 352c4762a1bSJed Brown if (user->reactions) { 353c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = MatZeroEntries(Pmat);CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 359c4762a1bSJed Brown 360c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 361c4762a1bSJed Brown ierr = PetscArraycpy(user->tchemwork,x[i],user->Nspec+1);CHKERRQ(ierr); 362c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */ 363c4762a1bSJed Brown ierr = TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1);CHKERRQ(ierr); 364c4762a1bSJed Brown 365c4762a1bSJed Brown for (j=0; j<M; j++) user->Jdense[j + 0*M] /= user->Tini; /* Non-dimensionalize first column */ 366c4762a1bSJed Brown for (j=0; j<M; j++) user->Jdense[0 + j*M] /= user->Tini; /* Non-dimensionalize first row */ 367c4762a1bSJed Brown for (j=0; j<M; j++) user->rows[j] = i*M+j; 368c4762a1bSJed Brown ierr = MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES);CHKERRQ(ierr); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 373c4762a1bSJed Brown } else { 374c4762a1bSJed Brown ierr = MatZeroEntries(Pmat);CHKERRQ(ierr); 375c4762a1bSJed Brown } 376c4762a1bSJed Brown if (user->diffusion) { 377c4762a1bSJed Brown ierr = FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr);CHKERRQ(ierr); 378c4762a1bSJed Brown } 379c4762a1bSJed Brown if (Amat != Pmat) { 380c4762a1bSJed Brown ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 382c4762a1bSJed Brown } 383c4762a1bSJed Brown PetscFunctionReturn(0); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 387c4762a1bSJed Brown { 388c4762a1bSJed Brown PetscScalar **x,*xc; 389c4762a1bSJed Brown PetscErrorCode ierr; 390c4762a1bSJed Brown struct {const char *name; PetscReal massfrac;} initial[] = { 391c4762a1bSJed Brown {"CH4", 0.0948178320887}, 392c4762a1bSJed Brown {"O2", 0.189635664177}, 393c4762a1bSJed Brown {"N2", 0.706766236705}, 394c4762a1bSJed Brown {"AR", 0.00878026702874} 395c4762a1bSJed Brown }; 396c4762a1bSJed Brown PetscInt i,j,xs,xm; 397c4762a1bSJed Brown DM dm; 398c4762a1bSJed Brown 399c4762a1bSJed Brown PetscFunctionBeginUser; 400c4762a1bSJed Brown ierr = VecZeroEntries(X);CHKERRQ(ierr); 401c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 402c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 403c4762a1bSJed Brown 404c4762a1bSJed Brown ierr = DMDAGetCoordinateArray(dm,&xc);CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,X,&x);CHKERRQ(ierr); 406c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 407c4762a1bSJed Brown x[i][0] = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[i]); /* Non-dimensionalized by user->Tini */ 408c4762a1bSJed Brown for (j=0; j<sizeof(initial)/sizeof(initial[0]); j++) { 409c4762a1bSJed Brown int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name)); 410c4762a1bSJed Brown if (ispec < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",initial[j].name); 411c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",j,initial[j].name,initial[j].massfrac);CHKERRQ(ierr); 412c4762a1bSJed Brown x[i][1+ispec] = initial[j].massfrac; 413c4762a1bSJed Brown } 414c4762a1bSJed Brown } 415c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,X,&x);CHKERRQ(ierr); 416c4762a1bSJed Brown ierr = DMDARestoreCoordinateArray(dm,&xc);CHKERRQ(ierr); 417c4762a1bSJed Brown PetscFunctionReturn(0); 418c4762a1bSJed Brown } 419c4762a1bSJed Brown 420c4762a1bSJed Brown /* 421c4762a1bSJed Brown Routines for displaying the solutions 422c4762a1bSJed Brown */ 423c4762a1bSJed Brown typedef struct { 424c4762a1bSJed Brown PetscInt cell; 425c4762a1bSJed Brown User user; 426c4762a1bSJed Brown } UserLGCtx; 427c4762a1bSJed Brown 428c4762a1bSJed Brown static PetscErrorCode FormMoleFraction(UserLGCtx *ctx,Vec massf,Vec *molef) 429c4762a1bSJed Brown { 430c4762a1bSJed Brown User user = ctx->user; 431c4762a1bSJed Brown PetscErrorCode ierr; 432c4762a1bSJed Brown PetscReal *M,tM=0; 433c4762a1bSJed Brown PetscInt i,n = user->Nspec+1; 434c4762a1bSJed Brown PetscScalar *mof; 435c4762a1bSJed Brown const PetscScalar **maf; 436c4762a1bSJed Brown 437c4762a1bSJed Brown PetscFunctionBegin; 438c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,molef);CHKERRQ(ierr); 439c4762a1bSJed Brown ierr = PetscMalloc1(user->Nspec,&M);CHKERRQ(ierr); 440c4762a1bSJed Brown TC_getSmass(user->Nspec, M); 441c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(user->dm,massf,&maf);CHKERRQ(ierr); 442c4762a1bSJed Brown ierr = VecGetArray(*molef,&mof);CHKERRQ(ierr); 443c4762a1bSJed Brown mof[0] = maf[ctx->cell][0]; /* copy over temperature */ 444c4762a1bSJed Brown for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1]; 445c4762a1bSJed Brown for (i=1; i<n; i++) { 446c4762a1bSJed Brown mof[i] = maf[ctx->cell][i]/(M[i-1]*tM); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf);CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = VecRestoreArray(*molef,&mof);CHKERRQ(ierr); 450c4762a1bSJed Brown ierr = PetscFree(M);CHKERRQ(ierr); 451c4762a1bSJed Brown PetscFunctionReturn(0); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454c4762a1bSJed Brown static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx) 455c4762a1bSJed Brown { 456c4762a1bSJed Brown PetscErrorCode ierr; 457c4762a1bSJed Brown 458c4762a1bSJed Brown PetscFunctionBegin; 459c4762a1bSJed Brown ierr = PetscFree(uctx);CHKERRQ(ierr); 460c4762a1bSJed Brown PetscFunctionReturn(0); 461c4762a1bSJed Brown } 462c4762a1bSJed Brown 463c4762a1bSJed Brown /* 464c4762a1bSJed Brown Use TSMonitorLG to monitor the reactions in a particular cell 465c4762a1bSJed Brown */ 466c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell) 467c4762a1bSJed Brown { 468c4762a1bSJed Brown PetscErrorCode ierr; 469c4762a1bSJed Brown TSMonitorLGCtx ctx; 470c4762a1bSJed Brown char **snames; 471c4762a1bSJed Brown UserLGCtx *uctx; 472c4762a1bSJed Brown char label[128]; 473c4762a1bSJed Brown PetscReal temp,*xc; 474c4762a1bSJed Brown PetscMPIInt rank; 475c4762a1bSJed Brown 476c4762a1bSJed Brown PetscFunctionBegin; 477c4762a1bSJed Brown ierr = DMDAGetCoordinateArray(user->dm,&xc);CHKERRQ(ierr); 478c4762a1bSJed Brown temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]); /* Non-dimensionalized by user->Tini */ 479c4762a1bSJed Brown ierr = DMDARestoreCoordinateArray(user->dm,&xc);CHKERRQ(ierr); 480ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 481c4762a1bSJed Brown ierr = PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank);CHKERRQ(ierr); 482c4762a1bSJed Brown ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx);CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = DMDAGetFieldNames(user->dm,(const char * const **)&snames);CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames);CHKERRQ(ierr); 485c4762a1bSJed Brown ierr = PetscNew(&uctx);CHKERRQ(ierr); 486c4762a1bSJed Brown uctx->cell = cell; 487c4762a1bSJed Brown uctx->user = user; 488c4762a1bSJed Brown ierr = TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx);CHKERRQ(ierr); 489c4762a1bSJed Brown ierr = TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr); 490c4762a1bSJed Brown PetscFunctionReturn(0); 491c4762a1bSJed Brown } 492