1e77315bcSPatrick Sanan 2e77315bcSPatrick Sanan static char help[] ="Nonlinear Radiative Transport PDE with multigrid in 3d.\n\ 3e77315bcSPatrick Sanan Uses 3-dimensional distributed arrays.\n\ 4e77315bcSPatrick Sanan A 3-dim simplified Radiative Transport test problem is used, with analytic Jacobian. \n\ 5e77315bcSPatrick Sanan \n\ 6e77315bcSPatrick Sanan Solves the linear systems via multilevel methods \n\ 7e77315bcSPatrick Sanan \n\ 8e77315bcSPatrick Sanan The command line\n\ 9e77315bcSPatrick Sanan options are:\n\ 10e77315bcSPatrick Sanan -tleft <tl>, where <tl> indicates the left Diriclet BC \n\ 11e77315bcSPatrick Sanan -tright <tr>, where <tr> indicates the right Diriclet BC \n\ 12e77315bcSPatrick Sanan -beta <beta>, where <beta> indicates the exponent in T \n\n"; 13e77315bcSPatrick Sanan 14e77315bcSPatrick Sanan /* 15e77315bcSPatrick Sanan 16e77315bcSPatrick Sanan This example models the partial differential equation 17e77315bcSPatrick Sanan 18e77315bcSPatrick Sanan - Div(alpha* T^beta (GRAD T)) = 0. 19e77315bcSPatrick Sanan 20e77315bcSPatrick Sanan where beta = 2.5 and alpha = 1.0 21e77315bcSPatrick Sanan 22e77315bcSPatrick Sanan BC: T_left = 1.0, T_right = 0.1, dT/dn_top = dTdn_bottom = dT/dn_up = dT/dn_down = 0. 23e77315bcSPatrick Sanan 24e77315bcSPatrick Sanan in the unit square, which is uniformly discretized in each of x and 25e77315bcSPatrick Sanan y in this simple encoding. The degrees of freedom are cell centered. 26e77315bcSPatrick Sanan 27e77315bcSPatrick Sanan A finite volume approximation with the usual 7-point stencil 28e77315bcSPatrick Sanan is used to discretize the boundary value problem to obtain a 29e77315bcSPatrick Sanan nonlinear system of equations. 30e77315bcSPatrick Sanan 31e77315bcSPatrick Sanan This code was contributed by Nickolas Jovanovic based on ex18.c 32e77315bcSPatrick Sanan 33e77315bcSPatrick Sanan */ 34e77315bcSPatrick Sanan 35e77315bcSPatrick Sanan #include <petscsnes.h> 36e77315bcSPatrick Sanan #include <petscdm.h> 37e77315bcSPatrick Sanan #include <petscdmda.h> 38e77315bcSPatrick Sanan 39e77315bcSPatrick Sanan /* User-defined application context */ 40e77315bcSPatrick Sanan 41e77315bcSPatrick Sanan typedef struct { 42e77315bcSPatrick Sanan PetscReal tleft,tright; /* Dirichlet boundary conditions */ 43e77315bcSPatrick Sanan PetscReal beta,bm1,coef; /* nonlinear diffusivity parameterizations */ 44e77315bcSPatrick Sanan } AppCtx; 45e77315bcSPatrick Sanan 46e77315bcSPatrick Sanan #define POWFLOP 5 /* assume a pow() takes five flops */ 47e77315bcSPatrick Sanan 48e77315bcSPatrick Sanan extern PetscErrorCode FormInitialGuess(SNES,Vec,void*); 49e77315bcSPatrick Sanan extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 50e77315bcSPatrick Sanan extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 51e77315bcSPatrick Sanan 52e77315bcSPatrick Sanan int main(int argc,char **argv) 53e77315bcSPatrick Sanan { 54e77315bcSPatrick Sanan SNES snes; 55e77315bcSPatrick Sanan AppCtx user; 56e77315bcSPatrick Sanan PetscInt its,lits; 57e77315bcSPatrick Sanan PetscReal litspit; 58e77315bcSPatrick Sanan DM da; 59e77315bcSPatrick Sanan 60*327415f7SBarry Smith PetscFunctionBeginUser; 619566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 62e77315bcSPatrick Sanan /* set problem parameters */ 63e77315bcSPatrick Sanan user.tleft = 1.0; 64e77315bcSPatrick Sanan user.tright = 0.1; 65e77315bcSPatrick Sanan user.beta = 2.5; 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-tleft",&user.tleft,NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-tright",&user.tright,NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-beta",&user.beta,NULL)); 69e77315bcSPatrick Sanan user.bm1 = user.beta - 1.0; 70e77315bcSPatrick Sanan user.coef = user.beta/2.0; 71e77315bcSPatrick Sanan 72e77315bcSPatrick Sanan /* 73e77315bcSPatrick Sanan Set the DMDA (grid structure) for the grids. 74e77315bcSPatrick Sanan */ 759566063dSJacob Faibussowitsch PetscCall(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,5,5,5,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,0,&da)); 769566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 779566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 789566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da,&user)); 79e77315bcSPatrick Sanan 80e77315bcSPatrick Sanan /* 81e77315bcSPatrick Sanan Create the nonlinear solver 82e77315bcSPatrick Sanan */ 839566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 849566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes,da)); 859566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,NULL,FormFunction,&user)); 869566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,NULL,NULL,FormJacobian,&user)); 879566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 889566063dSJacob Faibussowitsch PetscCall(SNESSetComputeInitialGuess(snes,FormInitialGuess,NULL)); 89e77315bcSPatrick Sanan 909566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,NULL)); 919566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 929566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes,&lits)); 93e77315bcSPatrick Sanan litspit = ((PetscReal)lits)/((PetscReal)its); 9463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its)); 9563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %" PetscInt_FMT "\n",lits)); 969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / SNES = %e\n",(double)litspit)); 97e77315bcSPatrick Sanan 989566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 101b122ec5aSJacob Faibussowitsch return 0; 102e77315bcSPatrick Sanan } 103e77315bcSPatrick Sanan /* -------------------- Form initial approximation ----------------- */ 104e77315bcSPatrick Sanan PetscErrorCode FormInitialGuess(SNES snes,Vec X,void *ctx) 105e77315bcSPatrick Sanan { 106e77315bcSPatrick Sanan AppCtx *user; 107e77315bcSPatrick Sanan PetscInt i,j,k,xs,ys,xm,ym,zs,zm; 108e77315bcSPatrick Sanan PetscScalar ***x; 109e77315bcSPatrick Sanan DM da; 110e77315bcSPatrick Sanan 111e77315bcSPatrick Sanan PetscFunctionBeginUser; 1129566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&da)); 1139566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da,&user)); 1149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 1159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,X,&x)); 116e77315bcSPatrick Sanan 117e77315bcSPatrick Sanan /* Compute initial guess */ 118e77315bcSPatrick Sanan for (k=zs; k<zs+zm; k++) { 119e77315bcSPatrick Sanan for (j=ys; j<ys+ym; j++) { 120e77315bcSPatrick Sanan for (i=xs; i<xs+xm; i++) { 121e77315bcSPatrick Sanan x[k][j][i] = user->tleft; 122e77315bcSPatrick Sanan } 123e77315bcSPatrick Sanan } 124e77315bcSPatrick Sanan } 1259566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,X,&x)); 126e77315bcSPatrick Sanan PetscFunctionReturn(0); 127e77315bcSPatrick Sanan } 128e77315bcSPatrick Sanan /* -------------------- Evaluate Function F(x) --------------------- */ 129e77315bcSPatrick Sanan PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr) 130e77315bcSPatrick Sanan { 131e77315bcSPatrick Sanan AppCtx *user = (AppCtx*)ptr; 132e77315bcSPatrick Sanan PetscInt i,j,k,mx,my,mz,xs,ys,zs,xm,ym,zm; 133e77315bcSPatrick Sanan PetscScalar zero = 0.0,one = 1.0; 134e77315bcSPatrick Sanan PetscScalar hx,hy,hz,hxhydhz,hyhzdhx,hzhxdhy; 135e77315bcSPatrick Sanan PetscScalar t0,tn,ts,te,tw,an,as,ae,aw,dn,ds,de,dw,fn = 0.0,fs = 0.0,fe =0.0,fw = 0.0; 136e77315bcSPatrick Sanan PetscScalar tleft,tright,beta,td,ad,dd,fd=0.0,tu,au,du=0.0,fu=0.0; 137e77315bcSPatrick Sanan PetscScalar ***x,***f; 138e77315bcSPatrick Sanan Vec localX; 139e77315bcSPatrick Sanan DM da; 140e77315bcSPatrick Sanan 141e77315bcSPatrick Sanan PetscFunctionBeginUser; 1429566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&da)); 1439566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localX)); 1449566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 145e77315bcSPatrick Sanan hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1); hz = one/(PetscReal)(mz-1); 146e77315bcSPatrick Sanan hxhydhz = hx*hy/hz; hyhzdhx = hy*hz/hx; hzhxdhy = hz*hx/hy; 147e77315bcSPatrick Sanan tleft = user->tleft; tright = user->tright; 148e77315bcSPatrick Sanan beta = user->beta; 149e77315bcSPatrick Sanan 150e77315bcSPatrick Sanan /* Get ghost points */ 1519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 1529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 1539566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 1549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,localX,&x)); 1559566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 156e77315bcSPatrick Sanan 157e77315bcSPatrick Sanan /* Evaluate function */ 158e77315bcSPatrick Sanan for (k=zs; k<zs+zm; k++) { 159e77315bcSPatrick Sanan for (j=ys; j<ys+ym; j++) { 160e77315bcSPatrick Sanan for (i=xs; i<xs+xm; i++) { 161e77315bcSPatrick Sanan t0 = x[k][j][i]; 162e77315bcSPatrick Sanan 163e77315bcSPatrick Sanan if (i > 0 && i < mx-1 && j > 0 && j < my-1 && k > 0 && k < mz-1) { 164e77315bcSPatrick Sanan 165e77315bcSPatrick Sanan /* general interior volume */ 166e77315bcSPatrick Sanan 167e77315bcSPatrick Sanan tw = x[k][j][i-1]; 168e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 169e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 170e77315bcSPatrick Sanan fw = dw*(t0 - tw); 171e77315bcSPatrick Sanan 172e77315bcSPatrick Sanan te = x[k][j][i+1]; 173e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 174e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 175e77315bcSPatrick Sanan fe = de*(te - t0); 176e77315bcSPatrick Sanan 177e77315bcSPatrick Sanan ts = x[k][j-1][i]; 178e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 179e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 180e77315bcSPatrick Sanan fs = ds*(t0 - ts); 181e77315bcSPatrick Sanan 182e77315bcSPatrick Sanan tn = x[k][j+1][i]; 183e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 184e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 185e77315bcSPatrick Sanan fn = dn*(tn - t0); 186e77315bcSPatrick Sanan 187e77315bcSPatrick Sanan td = x[k-1][j][i]; 188e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 189e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 190e77315bcSPatrick Sanan fd = dd*(t0 - td); 191e77315bcSPatrick Sanan 192e77315bcSPatrick Sanan tu = x[k+1][j][i]; 193e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 194e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 195e77315bcSPatrick Sanan fu = du*(tu - t0); 196e77315bcSPatrick Sanan 197e77315bcSPatrick Sanan } else if (i == 0) { 198e77315bcSPatrick Sanan 199e77315bcSPatrick Sanan /* left-hand (west) boundary */ 200e77315bcSPatrick Sanan tw = tleft; 201e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 202e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 203e77315bcSPatrick Sanan fw = dw*(t0 - tw); 204e77315bcSPatrick Sanan 205e77315bcSPatrick Sanan te = x[k][j][i+1]; 206e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 207e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 208e77315bcSPatrick Sanan fe = de*(te - t0); 209e77315bcSPatrick Sanan 210e77315bcSPatrick Sanan if (j > 0) { 211e77315bcSPatrick Sanan ts = x[k][j-1][i]; 212e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 213e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 214e77315bcSPatrick Sanan fs = ds*(t0 - ts); 215e77315bcSPatrick Sanan } else { 216e77315bcSPatrick Sanan fs = zero; 217e77315bcSPatrick Sanan } 218e77315bcSPatrick Sanan 219e77315bcSPatrick Sanan if (j < my-1) { 220e77315bcSPatrick Sanan tn = x[k][j+1][i]; 221e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 222e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 223e77315bcSPatrick Sanan fn = dn*(tn - t0); 224e77315bcSPatrick Sanan } else { 225e77315bcSPatrick Sanan fn = zero; 226e77315bcSPatrick Sanan } 227e77315bcSPatrick Sanan 228e77315bcSPatrick Sanan if (k > 0) { 229e77315bcSPatrick Sanan td = x[k-1][j][i]; 230e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 231e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 232e77315bcSPatrick Sanan fd = dd*(t0 - td); 233e77315bcSPatrick Sanan } else { 234e77315bcSPatrick Sanan fd = zero; 235e77315bcSPatrick Sanan } 236e77315bcSPatrick Sanan 237e77315bcSPatrick Sanan if (k < mz-1) { 238e77315bcSPatrick Sanan tu = x[k+1][j][i]; 239e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 240e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 241e77315bcSPatrick Sanan fu = du*(tu - t0); 242e77315bcSPatrick Sanan } else { 243e77315bcSPatrick Sanan fu = zero; 244e77315bcSPatrick Sanan } 245e77315bcSPatrick Sanan 246e77315bcSPatrick Sanan } else if (i == mx-1) { 247e77315bcSPatrick Sanan 248e77315bcSPatrick Sanan /* right-hand (east) boundary */ 249e77315bcSPatrick Sanan tw = x[k][j][i-1]; 250e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 251e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 252e77315bcSPatrick Sanan fw = dw*(t0 - tw); 253e77315bcSPatrick Sanan 254e77315bcSPatrick Sanan te = tright; 255e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 256e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 257e77315bcSPatrick Sanan fe = de*(te - t0); 258e77315bcSPatrick Sanan 259e77315bcSPatrick Sanan if (j > 0) { 260e77315bcSPatrick Sanan ts = x[k][j-1][i]; 261e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 262e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 263e77315bcSPatrick Sanan fs = ds*(t0 - ts); 264e77315bcSPatrick Sanan } else { 265e77315bcSPatrick Sanan fs = zero; 266e77315bcSPatrick Sanan } 267e77315bcSPatrick Sanan 268e77315bcSPatrick Sanan if (j < my-1) { 269e77315bcSPatrick Sanan tn = x[k][j+1][i]; 270e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 271e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 272e77315bcSPatrick Sanan fn = dn*(tn - t0); 273e77315bcSPatrick Sanan } else { 274e77315bcSPatrick Sanan fn = zero; 275e77315bcSPatrick Sanan } 276e77315bcSPatrick Sanan 277e77315bcSPatrick Sanan if (k > 0) { 278e77315bcSPatrick Sanan td = x[k-1][j][i]; 279e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 280e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 281e77315bcSPatrick Sanan fd = dd*(t0 - td); 282e77315bcSPatrick Sanan } else { 283e77315bcSPatrick Sanan fd = zero; 284e77315bcSPatrick Sanan } 285e77315bcSPatrick Sanan 286e77315bcSPatrick Sanan if (k < mz-1) { 287e77315bcSPatrick Sanan tu = x[k+1][j][i]; 288e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 289e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 290e77315bcSPatrick Sanan fu = du*(tu - t0); 291e77315bcSPatrick Sanan } else { 292e77315bcSPatrick Sanan fu = zero; 293e77315bcSPatrick Sanan } 294e77315bcSPatrick Sanan 295e77315bcSPatrick Sanan } else if (j == 0) { 296e77315bcSPatrick Sanan 297e77315bcSPatrick Sanan /* bottom (south) boundary, and i <> 0 or mx-1 */ 298e77315bcSPatrick Sanan tw = x[k][j][i-1]; 299e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 300e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 301e77315bcSPatrick Sanan fw = dw*(t0 - tw); 302e77315bcSPatrick Sanan 303e77315bcSPatrick Sanan te = x[k][j][i+1]; 304e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 305e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 306e77315bcSPatrick Sanan fe = de*(te - t0); 307e77315bcSPatrick Sanan 308e77315bcSPatrick Sanan fs = zero; 309e77315bcSPatrick Sanan 310e77315bcSPatrick Sanan tn = x[k][j+1][i]; 311e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 312e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 313e77315bcSPatrick Sanan fn = dn*(tn - t0); 314e77315bcSPatrick Sanan 315e77315bcSPatrick Sanan if (k > 0) { 316e77315bcSPatrick Sanan td = x[k-1][j][i]; 317e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 318e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 319e77315bcSPatrick Sanan fd = dd*(t0 - td); 320e77315bcSPatrick Sanan } else { 321e77315bcSPatrick Sanan fd = zero; 322e77315bcSPatrick Sanan } 323e77315bcSPatrick Sanan 324e77315bcSPatrick Sanan if (k < mz-1) { 325e77315bcSPatrick Sanan tu = x[k+1][j][i]; 326e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 327e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 328e77315bcSPatrick Sanan fu = du*(tu - t0); 329e77315bcSPatrick Sanan } else { 330e77315bcSPatrick Sanan fu = zero; 331e77315bcSPatrick Sanan } 332e77315bcSPatrick Sanan 333e77315bcSPatrick Sanan } else if (j == my-1) { 334e77315bcSPatrick Sanan 335e77315bcSPatrick Sanan /* top (north) boundary, and i <> 0 or mx-1 */ 336e77315bcSPatrick Sanan tw = x[k][j][i-1]; 337e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 338e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 339e77315bcSPatrick Sanan fw = dw*(t0 - tw); 340e77315bcSPatrick Sanan 341e77315bcSPatrick Sanan te = x[k][j][i+1]; 342e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 343e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 344e77315bcSPatrick Sanan fe = de*(te - t0); 345e77315bcSPatrick Sanan 346e77315bcSPatrick Sanan ts = x[k][j-1][i]; 347e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 348e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 349e77315bcSPatrick Sanan fs = ds*(t0 - ts); 350e77315bcSPatrick Sanan 351e77315bcSPatrick Sanan fn = zero; 352e77315bcSPatrick Sanan 353e77315bcSPatrick Sanan if (k > 0) { 354e77315bcSPatrick Sanan td = x[k-1][j][i]; 355e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 356e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 357e77315bcSPatrick Sanan fd = dd*(t0 - td); 358e77315bcSPatrick Sanan } else { 359e77315bcSPatrick Sanan fd = zero; 360e77315bcSPatrick Sanan } 361e77315bcSPatrick Sanan 362e77315bcSPatrick Sanan if (k < mz-1) { 363e77315bcSPatrick Sanan tu = x[k+1][j][i]; 364e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 365e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 366e77315bcSPatrick Sanan fu = du*(tu - t0); 367e77315bcSPatrick Sanan } else { 368e77315bcSPatrick Sanan fu = zero; 369e77315bcSPatrick Sanan } 370e77315bcSPatrick Sanan 371e77315bcSPatrick Sanan } else if (k == 0) { 372e77315bcSPatrick Sanan 373e77315bcSPatrick Sanan /* down boundary (interior only) */ 374e77315bcSPatrick Sanan tw = x[k][j][i-1]; 375e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 376e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 377e77315bcSPatrick Sanan fw = dw*(t0 - tw); 378e77315bcSPatrick Sanan 379e77315bcSPatrick Sanan te = x[k][j][i+1]; 380e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 381e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 382e77315bcSPatrick Sanan fe = de*(te - t0); 383e77315bcSPatrick Sanan 384e77315bcSPatrick Sanan ts = x[k][j-1][i]; 385e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 386e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 387e77315bcSPatrick Sanan fs = ds*(t0 - ts); 388e77315bcSPatrick Sanan 389e77315bcSPatrick Sanan tn = x[k][j+1][i]; 390e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 391e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 392e77315bcSPatrick Sanan fn = dn*(tn - t0); 393e77315bcSPatrick Sanan 394e77315bcSPatrick Sanan fd = zero; 395e77315bcSPatrick Sanan 396e77315bcSPatrick Sanan tu = x[k+1][j][i]; 397e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 398e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 399e77315bcSPatrick Sanan fu = du*(tu - t0); 400e77315bcSPatrick Sanan 401e77315bcSPatrick Sanan } else if (k == mz-1) { 402e77315bcSPatrick Sanan 403e77315bcSPatrick Sanan /* up boundary (interior only) */ 404e77315bcSPatrick Sanan tw = x[k][j][i-1]; 405e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 406e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 407e77315bcSPatrick Sanan fw = dw*(t0 - tw); 408e77315bcSPatrick Sanan 409e77315bcSPatrick Sanan te = x[k][j][i+1]; 410e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 411e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 412e77315bcSPatrick Sanan fe = de*(te - t0); 413e77315bcSPatrick Sanan 414e77315bcSPatrick Sanan ts = x[k][j-1][i]; 415e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 416e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 417e77315bcSPatrick Sanan fs = ds*(t0 - ts); 418e77315bcSPatrick Sanan 419e77315bcSPatrick Sanan tn = x[k][j+1][i]; 420e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 421e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 422e77315bcSPatrick Sanan fn = dn*(tn - t0); 423e77315bcSPatrick Sanan 424e77315bcSPatrick Sanan td = x[k-1][j][i]; 425e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 426e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 427e77315bcSPatrick Sanan fd = dd*(t0 - td); 428e77315bcSPatrick Sanan 429e77315bcSPatrick Sanan fu = zero; 430e77315bcSPatrick Sanan } 431e77315bcSPatrick Sanan 432e77315bcSPatrick Sanan f[k][j][i] = -hyhzdhx*(fe-fw) - hzhxdhy*(fn-fs) - hxhydhz*(fu-fd); 433e77315bcSPatrick Sanan } 434e77315bcSPatrick Sanan } 435e77315bcSPatrick Sanan } 4369566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,localX,&x)); 4379566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 4389566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localX)); 4399566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((22.0 + 4.0*POWFLOP)*ym*xm)); 440e77315bcSPatrick Sanan PetscFunctionReturn(0); 441e77315bcSPatrick Sanan } 442e77315bcSPatrick Sanan /* -------------------- Evaluate Jacobian F(x) --------------------- */ 443e77315bcSPatrick Sanan PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr) 444e77315bcSPatrick Sanan { 445e77315bcSPatrick Sanan AppCtx *user = (AppCtx*)ptr; 446e77315bcSPatrick Sanan PetscInt i,j,k,mx,my,mz,xs,ys,zs,xm,ym,zm; 447e77315bcSPatrick Sanan PetscScalar one = 1.0; 448e77315bcSPatrick Sanan PetscScalar hx,hy,hz,hxhydhz,hyhzdhx,hzhxdhy; 449e77315bcSPatrick Sanan PetscScalar t0,tn,ts,te,tw,an,as,ae,aw,dn,ds,de,dw; 450e77315bcSPatrick Sanan PetscScalar tleft,tright,beta,td,ad,dd,tu,au,du,v[7],bm1,coef; 451e77315bcSPatrick Sanan PetscScalar ***x,bn,bs,be,bw,bu,bd,gn,gs,ge,gw,gu,gd; 452e77315bcSPatrick Sanan Vec localX; 453e77315bcSPatrick Sanan MatStencil c[7],row; 454e77315bcSPatrick Sanan DM da; 455e77315bcSPatrick Sanan 456e77315bcSPatrick Sanan PetscFunctionBeginUser; 4579566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&da)); 4589566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localX)); 4599566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,NULL,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 460e77315bcSPatrick Sanan hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1); hz = one/(PetscReal)(mz-1); 461e77315bcSPatrick Sanan hxhydhz = hx*hy/hz; hyhzdhx = hy*hz/hx; hzhxdhy = hz*hx/hy; 462e77315bcSPatrick Sanan tleft = user->tleft; tright = user->tright; 463e77315bcSPatrick Sanan beta = user->beta; bm1 = user->bm1; coef = user->coef; 464e77315bcSPatrick Sanan 465e77315bcSPatrick Sanan /* Get ghost points */ 4669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 4679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 4689566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 4699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,localX,&x)); 470e77315bcSPatrick Sanan 471e77315bcSPatrick Sanan /* Evaluate Jacobian of function */ 472e77315bcSPatrick Sanan for (k=zs; k<zs+zm; k++) { 473e77315bcSPatrick Sanan for (j=ys; j<ys+ym; j++) { 474e77315bcSPatrick Sanan for (i=xs; i<xs+xm; i++) { 475e77315bcSPatrick Sanan t0 = x[k][j][i]; 476e77315bcSPatrick Sanan row.k = k; row.j = j; row.i = i; 477e77315bcSPatrick Sanan if (i > 0 && i < mx-1 && j > 0 && j < my-1 && k > 0 && k < mz-1) { 478e77315bcSPatrick Sanan 479e77315bcSPatrick Sanan /* general interior volume */ 480e77315bcSPatrick Sanan 481e77315bcSPatrick Sanan tw = x[k][j][i-1]; 482e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 483e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 484e77315bcSPatrick Sanan /* dw = bw * aw */ 485e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 486e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 487e77315bcSPatrick Sanan 488e77315bcSPatrick Sanan te = x[k][j][i+1]; 489e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 490e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 491e77315bcSPatrick Sanan /* de = be * ae; */ 492e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 493e77315bcSPatrick Sanan ge = coef*be*(te - t0); 494e77315bcSPatrick Sanan 495e77315bcSPatrick Sanan ts = x[k][j-1][i]; 496e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 497e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 498e77315bcSPatrick Sanan /* ds = bs * as; */ 499e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 500e77315bcSPatrick Sanan gs = coef*bs*(t0 - ts); 501e77315bcSPatrick Sanan 502e77315bcSPatrick Sanan tn = x[k][j+1][i]; 503e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 504e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 505e77315bcSPatrick Sanan /* dn = bn * an; */ 506e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 507e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 508e77315bcSPatrick Sanan 509e77315bcSPatrick Sanan td = x[k-1][j][i]; 510e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 511e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 512e77315bcSPatrick Sanan /* dd = bd * ad; */ 513e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 514e77315bcSPatrick Sanan gd = coef*bd*(t0 - td); 515e77315bcSPatrick Sanan 516e77315bcSPatrick Sanan tu = x[k+1][j][i]; 517e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 518e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 519e77315bcSPatrick Sanan /* du = bu * au; */ 520e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 521e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 522e77315bcSPatrick Sanan 523e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; v[0] = -hxhydhz*(dd - gd); 524e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 525e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 526e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 527e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 528e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 529e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 530e77315bcSPatrick Sanan c[4].k = k; c[4].j = j; c[4].i = i+1; 531e77315bcSPatrick Sanan v[4] = -hyhzdhx*(de + ge); 532e77315bcSPatrick Sanan c[5].k = k; c[5].j = j+1; c[5].i = i; 533e77315bcSPatrick Sanan v[5] = -hzhxdhy*(dn + gn); 534e77315bcSPatrick Sanan c[6].k = k+1; c[6].j = j; c[6].i = i; 535e77315bcSPatrick Sanan v[6] = -hxhydhz*(du + gu); 5369566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,7,c,v,INSERT_VALUES)); 537e77315bcSPatrick Sanan 538e77315bcSPatrick Sanan } else if (i == 0) { 539e77315bcSPatrick Sanan 540e77315bcSPatrick Sanan /* left-hand plane boundary */ 541e77315bcSPatrick Sanan tw = tleft; 542e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 543e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 544e77315bcSPatrick Sanan /* dw = bw * aw */ 545e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 546e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 547e77315bcSPatrick Sanan 548e77315bcSPatrick Sanan te = x[k][j][i+1]; 549e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 550e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 551e77315bcSPatrick Sanan /* de = be * ae; */ 552e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 553e77315bcSPatrick Sanan ge = coef*be*(te - t0); 554e77315bcSPatrick Sanan 555e77315bcSPatrick Sanan /* left-hand bottom edge */ 556e77315bcSPatrick Sanan if (j == 0) { 557e77315bcSPatrick Sanan 558e77315bcSPatrick Sanan tn = x[k][j+1][i]; 559e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 560e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 561e77315bcSPatrick Sanan /* dn = bn * an; */ 562e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 563e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 564e77315bcSPatrick Sanan 565e77315bcSPatrick Sanan /* left-hand bottom down corner */ 566e77315bcSPatrick Sanan if (k == 0) { 567e77315bcSPatrick Sanan 568e77315bcSPatrick Sanan tu = x[k+1][j][i]; 569e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 570e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 571e77315bcSPatrick Sanan /* du = bu * au; */ 572e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 573e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 574e77315bcSPatrick Sanan 575e77315bcSPatrick Sanan c[0].k = k; c[0].j = j; c[0].i = i; 576e77315bcSPatrick Sanan v[0] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 577e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i+1; 578e77315bcSPatrick Sanan v[1] = -hyhzdhx*(de + ge); 579e77315bcSPatrick Sanan c[2].k = k; c[2].j = j+1; c[2].i = i; 580e77315bcSPatrick Sanan v[2] = -hzhxdhy*(dn + gn); 581e77315bcSPatrick Sanan c[3].k = k+1; c[3].j = j; c[3].i = i; 582e77315bcSPatrick Sanan v[3] = -hxhydhz*(du + gu); 5839566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 584e77315bcSPatrick Sanan 585e77315bcSPatrick Sanan /* left-hand bottom interior edge */ 586e77315bcSPatrick Sanan } else if (k < mz-1) { 587e77315bcSPatrick Sanan 588e77315bcSPatrick Sanan tu = x[k+1][j][i]; 589e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 590e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 591e77315bcSPatrick Sanan /* du = bu * au; */ 592e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 593e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 594e77315bcSPatrick Sanan 595e77315bcSPatrick Sanan td = x[k-1][j][i]; 596e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 597e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 598e77315bcSPatrick Sanan /* dd = bd * ad; */ 599e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 600e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 601e77315bcSPatrick Sanan 602e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 603e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 604e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i; 605e77315bcSPatrick Sanan v[1] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 606e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i+1; 607e77315bcSPatrick Sanan v[2] = -hyhzdhx*(de + ge); 608e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 609e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 610e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 611e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 6129566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 613e77315bcSPatrick Sanan 614e77315bcSPatrick Sanan /* left-hand bottom up corner */ 615e77315bcSPatrick Sanan } else { 616e77315bcSPatrick Sanan 617e77315bcSPatrick Sanan td = x[k-1][j][i]; 618e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 619e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 620e77315bcSPatrick Sanan /* dd = bd * ad; */ 621e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 622e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 623e77315bcSPatrick Sanan 624e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 625e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 626e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i; 627e77315bcSPatrick Sanan v[1] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 628e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i+1; 629e77315bcSPatrick Sanan v[2] = -hyhzdhx*(de + ge); 630e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 631e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 6329566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 633e77315bcSPatrick Sanan } 634e77315bcSPatrick Sanan 635e77315bcSPatrick Sanan /* left-hand top edge */ 636e77315bcSPatrick Sanan } else if (j == my-1) { 637e77315bcSPatrick Sanan 638e77315bcSPatrick Sanan ts = x[k][j-1][i]; 639e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 640e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 641e77315bcSPatrick Sanan /* ds = bs * as; */ 642e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 643e77315bcSPatrick Sanan gs = coef*bs*(ts - t0); 644e77315bcSPatrick Sanan 645e77315bcSPatrick Sanan /* left-hand top down corner */ 646e77315bcSPatrick Sanan if (k == 0) { 647e77315bcSPatrick Sanan 648e77315bcSPatrick Sanan tu = x[k+1][j][i]; 649e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 650e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 651e77315bcSPatrick Sanan /* du = bu * au; */ 652e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 653e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 654e77315bcSPatrick Sanan 655e77315bcSPatrick Sanan c[0].k = k; c[0].j = j-1; c[0].i = i; 656e77315bcSPatrick Sanan v[0] = -hzhxdhy*(ds - gs); 657e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i; 658e77315bcSPatrick Sanan v[1] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 659e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i+1; 660e77315bcSPatrick Sanan v[2] = -hyhzdhx*(de + ge); 661e77315bcSPatrick Sanan c[3].k = k+1; c[3].j = j; c[3].i = i; 662e77315bcSPatrick Sanan v[3] = -hxhydhz*(du + gu); 6639566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 664e77315bcSPatrick Sanan 665e77315bcSPatrick Sanan /* left-hand top interior edge */ 666e77315bcSPatrick Sanan } else if (k < mz-1) { 667e77315bcSPatrick Sanan 668e77315bcSPatrick Sanan tu = x[k+1][j][i]; 669e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 670e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 671e77315bcSPatrick Sanan /* du = bu * au; */ 672e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 673e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 674e77315bcSPatrick Sanan 675e77315bcSPatrick Sanan td = x[k-1][j][i]; 676e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 677e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 678e77315bcSPatrick Sanan /* dd = bd * ad; */ 679e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 680e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 681e77315bcSPatrick Sanan 682e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 683e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 684e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 685e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 686e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 687e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 688e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 689e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 690e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 691e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 6929566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 693e77315bcSPatrick Sanan 694e77315bcSPatrick Sanan /* left-hand top up corner */ 695e77315bcSPatrick Sanan } else { 696e77315bcSPatrick Sanan 697e77315bcSPatrick Sanan td = x[k-1][j][i]; 698e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 699e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 700e77315bcSPatrick Sanan /* dd = bd * ad; */ 701e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 702e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 703e77315bcSPatrick Sanan 704e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 705e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 706e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 707e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 708e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 709e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 710e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 711e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 7129566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 713e77315bcSPatrick Sanan } 714e77315bcSPatrick Sanan 715e77315bcSPatrick Sanan } else { 716e77315bcSPatrick Sanan 717e77315bcSPatrick Sanan ts = x[k][j-1][i]; 718e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 719e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 720e77315bcSPatrick Sanan /* ds = bs * as; */ 721e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 722e77315bcSPatrick Sanan gs = coef*bs*(t0 - ts); 723e77315bcSPatrick Sanan 724e77315bcSPatrick Sanan tn = x[k][j+1][i]; 725e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 726e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 727e77315bcSPatrick Sanan /* dn = bn * an; */ 728e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 729e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 730e77315bcSPatrick Sanan 731e77315bcSPatrick Sanan /* left-hand down interior edge */ 732e77315bcSPatrick Sanan if (k == 0) { 733e77315bcSPatrick Sanan 734e77315bcSPatrick Sanan tu = x[k+1][j][i]; 735e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 736e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 737e77315bcSPatrick Sanan /* du = bu * au; */ 738e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 739e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 740e77315bcSPatrick Sanan 741e77315bcSPatrick Sanan c[0].k = k; c[0].j = j-1; c[0].i = i; 742e77315bcSPatrick Sanan v[0] = -hzhxdhy*(ds - gs); 743e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i; 744e77315bcSPatrick Sanan v[1] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 745e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i+1; 746e77315bcSPatrick Sanan v[2] = -hyhzdhx*(de + ge); 747e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 748e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 749e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 750e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 7519566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 752e77315bcSPatrick Sanan 753e77315bcSPatrick Sanan } else if (k == mz-1) { /* left-hand up interior edge */ 754e77315bcSPatrick Sanan 755e77315bcSPatrick Sanan td = x[k-1][j][i]; 756e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 757e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 758e77315bcSPatrick Sanan /* dd = bd * ad; */ 759e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 760e77315bcSPatrick Sanan gd = coef*bd*(t0 - td); 761e77315bcSPatrick Sanan 762e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 763e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 764e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 765e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 766e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 767e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 768e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 769e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 770e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 771e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 7729566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 773e77315bcSPatrick Sanan } else { /* left-hand interior plane */ 774e77315bcSPatrick Sanan 775e77315bcSPatrick Sanan td = x[k-1][j][i]; 776e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 777e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 778e77315bcSPatrick Sanan /* dd = bd * ad; */ 779e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 780e77315bcSPatrick Sanan gd = coef*bd*(t0 - td); 781e77315bcSPatrick Sanan 782e77315bcSPatrick Sanan tu = x[k+1][j][i]; 783e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 784e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 785e77315bcSPatrick Sanan /* du = bu * au; */ 786e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 787e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 788e77315bcSPatrick Sanan 789e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 790e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 791e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 792e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 793e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 794e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 795e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 796e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 797e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 798e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 799e77315bcSPatrick Sanan c[5].k = k+1; c[5].j = j; c[5].i = i; 800e77315bcSPatrick Sanan v[5] = -hxhydhz*(du + gu); 8019566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES)); 802e77315bcSPatrick Sanan } 803e77315bcSPatrick Sanan } 804e77315bcSPatrick Sanan 805e77315bcSPatrick Sanan } else if (i == mx-1) { 806e77315bcSPatrick Sanan 807e77315bcSPatrick Sanan /* right-hand plane boundary */ 808e77315bcSPatrick Sanan tw = x[k][j][i-1]; 809e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 810e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 811e77315bcSPatrick Sanan /* dw = bw * aw */ 812e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 813e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 814e77315bcSPatrick Sanan 815e77315bcSPatrick Sanan te = tright; 816e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 817e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 818e77315bcSPatrick Sanan /* de = be * ae; */ 819e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 820e77315bcSPatrick Sanan ge = coef*be*(te - t0); 821e77315bcSPatrick Sanan 822e77315bcSPatrick Sanan /* right-hand bottom edge */ 823e77315bcSPatrick Sanan if (j == 0) { 824e77315bcSPatrick Sanan 825e77315bcSPatrick Sanan tn = x[k][j+1][i]; 826e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 827e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 828e77315bcSPatrick Sanan /* dn = bn * an; */ 829e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 830e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 831e77315bcSPatrick Sanan 832e77315bcSPatrick Sanan /* right-hand bottom down corner */ 833e77315bcSPatrick Sanan if (k == 0) { 834e77315bcSPatrick Sanan 835e77315bcSPatrick Sanan tu = x[k+1][j][i]; 836e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 837e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 838e77315bcSPatrick Sanan /* du = bu * au; */ 839e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 840e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 841e77315bcSPatrick Sanan 842e77315bcSPatrick Sanan c[0].k = k; c[0].j = j; c[0].i = i-1; 843e77315bcSPatrick Sanan v[0] = -hyhzdhx*(dw - gw); 844e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i; 845e77315bcSPatrick Sanan v[1] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 846e77315bcSPatrick Sanan c[2].k = k; c[2].j = j+1; c[2].i = i; 847e77315bcSPatrick Sanan v[2] = -hzhxdhy*(dn + gn); 848e77315bcSPatrick Sanan c[3].k = k+1; c[3].j = j; c[3].i = i; 849e77315bcSPatrick Sanan v[3] = -hxhydhz*(du + gu); 8509566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 851e77315bcSPatrick Sanan 852e77315bcSPatrick Sanan /* right-hand bottom interior edge */ 853e77315bcSPatrick Sanan } else if (k < mz-1) { 854e77315bcSPatrick Sanan 855e77315bcSPatrick Sanan tu = x[k+1][j][i]; 856e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 857e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 858e77315bcSPatrick Sanan /* du = bu * au; */ 859e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 860e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 861e77315bcSPatrick Sanan 862e77315bcSPatrick Sanan td = x[k-1][j][i]; 863e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 864e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 865e77315bcSPatrick Sanan /* dd = bd * ad; */ 866e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 867e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 868e77315bcSPatrick Sanan 869e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 870e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 871e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 872e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 873e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 874e77315bcSPatrick Sanan v[2] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 875e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 876e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 877e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 878e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 8799566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 880e77315bcSPatrick Sanan 881e77315bcSPatrick Sanan /* right-hand bottom up corner */ 882e77315bcSPatrick Sanan } else { 883e77315bcSPatrick Sanan 884e77315bcSPatrick Sanan td = x[k-1][j][i]; 885e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 886e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 887e77315bcSPatrick Sanan /* dd = bd * ad; */ 888e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 889e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 890e77315bcSPatrick Sanan 891e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 892e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 893e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 894e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 895e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 896e77315bcSPatrick Sanan v[2] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 897e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 898e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 8999566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 900e77315bcSPatrick Sanan } 901e77315bcSPatrick Sanan 902e77315bcSPatrick Sanan /* right-hand top edge */ 903e77315bcSPatrick Sanan } else if (j == my-1) { 904e77315bcSPatrick Sanan 905e77315bcSPatrick Sanan ts = x[k][j-1][i]; 906e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 907e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 908e77315bcSPatrick Sanan /* ds = bs * as; */ 909e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 910e77315bcSPatrick Sanan gs = coef*bs*(ts - t0); 911e77315bcSPatrick Sanan 912e77315bcSPatrick Sanan /* right-hand top down corner */ 913e77315bcSPatrick Sanan if (k == 0) { 914e77315bcSPatrick Sanan 915e77315bcSPatrick Sanan tu = x[k+1][j][i]; 916e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 917e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 918e77315bcSPatrick Sanan /* du = bu * au; */ 919e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 920e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 921e77315bcSPatrick Sanan 922e77315bcSPatrick Sanan c[0].k = k; c[0].j = j-1; c[0].i = i; 923e77315bcSPatrick Sanan v[0] = -hzhxdhy*(ds - gs); 924e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 925e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 926e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 927e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 928e77315bcSPatrick Sanan c[3].k = k+1; c[3].j = j; c[3].i = i; 929e77315bcSPatrick Sanan v[3] = -hxhydhz*(du + gu); 9309566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 931e77315bcSPatrick Sanan 932e77315bcSPatrick Sanan /* right-hand top interior edge */ 933e77315bcSPatrick Sanan } else if (k < mz-1) { 934e77315bcSPatrick Sanan 935e77315bcSPatrick Sanan tu = x[k+1][j][i]; 936e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 937e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 938e77315bcSPatrick Sanan /* du = bu * au; */ 939e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 940e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 941e77315bcSPatrick Sanan 942e77315bcSPatrick Sanan td = x[k-1][j][i]; 943e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 944e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 945e77315bcSPatrick Sanan /* dd = bd * ad; */ 946e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 947e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 948e77315bcSPatrick Sanan 949e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 950e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 951e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 952e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 953e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 954e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 955e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 956e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 957e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 958e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 9599566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 960e77315bcSPatrick Sanan 961e77315bcSPatrick Sanan /* right-hand top up corner */ 962e77315bcSPatrick Sanan } else { 963e77315bcSPatrick Sanan 964e77315bcSPatrick Sanan td = x[k-1][j][i]; 965e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 966e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 967e77315bcSPatrick Sanan /* dd = bd * ad; */ 968e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 969e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 970e77315bcSPatrick Sanan 971e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 972e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 973e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 974e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 975e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 976e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 977e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 978e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 9799566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES)); 980e77315bcSPatrick Sanan } 981e77315bcSPatrick Sanan 982e77315bcSPatrick Sanan } else { 983e77315bcSPatrick Sanan 984e77315bcSPatrick Sanan ts = x[k][j-1][i]; 985e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 986e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 987e77315bcSPatrick Sanan /* ds = bs * as; */ 988e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 989e77315bcSPatrick Sanan gs = coef*bs*(t0 - ts); 990e77315bcSPatrick Sanan 991e77315bcSPatrick Sanan tn = x[k][j+1][i]; 992e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 993e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 994e77315bcSPatrick Sanan /* dn = bn * an; */ 995e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 996e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 997e77315bcSPatrick Sanan 998e77315bcSPatrick Sanan /* right-hand down interior edge */ 999e77315bcSPatrick Sanan if (k == 0) { 1000e77315bcSPatrick Sanan 1001e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1002e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1003e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1004e77315bcSPatrick Sanan /* du = bu * au; */ 1005e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1006e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1007e77315bcSPatrick Sanan 1008e77315bcSPatrick Sanan c[0].k = k; c[0].j = j-1; c[0].i = i; 1009e77315bcSPatrick Sanan v[0] = -hzhxdhy*(ds - gs); 1010e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 1011e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 1012e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 1013e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 1014e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 1015e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 1016e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 1017e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 10189566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 1019e77315bcSPatrick Sanan 1020e77315bcSPatrick Sanan } else if (k == mz-1) { /* right-hand up interior edge */ 1021e77315bcSPatrick Sanan 1022e77315bcSPatrick Sanan td = x[k-1][j][i]; 1023e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1024e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1025e77315bcSPatrick Sanan /* dd = bd * ad; */ 1026e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1027e77315bcSPatrick Sanan gd = coef*bd*(t0 - td); 1028e77315bcSPatrick Sanan 1029e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1030e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1031e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 1032e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 1033e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 1034e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 1035e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 1036e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 1037e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 1038e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 10399566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 1040e77315bcSPatrick Sanan 1041e77315bcSPatrick Sanan } else { /* right-hand interior plane */ 1042e77315bcSPatrick Sanan 1043e77315bcSPatrick Sanan td = x[k-1][j][i]; 1044e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1045e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1046e77315bcSPatrick Sanan /* dd = bd * ad; */ 1047e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1048e77315bcSPatrick Sanan gd = coef*bd*(t0 - td); 1049e77315bcSPatrick Sanan 1050e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1051e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1052e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1053e77315bcSPatrick Sanan /* du = bu * au; */ 1054e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1055e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1056e77315bcSPatrick Sanan 1057e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1058e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1059e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 1060e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 1061e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 1062e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 1063e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 1064e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 1065e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 1066e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 1067e77315bcSPatrick Sanan c[5].k = k+1; c[5].j = j; c[5].i = i; 1068e77315bcSPatrick Sanan v[5] = -hxhydhz*(du + gu); 10699566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES)); 1070e77315bcSPatrick Sanan } 1071e77315bcSPatrick Sanan } 1072e77315bcSPatrick Sanan 1073e77315bcSPatrick Sanan } else if (j == 0) { 1074e77315bcSPatrick Sanan 1075e77315bcSPatrick Sanan tw = x[k][j][i-1]; 1076e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 1077e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 1078e77315bcSPatrick Sanan /* dw = bw * aw */ 1079e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 1080e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 1081e77315bcSPatrick Sanan 1082e77315bcSPatrick Sanan te = x[k][j][i+1]; 1083e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 1084e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 1085e77315bcSPatrick Sanan /* de = be * ae; */ 1086e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 1087e77315bcSPatrick Sanan ge = coef*be*(te - t0); 1088e77315bcSPatrick Sanan 1089e77315bcSPatrick Sanan tn = x[k][j+1][i]; 1090e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 1091e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 1092e77315bcSPatrick Sanan /* dn = bn * an; */ 1093e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 1094e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 1095e77315bcSPatrick Sanan 1096e77315bcSPatrick Sanan /* bottom down interior edge */ 1097e77315bcSPatrick Sanan if (k == 0) { 1098e77315bcSPatrick Sanan 1099e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1100e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1101e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1102e77315bcSPatrick Sanan /* du = bu * au; */ 1103e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1104e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1105e77315bcSPatrick Sanan 1106e77315bcSPatrick Sanan c[0].k = k; c[0].j = j; c[0].i = i-1; 1107e77315bcSPatrick Sanan v[0] = -hyhzdhx*(dw - gw); 1108e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i; 1109e77315bcSPatrick Sanan v[1] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 1110e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i+1; 1111e77315bcSPatrick Sanan v[2] = -hyhzdhx*(de + ge); 1112e77315bcSPatrick Sanan c[3].k = k; c[3].j = j+1; c[3].i = i; 1113e77315bcSPatrick Sanan v[3] = -hzhxdhy*(dn + gn); 1114e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 1115e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 11169566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 1117e77315bcSPatrick Sanan 1118e77315bcSPatrick Sanan } else if (k == mz-1) { /* bottom up interior edge */ 1119e77315bcSPatrick Sanan 1120e77315bcSPatrick Sanan td = x[k-1][j][i]; 1121e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1122e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1123e77315bcSPatrick Sanan /* dd = bd * ad; */ 1124e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1125e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 1126e77315bcSPatrick Sanan 1127e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1128e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1129e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 1130e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 1131e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 1132e77315bcSPatrick Sanan v[2] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 1133e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 1134e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 1135e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 1136e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 11379566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 1138e77315bcSPatrick Sanan 1139e77315bcSPatrick Sanan } else { /* bottom interior plane */ 1140e77315bcSPatrick Sanan 1141e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1142e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1143e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1144e77315bcSPatrick Sanan /* du = bu * au; */ 1145e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1146e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1147e77315bcSPatrick Sanan 1148e77315bcSPatrick Sanan td = x[k-1][j][i]; 1149e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1150e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1151e77315bcSPatrick Sanan /* dd = bd * ad; */ 1152e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1153e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 1154e77315bcSPatrick Sanan 1155e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1156e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1157e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 1158e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 1159e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 1160e77315bcSPatrick Sanan v[2] = hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 1161e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 1162e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 1163e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 1164e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 1165e77315bcSPatrick Sanan c[5].k = k+1; c[5].j = j; c[5].i = i; 1166e77315bcSPatrick Sanan v[5] = -hxhydhz*(du + gu); 11679566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES)); 1168e77315bcSPatrick Sanan } 1169e77315bcSPatrick Sanan 1170e77315bcSPatrick Sanan } else if (j == my-1) { 1171e77315bcSPatrick Sanan 1172e77315bcSPatrick Sanan tw = x[k][j][i-1]; 1173e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 1174e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 1175e77315bcSPatrick Sanan /* dw = bw * aw */ 1176e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 1177e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 1178e77315bcSPatrick Sanan 1179e77315bcSPatrick Sanan te = x[k][j][i+1]; 1180e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 1181e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 1182e77315bcSPatrick Sanan /* de = be * ae; */ 1183e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 1184e77315bcSPatrick Sanan ge = coef*be*(te - t0); 1185e77315bcSPatrick Sanan 1186e77315bcSPatrick Sanan ts = x[k][j-1][i]; 1187e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 1188e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 1189e77315bcSPatrick Sanan /* ds = bs * as; */ 1190e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 1191e77315bcSPatrick Sanan gs = coef*bs*(t0 - ts); 1192e77315bcSPatrick Sanan 1193e77315bcSPatrick Sanan /* top down interior edge */ 1194e77315bcSPatrick Sanan if (k == 0) { 1195e77315bcSPatrick Sanan 1196e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1197e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1198e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1199e77315bcSPatrick Sanan /* du = bu * au; */ 1200e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1201e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1202e77315bcSPatrick Sanan 1203e77315bcSPatrick Sanan c[0].k = k; c[0].j = j-1; c[0].i = i; 1204e77315bcSPatrick Sanan v[0] = -hzhxdhy*(ds - gs); 1205e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 1206e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 1207e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 1208e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 1209e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 1210e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 1211e77315bcSPatrick Sanan c[4].k = k+1; c[4].j = j; c[4].i = i; 1212e77315bcSPatrick Sanan v[4] = -hxhydhz*(du + gu); 12139566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 1214e77315bcSPatrick Sanan 1215e77315bcSPatrick Sanan } else if (k == mz-1) { /* top up interior edge */ 1216e77315bcSPatrick Sanan 1217e77315bcSPatrick Sanan td = x[k-1][j][i]; 1218e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1219e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1220e77315bcSPatrick Sanan /* dd = bd * ad; */ 1221e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1222e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 1223e77315bcSPatrick Sanan 1224e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1225e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1226e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 1227e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 1228e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 1229e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 1230e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 1231e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 1232e77315bcSPatrick Sanan c[4].k = k; c[4].j = j; c[4].i = i+1; 1233e77315bcSPatrick Sanan v[4] = -hyhzdhx*(de + ge); 12349566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES)); 1235e77315bcSPatrick Sanan 1236e77315bcSPatrick Sanan } else { /* top interior plane */ 1237e77315bcSPatrick Sanan 1238e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1239e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1240e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1241e77315bcSPatrick Sanan /* du = bu * au; */ 1242e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1243e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1244e77315bcSPatrick Sanan 1245e77315bcSPatrick Sanan td = x[k-1][j][i]; 1246e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1247e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1248e77315bcSPatrick Sanan /* dd = bd * ad; */ 1249e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1250e77315bcSPatrick Sanan gd = coef*bd*(td - t0); 1251e77315bcSPatrick Sanan 1252e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1253e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1254e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 1255e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 1256e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 1257e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 1258e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 1259e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu); 1260e77315bcSPatrick Sanan c[4].k = k; c[4].j = j; c[4].i = i+1; 1261e77315bcSPatrick Sanan v[4] = -hyhzdhx*(de + ge); 1262e77315bcSPatrick Sanan c[5].k = k+1; c[5].j = j; c[5].i = i; 1263e77315bcSPatrick Sanan v[5] = -hxhydhz*(du + gu); 12649566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES)); 1265e77315bcSPatrick Sanan } 1266e77315bcSPatrick Sanan 1267e77315bcSPatrick Sanan } else if (k == 0) { 1268e77315bcSPatrick Sanan 1269e77315bcSPatrick Sanan /* down interior plane */ 1270e77315bcSPatrick Sanan 1271e77315bcSPatrick Sanan tw = x[k][j][i-1]; 1272e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 1273e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 1274e77315bcSPatrick Sanan /* dw = bw * aw */ 1275e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 1276e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 1277e77315bcSPatrick Sanan 1278e77315bcSPatrick Sanan te = x[k][j][i+1]; 1279e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 1280e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 1281e77315bcSPatrick Sanan /* de = be * ae; */ 1282e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 1283e77315bcSPatrick Sanan ge = coef*be*(te - t0); 1284e77315bcSPatrick Sanan 1285e77315bcSPatrick Sanan ts = x[k][j-1][i]; 1286e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 1287e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 1288e77315bcSPatrick Sanan /* ds = bs * as; */ 1289e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 1290e77315bcSPatrick Sanan gs = coef*bs*(t0 - ts); 1291e77315bcSPatrick Sanan 1292e77315bcSPatrick Sanan tn = x[k][j+1][i]; 1293e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 1294e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 1295e77315bcSPatrick Sanan /* dn = bn * an; */ 1296e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 1297e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 1298e77315bcSPatrick Sanan 1299e77315bcSPatrick Sanan tu = x[k+1][j][i]; 1300e77315bcSPatrick Sanan au = 0.5*(t0 + tu); 1301e77315bcSPatrick Sanan bu = PetscPowScalar(au,bm1); 1302e77315bcSPatrick Sanan /* du = bu * au; */ 1303e77315bcSPatrick Sanan du = PetscPowScalar(au,beta); 1304e77315bcSPatrick Sanan gu = coef*bu*(tu - t0); 1305e77315bcSPatrick Sanan 1306e77315bcSPatrick Sanan c[0].k = k; c[0].j = j-1; c[0].i = i; 1307e77315bcSPatrick Sanan v[0] = -hzhxdhy*(ds - gs); 1308e77315bcSPatrick Sanan c[1].k = k; c[1].j = j; c[1].i = i-1; 1309e77315bcSPatrick Sanan v[1] = -hyhzdhx*(dw - gw); 1310e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i; 1311e77315bcSPatrick Sanan v[2] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu); 1312e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i+1; 1313e77315bcSPatrick Sanan v[3] = -hyhzdhx*(de + ge); 1314e77315bcSPatrick Sanan c[4].k = k; c[4].j = j+1; c[4].i = i; 1315e77315bcSPatrick Sanan v[4] = -hzhxdhy*(dn + gn); 1316e77315bcSPatrick Sanan c[5].k = k+1; c[5].j = j; c[5].i = i; 1317e77315bcSPatrick Sanan v[5] = -hxhydhz*(du + gu); 13189566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES)); 1319e77315bcSPatrick Sanan 1320e77315bcSPatrick Sanan } else if (k == mz-1) { 1321e77315bcSPatrick Sanan 1322e77315bcSPatrick Sanan /* up interior plane */ 1323e77315bcSPatrick Sanan 1324e77315bcSPatrick Sanan tw = x[k][j][i-1]; 1325e77315bcSPatrick Sanan aw = 0.5*(t0 + tw); 1326e77315bcSPatrick Sanan bw = PetscPowScalar(aw,bm1); 1327e77315bcSPatrick Sanan /* dw = bw * aw */ 1328e77315bcSPatrick Sanan dw = PetscPowScalar(aw,beta); 1329e77315bcSPatrick Sanan gw = coef*bw*(t0 - tw); 1330e77315bcSPatrick Sanan 1331e77315bcSPatrick Sanan te = x[k][j][i+1]; 1332e77315bcSPatrick Sanan ae = 0.5*(t0 + te); 1333e77315bcSPatrick Sanan be = PetscPowScalar(ae,bm1); 1334e77315bcSPatrick Sanan /* de = be * ae; */ 1335e77315bcSPatrick Sanan de = PetscPowScalar(ae,beta); 1336e77315bcSPatrick Sanan ge = coef*be*(te - t0); 1337e77315bcSPatrick Sanan 1338e77315bcSPatrick Sanan ts = x[k][j-1][i]; 1339e77315bcSPatrick Sanan as = 0.5*(t0 + ts); 1340e77315bcSPatrick Sanan bs = PetscPowScalar(as,bm1); 1341e77315bcSPatrick Sanan /* ds = bs * as; */ 1342e77315bcSPatrick Sanan ds = PetscPowScalar(as,beta); 1343e77315bcSPatrick Sanan gs = coef*bs*(t0 - ts); 1344e77315bcSPatrick Sanan 1345e77315bcSPatrick Sanan tn = x[k][j+1][i]; 1346e77315bcSPatrick Sanan an = 0.5*(t0 + tn); 1347e77315bcSPatrick Sanan bn = PetscPowScalar(an,bm1); 1348e77315bcSPatrick Sanan /* dn = bn * an; */ 1349e77315bcSPatrick Sanan dn = PetscPowScalar(an,beta); 1350e77315bcSPatrick Sanan gn = coef*bn*(tn - t0); 1351e77315bcSPatrick Sanan 1352e77315bcSPatrick Sanan td = x[k-1][j][i]; 1353e77315bcSPatrick Sanan ad = 0.5*(t0 + td); 1354e77315bcSPatrick Sanan bd = PetscPowScalar(ad,bm1); 1355e77315bcSPatrick Sanan /* dd = bd * ad; */ 1356e77315bcSPatrick Sanan dd = PetscPowScalar(ad,beta); 1357e77315bcSPatrick Sanan gd = coef*bd*(t0 - td); 1358e77315bcSPatrick Sanan 1359e77315bcSPatrick Sanan c[0].k = k-1; c[0].j = j; c[0].i = i; 1360e77315bcSPatrick Sanan v[0] = -hxhydhz*(dd - gd); 1361e77315bcSPatrick Sanan c[1].k = k; c[1].j = j-1; c[1].i = i; 1362e77315bcSPatrick Sanan v[1] = -hzhxdhy*(ds - gs); 1363e77315bcSPatrick Sanan c[2].k = k; c[2].j = j; c[2].i = i-1; 1364e77315bcSPatrick Sanan v[2] = -hyhzdhx*(dw - gw); 1365e77315bcSPatrick Sanan c[3].k = k; c[3].j = j; c[3].i = i; 1366e77315bcSPatrick Sanan v[3] = hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd); 1367e77315bcSPatrick Sanan c[4].k = k; c[4].j = j; c[4].i = i+1; 1368e77315bcSPatrick Sanan v[4] = -hyhzdhx*(de + ge); 1369e77315bcSPatrick Sanan c[5].k = k; c[5].j = j+1; c[5].i = i; 1370e77315bcSPatrick Sanan v[5] = -hzhxdhy*(dn + gn); 13719566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES)); 1372e77315bcSPatrick Sanan } 1373e77315bcSPatrick Sanan } 1374e77315bcSPatrick Sanan } 1375e77315bcSPatrick Sanan } 13769566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 13779566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 1378e77315bcSPatrick Sanan if (jac != J) { 13799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 13809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 1381e77315bcSPatrick Sanan } 13829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,localX,&x)); 13839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localX)); 1384e77315bcSPatrick Sanan 13859566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((41.0 + 8.0*POWFLOP)*xm*ym)); 1386e77315bcSPatrick Sanan PetscFunctionReturn(0); 1387e77315bcSPatrick Sanan } 1388e77315bcSPatrick Sanan 1389e77315bcSPatrick Sanan /*TEST 1390e77315bcSPatrick Sanan 1391e77315bcSPatrick Sanan test: 1392e77315bcSPatrick Sanan nsize: 4 1393e77315bcSPatrick Sanan args: -snes_monitor_short -pc_mg_type full -ksp_type fgmres -pc_type mg -snes_view -pc_mg_levels 2 -pc_mg_galerkin pmat 1394e77315bcSPatrick Sanan requires: !single 1395e77315bcSPatrick Sanan 1396e77315bcSPatrick Sanan TEST*/ 1397