1c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Solves the equation 5c4762a1bSJed Brown 6c4762a1bSJed Brown u_t = - kappa \Delta \Delta u 7c4762a1bSJed Brown Periodic boundary conditions 8c4762a1bSJed Brown 9c4762a1bSJed Brown Evolve the biharmonic heat equation: 10c4762a1bSJed Brown --------------- 11c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -mymonitor 12c4762a1bSJed Brown 13c4762a1bSJed Brown Evolve with the restriction that -1 <= u <= 1; i.e. as a variational inequality 14c4762a1bSJed Brown --------------- 15c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -mymonitor 16c4762a1bSJed Brown 17c4762a1bSJed Brown u_t = kappa \Delta \Delta u + 6.*u*(u_x)^2 + (3*u^2 - 12) \Delta u 18c4762a1bSJed Brown -1 <= u <= 1 19c4762a1bSJed Brown Periodic boundary conditions 20c4762a1bSJed Brown 21c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double well Initial hump shrinks then grows 22c4762a1bSJed Brown --------------- 23c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 6 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -ts_monitor_draw_solution --mymonitor 24c4762a1bSJed Brown 25c4762a1bSJed Brown Initial hump neither shrinks nor grows when degenerate (otherwise similar solution) 26c4762a1bSJed Brown 27c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 6 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -degenerate -ts_monitor_draw_solution --mymonitor 28c4762a1bSJed Brown 29c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 6 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -snes_vi_ignore_function_sign -ts_monitor_draw_solution --mymonitor 30c4762a1bSJed Brown 31c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double obstacle 32c4762a1bSJed Brown --------------- 33c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -energy 2 -snes_linesearch_monitor -ts_monitor_draw_solution --mymonitor 34c4762a1bSJed Brown 35c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic + double well (never shrinks and then grows) 36c4762a1bSJed Brown --------------- 37c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu --snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .0001 -ts_dt 5.96046e-06 -cahn-hillard -energy 3 -snes_linesearch_monitor -theta .00000001 -ts_monitor_draw_solution --ts_max_time 1. -mymonitor 38c4762a1bSJed Brown 39c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu --snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .0001 -ts_dt 5.96046e-06 -cahn-hillard -energy 3 -snes_linesearch_monitor -theta .00000001 -ts_monitor_draw_solution --ts_max_time 1. -degenerate -mymonitor 40c4762a1bSJed Brown 41c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic + double obstacle (never shrinks, never grows) 42c4762a1bSJed Brown --------------- 43c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu --snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -energy 4 -snes_linesearch_monitor -theta .00000001 -ts_monitor_draw_solution --mymonitor 44c4762a1bSJed Brown 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown #include <petscdm.h> 47c4762a1bSJed Brown #include <petscdmda.h> 48c4762a1bSJed Brown #include <petscts.h> 49c4762a1bSJed Brown #include <petscdraw.h> 50c4762a1bSJed Brown 51c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, void *), FormInitialSolution(DM, Vec), MyMonitor(TS, PetscInt, PetscReal, Vec, void *), MyDestroy(void **), FormJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 529371c9d4SSatish Balay typedef struct { 539371c9d4SSatish Balay PetscBool cahnhillard; 549371c9d4SSatish Balay PetscBool degenerate; 559371c9d4SSatish Balay PetscReal kappa; 569371c9d4SSatish Balay PetscInt energy; 579371c9d4SSatish Balay PetscReal tol; 589371c9d4SSatish Balay PetscReal theta, theta_c; 599371c9d4SSatish Balay PetscInt truncation; 609371c9d4SSatish Balay PetscBool netforce; 619371c9d4SSatish Balay PetscDrawViewPorts *ports; 629371c9d4SSatish Balay } UserCtx; 63c4762a1bSJed Brown 64d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 65d71ae5a4SJacob Faibussowitsch { 66c4762a1bSJed Brown TS ts; /* nonlinear solver */ 67c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 68c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 69c4762a1bSJed Brown PetscInt steps, Mx; 70c4762a1bSJed Brown DM da; 71c4762a1bSJed Brown PetscReal dt; 72c4762a1bSJed Brown PetscBool mymonitor; 73c4762a1bSJed Brown UserCtx ctx; 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 76c4762a1bSJed Brown Initialize program 77c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 78327415f7SBarry Smith PetscFunctionBeginUser; 79c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 80c4762a1bSJed Brown ctx.kappa = 1.0; 819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &ctx.kappa, NULL)); 82c4762a1bSJed Brown ctx.degenerate = PETSC_FALSE; 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-degenerate", &ctx.degenerate, NULL)); 84c4762a1bSJed Brown ctx.cahnhillard = PETSC_FALSE; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-cahn-hillard", &ctx.cahnhillard, NULL)); 86c4762a1bSJed Brown ctx.netforce = PETSC_FALSE; 879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-netforce", &ctx.netforce, NULL)); 88c4762a1bSJed Brown ctx.energy = 1; 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-energy", &ctx.energy, NULL)); 90c4762a1bSJed Brown ctx.tol = 1.0e-8; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &ctx.tol, NULL)); 92c4762a1bSJed Brown ctx.theta = .001; 93c4762a1bSJed Brown ctx.theta_c = 1.0; 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta", &ctx.theta, NULL)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-theta_c", &ctx.theta_c, NULL)); 96c4762a1bSJed Brown ctx.truncation = 1; 979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-truncation", &ctx.truncation, NULL)); 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-mymonitor", &mymonitor)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 102c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1039566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, 1, 2, NULL, &da)); 1049566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1059566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1069566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "Biharmonic heat equation: u")); 1079566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 108c4762a1bSJed Brown dt = 1.0 / (10. * ctx.kappa * Mx * Mx * Mx * Mx); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 111c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 112c4762a1bSJed Brown vectors that are the same types 113c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1149566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 1159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Create timestepping solver context 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1209566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1219566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 1229566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1239566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &ctx)); 1249566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ)); 1259566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J)); 1269566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, J, J, FormJacobian, &ctx)); 1279566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, .02)); 1289566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_INTERPOLATE)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 132c4762a1bSJed Brown 133c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 134c4762a1bSJed Brown routine. User can override with: 135c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 136c4762a1bSJed Brown (unless user explicitly sets preconditioner) 137*7addb90fSBarry Smith -snes_mf_operator : form matrix used to construct the preconditioner as set by the user, 138c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 139c4762a1bSJed Brown products within Newton-Krylov method 140c4762a1bSJed Brown 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Customize nonlinear solver 144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1459566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148c4762a1bSJed Brown Set initial conditions 149c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1509566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x)); 1519566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 1529566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown if (mymonitor) { 155c4762a1bSJed Brown ctx.ports = NULL; 1569566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyMonitor, &ctx, MyDestroy)); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown Set runtime options 161c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1629566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 165c4762a1bSJed Brown Solve nonlinear system 166c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1679566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1689566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 1699566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_BINARY_WORLD)); 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 173c4762a1bSJed Brown are no longer needed. 174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1789566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 180c4762a1bSJed Brown 1819566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 182b122ec5aSJacob Faibussowitsch return 0; 183c4762a1bSJed Brown } 184c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 187c4762a1bSJed Brown 188c4762a1bSJed Brown Input Parameters: 189c4762a1bSJed Brown . ts - the TS context 190c4762a1bSJed Brown . X - input vector 191c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 192c4762a1bSJed Brown 193c4762a1bSJed Brown Output Parameter: 194c4762a1bSJed Brown . F - function vector 195c4762a1bSJed Brown */ 196d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ptr) 197d71ae5a4SJacob Faibussowitsch { 198c4762a1bSJed Brown DM da; 199c4762a1bSJed Brown PetscInt i, Mx, xs, xm; 200c4762a1bSJed Brown PetscReal hx, sx; 201c4762a1bSJed Brown PetscScalar *x, *f, c, r, l; 202c4762a1bSJed Brown Vec localX; 203c4762a1bSJed Brown UserCtx *ctx = (UserCtx *)ptr; 204c4762a1bSJed Brown PetscReal tol = ctx->tol, theta = ctx->theta, theta_c = ctx->theta_c, a, b; /* a and b are used in the cubic truncation of the log function */ 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2089566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 2099566063dSJacob 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)); 210c4762a1bSJed Brown 2119371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx; 2129371c9d4SSatish Balay sx = 1.0 / (hx * hx); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 216c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 217c4762a1bSJed Brown By placing code between these two statements, computations can be 218c4762a1bSJed Brown done while messages are in transition. 219c4762a1bSJed Brown */ 2209566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 2219566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* 224c4762a1bSJed Brown Get pointers to vector data 225c4762a1bSJed Brown */ 2269566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 2279566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* 230c4762a1bSJed Brown Get local grid boundaries 231c4762a1bSJed Brown */ 2329566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* 235c4762a1bSJed Brown Compute function over the locally owned part of the grid 236c4762a1bSJed Brown */ 237c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 238c4762a1bSJed Brown if (ctx->degenerate) { 239c4762a1bSJed Brown c = (1. - x[i] * x[i]) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 240c4762a1bSJed Brown r = (1. - x[i + 1] * x[i + 1]) * (x[i] + x[i + 2] - 2.0 * x[i + 1]) * sx; 241c4762a1bSJed Brown l = (1. - x[i - 1] * x[i - 1]) * (x[i - 2] + x[i] - 2.0 * x[i - 1]) * sx; 242c4762a1bSJed Brown } else { 243c4762a1bSJed Brown c = (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 244c4762a1bSJed Brown r = (x[i] + x[i + 2] - 2.0 * x[i + 1]) * sx; 245c4762a1bSJed Brown l = (x[i - 2] + x[i] - 2.0 * x[i - 1]) * sx; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown f[i] = -ctx->kappa * (l + r - 2.0 * c) * sx; 248c4762a1bSJed Brown if (ctx->cahnhillard) { 249c4762a1bSJed Brown switch (ctx->energy) { 250d71ae5a4SJacob Faibussowitsch case 1: /* double well */ 251d71ae5a4SJacob Faibussowitsch f[i] += 6. * .25 * x[i] * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (3. * x[i] * x[i] - 1.) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 252d71ae5a4SJacob Faibussowitsch break; 253d71ae5a4SJacob Faibussowitsch case 2: /* double obstacle */ 254d71ae5a4SJacob Faibussowitsch f[i] += -(x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 255d71ae5a4SJacob Faibussowitsch break; 256c4762a1bSJed Brown case 3: /* logarithmic + double well */ 257c4762a1bSJed Brown f[i] += 6. * .25 * x[i] * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (3. * x[i] * x[i] - 1.) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 258c4762a1bSJed Brown if (ctx->truncation == 2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */ 259c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0 * tol) f[i] += (.25 * theta / (tol - tol * tol)) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 260c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0 * tol) f[i] += (.25 * theta / (tol - tol * tol)) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 261c4762a1bSJed Brown else f[i] += 2.0 * theta * x[i] / ((1.0 - x[i] * x[i]) * (1.0 - x[i] * x[i])) * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (theta / (1.0 - x[i] * x[i])) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 262c4762a1bSJed Brown } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */ 263c4762a1bSJed Brown a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol)); 264c4762a1bSJed Brown b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol); 265c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0 * tol) f[i] += -1.0 * a * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (-1.0 * a * x[i] + b) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 266c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0 * tol) f[i] += 1.0 * a * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (a * x[i] + b) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 267c4762a1bSJed Brown else f[i] += 2.0 * theta * x[i] / ((1.0 - x[i] * x[i]) * (1.0 - x[i] * x[i])) * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (theta / (1.0 - x[i] * x[i])) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown break; 270c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 271c4762a1bSJed Brown f[i] += -theta_c * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 272c4762a1bSJed Brown if (ctx->truncation == 2) { /* quadratic */ 273c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0 * tol) f[i] += (.25 * theta / (tol - tol * tol)) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 274c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0 * tol) f[i] += (.25 * theta / (tol - tol * tol)) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 275c4762a1bSJed Brown else f[i] += 2.0 * theta * x[i] / ((1.0 - x[i] * x[i]) * (1.0 - x[i] * x[i])) * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (theta / (1.0 - x[i] * x[i])) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 276c4762a1bSJed Brown } else { /* cubic */ 277c4762a1bSJed Brown a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol)); 278c4762a1bSJed Brown b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol); 279c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0 * tol) f[i] += -1.0 * a * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (-1.0 * a * x[i] + b) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 280c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0 * tol) f[i] += 1.0 * a * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (a * x[i] + b) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 281c4762a1bSJed Brown else f[i] += 2.0 * theta * x[i] / ((1.0 - x[i] * x[i]) * (1.0 - x[i] * x[i])) * .25 * (x[i + 1] - x[i - 1]) * (x[i + 1] - x[i - 1]) * sx + (theta / (1.0 - x[i] * x[i])) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 282c4762a1bSJed Brown } 283c4762a1bSJed Brown break; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* 289c4762a1bSJed Brown Restore vectors 290c4762a1bSJed Brown */ 2919566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 2929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 298c4762a1bSJed Brown /* 299c4762a1bSJed Brown FormJacobian - Evaluates nonlinear function's Jacobian 300c4762a1bSJed Brown 301c4762a1bSJed Brown */ 302d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(TS ts, PetscReal ftime, Vec X, Mat A, Mat B, void *ptr) 303d71ae5a4SJacob Faibussowitsch { 304c4762a1bSJed Brown DM da; 305c4762a1bSJed Brown PetscInt i, Mx, xs, xm; 306c4762a1bSJed Brown MatStencil row, cols[5]; 307c4762a1bSJed Brown PetscReal hx, sx; 308c4762a1bSJed Brown PetscScalar *x, vals[5]; 309c4762a1bSJed Brown Vec localX; 310c4762a1bSJed Brown UserCtx *ctx = (UserCtx *)ptr; 311c4762a1bSJed Brown 312c4762a1bSJed Brown PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 3149566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 3159566063dSJacob 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)); 316c4762a1bSJed Brown 3179371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx; 3189371c9d4SSatish Balay sx = 1.0 / (hx * hx); 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* 321c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 322c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 323c4762a1bSJed Brown By placing code between these two statements, computations can be 324c4762a1bSJed Brown done while messages are in transition. 325c4762a1bSJed Brown */ 3269566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 3279566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* 330c4762a1bSJed Brown Get pointers to vector data 331c4762a1bSJed Brown */ 3329566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 333c4762a1bSJed Brown 334c4762a1bSJed Brown /* 335c4762a1bSJed Brown Get local grid boundaries 336c4762a1bSJed Brown */ 3379566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 338c4762a1bSJed Brown 339c4762a1bSJed Brown /* 340c4762a1bSJed Brown Compute function over the locally owned part of the grid 341c4762a1bSJed Brown */ 342c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 343c4762a1bSJed Brown row.i = i; 344c4762a1bSJed Brown if (ctx->degenerate) { 345c4762a1bSJed Brown /*PetscScalar c,r,l; 346c4762a1bSJed Brown c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 347c4762a1bSJed Brown r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx; 348c4762a1bSJed Brown l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */ 349c4762a1bSJed Brown } else { 3509371c9d4SSatish Balay cols[0].i = i - 2; 3519371c9d4SSatish Balay vals[0] = -ctx->kappa * sx * sx; 3529371c9d4SSatish Balay cols[1].i = i - 1; 3539371c9d4SSatish Balay vals[1] = 4.0 * ctx->kappa * sx * sx; 3549371c9d4SSatish Balay cols[2].i = i; 3559371c9d4SSatish Balay vals[2] = -6.0 * ctx->kappa * sx * sx; 3569371c9d4SSatish Balay cols[3].i = i + 1; 3579371c9d4SSatish Balay vals[3] = 4.0 * ctx->kappa * sx * sx; 3589371c9d4SSatish Balay cols[4].i = i + 2; 3599371c9d4SSatish Balay vals[4] = -ctx->kappa * sx * sx; 360c4762a1bSJed Brown } 3619566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(B, 1, &row, 5, cols, vals, INSERT_VALUES)); 362c4762a1bSJed Brown 363c4762a1bSJed Brown if (ctx->cahnhillard) { 364c4762a1bSJed Brown switch (ctx->energy) { 365c4762a1bSJed Brown case 1: /* double well */ 366c4762a1bSJed Brown /* f[i] += 6.*.25*x[i]*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (3.*x[i]*x[i] - 1.)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; */ 367c4762a1bSJed Brown break; 368c4762a1bSJed Brown case 2: /* double obstacle */ 369c4762a1bSJed Brown /* f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */ 370c4762a1bSJed Brown break; 371d71ae5a4SJacob Faibussowitsch case 3: /* logarithmic + double well */ 372d71ae5a4SJacob Faibussowitsch break; 373d71ae5a4SJacob Faibussowitsch case 4: /* logarithmic + double obstacle */ 374d71ae5a4SJacob Faibussowitsch break; 375c4762a1bSJed Brown } 376c4762a1bSJed Brown } 377c4762a1bSJed Brown } 378c4762a1bSJed Brown 379c4762a1bSJed Brown /* 380c4762a1bSJed Brown Restore vectors 381c4762a1bSJed Brown */ 3829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 3839566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 3849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 3859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 386c4762a1bSJed Brown if (A != B) { 3879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 389c4762a1bSJed Brown } 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 393d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec U) 394d71ae5a4SJacob Faibussowitsch { 395c4762a1bSJed Brown PetscInt i, xs, xm, Mx, N, scale; 396c4762a1bSJed Brown PetscScalar *u; 397c4762a1bSJed Brown PetscReal r, hx, x; 398c4762a1bSJed Brown const PetscScalar *f; 399c4762a1bSJed Brown Vec finesolution; 400c4762a1bSJed Brown PetscViewer viewer; 401c4762a1bSJed Brown 402c4762a1bSJed Brown PetscFunctionBegin; 4039566063dSJacob 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)); 404c4762a1bSJed Brown 405c4762a1bSJed Brown hx = 1.0 / (PetscReal)Mx; 406c4762a1bSJed Brown 407c4762a1bSJed Brown /* 408c4762a1bSJed Brown Get pointers to vector data 409c4762a1bSJed Brown */ 4109566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 411c4762a1bSJed Brown 412c4762a1bSJed Brown /* 413c4762a1bSJed Brown Get local grid boundaries 414c4762a1bSJed Brown */ 4159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 416c4762a1bSJed Brown 417c4762a1bSJed Brown /* 418c4762a1bSJed Brown Seee heat.c for how to generate InitialSolution.heat 419c4762a1bSJed Brown */ 4209566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "InitialSolution.heat", FILE_MODE_READ, &viewer)); 4219566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &finesolution)); 4229566063dSJacob Faibussowitsch PetscCall(VecLoad(finesolution, viewer)); 4239566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 4249566063dSJacob Faibussowitsch PetscCall(VecGetSize(finesolution, &N)); 425c4762a1bSJed Brown scale = N / Mx; 4269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(finesolution, &f)); 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* 429c4762a1bSJed Brown Compute function over the locally owned part of the grid 430c4762a1bSJed Brown */ 431c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 432c4762a1bSJed Brown x = i * hx; 433c4762a1bSJed Brown r = PetscSqrtReal((x - .5) * (x - .5)); 434c4762a1bSJed Brown if (r < .125) u[i] = 1.0; 435c4762a1bSJed Brown else u[i] = -.5; 436c4762a1bSJed Brown 437c4762a1bSJed Brown /* With the initial condition above the method is first order in space */ 438c4762a1bSJed Brown /* this is a smooth initial condition so the method becomes second order in space */ 439c4762a1bSJed Brown /*u[i] = PetscSinScalar(2*PETSC_PI*x); */ 440c4762a1bSJed Brown u[i] = f[scale * i]; 441c4762a1bSJed Brown } 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(finesolution, &f)); 4439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&finesolution)); 444c4762a1bSJed Brown 445c4762a1bSJed Brown /* 446c4762a1bSJed Brown Restore vectors 447c4762a1bSJed Brown */ 4489566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 4493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 450c4762a1bSJed Brown } 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* 453c4762a1bSJed Brown This routine is not parallel 454c4762a1bSJed Brown */ 455d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMonitor(TS ts, PetscInt step, PetscReal time, Vec U, void *ptr) 456d71ae5a4SJacob Faibussowitsch { 457c4762a1bSJed Brown UserCtx *ctx = (UserCtx *)ptr; 458c4762a1bSJed Brown PetscDrawLG lg; 459c4762a1bSJed Brown PetscScalar *u, l, r, c; 460c4762a1bSJed Brown PetscInt Mx, i, xs, xm, cnt; 461c4762a1bSJed Brown PetscReal x, y, hx, pause, sx, len, max, xx[4], yy[4], xx_netforce, yy_netforce, yup, ydown, y2, len2; 462c4762a1bSJed Brown PetscDraw draw; 463c4762a1bSJed Brown Vec localU; 464c4762a1bSJed Brown DM da; 465c4762a1bSJed Brown int colors[] = {PETSC_DRAW_YELLOW, PETSC_DRAW_RED, PETSC_DRAW_BLUE, PETSC_DRAW_PLUM, PETSC_DRAW_BLACK}; 466c4762a1bSJed Brown /* 467c4762a1bSJed Brown const char *const legend[3][3] = {{"-kappa (\\grad u,\\grad u)","(1 - u^2)^2"},{"-kappa (\\grad u,\\grad u)","(1 - u^2)"},{"-kappa (\\grad u,\\grad u)","logarithmic"}}; 468c4762a1bSJed Brown */ 469c4762a1bSJed Brown PetscDrawAxis axis; 470c4762a1bSJed Brown PetscDrawViewPorts *ports; 471c4762a1bSJed Brown PetscReal tol = ctx->tol, theta = ctx->theta, theta_c = ctx->theta_c, a, b; /* a and b are used in the cubic truncation of the log function */ 472c4762a1bSJed Brown PetscReal vbounds[] = {-1.1, 1.1}; 473c4762a1bSJed Brown 474c4762a1bSJed Brown PetscFunctionBegin; 4759566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, vbounds)); 4769566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 800, 600)); 4779566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 4789566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 4799371c9d4SSatish Balay 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)); 4809566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 4819371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx; 4829371c9d4SSatish Balay sx = 1.0 / (hx * hx); 4839566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 4849566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 4859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 486c4762a1bSJed Brown 4879566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, &lg)); 4889566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetDraw(lg, &draw)); 4899566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 49048a46eb9SPierre Jolivet if (!ctx->ports) PetscCall(PetscDrawViewPortsCreateRect(draw, 1, 3, &ctx->ports)); 491c4762a1bSJed Brown ports = ctx->ports; 4929566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lg, &axis)); 4939566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 494c4762a1bSJed Brown 4959371c9d4SSatish Balay xx[0] = 0.0; 4969371c9d4SSatish Balay xx[1] = 1.0; 4979371c9d4SSatish Balay cnt = 2; 4989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-zoom", xx, &cnt, NULL)); 4999371c9d4SSatish Balay xs = xx[0] / hx; 5009371c9d4SSatish Balay xm = (xx[1] - xx[0]) / hx; 501c4762a1bSJed Brown 502c4762a1bSJed Brown /* 503c4762a1bSJed Brown Plot the energies 504c4762a1bSJed Brown */ 5059566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3))); 5069566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors + 1)); 5079566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 2)); 508c4762a1bSJed Brown x = hx * xs; 509c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 510c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = x; 511c4762a1bSJed Brown if (ctx->degenerate) yy[0] = PetscRealPart(.25 * (1. - u[i] * u[i]) * ctx->kappa * (u[i - 1] - u[i + 1]) * (u[i - 1] - u[i + 1]) * sx); 512c4762a1bSJed Brown else yy[0] = PetscRealPart(.25 * ctx->kappa * (u[i - 1] - u[i + 1]) * (u[i - 1] - u[i + 1]) * sx); 513c4762a1bSJed Brown 514c4762a1bSJed Brown if (ctx->cahnhillard) { 515c4762a1bSJed Brown switch (ctx->energy) { 516d71ae5a4SJacob Faibussowitsch case 1: /* double well */ 517d71ae5a4SJacob Faibussowitsch yy[1] = .25 * PetscRealPart((1. - u[i] * u[i]) * (1. - u[i] * u[i])); 518d71ae5a4SJacob Faibussowitsch break; 519d71ae5a4SJacob Faibussowitsch case 2: /* double obstacle */ 520d71ae5a4SJacob Faibussowitsch yy[1] = .5 * PetscRealPart(1. - u[i] * u[i]); 521d71ae5a4SJacob Faibussowitsch break; 522c4762a1bSJed Brown case 3: /* logarithm + double well */ 523c4762a1bSJed Brown yy[1] = .25 * PetscRealPart((1. - u[i] * u[i]) * (1. - u[i] * u[i])); 524c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) yy[2] = .5 * theta * (2.0 * tol * PetscLogReal(tol) + PetscRealPart(1.0 - u[i]) * PetscLogReal(PetscRealPart(1. - u[i]) / 2.0)); 525c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) yy[2] = .5 * theta * (PetscRealPart(1.0 + u[i]) * PetscLogReal(PetscRealPart(1.0 + u[i]) / 2.0) + 2.0 * tol * PetscLogReal(tol)); 526c4762a1bSJed Brown else yy[2] = .5 * theta * (PetscRealPart(1.0 + u[i]) * PetscLogReal(PetscRealPart(1.0 + u[i]) / 2.0) + PetscRealPart(1.0 - u[i]) * PetscLogReal(PetscRealPart(1.0 - u[i]) / 2.0)); 527c4762a1bSJed Brown break; 528c4762a1bSJed Brown case 4: /* logarithm + double obstacle */ 529c4762a1bSJed Brown yy[1] = .5 * theta_c * PetscRealPart(1.0 - u[i] * u[i]); 530c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) yy[2] = .5 * theta * (2.0 * tol * PetscLogReal(tol) + PetscRealPart(1.0 - u[i]) * PetscLogReal(PetscRealPart(1. - u[i]) / 2.0)); 531c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) yy[2] = .5 * theta * (PetscRealPart(1.0 + u[i]) * PetscLogReal(PetscRealPart(1.0 + u[i]) / 2.0) + 2.0 * tol * PetscLogReal(tol)); 532c4762a1bSJed Brown else yy[2] = .5 * theta * (PetscRealPart(1.0 + u[i]) * PetscLogReal(PetscRealPart(1.0 + u[i]) / 2.0) + PetscRealPart(1.0 - u[i]) * PetscLogReal(PetscRealPart(1.0 - u[i]) / 2.0)); 533c4762a1bSJed Brown break; 534d71ae5a4SJacob Faibussowitsch default: 535d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "It will always be one of the values"); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown } 5389566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 539c4762a1bSJed Brown x += hx; 540c4762a1bSJed Brown } 5419566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause)); 5429566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 5439566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Energy", "", "")); 5449566063dSJacob Faibussowitsch /* PetscCall(PetscDrawLGSetLegend(lg,legend[ctx->energy-1])); */ 5459566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 546c4762a1bSJed Brown 547c4762a1bSJed Brown /* 548c4762a1bSJed Brown Plot the forces 549c4762a1bSJed Brown */ 5509566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3))); 5519566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors + 1)); 5529566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 1)); 5539566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 554c4762a1bSJed Brown x = xs * hx; 555c4762a1bSJed Brown max = 0.; 556c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 557c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = xx[3] = x; 558c4762a1bSJed Brown xx_netforce = x; 559c4762a1bSJed Brown if (ctx->degenerate) { 560c4762a1bSJed Brown c = (1. - u[i] * u[i]) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 561c4762a1bSJed Brown r = (1. - u[i + 1] * u[i + 1]) * (u[i] + u[i + 2] - 2.0 * u[i + 1]) * sx; 562c4762a1bSJed Brown l = (1. - u[i - 1] * u[i - 1]) * (u[i - 2] + u[i] - 2.0 * u[i - 1]) * sx; 563c4762a1bSJed Brown } else { 564c4762a1bSJed Brown c = (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 565c4762a1bSJed Brown r = (u[i] + u[i + 2] - 2.0 * u[i + 1]) * sx; 566c4762a1bSJed Brown l = (u[i - 2] + u[i] - 2.0 * u[i - 1]) * sx; 567c4762a1bSJed Brown } 568c4762a1bSJed Brown yy[0] = PetscRealPart(-ctx->kappa * (l + r - 2.0 * c) * sx); 569c4762a1bSJed Brown yy_netforce = yy[0]; 570c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[0])); 571c4762a1bSJed Brown if (ctx->cahnhillard) { 572c4762a1bSJed Brown switch (ctx->energy) { 573d71ae5a4SJacob Faibussowitsch case 1: /* double well */ 574d71ae5a4SJacob Faibussowitsch yy[1] = PetscRealPart(6. * .25 * u[i] * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (3. * u[i] * u[i] - 1.) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 575d71ae5a4SJacob Faibussowitsch break; 576d71ae5a4SJacob Faibussowitsch case 2: /* double obstacle */ 577d71ae5a4SJacob Faibussowitsch yy[1] = -PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 578d71ae5a4SJacob Faibussowitsch break; 579c4762a1bSJed Brown case 3: /* logarithmic + double well */ 580c4762a1bSJed Brown yy[1] = PetscRealPart(6. * .25 * u[i] * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (3. * u[i] * u[i] - 1.) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 581c4762a1bSJed Brown if (ctx->truncation == 2) { /* quadratic */ 582c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) yy[2] = (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 583c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) yy[2] = (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 584c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 585c4762a1bSJed Brown } else { /* cubic */ 586c4762a1bSJed Brown a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol)); 587c4762a1bSJed Brown b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol); 588c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) yy[2] = PetscRealPart(-1.0 * a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (-1.0 * a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 589c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) yy[2] = PetscRealPart(1.0 * a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 590c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 591c4762a1bSJed Brown } 592c4762a1bSJed Brown break; 593c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 594c4762a1bSJed Brown yy[1] = theta_c * PetscRealPart(-(u[i - 1] + u[i + 1] - 2.0 * u[i])) * sx; 595c4762a1bSJed Brown if (ctx->truncation == 2) { 596c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) yy[2] = (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 597c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) yy[2] = (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 598c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 5999371c9d4SSatish Balay } else { 600c4762a1bSJed Brown a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol)); 601c4762a1bSJed Brown b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol); 602c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) yy[2] = PetscRealPart(-1.0 * a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (-1.0 * a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 603c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) yy[2] = PetscRealPart(1.0 * a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 604c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 605c4762a1bSJed Brown } 606c4762a1bSJed Brown break; 607d71ae5a4SJacob Faibussowitsch default: 608d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "It will always be one of the values"); 609c4762a1bSJed Brown } 610c4762a1bSJed Brown if (ctx->energy < 3) { 611c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[1])); 612c4762a1bSJed Brown yy[2] = yy[0] + yy[1]; 613c4762a1bSJed Brown yy_netforce = yy[2]; 614c4762a1bSJed Brown } else { 615c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[1] + yy[2])); 616c4762a1bSJed Brown yy[3] = yy[0] + yy[1] + yy[2]; 617c4762a1bSJed Brown yy_netforce = yy[3]; 618c4762a1bSJed Brown } 619c4762a1bSJed Brown } 620c4762a1bSJed Brown if (ctx->netforce) { 6219566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, &xx_netforce, &yy_netforce)); 622c4762a1bSJed Brown } else { 6239566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 624c4762a1bSJed Brown } 625c4762a1bSJed Brown x += hx; 626c4762a1bSJed Brown /*if (max > 7200150000.0) */ 627c4762a1bSJed Brown /* printf("max very big when i = %d\n",i); */ 628c4762a1bSJed Brown } 6299566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Right hand side", "", "")); 6309566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, NULL)); 6319566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 632c4762a1bSJed Brown 633c4762a1bSJed Brown /* 634c4762a1bSJed Brown Plot the solution 635c4762a1bSJed Brown */ 6369566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 1)); 6379566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 0)); 6389566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 639c4762a1bSJed Brown x = hx * xs; 6409566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLimits(lg, x, x + (xm - 1) * hx, -1.1, 1.1)); 6419566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors)); 642c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 643c4762a1bSJed Brown xx[0] = x; 644c4762a1bSJed Brown yy[0] = PetscRealPart(u[i]); 6459566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 646c4762a1bSJed Brown x += hx; 647c4762a1bSJed Brown } 6489566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution", "", "")); 6499566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 650c4762a1bSJed Brown 651c4762a1bSJed Brown /* 652c4762a1bSJed Brown Print the forces as arrows on the solution 653c4762a1bSJed Brown */ 654c4762a1bSJed Brown x = hx * xs; 655c4762a1bSJed Brown cnt = xm / 60; 656c4762a1bSJed Brown cnt = (!cnt) ? 1 : cnt; 657c4762a1bSJed Brown 658c4762a1bSJed Brown for (i = xs; i < xs + xm; i += cnt) { 659c4762a1bSJed Brown y = yup = ydown = PetscRealPart(u[i]); 660c4762a1bSJed Brown c = (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; 661c4762a1bSJed Brown r = (u[i] + u[i + 2] - 2.0 * u[i + 1]) * sx; 662c4762a1bSJed Brown l = (u[i - 2] + u[i] - 2.0 * u[i - 1]) * sx; 663c4762a1bSJed Brown len = -.5 * PetscRealPart(ctx->kappa * (l + r - 2.0 * c) * sx) / max; 6649566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_RED)); 665c4762a1bSJed Brown if (ctx->cahnhillard) { 666c4762a1bSJed Brown if (len < 0.) ydown += len; 667c4762a1bSJed Brown else yup += len; 668c4762a1bSJed Brown 669c4762a1bSJed Brown switch (ctx->energy) { 670d71ae5a4SJacob Faibussowitsch case 1: /* double well */ 671d71ae5a4SJacob Faibussowitsch len = .5 * PetscRealPart(6. * .25 * u[i] * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (3. * u[i] * u[i] - 1.) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max; 672d71ae5a4SJacob Faibussowitsch break; 673d71ae5a4SJacob Faibussowitsch case 2: /* double obstacle */ 674d71ae5a4SJacob Faibussowitsch len = -.5 * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max; 675d71ae5a4SJacob Faibussowitsch break; 676c4762a1bSJed Brown case 3: /* logarithmic + double well */ 677c4762a1bSJed Brown len = .5 * PetscRealPart(6. * .25 * u[i] * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (3. * u[i] * u[i] - 1.) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max; 678c4762a1bSJed Brown if (len < 0.) ydown += len; 679c4762a1bSJed Brown else yup += len; 680c4762a1bSJed Brown 681c4762a1bSJed Brown if (ctx->truncation == 2) { /* quadratic */ 682c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) len2 = .5 * (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max; 683c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) len2 = .5 * (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max; 684c4762a1bSJed Brown else len2 = PetscRealPart(.5 * (2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max); 685c4762a1bSJed Brown } else { /* cubic */ 686c4762a1bSJed Brown a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol)); 687c4762a1bSJed Brown b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol); 688c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) len2 = PetscRealPart(.5 * (-1.0 * a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (-1.0 * a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max); 689c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) len2 = PetscRealPart(.5 * (a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max); 690c4762a1bSJed Brown else len2 = PetscRealPart(.5 * (2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max); 691c4762a1bSJed Brown } 692c4762a1bSJed Brown y2 = len < 0 ? ydown : yup; 6939566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y2, x, y2 + len2, PETSC_DRAW_PLUM)); 694c4762a1bSJed Brown break; 695c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 696c4762a1bSJed Brown len = -.5 * theta_c * PetscRealPart(-(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max); 697c4762a1bSJed Brown if (len < 0.) ydown += len; 698c4762a1bSJed Brown else yup += len; 699c4762a1bSJed Brown 700c4762a1bSJed Brown if (ctx->truncation == 2) { /* quadratic */ 701c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) len2 = .5 * (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max; 702c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) len2 = .5 * (.25 * theta / (tol - tol * tol)) * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max; 703c4762a1bSJed Brown else len2 = PetscRealPart(.5 * (2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max); 704c4762a1bSJed Brown } else { /* cubic */ 705c4762a1bSJed Brown a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol)); 706c4762a1bSJed Brown b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol); 707c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0 * tol) len2 = .5 * PetscRealPart(-1.0 * a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (-1.0 * a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max; 708c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0 * tol) len2 = .5 * PetscRealPart(a * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (a * u[i] + b) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max; 709c4762a1bSJed Brown else len2 = .5 * PetscRealPart(2.0 * theta * u[i] / ((1.0 - u[i] * u[i]) * (1.0 - u[i] * u[i])) * .25 * (u[i + 1] - u[i - 1]) * (u[i + 1] - u[i - 1]) * sx + (theta / (1.0 - u[i] * u[i])) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max; 710c4762a1bSJed Brown } 711c4762a1bSJed Brown y2 = len < 0 ? ydown : yup; 7129566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y2, x, y2 + len2, PETSC_DRAW_PLUM)); 713c4762a1bSJed Brown break; 714c4762a1bSJed Brown } 7159566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_BLUE)); 716c4762a1bSJed Brown } 717c4762a1bSJed Brown x += cnt * hx; 718c4762a1bSJed Brown } 7199566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &x)); 7209566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 7219566063dSJacob Faibussowitsch PetscCall(PetscDrawStringSetSize(draw, .2, .2)); 7229566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 7239566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause)); 7249566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 7253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 726c4762a1bSJed Brown } 727c4762a1bSJed Brown 728d71ae5a4SJacob Faibussowitsch PetscErrorCode MyDestroy(void **ptr) 729d71ae5a4SJacob Faibussowitsch { 730c4762a1bSJed Brown UserCtx *ctx = *(UserCtx **)ptr; 731c4762a1bSJed Brown 732c4762a1bSJed Brown PetscFunctionBegin; 7339566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ctx->ports)); 7343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 735c4762a1bSJed Brown } 736c4762a1bSJed Brown 737c4762a1bSJed Brown /*TEST 738c4762a1bSJed Brown 739c4762a1bSJed Brown test: 740c4762a1bSJed Brown TODO: currently requires initial condition file generated by heat 741c4762a1bSJed Brown 742c4762a1bSJed Brown TEST*/ 743