1*c4762a1bSJed Brown static char help[] = "Time-Dependent Reactive Flow example in 2D with Darcy Flow"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /*T 4*c4762a1bSJed Brown Concepts: Solving a multicomponent time-dependent reactive flow system 5*c4762a1bSJed Brown Concepts: DMDA with timestepping 6*c4762a1bSJed Brown Processors: n 7*c4762a1bSJed Brown T*/ 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown /* 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown This example solves the elementary chemical reaction: 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown SP_A + 2SP_B = SP_C 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown Subject to predetermined flow modeled as though it were in a porous media. 16*c4762a1bSJed Brown This flow is modeled as advection and diffusion of the three species as 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown v = porosity*saturation*grad(q) 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown and the time-dependent equation solved for each species as 21*c4762a1bSJed Brown advection + diffusion + reaction: 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown v dot grad SP_i + dSP_i / dt + dispersivity*div*grad(SP_i) + R(SP_i) = 0 24*c4762a1bSJed Brown 25*c4762a1bSJed Brown The following options are available: 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown -length_x - The length of the domain in the direction of the flow 28*c4762a1bSJed Brown -length_y - The length of the domain in the direction orthogonal to the flow 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown -gradq_inflow - The inflow speed as if the saturation and porosity were 1. 31*c4762a1bSJed Brown -saturation - The saturation of the porous medium. 32*c4762a1bSJed Brown -porosity - The porosity of the medium. 33*c4762a1bSJed Brown -dispersivity - The dispersivity of the flow. 34*c4762a1bSJed Brown -rate_constant - The rate constant for the chemical reaction 35*c4762a1bSJed Brown -stoich - The stoichiometry matrix for the reaction 36*c4762a1bSJed Brown -sp_inflow - The species concentrations at the inflow 37*c4762a1bSJed Brown -sp_0 - The species concentrations at time 0 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown */ 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown #include <petscdm.h> 42*c4762a1bSJed Brown #include <petscdmda.h> 43*c4762a1bSJed Brown #include <petscsnes.h> 44*c4762a1bSJed Brown #include <petscts.h> 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown #define N_SPECIES 3 47*c4762a1bSJed Brown #define N_REACTIONS 1 48*c4762a1bSJed Brown #define DIM 2 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown #define stoich(i, j) ctx->stoichiometry[N_SPECIES*i + j] 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown typedef struct { 53*c4762a1bSJed Brown PetscScalar sp[N_SPECIES]; 54*c4762a1bSJed Brown } Field; 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown typedef struct { 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown Field x_inflow; 60*c4762a1bSJed Brown Field x_0; 61*c4762a1bSJed Brown PetscReal stoichiometry[N_SPECIES*N_REACTIONS]; 62*c4762a1bSJed Brown PetscReal porosity; 63*c4762a1bSJed Brown PetscReal dispersivity; 64*c4762a1bSJed Brown PetscReal saturation; 65*c4762a1bSJed Brown PetscReal rate_constant[N_REACTIONS]; 66*c4762a1bSJed Brown PetscReal gradq_inflow; 67*c4762a1bSJed Brown PetscReal length[DIM]; 68*c4762a1bSJed Brown } AppCtx; 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(DM da,AppCtx *ctx,Vec X); 71*c4762a1bSJed Brown extern PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,AppCtx*); 72*c4762a1bSJed Brown extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 73*c4762a1bSJed Brown extern PetscErrorCode ReactingFlowPostCheck(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown PetscErrorCode SetFromOptions(AppCtx * ctx) 77*c4762a1bSJed Brown { 78*c4762a1bSJed Brown PetscErrorCode ierr; 79*c4762a1bSJed Brown PetscInt i,j; 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown PetscFunctionBeginUser; 82*c4762a1bSJed Brown ctx->dispersivity = 0.5; 83*c4762a1bSJed Brown ctx->porosity = 0.25; 84*c4762a1bSJed Brown ctx->saturation = 1.0; 85*c4762a1bSJed Brown ctx->gradq_inflow = 1.0; 86*c4762a1bSJed Brown ctx->rate_constant[0] = 0.5; 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown ctx->length[0] = 100.0; 89*c4762a1bSJed Brown ctx->length[1] = 100.0; 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ctx->x_0.sp[0] = 0.0; 92*c4762a1bSJed Brown ctx->x_0.sp[1] = 0.0; 93*c4762a1bSJed Brown ctx->x_0.sp[2] = 0.0; 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown ctx->x_inflow.sp[0] = 0.05; 96*c4762a1bSJed Brown ctx->x_inflow.sp[1] = 0.05; 97*c4762a1bSJed Brown ctx->x_inflow.sp[2] = 0.0; 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown for (i = 0; i < N_REACTIONS; i++) { 100*c4762a1bSJed Brown for (j = 0; j < N_SPECIES; j++) stoich(i, j) = 0.; 101*c4762a1bSJed Brown } 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown /* set up a pretty easy example */ 104*c4762a1bSJed Brown stoich(0, 0) = -1.; 105*c4762a1bSJed Brown stoich(0, 1) = -2.; 106*c4762a1bSJed Brown stoich(0, 2) = 1.; 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown PetscInt as = N_SPECIES; 109*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-length_x",&ctx->length[0],NULL);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-length_y",&ctx->length[1],NULL);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-porosity",&ctx->porosity,NULL);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-saturation",&ctx->saturation,NULL);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-dispersivity",&ctx->dispersivity,NULL);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-gradq_inflow",&ctx->gradq_inflow,NULL);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-rate_constant",&ctx->rate_constant[0],NULL);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = PetscOptionsGetRealArray(NULL,NULL,"-sp_inflow",ctx->x_inflow.sp,&as,NULL);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = PetscOptionsGetRealArray(NULL,NULL,"-sp_0",ctx->x_0.sp,&as,NULL);CHKERRQ(ierr); 118*c4762a1bSJed Brown as = N_SPECIES; 119*c4762a1bSJed Brown ierr = PetscOptionsGetRealArray(NULL,NULL,"-stoich",ctx->stoichiometry,&as,NULL);CHKERRQ(ierr); 120*c4762a1bSJed Brown PetscFunctionReturn(0); 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown int main(int argc,char **argv) 124*c4762a1bSJed Brown { 125*c4762a1bSJed Brown TS ts; 126*c4762a1bSJed Brown SNES snes; 127*c4762a1bSJed Brown SNESLineSearch linesearch; 128*c4762a1bSJed Brown Vec x; 129*c4762a1bSJed Brown AppCtx ctx; 130*c4762a1bSJed Brown PetscErrorCode ierr; 131*c4762a1bSJed Brown DM da; 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 134*c4762a1bSJed Brown ierr = SetFromOptions(&ctx);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = TSSetIFunction(ts, NULL, FormIFunction, &ctx);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,N_SPECIES,1,NULL,NULL,&da);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 141*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"species A");CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"species B");CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = DMDASetFieldName(da,2,"species C");CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = DMSetApplicationContext(da,&ctx);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = FormInitialGuess(da, &ctx, x);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown ierr = TSSetDM(ts, da);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,1000.0);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,1.0);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = SNESLineSearchSetPostCheck(linesearch, ReactingFlowPostCheck, (void*)&ctx);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = PetscFinalize(); 167*c4762a1bSJed Brown PetscFunctionReturn(0); 168*c4762a1bSJed Brown } 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(DM da,AppCtx *ctx,Vec X) 173*c4762a1bSJed Brown { 174*c4762a1bSJed Brown PetscInt i,j,l,Mx,My,xs,ys,xm,ym; 175*c4762a1bSJed Brown PetscErrorCode ierr; 176*c4762a1bSJed Brown Field **x; 177*c4762a1bSJed Brown 178*c4762a1bSJed Brown PetscFunctionBeginUser; 179*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE, 180*c4762a1bSJed Brown PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE); 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 186*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 187*c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 188*c4762a1bSJed Brown if (i == 0) { 189*c4762a1bSJed Brown if (l == 0) x[j][i].sp[l] = (ctx->x_inflow.sp[l]*((PetscScalar)j) / (My - 1)); 190*c4762a1bSJed Brown else if (l == 1) x[j][i].sp[l] = (ctx->x_inflow.sp[l]*(1. - ((PetscScalar)j) / (My - 1))); 191*c4762a1bSJed Brown else x[j][i].sp[l] = ctx->x_0.sp[l]; 192*c4762a1bSJed Brown } 193*c4762a1bSJed Brown } 194*c4762a1bSJed Brown } 195*c4762a1bSJed Brown } 196*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr); 197*c4762a1bSJed Brown PetscFunctionReturn(0); 198*c4762a1bSJed Brown 199*c4762a1bSJed Brown } 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscScalar ptime,Field **x,Field **xt,Field **f,AppCtx *ctx) 202*c4762a1bSJed Brown { 203*c4762a1bSJed Brown PetscErrorCode ierr; 204*c4762a1bSJed Brown PetscInt i,j,l,m; 205*c4762a1bSJed Brown PetscReal hx,hy,dhx,dhy,hxdhy,hydhx,scale; 206*c4762a1bSJed Brown PetscScalar u,uxx,uyy; 207*c4762a1bSJed Brown PetscScalar vx, vy,sxp,syp,sxm,sym,avx,vxp,vxm,avy,vyp,vym,f_advect; 208*c4762a1bSJed Brown PetscScalar rate; 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown PetscFunctionBeginUser; 211*c4762a1bSJed Brown hx = ctx->length[0]/((PetscReal)(info->mx-1)); 212*c4762a1bSJed Brown hy = ctx->length[1]/((PetscReal)(info->my-1)); 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown dhx = 1. / hx; 215*c4762a1bSJed Brown dhy = 1. / hy; 216*c4762a1bSJed Brown hxdhy = hx/hy; 217*c4762a1bSJed Brown hydhx = hy/hx; 218*c4762a1bSJed Brown scale = hx*hy; 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown for (j=info->ys; j<info->ys+info->ym; j++) { 221*c4762a1bSJed Brown for (i=info->xs; i<info->xs+info->xm; i++) { 222*c4762a1bSJed Brown vx = ctx->gradq_inflow*ctx->porosity*ctx->saturation; 223*c4762a1bSJed Brown vy = 0.0*dhy; 224*c4762a1bSJed Brown avx = PetscAbsScalar(vx); 225*c4762a1bSJed Brown vxp = .5*(vx+avx); vxm = .5*(vx-avx); 226*c4762a1bSJed Brown avy = PetscAbsScalar(vy); 227*c4762a1bSJed Brown vyp = .5*(vy+avy); vym = .5*(vy-avy); 228*c4762a1bSJed Brown /* chemical reactions */ 229*c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 230*c4762a1bSJed Brown /* determine the velocites as the gradients of the pressure */ 231*c4762a1bSJed Brown if (i == 0) { 232*c4762a1bSJed Brown sxp = (x[j][i+1].sp[l] - x[j][i].sp[l])*dhx; 233*c4762a1bSJed Brown sxm = sxp; 234*c4762a1bSJed Brown } else if (i == info->mx - 1) { 235*c4762a1bSJed Brown sxp = (x[j][i].sp[l] - x[j][i-1].sp[l])*dhx; 236*c4762a1bSJed Brown sxm = sxp; 237*c4762a1bSJed Brown } else { 238*c4762a1bSJed Brown sxm = (x[j][i+1].sp[l] - x[j][i].sp[l])*dhx; 239*c4762a1bSJed Brown sxp = (x[j][i].sp[l] - x[j][i-1].sp[l])*dhx; 240*c4762a1bSJed Brown } 241*c4762a1bSJed Brown if (j == 0) { 242*c4762a1bSJed Brown syp = (x[j+1][i].sp[l] - x[j][i].sp[l])*dhy; 243*c4762a1bSJed Brown sym = syp; 244*c4762a1bSJed Brown } else if (j == info->my - 1) { 245*c4762a1bSJed Brown syp = (x[j][i].sp[l] - x[j-1][i].sp[l])*dhy; 246*c4762a1bSJed Brown sym = syp; 247*c4762a1bSJed Brown } else { 248*c4762a1bSJed Brown sym = (x[j+1][i].sp[l] - x[j][i].sp[l])*dhy; 249*c4762a1bSJed Brown syp = (x[j][i].sp[l] - x[j-1][i].sp[l])*dhy; 250*c4762a1bSJed Brown } /* 4 flops per species*point */ 251*c4762a1bSJed Brown 252*c4762a1bSJed Brown if (i == 0) { 253*c4762a1bSJed Brown if (l == 0) f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l]*((PetscScalar)j) / (info->my - 1)); 254*c4762a1bSJed Brown else if (l == 1) f[j][i].sp[l] = (x[j][i].sp[l] - ctx->x_inflow.sp[l]*(1. - ((PetscScalar)j) / (info->my - 1))); 255*c4762a1bSJed Brown else f[j][i].sp[l] = x[j][i].sp[l]; 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown } else { 258*c4762a1bSJed Brown f[j][i].sp[l] = xt[j][i].sp[l]*scale; 259*c4762a1bSJed Brown u = x[j][i].sp[l]; 260*c4762a1bSJed Brown if (j == 0) uyy = u - x[j+1][i].sp[l]; 261*c4762a1bSJed Brown else if (j == info->my - 1) uyy = u - x[j-1][i].sp[l]; 262*c4762a1bSJed Brown else uyy = (2.0*u - x[j-1][i].sp[l] - x[j+1][i].sp[l])*hxdhy; 263*c4762a1bSJed Brown 264*c4762a1bSJed Brown if (i != info->mx - 1) uxx = (2.0*u - x[j][i-1].sp[l] - x[j][i+1].sp[l])*hydhx; 265*c4762a1bSJed Brown else uxx = u - x[j][i-1].sp[l]; 266*c4762a1bSJed Brown /* 10 flops per species*point */ 267*c4762a1bSJed Brown 268*c4762a1bSJed Brown f_advect = 0.; 269*c4762a1bSJed Brown f_advect += scale*(vxp*sxp + vxm*sxm); 270*c4762a1bSJed Brown f_advect += scale*(vyp*syp + vym*sym); 271*c4762a1bSJed Brown f[j][i].sp[l] += f_advect + ctx->dispersivity*(uxx + uyy); 272*c4762a1bSJed Brown /* 14 flops per species*point */ 273*c4762a1bSJed Brown } 274*c4762a1bSJed Brown } 275*c4762a1bSJed Brown /* reaction */ 276*c4762a1bSJed Brown if (i != 0) { 277*c4762a1bSJed Brown for (m = 0; m < N_REACTIONS; m++) { 278*c4762a1bSJed Brown rate = ctx->rate_constant[m]; 279*c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 280*c4762a1bSJed Brown if (stoich(m, l) < 0) { 281*c4762a1bSJed Brown /* assume an elementary reaction */ 282*c4762a1bSJed Brown rate *= PetscPowScalar(x[j][i].sp[l], PetscAbsScalar(stoich(m, l))); 283*c4762a1bSJed Brown /* ~10 flops per reaction per species per point */ 284*c4762a1bSJed Brown } 285*c4762a1bSJed Brown } 286*c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 287*c4762a1bSJed Brown f[j][i].sp[l] += -scale*stoich(m, l)*rate; /* Reaction term */ 288*c4762a1bSJed Brown /* 3 flops per reaction per species per point */ 289*c4762a1bSJed Brown } 290*c4762a1bSJed Brown } 291*c4762a1bSJed Brown } 292*c4762a1bSJed Brown } 293*c4762a1bSJed Brown } 294*c4762a1bSJed Brown ierr = PetscLogFlops((N_SPECIES*(28 + 13*N_REACTIONS))*info->ym*info->xm);CHKERRQ(ierr); 295*c4762a1bSJed Brown PetscFunctionReturn(0); 296*c4762a1bSJed Brown } 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown 299*c4762a1bSJed Brown PetscErrorCode ReactingFlowPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_y, PetscBool *changed_w, void *vctx) 300*c4762a1bSJed Brown { 301*c4762a1bSJed Brown PetscInt i,j,l,Mx,My,xs,ys,xm,ym; 302*c4762a1bSJed Brown PetscErrorCode ierr; 303*c4762a1bSJed Brown Field **x; 304*c4762a1bSJed Brown SNES snes; 305*c4762a1bSJed Brown DM da; 306*c4762a1bSJed Brown PetscScalar min; 307*c4762a1bSJed Brown 308*c4762a1bSJed Brown PetscFunctionBeginUser; 309*c4762a1bSJed Brown *changed_w = PETSC_FALSE; 310*c4762a1bSJed Brown ierr = VecMin(X,NULL,&min);CHKERRQ(ierr); 311*c4762a1bSJed Brown if (min >= 0.) PetscFunctionReturn(0); 312*c4762a1bSJed Brown 313*c4762a1bSJed Brown *changed_w = PETSC_TRUE; 314*c4762a1bSJed Brown ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 315*c4762a1bSJed Brown ierr = SNESGetDM(snes,&da);CHKERRQ(ierr); 316*c4762a1bSJed Brown ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE, 317*c4762a1bSJed Brown PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE); 318*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,W,&x);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 320*c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 321*c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 322*c4762a1bSJed Brown for (l = 0; l < N_SPECIES; l++) { 323*c4762a1bSJed Brown if (x[j][i].sp[l] < 0.) x[j][i].sp[l] = 0.; 324*c4762a1bSJed Brown } 325*c4762a1bSJed Brown } 326*c4762a1bSJed Brown } 327*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,W,&x);CHKERRQ(ierr); 328*c4762a1bSJed Brown PetscFunctionReturn(0); 329*c4762a1bSJed Brown } 330*c4762a1bSJed Brown 331*c4762a1bSJed Brown 332*c4762a1bSJed Brown PetscErrorCode FormIFunction(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *user) 333*c4762a1bSJed Brown { 334*c4762a1bSJed Brown DMDALocalInfo info; 335*c4762a1bSJed Brown Field **u,**udot,**fu; 336*c4762a1bSJed Brown PetscErrorCode ierr; 337*c4762a1bSJed Brown Vec localX; 338*c4762a1bSJed Brown DM da; 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown PetscFunctionBeginUser; 341*c4762a1bSJed Brown ierr = TSGetDM(ts,(DM*)&da);CHKERRQ(ierr); 342*c4762a1bSJed Brown ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr); 343*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 344*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 345*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 346*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,localX,&u);CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da,Xdot,&udot);CHKERRQ(ierr); 348*c4762a1bSJed Brown ierr = DMDAVecGetArray(da,F,&fu);CHKERRQ(ierr); 349*c4762a1bSJed Brown ierr = FormIFunctionLocal(&info,ptime,u,udot,fu,(AppCtx*)user);CHKERRQ(ierr); 350*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,localX,&u);CHKERRQ(ierr); 351*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da,Xdot,&udot);CHKERRQ(ierr); 352*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da,F,&fu);CHKERRQ(ierr); 353*c4762a1bSJed Brown ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr); 354*c4762a1bSJed Brown PetscFunctionReturn(0); 355*c4762a1bSJed Brown } 356