1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d. Modified from ex13.c for illustrating how to solve DAEs. \n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown u_t = uxx + uyy 5c4762a1bSJed Brown 0 < x < 1, 0 < y < 1; 6c4762a1bSJed Brown At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125 7c4762a1bSJed Brown u(x,y) = 0.0 if r >= .125 8c4762a1bSJed Brown 9c4762a1bSJed Brown Boundary conditions: 10c4762a1bSJed Brown Drichlet BC: 11c4762a1bSJed Brown At x=0, x=1, y=0, y=1: u = 0.0 12c4762a1bSJed Brown 13c4762a1bSJed Brown Neumann BC: 14c4762a1bSJed Brown At x=0, x=1: du(x,y,t)/dx = 0 15c4762a1bSJed Brown At y=0, y=1: du(x,y,t)/dy = 0 16c4762a1bSJed Brown 17c4762a1bSJed Brown mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor 18c4762a1bSJed Brown ./ex15 -da_grid_x 40 -da_grid_y 40 -draw_pause .1 -boundary 1 -ts_monitor_draw_solution 19c4762a1bSJed Brown ./ex15 -da_grid_x 40 -da_grid_y 40 -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9 20c4762a1bSJed Brown 21c4762a1bSJed Brown */ 22c4762a1bSJed Brown 23c4762a1bSJed Brown #include <petscdm.h> 24c4762a1bSJed Brown #include <petscdmda.h> 25c4762a1bSJed Brown #include <petscts.h> 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown User-defined data structures and routines 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* AppCtx: used by FormIFunction() and FormIJacobian() */ 32c4762a1bSJed Brown typedef struct { 33c4762a1bSJed Brown DM da; 34c4762a1bSJed Brown PetscInt nstencilpts; /* number of stencil points: 5 or 9 */ 35c4762a1bSJed Brown PetscReal c; 36c4762a1bSJed Brown PetscInt boundary; /* Type of boundary condition */ 37c4762a1bSJed Brown PetscBool viewJacobian; 38c4762a1bSJed Brown } AppCtx; 39c4762a1bSJed Brown 40c4762a1bSJed Brown extern PetscErrorCode FormIFunction(TS, PetscReal, Vec, Vec, Vec, void *); 41c4762a1bSJed Brown extern PetscErrorCode FormIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 42c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(Vec, void *); 43c4762a1bSJed Brown 44*9371c9d4SSatish Balay int main(int argc, char **argv) { 45c4762a1bSJed Brown TS ts; /* nonlinear solver */ 46c4762a1bSJed Brown Vec u, r; /* solution, residual vectors */ 47c4762a1bSJed Brown Mat J, Jmf = NULL; /* Jacobian matrices */ 48c4762a1bSJed Brown DM da; 49c4762a1bSJed Brown PetscReal dt; 50c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 51c4762a1bSJed Brown SNES snes; 52c4762a1bSJed Brown PetscInt Jtype; /* Jacobian type 53c4762a1bSJed Brown 0: user provide Jacobian; 54c4762a1bSJed Brown 1: slow finite difference; 55c4762a1bSJed Brown 2: fd with coloring; */ 56c4762a1bSJed Brown 57327415f7SBarry Smith PetscFunctionBeginUser; 589566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 59c4762a1bSJed Brown /* Initialize user application context */ 60c4762a1bSJed Brown user.da = NULL; 61c4762a1bSJed Brown user.nstencilpts = 5; 62c4762a1bSJed Brown user.c = -30.0; 63c4762a1bSJed Brown user.boundary = 0; /* 0: Drichlet BC; 1: Neumann BC */ 64c4762a1bSJed Brown user.viewJacobian = PETSC_FALSE; 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nstencilpts", &user.nstencilpts, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-boundary", &user.boundary, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-viewJacobian", &user.viewJacobian)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 71c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 72c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 73c4762a1bSJed Brown if (user.nstencilpts == 5) { 749566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 11, 11, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 75c4762a1bSJed Brown } else if (user.nstencilpts == 9) { 769566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 11, 11, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 7763a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "nstencilpts %" PetscInt_FMT " is not supported", user.nstencilpts); 789566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 799566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 80c4762a1bSJed Brown user.da = da; 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 83c4762a1bSJed Brown Extract global vectors from DMDA; 84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 859566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &u)); 869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89c4762a1bSJed Brown Create timestepping solver context 90c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 919566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 929566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 939566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 949566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 959566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, r, FormIFunction, &user)); 969566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0)); 979566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Set initial conditions 101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1029566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(u, &user)); 1039566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, u)); 104c4762a1bSJed Brown dt = .01; 1059566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Set Jacobian evaluation routine 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1109566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ)); 1119566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J)); 112c4762a1bSJed Brown Jtype = 0; 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Jtype", &Jtype, NULL)); 114c4762a1bSJed Brown if (Jtype == 0) { /* use user provided Jacobian evaluation routine */ 11563a3b9bcSJacob Faibussowitsch PetscCheck(user.nstencilpts == 5, PETSC_COMM_WORLD, PETSC_ERR_SUP, "user Jacobian routine FormIJacobian() does not support nstencilpts=%" PetscInt_FMT, user.nstencilpts); 1169566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user)); 117c4762a1bSJed Brown } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */ 1189566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1199566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(snes, &Jmf)); 120c4762a1bSJed Brown if (Jtype == 1) { /* slow finite difference J; */ 1219566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Jmf, J, SNESComputeJacobianDefault, NULL)); 122c4762a1bSJed Brown } else if (Jtype == 2) { /* Use coloring to compute finite difference J efficiently */ 1239566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Jmf, J, SNESComputeJacobianDefaultColor, 0)); 124c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Jtype is not supported"); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128c4762a1bSJed Brown Sets various TS parameters from user options 129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1309566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 133c4762a1bSJed Brown Solve nonlinear system 134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1359566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Free work space. 139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jmf)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1449566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 148b122ec5aSJacob Faibussowitsch return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 152c4762a1bSJed Brown /* 153c4762a1bSJed Brown FormIFunction = Udot - RHSFunction 154c4762a1bSJed Brown */ 155*9371c9d4SSatish Balay PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) { 156c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 157c4762a1bSJed Brown DM da = (DM)user->da; 158c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 159c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 160c4762a1bSJed Brown PetscScalar u, uxx, uyy, **uarray, **f, **udot; 161c4762a1bSJed Brown Vec localU; 162c4762a1bSJed Brown 163c4762a1bSJed Brown PetscFunctionBeginUser; 1649566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 1659566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 166c4762a1bSJed Brown 167*9371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1); 168*9371c9d4SSatish Balay sx = 1.0 / (hx * hx); 169*9371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1); 170*9371c9d4SSatish Balay sy = 1.0 / (hy * hy); 1713c633725SBarry Smith PetscCheck(user->nstencilpts != 9 || hx == hy, PETSC_COMM_WORLD, PETSC_ERR_SUP, "hx must equal hy when nstencilpts = 9 for this example"); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* 174c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 175c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 176c4762a1bSJed Brown By placing code between these two statements, computations can be 177c4762a1bSJed Brown done while messages are in transition. 178c4762a1bSJed Brown */ 1799566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 1809566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* Get pointers to vector data */ 1839566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &uarray)); 1849566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 1859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Udot, &udot)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Get local grid boundaries */ 1889566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 191c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 192c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 193c4762a1bSJed Brown /* Boundary conditions */ 194c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 195c4762a1bSJed Brown if (user->boundary == 0) { /* Drichlet BC */ 196c4762a1bSJed Brown f[j][i] = uarray[j][i]; /* F = U */ 197c4762a1bSJed Brown } else { /* Neumann BC */ 198c4762a1bSJed Brown if (i == 0 && j == 0) { /* SW corner */ 199c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i + 1]; 200c4762a1bSJed Brown } else if (i == Mx - 1 && j == 0) { /* SE corner */ 201c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i - 1]; 202c4762a1bSJed Brown } else if (i == 0 && j == My - 1) { /* NW corner */ 203c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i + 1]; 204c4762a1bSJed Brown } else if (i == Mx - 1 && j == My - 1) { /* NE corner */ 205c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i - 1]; 206c4762a1bSJed Brown } else if (i == 0) { /* Left */ 207c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i + 1]; 208c4762a1bSJed Brown } else if (i == Mx - 1) { /* Right */ 209c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i - 1]; 210c4762a1bSJed Brown } else if (j == 0) { /* Bottom */ 211c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i]; 212c4762a1bSJed Brown } else if (j == My - 1) { /* Top */ 213c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i]; 214c4762a1bSJed Brown } 215c4762a1bSJed Brown } 216c4762a1bSJed Brown } else { /* Interior */ 217c4762a1bSJed Brown u = uarray[j][i]; 218c4762a1bSJed Brown /* 5-point stencil */ 219c4762a1bSJed Brown uxx = (-2.0 * u + uarray[j][i - 1] + uarray[j][i + 1]); 220c4762a1bSJed Brown uyy = (-2.0 * u + uarray[j - 1][i] + uarray[j + 1][i]); 221c4762a1bSJed Brown if (user->nstencilpts == 9) { 222c4762a1bSJed Brown /* 9-point stencil: assume hx=hy */ 223c4762a1bSJed Brown uxx = 2.0 * uxx / 3.0 + (0.5 * (uarray[j - 1][i - 1] + uarray[j - 1][i + 1] + uarray[j + 1][i - 1] + uarray[j + 1][i + 1]) - 2.0 * u) / 6.0; 224c4762a1bSJed Brown uyy = 2.0 * uyy / 3.0 + (0.5 * (uarray[j - 1][i - 1] + uarray[j - 1][i + 1] + uarray[j + 1][i - 1] + uarray[j + 1][i + 1]) - 2.0 * u) / 6.0; 225c4762a1bSJed Brown } 226c4762a1bSJed Brown f[j][i] = udot[j][i] - (uxx * sx + uyy * sy); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* Restore vectors */ 2329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &uarray)); 2339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Udot, &udot)); 2359566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 2369566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm)); 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 241c4762a1bSJed Brown /* 242c4762a1bSJed Brown FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot 243c4762a1bSJed Brown This routine is not used with option '-use_coloring' 244c4762a1bSJed Brown */ 245*9371c9d4SSatish Balay PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, void *ctx) { 246c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym, nc; 247c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 248c4762a1bSJed Brown DM da = (DM)user->da; 249c4762a1bSJed Brown MatStencil col[5], row; 250c4762a1bSJed Brown PetscScalar vals[5], hx, hy, sx, sy; 251c4762a1bSJed Brown 252c4762a1bSJed Brown PetscFunctionBeginUser; 2539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 2549566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 255c4762a1bSJed Brown 256*9371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1); 257*9371c9d4SSatish Balay sx = 1.0 / (hx * hx); 258*9371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1); 259*9371c9d4SSatish Balay sy = 1.0 / (hy * hy); 260c4762a1bSJed Brown 261c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 262c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 263c4762a1bSJed Brown nc = 0; 264*9371c9d4SSatish Balay row.j = j; 265*9371c9d4SSatish Balay row.i = i; 266c4762a1bSJed Brown if (user->boundary == 0 && (i == 0 || i == Mx - 1 || j == 0 || j == My - 1)) { 267*9371c9d4SSatish Balay col[nc].j = j; 268*9371c9d4SSatish Balay col[nc].i = i; 269*9371c9d4SSatish Balay vals[nc++] = 1.0; 270c4762a1bSJed Brown 271c4762a1bSJed Brown } else if (user->boundary > 0 && i == 0) { /* Left Neumann */ 272*9371c9d4SSatish Balay col[nc].j = j; 273*9371c9d4SSatish Balay col[nc].i = i; 274*9371c9d4SSatish Balay vals[nc++] = 1.0; 275*9371c9d4SSatish Balay col[nc].j = j; 276*9371c9d4SSatish Balay col[nc].i = i + 1; 277*9371c9d4SSatish Balay vals[nc++] = -1.0; 278c4762a1bSJed Brown } else if (user->boundary > 0 && i == Mx - 1) { /* Right Neumann */ 279*9371c9d4SSatish Balay col[nc].j = j; 280*9371c9d4SSatish Balay col[nc].i = i; 281*9371c9d4SSatish Balay vals[nc++] = 1.0; 282*9371c9d4SSatish Balay col[nc].j = j; 283*9371c9d4SSatish Balay col[nc].i = i - 1; 284*9371c9d4SSatish Balay vals[nc++] = -1.0; 285c4762a1bSJed Brown } else if (user->boundary > 0 && j == 0) { /* Bottom Neumann */ 286*9371c9d4SSatish Balay col[nc].j = j; 287*9371c9d4SSatish Balay col[nc].i = i; 288*9371c9d4SSatish Balay vals[nc++] = 1.0; 289*9371c9d4SSatish Balay col[nc].j = j + 1; 290*9371c9d4SSatish Balay col[nc].i = i; 291*9371c9d4SSatish Balay vals[nc++] = -1.0; 292c4762a1bSJed Brown } else if (user->boundary > 0 && j == My - 1) { /* Top Neumann */ 293*9371c9d4SSatish Balay col[nc].j = j; 294*9371c9d4SSatish Balay col[nc].i = i; 295*9371c9d4SSatish Balay vals[nc++] = 1.0; 296*9371c9d4SSatish Balay col[nc].j = j - 1; 297*9371c9d4SSatish Balay col[nc].i = i; 298*9371c9d4SSatish Balay vals[nc++] = -1.0; 299c4762a1bSJed Brown } else { /* Interior */ 300*9371c9d4SSatish Balay col[nc].j = j - 1; 301*9371c9d4SSatish Balay col[nc].i = i; 302*9371c9d4SSatish Balay vals[nc++] = -sy; 303*9371c9d4SSatish Balay col[nc].j = j; 304*9371c9d4SSatish Balay col[nc].i = i - 1; 305*9371c9d4SSatish Balay vals[nc++] = -sx; 306*9371c9d4SSatish Balay col[nc].j = j; 307*9371c9d4SSatish Balay col[nc].i = i; 308*9371c9d4SSatish Balay vals[nc++] = 2.0 * (sx + sy) + a; 309*9371c9d4SSatish Balay col[nc].j = j; 310*9371c9d4SSatish Balay col[nc].i = i + 1; 311*9371c9d4SSatish Balay vals[nc++] = -sx; 312*9371c9d4SSatish Balay col[nc].j = j + 1; 313*9371c9d4SSatish Balay col[nc].i = i; 314*9371c9d4SSatish Balay vals[nc++] = -sy; 315c4762a1bSJed Brown } 3169566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, vals, INSERT_VALUES)); 317c4762a1bSJed Brown } 318c4762a1bSJed Brown } 3199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 3209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 321c4762a1bSJed Brown if (J != Jpre) { 3229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 326c4762a1bSJed Brown if (user->viewJacobian) { 3279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)Jpre), "Jpre:\n")); 3289566063dSJacob Faibussowitsch PetscCall(MatView(Jpre, PETSC_VIEWER_STDOUT_WORLD)); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown PetscFunctionReturn(0); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 334*9371c9d4SSatish Balay PetscErrorCode FormInitialSolution(Vec U, void *ptr) { 335c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 336c4762a1bSJed Brown DM da = user->da; 337c4762a1bSJed Brown PetscReal c = user->c; 338c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 339c4762a1bSJed Brown PetscScalar **u; 340c4762a1bSJed Brown PetscReal hx, hy, x, y, r; 341c4762a1bSJed Brown 342c4762a1bSJed Brown PetscFunctionBeginUser; 3439566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 344c4762a1bSJed Brown 345c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 346c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* Get pointers to vector data */ 3499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* Get local grid boundaries */ 3529566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 353c4762a1bSJed Brown 354c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 355c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 356c4762a1bSJed Brown y = j * hy; 357c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 358c4762a1bSJed Brown x = i * hx; 359c4762a1bSJed Brown r = PetscSqrtReal((x - .5) * (x - .5) + (y - .5) * (y - .5)); 360c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(c * r * r * r); 361c4762a1bSJed Brown else u[j][i] = 0.0; 362c4762a1bSJed Brown } 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365c4762a1bSJed Brown /* Restore vectors */ 3669566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 367c4762a1bSJed Brown PetscFunctionReturn(0); 368c4762a1bSJed Brown } 369c4762a1bSJed Brown 370c4762a1bSJed Brown /*TEST 371c4762a1bSJed Brown 372c4762a1bSJed Brown test: 373c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor 374c4762a1bSJed Brown 375c4762a1bSJed Brown test: 376c4762a1bSJed Brown suffix: 2 377c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor 378c4762a1bSJed Brown 379c4762a1bSJed Brown test: 380c4762a1bSJed Brown suffix: 3 381c4762a1bSJed Brown requires: !single 382c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 383c4762a1bSJed Brown 384c4762a1bSJed Brown test: 385c4762a1bSJed Brown suffix: 4 386c4762a1bSJed Brown requires: !single 387c4762a1bSJed Brown nsize: 2 388c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 389c4762a1bSJed Brown 390c4762a1bSJed Brown test: 391c4762a1bSJed Brown suffix: 5 392c4762a1bSJed Brown nsize: 1 393c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor 394c4762a1bSJed Brown 395c4762a1bSJed Brown TEST*/ 396