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