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 44d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 45d71ae5a4SJacob Faibussowitsch { 46c4762a1bSJed Brown TS ts; /* nonlinear solver */ 47c4762a1bSJed Brown Vec u, r; /* solution, residual vectors */ 48c4762a1bSJed Brown Mat J, Jmf = NULL; /* Jacobian matrices */ 49c4762a1bSJed Brown DM da; 50c4762a1bSJed Brown PetscReal dt; 51c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 52c4762a1bSJed Brown SNES snes; 53c4762a1bSJed Brown PetscInt Jtype; /* Jacobian type 54c4762a1bSJed Brown 0: user provide Jacobian; 55c4762a1bSJed Brown 1: slow finite difference; 56c4762a1bSJed Brown 2: fd with coloring; */ 57c4762a1bSJed Brown 58327415f7SBarry Smith PetscFunctionBeginUser; 599566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 60c4762a1bSJed Brown /* Initialize user application context */ 61c4762a1bSJed Brown user.da = NULL; 62c4762a1bSJed Brown user.nstencilpts = 5; 63c4762a1bSJed Brown user.c = -30.0; 64c4762a1bSJed Brown user.boundary = 0; /* 0: Drichlet BC; 1: Neumann BC */ 65c4762a1bSJed Brown user.viewJacobian = PETSC_FALSE; 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nstencilpts", &user.nstencilpts, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-boundary", &user.boundary, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-viewJacobian", &user.viewJacobian)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74c4762a1bSJed Brown if (user.nstencilpts == 5) { 759566063dSJacob 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)); 76c4762a1bSJed Brown } else if (user.nstencilpts == 9) { 779566063dSJacob 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)); 7863a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "nstencilpts %" PetscInt_FMT " is not supported", user.nstencilpts); 799566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 809566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 81c4762a1bSJed Brown user.da = da; 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Extract global vectors from DMDA; 85c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 869566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &u)); 879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 90c4762a1bSJed Brown Create timestepping solver context 91c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 929566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 939566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 949566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 959566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 969566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, r, FormIFunction, &user)); 979566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0)); 989566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Set initial conditions 102c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1039566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(u, &user)); 1049566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, u)); 105c4762a1bSJed Brown dt = .01; 1069566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Set Jacobian evaluation routine 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1119566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ)); 1129566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J)); 113c4762a1bSJed Brown Jtype = 0; 1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-Jtype", &Jtype, NULL)); 115c4762a1bSJed Brown if (Jtype == 0) { /* use user provided Jacobian evaluation routine */ 11663a3b9bcSJacob Faibussowitsch PetscCheck(user.nstencilpts == 5, PETSC_COMM_WORLD, PETSC_ERR_SUP, "user Jacobian routine FormIJacobian() does not support nstencilpts=%" PetscInt_FMT, user.nstencilpts); 1179566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user)); 118c4762a1bSJed Brown } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */ 1199566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1209566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(snes, &Jmf)); 121c4762a1bSJed Brown if (Jtype == 1) { /* slow finite difference J; */ 1229566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Jmf, J, SNESComputeJacobianDefault, NULL)); 123c4762a1bSJed Brown } else if (Jtype == 2) { /* Use coloring to compute finite difference J efficiently */ 1249566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Jmf, J, SNESComputeJacobianDefaultColor, 0)); 125c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Jtype is not supported"); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Sets various TS parameters from user options 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1319566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 134c4762a1bSJed Brown Solve nonlinear system 135c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1369566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Free work space. 140c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jmf)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1459566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 149b122ec5aSJacob Faibussowitsch return 0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 153c4762a1bSJed Brown /* 154c4762a1bSJed Brown FormIFunction = Udot - RHSFunction 155c4762a1bSJed Brown */ 156d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 157d71ae5a4SJacob Faibussowitsch { 158c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 159c4762a1bSJed Brown DM da = (DM)user->da; 160c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 161c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 162c4762a1bSJed Brown PetscScalar u, uxx, uyy, **uarray, **f, **udot; 163c4762a1bSJed Brown Vec localU; 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscFunctionBeginUser; 1669566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 1679566063dSJacob 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)); 168c4762a1bSJed Brown 1699371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1); 1709371c9d4SSatish Balay sx = 1.0 / (hx * hx); 1719371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1); 1729371c9d4SSatish Balay sy = 1.0 / (hy * hy); 1733c633725SBarry Smith PetscCheck(user->nstencilpts != 9 || hx == hy, PETSC_COMM_WORLD, PETSC_ERR_SUP, "hx must equal hy when nstencilpts = 9 for this example"); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 177c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 178c4762a1bSJed Brown By placing code between these two statements, computations can be 179c4762a1bSJed Brown done while messages are in transition. 180c4762a1bSJed Brown */ 1819566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 1829566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* Get pointers to vector data */ 1859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &uarray)); 1869566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 1879566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Udot, &udot)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Get local grid boundaries */ 1909566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 193c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 194c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 195c4762a1bSJed Brown /* Boundary conditions */ 196c4762a1bSJed Brown if (i == 0 || j == 0 || i == Mx - 1 || j == My - 1) { 197c4762a1bSJed Brown if (user->boundary == 0) { /* Drichlet BC */ 198c4762a1bSJed Brown f[j][i] = uarray[j][i]; /* F = U */ 199c4762a1bSJed Brown } else { /* Neumann BC */ 200c4762a1bSJed Brown if (i == 0 && j == 0) { /* SW corner */ 201c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i + 1]; 202c4762a1bSJed Brown } else if (i == Mx - 1 && j == 0) { /* SE corner */ 203c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i - 1]; 204c4762a1bSJed Brown } else if (i == 0 && j == My - 1) { /* NW corner */ 205c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i + 1]; 206c4762a1bSJed Brown } else if (i == Mx - 1 && j == My - 1) { /* NE corner */ 207c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i - 1]; 208c4762a1bSJed Brown } else if (i == 0) { /* Left */ 209c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i + 1]; 210c4762a1bSJed Brown } else if (i == Mx - 1) { /* Right */ 211c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j][i - 1]; 212c4762a1bSJed Brown } else if (j == 0) { /* Bottom */ 213c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j + 1][i]; 214c4762a1bSJed Brown } else if (j == My - 1) { /* Top */ 215c4762a1bSJed Brown f[j][i] = uarray[j][i] - uarray[j - 1][i]; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown } 218c4762a1bSJed Brown } else { /* Interior */ 219c4762a1bSJed Brown u = uarray[j][i]; 220c4762a1bSJed Brown /* 5-point stencil */ 221c4762a1bSJed Brown uxx = (-2.0 * u + uarray[j][i - 1] + uarray[j][i + 1]); 222c4762a1bSJed Brown uyy = (-2.0 * u + uarray[j - 1][i] + uarray[j + 1][i]); 223c4762a1bSJed Brown if (user->nstencilpts == 9) { 224c4762a1bSJed Brown /* 9-point stencil: assume hx=hy */ 225c4762a1bSJed 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; 226c4762a1bSJed 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; 227c4762a1bSJed Brown } 228c4762a1bSJed Brown f[j][i] = udot[j][i] - (uxx * sx + uyy * sy); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Restore vectors */ 2349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &uarray)); 2359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2369566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Udot, &udot)); 2379566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 2389566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(11.0 * ym * xm)); 239*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 243c4762a1bSJed Brown /* 244c4762a1bSJed Brown FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot 245c4762a1bSJed Brown This routine is not used with option '-use_coloring' 246c4762a1bSJed Brown */ 247d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat J, Mat Jpre, void *ctx) 248d71ae5a4SJacob Faibussowitsch { 249c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym, nc; 250c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 251c4762a1bSJed Brown DM da = (DM)user->da; 252c4762a1bSJed Brown MatStencil col[5], row; 253c4762a1bSJed Brown PetscScalar vals[5], hx, hy, sx, sy; 254c4762a1bSJed Brown 255c4762a1bSJed Brown PetscFunctionBeginUser; 2569566063dSJacob 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)); 2579566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 258c4762a1bSJed Brown 2599371c9d4SSatish Balay hx = 1.0 / (PetscReal)(Mx - 1); 2609371c9d4SSatish Balay sx = 1.0 / (hx * hx); 2619371c9d4SSatish Balay hy = 1.0 / (PetscReal)(My - 1); 2629371c9d4SSatish Balay sy = 1.0 / (hy * hy); 263c4762a1bSJed Brown 264c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 265c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 266c4762a1bSJed Brown nc = 0; 2679371c9d4SSatish Balay row.j = j; 2689371c9d4SSatish Balay row.i = i; 269c4762a1bSJed Brown if (user->boundary == 0 && (i == 0 || i == Mx - 1 || j == 0 || j == My - 1)) { 2709371c9d4SSatish Balay col[nc].j = j; 2719371c9d4SSatish Balay col[nc].i = i; 2729371c9d4SSatish Balay vals[nc++] = 1.0; 273c4762a1bSJed Brown 274c4762a1bSJed Brown } else if (user->boundary > 0 && i == 0) { /* Left Neumann */ 2759371c9d4SSatish Balay col[nc].j = j; 2769371c9d4SSatish Balay col[nc].i = i; 2779371c9d4SSatish Balay vals[nc++] = 1.0; 2789371c9d4SSatish Balay col[nc].j = j; 2799371c9d4SSatish Balay col[nc].i = i + 1; 2809371c9d4SSatish Balay vals[nc++] = -1.0; 281c4762a1bSJed Brown } else if (user->boundary > 0 && i == Mx - 1) { /* Right Neumann */ 2829371c9d4SSatish Balay col[nc].j = j; 2839371c9d4SSatish Balay col[nc].i = i; 2849371c9d4SSatish Balay vals[nc++] = 1.0; 2859371c9d4SSatish Balay col[nc].j = j; 2869371c9d4SSatish Balay col[nc].i = i - 1; 2879371c9d4SSatish Balay vals[nc++] = -1.0; 288c4762a1bSJed Brown } else if (user->boundary > 0 && j == 0) { /* Bottom Neumann */ 2899371c9d4SSatish Balay col[nc].j = j; 2909371c9d4SSatish Balay col[nc].i = i; 2919371c9d4SSatish Balay vals[nc++] = 1.0; 2929371c9d4SSatish Balay col[nc].j = j + 1; 2939371c9d4SSatish Balay col[nc].i = i; 2949371c9d4SSatish Balay vals[nc++] = -1.0; 295c4762a1bSJed Brown } else if (user->boundary > 0 && j == My - 1) { /* Top Neumann */ 2969371c9d4SSatish Balay col[nc].j = j; 2979371c9d4SSatish Balay col[nc].i = i; 2989371c9d4SSatish Balay vals[nc++] = 1.0; 2999371c9d4SSatish Balay col[nc].j = j - 1; 3009371c9d4SSatish Balay col[nc].i = i; 3019371c9d4SSatish Balay vals[nc++] = -1.0; 302c4762a1bSJed Brown } else { /* Interior */ 3039371c9d4SSatish Balay col[nc].j = j - 1; 3049371c9d4SSatish Balay col[nc].i = i; 3059371c9d4SSatish Balay vals[nc++] = -sy; 3069371c9d4SSatish Balay col[nc].j = j; 3079371c9d4SSatish Balay col[nc].i = i - 1; 3089371c9d4SSatish Balay vals[nc++] = -sx; 3099371c9d4SSatish Balay col[nc].j = j; 3109371c9d4SSatish Balay col[nc].i = i; 3119371c9d4SSatish Balay vals[nc++] = 2.0 * (sx + sy) + a; 3129371c9d4SSatish Balay col[nc].j = j; 3139371c9d4SSatish Balay col[nc].i = i + 1; 3149371c9d4SSatish Balay vals[nc++] = -sx; 3159371c9d4SSatish Balay col[nc].j = j + 1; 3169371c9d4SSatish Balay col[nc].i = i; 3179371c9d4SSatish Balay vals[nc++] = -sy; 318c4762a1bSJed Brown } 3199566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, vals, INSERT_VALUES)); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown } 3229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 3239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 324c4762a1bSJed Brown if (J != Jpre) { 3259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown 329c4762a1bSJed Brown if (user->viewJacobian) { 3309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)Jpre), "Jpre:\n")); 3319566063dSJacob Faibussowitsch PetscCall(MatView(Jpre, PETSC_VIEWER_STDOUT_WORLD)); 332c4762a1bSJed Brown } 333*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 337d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(Vec U, void *ptr) 338d71ae5a4SJacob Faibussowitsch { 339c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 340c4762a1bSJed Brown DM da = user->da; 341c4762a1bSJed Brown PetscReal c = user->c; 342c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 343c4762a1bSJed Brown PetscScalar **u; 344c4762a1bSJed Brown PetscReal hx, hy, x, y, r; 345c4762a1bSJed Brown 346c4762a1bSJed Brown PetscFunctionBeginUser; 3479566063dSJacob 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)); 348c4762a1bSJed Brown 349c4762a1bSJed Brown hx = 1.0 / (PetscReal)(Mx - 1); 350c4762a1bSJed Brown hy = 1.0 / (PetscReal)(My - 1); 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* Get pointers to vector data */ 3539566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* Get local grid boundaries */ 3569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 359c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 360c4762a1bSJed Brown y = j * hy; 361c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 362c4762a1bSJed Brown x = i * hx; 363c4762a1bSJed Brown r = PetscSqrtReal((x - .5) * (x - .5) + (y - .5) * (y - .5)); 364c4762a1bSJed Brown if (r < .125) u[j][i] = PetscExpReal(c * r * r * r); 365c4762a1bSJed Brown else u[j][i] = 0.0; 366c4762a1bSJed Brown } 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* Restore vectors */ 3709566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 371*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 374c4762a1bSJed Brown /*TEST 375c4762a1bSJed Brown 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor 378c4762a1bSJed Brown 379c4762a1bSJed Brown test: 380c4762a1bSJed Brown suffix: 2 381c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor 382c4762a1bSJed Brown 383c4762a1bSJed Brown test: 384c4762a1bSJed Brown suffix: 3 385c4762a1bSJed Brown requires: !single 386c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 387c4762a1bSJed Brown 388c4762a1bSJed Brown test: 389c4762a1bSJed Brown suffix: 4 390c4762a1bSJed Brown requires: !single 391c4762a1bSJed Brown nsize: 2 392c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor 393c4762a1bSJed Brown 394c4762a1bSJed Brown test: 395c4762a1bSJed Brown suffix: 5 396c4762a1bSJed Brown nsize: 1 397c4762a1bSJed Brown args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor 398c4762a1bSJed Brown 399c4762a1bSJed Brown TEST*/ 400