1*c4762a1bSJed Brown static const char help[] = "Integrate chemistry using TChem.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscts.h> 4*c4762a1bSJed Brown #include <petscdmda.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown #if defined(PETSC_HAVE_TCHEM) 7*c4762a1bSJed Brown #if defined(MAX) 8*c4762a1bSJed Brown #undef MAX 9*c4762a1bSJed Brown #endif 10*c4762a1bSJed Brown #if defined(MIN) 11*c4762a1bSJed Brown #undef MIN 12*c4762a1bSJed Brown #endif 13*c4762a1bSJed Brown # include <TC_params.h> 14*c4762a1bSJed Brown # include <TC_interface.h> 15*c4762a1bSJed Brown #else 16*c4762a1bSJed Brown # error TChem is required for this example. Reconfigure PETSc using --download-tchem. 17*c4762a1bSJed Brown #endif 18*c4762a1bSJed Brown /* 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown This is an extension of extchem.c to solve the reaction equations independently in each cell of a one dimensional field 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown Obtain the three files into this directory 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown curl http://combustion.berkeley.edu/gri_mech/version30/files30/grimech30.dat > chem.inp 25*c4762a1bSJed Brown curl http://combustion.berkeley.edu/gri_mech/version30/files30/thermo30.dat > therm.dat 26*c4762a1bSJed Brown cp $PETSC_DIR/$PETSC_ARCH/externalpackages/tchem/data/periodictable.dat . 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown Run with 29*c4762a1bSJed 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 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown Options for visualizing the solution: 32*c4762a1bSJed Brown Watch certain variables in each cell evolve with time 33*c4762a1bSJed 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 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown Watch certain variables in all cells evolve with time 36*c4762a1bSJed Brown -da_refine 4 -ts_monitor_draw_solution -draw_fields_by_name Temp,H2 -draw_vec_mark_points -draw_pause -2 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown Keep the initial temperature distribution as one monitors the current temperature distribution 39*c4762a1bSJed Brown -ts_monitor_draw_solution_initial -draw_bounds .9,1.7 -draw_fields_by_name Temp 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown Save the images in a .gif (movie) file 42*c4762a1bSJed Brown -draw_save -draw_save_single_file 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown Compute the sensitivies of the solution of the first tempature on the initial conditions 45*c4762a1bSJed Brown -ts_adjoint_solve -ts_dt 1.e-5 -ts_type cn -ts_adjoint_view_solution draw 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown Turn off diffusion 48*c4762a1bSJed Brown -diffusion no 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown Turn off reactions 51*c4762a1bSJed Brown -reactions no 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown The solution for component i = 0 is the temperature. 55*c4762a1bSJed Brown 56*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 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown The mole fraction molef[i], i > 0, is the number of moles of a species/ total number of moles of all species 59*c4762a1bSJed Brown Define M[i] = mass per mole of species i then 60*c4762a1bSJed Brown molef[i] = massf[i]/(M[i]*(sum_j massf[j]/M[j])) 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown FormMoleFraction(User,massf,molef) converts the mass fraction solution of each species to the mole fraction of each species. 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown */ 65*c4762a1bSJed Brown typedef struct _User *User; 66*c4762a1bSJed Brown struct _User { 67*c4762a1bSJed Brown PetscReal pressure; 68*c4762a1bSJed Brown int Nspec; 69*c4762a1bSJed Brown int Nreac; 70*c4762a1bSJed Brown PetscReal Tini,dx; 71*c4762a1bSJed Brown PetscReal diffus; 72*c4762a1bSJed Brown DM dm; 73*c4762a1bSJed Brown PetscBool diffusion,reactions; 74*c4762a1bSJed Brown double *tchemwork; 75*c4762a1bSJed Brown double *Jdense; /* Dense array workspace where Tchem computes the Jacobian */ 76*c4762a1bSJed Brown PetscInt *rows; 77*c4762a1bSJed Brown }; 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS,User,PetscInt); 80*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*); 81*c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 82*c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 83*c4762a1bSJed Brown 84*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) 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown int main(int argc,char **argv) 87*c4762a1bSJed Brown { 88*c4762a1bSJed Brown TS ts; /* time integrator */ 89*c4762a1bSJed Brown TSAdapt adapt; 90*c4762a1bSJed Brown Vec X; /* solution vector */ 91*c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 92*c4762a1bSJed Brown PetscInt steps,ncells,xs,xm,i; 93*c4762a1bSJed Brown PetscErrorCode ierr; 94*c4762a1bSJed Brown PetscReal ftime,dt; 95*c4762a1bSJed Brown char chemfile[PETSC_MAX_PATH_LEN] = "chem.inp",thermofile[PETSC_MAX_PATH_LEN] = "therm.dat"; 96*c4762a1bSJed Brown struct _User user; 97*c4762a1bSJed Brown TSConvergedReason reason; 98*c4762a1bSJed Brown PetscBool showsolutions = PETSC_FALSE; 99*c4762a1bSJed Brown char **snames,*names; 100*c4762a1bSJed Brown Vec lambda; /* used with TSAdjoint for sensitivities */ 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 103*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Chemistry solver options","");CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscOptionsString("-chem","CHEMKIN input file","",chemfile,chemfile,sizeof(chemfile),NULL);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = PetscOptionsString("-thermo","NASA thermo input file","",thermofile,thermofile,sizeof(thermofile),NULL);CHKERRQ(ierr); 106*c4762a1bSJed Brown user.pressure = 1.01325e5; /* Pascal */ 107*c4762a1bSJed Brown ierr = PetscOptionsReal("-pressure","Pressure of reaction [Pa]","",user.pressure,&user.pressure,NULL);CHKERRQ(ierr); 108*c4762a1bSJed Brown user.Tini = 1550; 109*c4762a1bSJed Brown ierr = PetscOptionsReal("-Tini","Initial temperature [K]","",user.Tini,&user.Tini,NULL);CHKERRQ(ierr); 110*c4762a1bSJed Brown user.diffus = 100; 111*c4762a1bSJed Brown ierr = PetscOptionsReal("-diffus","Diffusion constant","",user.diffus,&user.diffus,NULL);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = PetscOptionsBool("-draw_solution","Plot the solution for each cell","",showsolutions,&showsolutions,NULL);CHKERRQ(ierr); 113*c4762a1bSJed Brown user.diffusion = PETSC_TRUE; 114*c4762a1bSJed Brown ierr = PetscOptionsBool("-diffusion","Have diffusion","",user.diffusion,&user.diffusion,NULL);CHKERRQ(ierr); 115*c4762a1bSJed Brown user.reactions = PETSC_TRUE; 116*c4762a1bSJed Brown ierr = PetscOptionsBool("-reactions","Have reactions","",user.reactions,&user.reactions,NULL);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown ierr = TC_initChem(chemfile, thermofile, 0, 1.0);TCCHKERRQ(ierr); 120*c4762a1bSJed Brown user.Nspec = TC_getNspec(); 121*c4762a1bSJed Brown user.Nreac = TC_getNreac(); 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,10,user.Nspec+1,1,NULL,&user.dm);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = DMSetUp(user.dm);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = DMDAGetInfo(user.dm,NULL,&ncells,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 127*c4762a1bSJed Brown user.dx = 1.0/ncells; /* Set the coordinates of the cell centers; note final ghost cell is at x coordinate 1.0 */ 128*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(user.dm,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown /* set the names of each field in the DMDA based on the species name */ 131*c4762a1bSJed Brown ierr = PetscMalloc1((user.Nspec+1)*LENGTHOFSPECNAME,&names);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = PetscStrcpy(names,"Temp");CHKERRQ(ierr); 133*c4762a1bSJed Brown TC_getSnames(user.Nspec,names+LENGTHOFSPECNAME);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = PetscMalloc1((user.Nspec+2),&snames);CHKERRQ(ierr); 135*c4762a1bSJed Brown for (i=0; i<user.Nspec+1; i++) snames[i] = names+i*LENGTHOFSPECNAME; 136*c4762a1bSJed Brown snames[user.Nspec+1] = NULL; 137*c4762a1bSJed Brown ierr = DMDASetFieldNames(user.dm,(const char * const *)snames);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = PetscFree(snames);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = PetscFree(names);CHKERRQ(ierr); 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown ierr = DMCreateMatrix(user.dm,&J);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.dm,&X);CHKERRQ(ierr); 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown ierr = PetscMalloc3(user.Nspec+1,&user.tchemwork,PetscSqr(user.Nspec+1),&user.Jdense,user.Nspec+1,&user.rows);CHKERRQ(ierr); 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148*c4762a1bSJed Brown Create timestepping solver context 149*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 150*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = TSSetDM(ts,user.dm);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = TSARKIMEXSetType(ts,TSARKIMEX4);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,J,J,FormRHSJacobian,&user);CHKERRQ(ierr); 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown ftime = 1.0; 159*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 161*c4762a1bSJed Brown 162*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 163*c4762a1bSJed Brown Set initial conditions 164*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165*c4762a1bSJed Brown ierr = FormInitialSolution(ts,X,&user);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 167*c4762a1bSJed Brown dt = 1e-10; /* Initial time step */ 168*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,1e-12,1e-4);CHKERRQ(ierr); /* Also available with -ts_adapt_dt_min/-ts_adapt_dt_max */ 171*c4762a1bSJed Brown ierr = TSSetMaxSNESFailures(ts,-1);CHKERRQ(ierr); /* Retry step an unlimited number of times */ 172*c4762a1bSJed Brown 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175*c4762a1bSJed Brown Pass information to graphical monitoring routine 176*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 177*c4762a1bSJed Brown if (showsolutions) { 178*c4762a1bSJed Brown ierr = DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 179*c4762a1bSJed Brown for (i=xs;i<xs+xm;i++) { 180*c4762a1bSJed Brown ierr = MonitorCell(ts,&user,i);CHKERRQ(ierr); 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown } 183*c4762a1bSJed Brown 184*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185*c4762a1bSJed Brown Set runtime options 186*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 187*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190*c4762a1bSJed Brown Set final conditions for sensitivities 191*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 192*c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.dm,&lambda);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,&lambda,NULL);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = VecSetValue(lambda,0,1.0,INSERT_VALUES);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = VecAssemblyBegin(lambda);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = VecAssemblyEnd(lambda);CHKERRQ(ierr); 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199*c4762a1bSJed Brown Solve ODE 200*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 201*c4762a1bSJed Brown ierr = TSSolve(ts,X);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 204*c4762a1bSJed Brown ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 205*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr); 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown { 208*c4762a1bSJed Brown Vec max; 209*c4762a1bSJed Brown const char * const *names; 210*c4762a1bSJed Brown PetscInt i; 211*c4762a1bSJed Brown const PetscReal *bmax; 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown ierr = TSMonitorEnvelopeGetBounds(ts,&max,NULL);CHKERRQ(ierr); 214*c4762a1bSJed Brown if (max) { 215*c4762a1bSJed Brown ierr = TSMonitorLGGetVariableNames(ts,&names);CHKERRQ(ierr); 216*c4762a1bSJed Brown if (names) { 217*c4762a1bSJed Brown ierr = VecGetArrayRead(max,&bmax);CHKERRQ(ierr); 218*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Species - maximum mass fraction\n");CHKERRQ(ierr); 219*c4762a1bSJed Brown for (i=1; i<user.Nspec; i++) { 220*c4762a1bSJed Brown if (bmax[i] > .01) {ierr = PetscPrintf(PETSC_COMM_SELF,"%s %g\n",names[i],bmax[i]);CHKERRQ(ierr);} 221*c4762a1bSJed Brown } 222*c4762a1bSJed Brown ierr = VecRestoreArrayRead(max,&bmax);CHKERRQ(ierr); 223*c4762a1bSJed Brown } 224*c4762a1bSJed Brown } 225*c4762a1bSJed Brown } 226*c4762a1bSJed Brown 227*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 228*c4762a1bSJed Brown Free work space. 229*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 230*c4762a1bSJed Brown TC_reset(); 231*c4762a1bSJed Brown ierr = DMDestroy(&user.dm);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 234*c4762a1bSJed Brown ierr = VecDestroy(&lambda);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = PetscFree3(user.tchemwork,user.Jdense,user.rows);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = PetscFinalize(); 238*c4762a1bSJed Brown return ierr; 239*c4762a1bSJed Brown } 240*c4762a1bSJed Brown 241*c4762a1bSJed Brown /* 242*c4762a1bSJed Brown Applies the second order centered difference diffusion operator on a one dimensional periodic domain 243*c4762a1bSJed Brown */ 244*c4762a1bSJed Brown static PetscErrorCode FormDiffusionFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 245*c4762a1bSJed Brown { 246*c4762a1bSJed Brown User user = (User)ptr; 247*c4762a1bSJed Brown PetscErrorCode ierr; 248*c4762a1bSJed Brown PetscScalar **f; 249*c4762a1bSJed Brown const PetscScalar **x; 250*c4762a1bSJed Brown DM dm; 251*c4762a1bSJed Brown PetscInt i,xs,xm,j,dof; 252*c4762a1bSJed Brown Vec Xlocal; 253*c4762a1bSJed Brown PetscReal idx; 254*c4762a1bSJed Brown 255*c4762a1bSJed Brown PetscFunctionBeginUser; 256*c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 257*c4762a1bSJed Brown ierr = DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 258*c4762a1bSJed Brown ierr = DMGetLocalVector(dm,&Xlocal);CHKERRQ(ierr); 259*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xlocal);CHKERRQ(ierr); 260*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xlocal);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(dm,Xlocal,&x);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,F,&f);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 264*c4762a1bSJed Brown 265*c4762a1bSJed Brown idx = 1.0*user->diffus/user->dx; 266*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 267*c4762a1bSJed Brown for (j=0; j<dof; j++) { 268*c4762a1bSJed Brown f[i][j] += idx*(x[i+1][j] - 2.0*x[i][j] + x[i-1][j]); 269*c4762a1bSJed Brown } 270*c4762a1bSJed Brown } 271*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(dm,Xlocal,&x);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,F,&f);CHKERRQ(ierr); 273*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm,&Xlocal);CHKERRQ(ierr); 274*c4762a1bSJed Brown PetscFunctionReturn(0); 275*c4762a1bSJed Brown } 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown /* 278*c4762a1bSJed Brown Produces the second order centered difference diffusion operator on a one dimensional periodic domain 279*c4762a1bSJed Brown */ 280*c4762a1bSJed Brown static PetscErrorCode FormDiffusionJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 281*c4762a1bSJed Brown { 282*c4762a1bSJed Brown User user = (User)ptr; 283*c4762a1bSJed Brown PetscErrorCode ierr; 284*c4762a1bSJed Brown DM dm; 285*c4762a1bSJed Brown PetscInt i,xs,xm,j,dof; 286*c4762a1bSJed Brown PetscReal idx,values[3]; 287*c4762a1bSJed Brown MatStencil row,col[3]; 288*c4762a1bSJed Brown 289*c4762a1bSJed Brown PetscFunctionBeginUser; 290*c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 291*c4762a1bSJed Brown ierr = DMDAGetInfo(dm,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 292*c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 293*c4762a1bSJed Brown 294*c4762a1bSJed Brown idx = 1.0*user->diffus/user->dx; 295*c4762a1bSJed Brown values[0] = idx; 296*c4762a1bSJed Brown values[1] = -2.0*idx; 297*c4762a1bSJed Brown values[2] = idx; 298*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 299*c4762a1bSJed Brown for (j=0; j<dof; j++) { 300*c4762a1bSJed Brown row.i = i; row.c = j; 301*c4762a1bSJed Brown col[0].i = i-1; col[0].c = j; 302*c4762a1bSJed Brown col[1].i = i; col[1].c = j; 303*c4762a1bSJed Brown col[2].i = i+1; col[2].c = j; 304*c4762a1bSJed Brown ierr = MatSetValuesStencil(Pmat,1,&row,3,col,values,ADD_VALUES);CHKERRQ(ierr); 305*c4762a1bSJed Brown } 306*c4762a1bSJed Brown } 307*c4762a1bSJed Brown ierr = MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 309*c4762a1bSJed Brown PetscFunctionReturn(0); 310*c4762a1bSJed Brown } 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 313*c4762a1bSJed Brown { 314*c4762a1bSJed Brown User user = (User)ptr; 315*c4762a1bSJed Brown PetscErrorCode ierr; 316*c4762a1bSJed Brown PetscScalar **f; 317*c4762a1bSJed Brown const PetscScalar **x; 318*c4762a1bSJed Brown DM dm; 319*c4762a1bSJed Brown PetscInt i,xs,xm; 320*c4762a1bSJed Brown 321*c4762a1bSJed Brown PetscFunctionBeginUser; 322*c4762a1bSJed Brown if (user->reactions) { 323*c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 324*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 325*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,F,&f);CHKERRQ(ierr); 326*c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 327*c4762a1bSJed Brown 328*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 329*c4762a1bSJed Brown ierr = PetscArraycpy(user->tchemwork,x[i],user->Nspec+1);CHKERRQ(ierr); 330*c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize */ 331*c4762a1bSJed Brown ierr = TC_getSrc(user->tchemwork,user->Nspec+1,f[i]);TCCHKERRQ(ierr); 332*c4762a1bSJed Brown f[i][0] /= user->Tini; /* Non-dimensionalize */ 333*c4762a1bSJed Brown } 334*c4762a1bSJed Brown 335*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 336*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,F,&f);CHKERRQ(ierr); 337*c4762a1bSJed Brown } else { 338*c4762a1bSJed Brown ierr = VecZeroEntries(F);CHKERRQ(ierr); 339*c4762a1bSJed Brown } 340*c4762a1bSJed Brown if (user->diffusion) { 341*c4762a1bSJed Brown ierr = FormDiffusionFunction(ts,t,X,F,ptr);CHKERRQ(ierr); 342*c4762a1bSJed Brown } 343*c4762a1bSJed Brown PetscFunctionReturn(0); 344*c4762a1bSJed Brown } 345*c4762a1bSJed Brown 346*c4762a1bSJed Brown static PetscErrorCode FormRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat,Mat Pmat,void *ptr) 347*c4762a1bSJed Brown { 348*c4762a1bSJed Brown User user = (User)ptr; 349*c4762a1bSJed Brown PetscErrorCode ierr; 350*c4762a1bSJed Brown const PetscScalar **x; 351*c4762a1bSJed Brown PetscInt M = user->Nspec+1,i,j,xs,xm; 352*c4762a1bSJed Brown DM dm; 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown PetscFunctionBeginUser; 355*c4762a1bSJed Brown if (user->reactions) { 356*c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 357*c4762a1bSJed Brown ierr = MatZeroEntries(Pmat);CHKERRQ(ierr); 358*c4762a1bSJed Brown ierr = MatSetOption(Pmat,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 359*c4762a1bSJed Brown ierr = MatSetOption(Pmat,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 362*c4762a1bSJed Brown 363*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 364*c4762a1bSJed Brown ierr = PetscArraycpy(user->tchemwork,x[i],user->Nspec+1);CHKERRQ(ierr); 365*c4762a1bSJed Brown user->tchemwork[0] *= user->Tini; /* Dimensionalize temperature (first row) because that is what Tchem wants */ 366*c4762a1bSJed Brown ierr = TC_getJacTYN(user->tchemwork,user->Nspec,user->Jdense,1);CHKERRQ(ierr); 367*c4762a1bSJed Brown 368*c4762a1bSJed Brown for (j=0; j<M; j++) user->Jdense[j + 0*M] /= user->Tini; /* Non-dimensionalize first column */ 369*c4762a1bSJed Brown for (j=0; j<M; j++) user->Jdense[0 + j*M] /= user->Tini; /* Non-dimensionalize first row */ 370*c4762a1bSJed Brown for (j=0; j<M; j++) user->rows[j] = i*M+j; 371*c4762a1bSJed Brown ierr = MatSetValues(Pmat,M,user->rows,M,user->rows,user->Jdense,INSERT_VALUES);CHKERRQ(ierr); 372*c4762a1bSJed Brown } 373*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(dm,X,&x);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 376*c4762a1bSJed Brown } else { 377*c4762a1bSJed Brown ierr = MatZeroEntries(Pmat);CHKERRQ(ierr); 378*c4762a1bSJed Brown } 379*c4762a1bSJed Brown if (user->diffusion) { 380*c4762a1bSJed Brown ierr = FormDiffusionJacobian(ts,t,X,Amat,Pmat,ptr);CHKERRQ(ierr); 381*c4762a1bSJed Brown } 382*c4762a1bSJed Brown if (Amat != Pmat) { 383*c4762a1bSJed Brown ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 384*c4762a1bSJed Brown ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 385*c4762a1bSJed Brown } 386*c4762a1bSJed Brown PetscFunctionReturn(0); 387*c4762a1bSJed Brown } 388*c4762a1bSJed Brown 389*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx) 390*c4762a1bSJed Brown { 391*c4762a1bSJed Brown PetscScalar **x,*xc; 392*c4762a1bSJed Brown PetscErrorCode ierr; 393*c4762a1bSJed Brown struct {const char *name; PetscReal massfrac;} initial[] = { 394*c4762a1bSJed Brown {"CH4", 0.0948178320887}, 395*c4762a1bSJed Brown {"O2", 0.189635664177}, 396*c4762a1bSJed Brown {"N2", 0.706766236705}, 397*c4762a1bSJed Brown {"AR", 0.00878026702874} 398*c4762a1bSJed Brown }; 399*c4762a1bSJed Brown PetscInt i,j,xs,xm; 400*c4762a1bSJed Brown DM dm; 401*c4762a1bSJed Brown 402*c4762a1bSJed Brown PetscFunctionBeginUser; 403*c4762a1bSJed Brown ierr = VecZeroEntries(X);CHKERRQ(ierr); 404*c4762a1bSJed Brown ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 405*c4762a1bSJed Brown ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr); 406*c4762a1bSJed Brown 407*c4762a1bSJed Brown ierr = DMDAGetCoordinateArray(dm,&xc);CHKERRQ(ierr); 408*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOF(dm,X,&x);CHKERRQ(ierr); 409*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 410*c4762a1bSJed Brown x[i][0] = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[i]); /* Non-dimensionalized by user->Tini */ 411*c4762a1bSJed Brown for (j=0; j<sizeof(initial)/sizeof(initial[0]); j++) { 412*c4762a1bSJed Brown int ispec = TC_getSpos(initial[j].name, strlen(initial[j].name)); 413*c4762a1bSJed Brown if (ispec < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"Could not find species %s",initial[j].name); 414*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Species %d: %s %g\n",j,initial[j].name,initial[j].massfrac);CHKERRQ(ierr); 415*c4762a1bSJed Brown x[i][1+ispec] = initial[j].massfrac; 416*c4762a1bSJed Brown } 417*c4762a1bSJed Brown } 418*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOF(dm,X,&x);CHKERRQ(ierr); 419*c4762a1bSJed Brown ierr = DMDARestoreCoordinateArray(dm,&xc);CHKERRQ(ierr); 420*c4762a1bSJed Brown PetscFunctionReturn(0); 421*c4762a1bSJed Brown } 422*c4762a1bSJed Brown 423*c4762a1bSJed Brown /* 424*c4762a1bSJed Brown Routines for displaying the solutions 425*c4762a1bSJed Brown */ 426*c4762a1bSJed Brown typedef struct { 427*c4762a1bSJed Brown PetscInt cell; 428*c4762a1bSJed Brown User user; 429*c4762a1bSJed Brown } UserLGCtx; 430*c4762a1bSJed Brown 431*c4762a1bSJed Brown static PetscErrorCode FormMoleFraction(UserLGCtx *ctx,Vec massf,Vec *molef) 432*c4762a1bSJed Brown { 433*c4762a1bSJed Brown User user = ctx->user; 434*c4762a1bSJed Brown PetscErrorCode ierr; 435*c4762a1bSJed Brown PetscReal *M,tM=0; 436*c4762a1bSJed Brown PetscInt i,n = user->Nspec+1; 437*c4762a1bSJed Brown PetscScalar *mof; 438*c4762a1bSJed Brown const PetscScalar **maf; 439*c4762a1bSJed Brown 440*c4762a1bSJed Brown PetscFunctionBegin; 441*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,n,molef);CHKERRQ(ierr); 442*c4762a1bSJed Brown ierr = PetscMalloc1(user->Nspec,&M);CHKERRQ(ierr); 443*c4762a1bSJed Brown TC_getSmass(user->Nspec, M); 444*c4762a1bSJed Brown ierr = DMDAVecGetArrayDOFRead(user->dm,massf,&maf);CHKERRQ(ierr); 445*c4762a1bSJed Brown ierr = VecGetArray(*molef,&mof);CHKERRQ(ierr); 446*c4762a1bSJed Brown mof[0] = maf[ctx->cell][0]; /* copy over temperature */ 447*c4762a1bSJed Brown for (i=1; i<n; i++) tM += maf[ctx->cell][i]/M[i-1]; 448*c4762a1bSJed Brown for (i=1; i<n; i++) { 449*c4762a1bSJed Brown mof[i] = maf[ctx->cell][i]/(M[i-1]*tM); 450*c4762a1bSJed Brown } 451*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayDOFRead(user->dm,massf,&maf);CHKERRQ(ierr); 452*c4762a1bSJed Brown ierr = VecRestoreArray(*molef,&mof);CHKERRQ(ierr); 453*c4762a1bSJed Brown ierr = PetscFree(M);CHKERRQ(ierr); 454*c4762a1bSJed Brown PetscFunctionReturn(0); 455*c4762a1bSJed Brown } 456*c4762a1bSJed Brown 457*c4762a1bSJed Brown static PetscErrorCode MonitorCellDestroy(UserLGCtx *uctx) 458*c4762a1bSJed Brown { 459*c4762a1bSJed Brown PetscErrorCode ierr; 460*c4762a1bSJed Brown 461*c4762a1bSJed Brown PetscFunctionBegin; 462*c4762a1bSJed Brown ierr = PetscFree(uctx);CHKERRQ(ierr); 463*c4762a1bSJed Brown PetscFunctionReturn(0); 464*c4762a1bSJed Brown } 465*c4762a1bSJed Brown 466*c4762a1bSJed Brown /* 467*c4762a1bSJed Brown Use TSMonitorLG to monitor the reactions in a particular cell 468*c4762a1bSJed Brown */ 469*c4762a1bSJed Brown static PetscErrorCode MonitorCell(TS ts,User user,PetscInt cell) 470*c4762a1bSJed Brown { 471*c4762a1bSJed Brown PetscErrorCode ierr; 472*c4762a1bSJed Brown TSMonitorLGCtx ctx; 473*c4762a1bSJed Brown char **snames; 474*c4762a1bSJed Brown UserLGCtx *uctx; 475*c4762a1bSJed Brown char label[128]; 476*c4762a1bSJed Brown PetscReal temp,*xc; 477*c4762a1bSJed Brown PetscMPIInt rank; 478*c4762a1bSJed Brown 479*c4762a1bSJed Brown PetscFunctionBegin; 480*c4762a1bSJed Brown ierr = DMDAGetCoordinateArray(user->dm,&xc);CHKERRQ(ierr); 481*c4762a1bSJed Brown temp = 1.0 + .05*PetscSinScalar(2.*PETSC_PI*xc[cell]); /* Non-dimensionalized by user->Tini */ 482*c4762a1bSJed Brown ierr = DMDARestoreCoordinateArray(user->dm,&xc);CHKERRQ(ierr); 483*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 484*c4762a1bSJed Brown ierr = PetscSNPrintf(label,sizeof(label),"Initial Temperature %g Cell %d Rank %d",(double)user->Tini*temp,(int)cell,rank);CHKERRQ(ierr); 485*c4762a1bSJed Brown ierr = TSMonitorLGCtxCreate(PETSC_COMM_SELF,NULL,label,PETSC_DECIDE,PETSC_DECIDE,600,400,1,&ctx);CHKERRQ(ierr); 486*c4762a1bSJed Brown ierr = DMDAGetFieldNames(user->dm,(const char * const **)&snames);CHKERRQ(ierr); 487*c4762a1bSJed Brown ierr = TSMonitorLGCtxSetVariableNames(ctx,(const char * const *)snames);CHKERRQ(ierr); 488*c4762a1bSJed Brown ierr = PetscNew(&uctx);CHKERRQ(ierr); 489*c4762a1bSJed Brown uctx->cell = cell; 490*c4762a1bSJed Brown uctx->user = user; 491*c4762a1bSJed Brown ierr = TSMonitorLGCtxSetTransform(ctx,(PetscErrorCode (*)(void*,Vec,Vec*))FormMoleFraction,(PetscErrorCode (*)(void*))MonitorCellDestroy,uctx);CHKERRQ(ierr); 492*c4762a1bSJed Brown ierr = TSMonitorSet(ts,TSMonitorLGSolution,ctx,(PetscErrorCode (*)(void**))TSMonitorLGCtxDestroy);CHKERRQ(ierr); 493*c4762a1bSJed Brown PetscFunctionReturn(0); 494*c4762a1bSJed Brown } 495