1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Solves the equation biharmonic equation in split form 6c4762a1bSJed Brown 7c4762a1bSJed Brown w = -kappa \Delta u 8c4762a1bSJed Brown u_t = \Delta w 9c4762a1bSJed Brown -1 <= u <= 1 10c4762a1bSJed Brown Periodic boundary conditions 11c4762a1bSJed Brown 12c4762a1bSJed Brown Evolve the biharmonic heat equation with bounds: (same as biharmonic) 13c4762a1bSJed Brown --------------- 14c4762a1bSJed Brown ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type beuler -da_refine 5 -draw_fields 1 -ts_dt 9.53674e-9 15c4762a1bSJed Brown 16c4762a1bSJed Brown w = -kappa \Delta u + u^3 - u 17c4762a1bSJed Brown u_t = \Delta w 18c4762a1bSJed Brown -1 <= u <= 1 19c4762a1bSJed Brown Periodic boundary conditions 20c4762a1bSJed Brown 21c4762a1bSJed Brown Evolve the Cahn-Hillard equations: (this fails after a few timesteps 12/17/2017) 22c4762a1bSJed Brown --------------- 23c4762a1bSJed Brown ./biharmonic2 -ts_monitor -snes_monitor -ts_monitor_draw_solution -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type beuler -da_refine 6 -draw_fields 1 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard 24c4762a1bSJed Brown 25c4762a1bSJed Brown */ 26c4762a1bSJed Brown #include <petscdm.h> 27c4762a1bSJed Brown #include <petscdmda.h> 28c4762a1bSJed Brown #include <petscts.h> 29c4762a1bSJed Brown #include <petscdraw.h> 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* 32c4762a1bSJed Brown User-defined routines 33c4762a1bSJed Brown */ 34c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, Vec, void *), FormInitialSolution(DM, Vec, PetscReal); 359371c9d4SSatish Balay typedef struct { 369371c9d4SSatish Balay PetscBool cahnhillard; 379371c9d4SSatish Balay PetscReal kappa; 389371c9d4SSatish Balay PetscInt energy; 399371c9d4SSatish Balay PetscReal tol; 409371c9d4SSatish Balay PetscReal theta; 419371c9d4SSatish Balay PetscReal theta_c; 429371c9d4SSatish Balay } UserCtx; 43c4762a1bSJed Brown 44d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 45d71ae5a4SJacob Faibussowitsch { 46c4762a1bSJed Brown TS ts; /* nonlinear solver */ 47c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 48c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 49c4762a1bSJed Brown PetscInt steps, Mx; 50c4762a1bSJed Brown DM da; 51c4762a1bSJed Brown MatFDColoring matfdcoloring; 52c4762a1bSJed Brown ISColoring iscoloring; 53c4762a1bSJed Brown PetscReal dt; 54c4762a1bSJed Brown PetscReal vbounds[] = {-100000, 100000, -1.1, 1.1}; 55c4762a1bSJed Brown SNES snes; 56c4762a1bSJed Brown UserCtx ctx; 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 59c4762a1bSJed Brown Initialize program 60c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 61327415f7SBarry Smith PetscFunctionBeginUser; 629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 63c4762a1bSJed Brown ctx.kappa = 1.0; 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &ctx.kappa, NULL)); 65c4762a1bSJed Brown ctx.cahnhillard = PETSC_FALSE; 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-cahn-hillard", &ctx.cahnhillard, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 2, vbounds)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 600, 600)); 70c4762a1bSJed Brown ctx.energy = 1; 719566063dSJacob Faibussowitsch /*PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL));*/ 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-energy", &ctx.energy, NULL)); 73c4762a1bSJed Brown ctx.tol = 1.0e-8; 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &ctx.tol, NULL)); 75c4762a1bSJed Brown ctx.theta = .001; 76c4762a1bSJed Brown ctx.theta_c = 1.0; 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta", &ctx.theta, NULL)); 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta_c", &ctx.theta_c, NULL)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 82c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 839566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, 2, 2, NULL, &da)); 849566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 859566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 869566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "Biharmonic heat equation: w = -kappa*u_xx")); 879566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "Biharmonic heat equation: u")); 889566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 89c4762a1bSJed Brown dt = 1.0 / (10. * ctx.kappa * Mx * Mx * Mx * Mx); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 93c4762a1bSJed Brown vectors that are the same types 94c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 959566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99c4762a1bSJed Brown Create timestepping solver context 100c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1019566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1029566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 1039566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1049566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormFunction, &ctx)); 1059566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, .02)); 1069566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_INTERPOLATE)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 110c4762a1bSJed Brown 111c4762a1bSJed Brown < Set Jacobian matrix data structure and default Jacobian evaluation 112c4762a1bSJed Brown routine. User can override with: 113c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 114c4762a1bSJed Brown (unless user explicitly sets preconditioner) 115c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 116c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 117c4762a1bSJed Brown products within Newton-Krylov method 118c4762a1bSJed Brown 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1209566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1219566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(da, IS_COLORING_GLOBAL, &iscoloring)); 1229566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ)); 1239566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J)); 1249566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(J, iscoloring, &matfdcoloring)); 1259566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(matfdcoloring, (PetscErrorCode(*)(void))SNESTSFormFunction, ts)); 1269566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(matfdcoloring)); 1279566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(J, iscoloring, matfdcoloring)); 1289566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 1299566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, SNESComputeJacobianDefaultColor, matfdcoloring)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132c4762a1bSJed Brown Customize nonlinear solver 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1349566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137c4762a1bSJed Brown Set initial conditions 138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1399566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x, ctx.kappa)); 1409566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 1419566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Set runtime options 145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1469566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 149c4762a1bSJed Brown Solve nonlinear system 150c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1519566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1529566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 156c4762a1bSJed Brown are no longer needed. 157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1599566063dSJacob Faibussowitsch PetscCall(MatFDColoringDestroy(&matfdcoloring)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1629566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 166b122ec5aSJacob Faibussowitsch return 0; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 1699371c9d4SSatish Balay typedef struct { 1709371c9d4SSatish Balay PetscScalar w, u; 1719371c9d4SSatish Balay } Field; 172c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 173c4762a1bSJed Brown /* 174c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 175c4762a1bSJed Brown 176c4762a1bSJed Brown Input Parameters: 177c4762a1bSJed Brown . ts - the TS context 178c4762a1bSJed Brown . X - input vector 179c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 180c4762a1bSJed Brown 181c4762a1bSJed Brown Output Parameter: 182c4762a1bSJed Brown . F - function vector 183c4762a1bSJed Brown */ 184d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec Xdot, Vec F, void *ptr) 185d71ae5a4SJacob Faibussowitsch { 186c4762a1bSJed Brown DM da; 187c4762a1bSJed Brown PetscInt i, Mx, xs, xm; 188c4762a1bSJed Brown PetscReal hx, sx; 189c4762a1bSJed Brown Field *x, *xdot, *f; 190c4762a1bSJed Brown Vec localX, localXdot; 191c4762a1bSJed Brown UserCtx *ctx = (UserCtx *)ptr; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBegin; 1949566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1959566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 1969566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localXdot)); 1979566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 198c4762a1bSJed Brown 1999371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx; 2009371c9d4SSatish Balay sx = 1.0 / (hx * hx); 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* 203c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 204c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 205c4762a1bSJed Brown By placing code between these two statements, computations can be 206c4762a1bSJed Brown done while messages are in transition. 207c4762a1bSJed Brown */ 2089566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 2099566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 2109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, Xdot, INSERT_VALUES, localXdot)); 2119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, Xdot, INSERT_VALUES, localXdot)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* 214c4762a1bSJed Brown Get pointers to vector data 215c4762a1bSJed Brown */ 2169566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 2179566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localXdot, &xdot)); 2189566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* 221c4762a1bSJed Brown Get local grid boundaries 222c4762a1bSJed Brown */ 2239566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* 226c4762a1bSJed Brown Compute function over the locally owned part of the grid 227c4762a1bSJed Brown */ 228c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 229c4762a1bSJed Brown f[i].w = x[i].w + ctx->kappa * (x[i - 1].u + x[i + 1].u - 2.0 * x[i].u) * sx; 230c4762a1bSJed Brown if (ctx->cahnhillard) { 231c4762a1bSJed Brown switch (ctx->energy) { 232d71ae5a4SJacob Faibussowitsch case 1: /* double well */ 233d71ae5a4SJacob Faibussowitsch f[i].w += -x[i].u * x[i].u * x[i].u + x[i].u; 234d71ae5a4SJacob Faibussowitsch break; 235d71ae5a4SJacob Faibussowitsch case 2: /* double obstacle */ 236d71ae5a4SJacob Faibussowitsch f[i].w += x[i].u; 237d71ae5a4SJacob Faibussowitsch break; 238c4762a1bSJed Brown case 3: /* logarithmic */ 239c4762a1bSJed Brown if (PetscRealPart(x[i].u) < -1.0 + 2.0 * ctx->tol) f[i].w += .5 * ctx->theta * (-PetscLogReal(ctx->tol) + PetscLogScalar((1.0 - x[i].u) / 2.0)) + ctx->theta_c * x[i].u; 240c4762a1bSJed Brown else if (PetscRealPart(x[i].u) > 1.0 - 2.0 * ctx->tol) f[i].w += .5 * ctx->theta * (-PetscLogScalar((1.0 + x[i].u) / 2.0) + PetscLogReal(ctx->tol)) + ctx->theta_c * x[i].u; 241c4762a1bSJed Brown else f[i].w += .5 * ctx->theta * (-PetscLogScalar((1.0 + x[i].u) / 2.0) + PetscLogScalar((1.0 - x[i].u) / 2.0)) + ctx->theta_c * x[i].u; 242c4762a1bSJed Brown break; 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 245c4762a1bSJed Brown f[i].u = xdot[i].u - (x[i - 1].w + x[i + 1].w - 2.0 * x[i].w) * sx; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* 249c4762a1bSJed Brown Restore vectors 250c4762a1bSJed Brown */ 2519566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localXdot, &xdot)); 2529566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 2539566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2549566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 2559566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localXdot)); 256*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 260d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec X, PetscReal kappa) 261d71ae5a4SJacob Faibussowitsch { 262c4762a1bSJed Brown PetscInt i, xs, xm, Mx, xgs, xgm; 263c4762a1bSJed Brown Field *x; 264c4762a1bSJed Brown PetscReal hx, xx, r, sx; 265c4762a1bSJed Brown Vec Xg; 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscFunctionBegin; 2689566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 269c4762a1bSJed Brown 270c4762a1bSJed Brown hx = 1.0 / (PetscReal)Mx; 271c4762a1bSJed Brown sx = 1.0 / (hx * hx); 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* 274c4762a1bSJed Brown Get pointers to vector data 275c4762a1bSJed Brown */ 2769566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(da, &Xg)); 2779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xg, &x)); 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown Get local grid boundaries 281c4762a1bSJed Brown */ 2829566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 2839566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &xgs, NULL, NULL, &xgm, NULL, NULL)); 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* 286c4762a1bSJed Brown Compute u function over the locally owned part of the grid including ghost points 287c4762a1bSJed Brown */ 288c4762a1bSJed Brown for (i = xgs; i < xgs + xgm; i++) { 289c4762a1bSJed Brown xx = i * hx; 290c4762a1bSJed Brown r = PetscSqrtReal((xx - .5) * (xx - .5)); 291c4762a1bSJed Brown if (r < .125) x[i].u = 1.0; 292c4762a1bSJed Brown else x[i].u = -.50; 293c4762a1bSJed Brown /* fill in x[i].w so that valgrind doesn't detect use of uninitialized memory */ 294c4762a1bSJed Brown x[i].w = 0; 295c4762a1bSJed Brown } 296c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) x[i].w = -kappa * (x[i - 1].u + x[i + 1].u - 2.0 * x[i].u) * sx; 297c4762a1bSJed Brown 298c4762a1bSJed Brown /* 299c4762a1bSJed Brown Restore vectors 300c4762a1bSJed Brown */ 3019566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xg, &x)); 302c4762a1bSJed Brown 303c4762a1bSJed Brown /* Grab only the global part of the vector */ 3049566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0)); 3059566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, Xg, ADD_VALUES, X)); 3069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, Xg, ADD_VALUES, X)); 3079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xg)); 308*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown /*TEST 312c4762a1bSJed Brown 313c4762a1bSJed Brown build: 314c4762a1bSJed Brown requires: !complex !single 315c4762a1bSJed Brown 316c4762a1bSJed Brown test: 317c4762a1bSJed Brown args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type beuler -da_refine 5 -ts_dt 9.53674e-9 -ts_max_steps 50 318c4762a1bSJed Brown requires: x 319c4762a1bSJed Brown 320c4762a1bSJed Brown TEST*/ 321