xref: /petsc/src/ts/tutorials/phasefield/biharmonic.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown   Solves the equation
6c4762a1bSJed Brown 
7c4762a1bSJed Brown     u_t = - kappa  \Delta \Delta u
8c4762a1bSJed Brown     Periodic boundary conditions
9c4762a1bSJed Brown 
10c4762a1bSJed Brown Evolve the biharmonic heat equation:
11c4762a1bSJed Brown ---------------
12c4762a1bSJed 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
13c4762a1bSJed Brown 
14c4762a1bSJed Brown Evolve with the restriction that -1 <= u <= 1; i.e. as a variational inequality
15c4762a1bSJed Brown ---------------
16c4762a1bSJed 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
17c4762a1bSJed Brown 
18c4762a1bSJed Brown    u_t =  kappa \Delta \Delta u +   6.*u*(u_x)^2 + (3*u^2 - 12) \Delta u
19c4762a1bSJed Brown     -1 <= u <= 1
20c4762a1bSJed Brown     Periodic boundary conditions
21c4762a1bSJed Brown 
22c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double well Initial hump shrinks then grows
23c4762a1bSJed Brown ---------------
24c4762a1bSJed 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
25c4762a1bSJed Brown 
26c4762a1bSJed Brown Initial hump neither shrinks nor grows when degenerate (otherwise similar solution)
27c4762a1bSJed Brown 
28c4762a1bSJed 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
29c4762a1bSJed Brown 
30c4762a1bSJed 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
31c4762a1bSJed Brown 
32c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double obstacle
33c4762a1bSJed Brown ---------------
34c4762a1bSJed 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
35c4762a1bSJed Brown 
36c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic + double well (never shrinks and then grows)
37c4762a1bSJed Brown ---------------
38c4762a1bSJed 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
39c4762a1bSJed Brown 
40c4762a1bSJed 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
41c4762a1bSJed Brown 
42c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic +  double obstacle (never shrinks, never grows)
43c4762a1bSJed Brown ---------------
44c4762a1bSJed 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
45c4762a1bSJed Brown 
46c4762a1bSJed Brown */
47c4762a1bSJed Brown #include <petscdm.h>
48c4762a1bSJed Brown #include <petscdmda.h>
49c4762a1bSJed Brown #include <petscts.h>
50c4762a1bSJed Brown #include <petscdraw.h>
51c4762a1bSJed Brown 
52c4762a1bSJed 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 *);
539371c9d4SSatish Balay typedef struct {
549371c9d4SSatish Balay   PetscBool           cahnhillard;
559371c9d4SSatish Balay   PetscBool           degenerate;
569371c9d4SSatish Balay   PetscReal           kappa;
579371c9d4SSatish Balay   PetscInt            energy;
589371c9d4SSatish Balay   PetscReal           tol;
599371c9d4SSatish Balay   PetscReal           theta, theta_c;
609371c9d4SSatish Balay   PetscInt            truncation;
619371c9d4SSatish Balay   PetscBool           netforce;
629371c9d4SSatish Balay   PetscDrawViewPorts *ports;
639371c9d4SSatish Balay } UserCtx;
64c4762a1bSJed Brown 
659371c9d4SSatish Balay int main(int argc, char **argv) {
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;
799566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, 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)
137c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix 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  */
1969371c9d4SSatish Balay PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ptr) {
197c4762a1bSJed Brown   DM           da;
198c4762a1bSJed Brown   PetscInt     i, Mx, xs, xm;
199c4762a1bSJed Brown   PetscReal    hx, sx;
200c4762a1bSJed Brown   PetscScalar *x, *f, c, r, l;
201c4762a1bSJed Brown   Vec          localX;
202c4762a1bSJed Brown   UserCtx     *ctx = (UserCtx *)ptr;
203c4762a1bSJed 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 */
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2079566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
2089566063dSJacob 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));
209c4762a1bSJed Brown 
2109371c9d4SSatish Balay   hx = 1.0 / (PetscReal)Mx;
2119371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /*
214c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
215c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
216c4762a1bSJed Brown      By placing code between these two statements, computations can be
217c4762a1bSJed Brown      done while messages are in transition.
218c4762a1bSJed Brown   */
2199566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
2209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /*
223c4762a1bSJed Brown      Get pointers to vector data
224c4762a1bSJed Brown   */
2259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localX, &x));
2269566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /*
229c4762a1bSJed Brown      Get local grid boundaries
230c4762a1bSJed Brown   */
2319566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /*
234c4762a1bSJed Brown      Compute function over the locally owned part of the grid
235c4762a1bSJed Brown   */
236c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
237c4762a1bSJed Brown     if (ctx->degenerate) {
238c4762a1bSJed Brown       c = (1. - x[i] * x[i]) * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx;
239c4762a1bSJed Brown       r = (1. - x[i + 1] * x[i + 1]) * (x[i] + x[i + 2] - 2.0 * x[i + 1]) * sx;
240c4762a1bSJed Brown       l = (1. - x[i - 1] * x[i - 1]) * (x[i - 2] + x[i] - 2.0 * x[i - 1]) * sx;
241c4762a1bSJed Brown     } else {
242c4762a1bSJed Brown       c = (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx;
243c4762a1bSJed Brown       r = (x[i] + x[i + 2] - 2.0 * x[i + 1]) * sx;
244c4762a1bSJed Brown       l = (x[i - 2] + x[i] - 2.0 * x[i - 1]) * sx;
245c4762a1bSJed Brown     }
246c4762a1bSJed Brown     f[i] = -ctx->kappa * (l + r - 2.0 * c) * sx;
247c4762a1bSJed Brown     if (ctx->cahnhillard) {
248c4762a1bSJed Brown       switch (ctx->energy) {
2499371c9d4SSatish Balay       case 1: /*  double well */ 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; break;
2509371c9d4SSatish Balay       case 2: /* double obstacle */ f[i] += -(x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; break;
251c4762a1bSJed Brown       case 3: /* logarithmic + double well */
252c4762a1bSJed 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;
253c4762a1bSJed Brown         if (ctx->truncation == 2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */
254c4762a1bSJed 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;
255c4762a1bSJed 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;
256c4762a1bSJed 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;
257c4762a1bSJed Brown         } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */
258c4762a1bSJed Brown           a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol));
259c4762a1bSJed Brown           b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol);
260c4762a1bSJed 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;
261c4762a1bSJed 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;
262c4762a1bSJed 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;
263c4762a1bSJed Brown         }
264c4762a1bSJed Brown         break;
265c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
266c4762a1bSJed Brown         f[i] += -theta_c * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx;
267c4762a1bSJed Brown         if (ctx->truncation == 2) { /* quadratic */
268c4762a1bSJed 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;
269c4762a1bSJed 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;
270c4762a1bSJed 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;
271c4762a1bSJed Brown         } else { /* cubic */
272c4762a1bSJed Brown           a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol));
273c4762a1bSJed Brown           b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol);
274c4762a1bSJed 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;
275c4762a1bSJed 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;
276c4762a1bSJed 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;
277c4762a1bSJed Brown         }
278c4762a1bSJed Brown         break;
279c4762a1bSJed Brown       }
280c4762a1bSJed Brown     }
281c4762a1bSJed Brown   }
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /*
284c4762a1bSJed Brown      Restore vectors
285c4762a1bSJed Brown   */
2869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
2879566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2889566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
289c4762a1bSJed Brown   PetscFunctionReturn(0);
290c4762a1bSJed Brown }
291c4762a1bSJed Brown 
292c4762a1bSJed Brown /* ------------------------------------------------------------------- */
293c4762a1bSJed Brown /*
294c4762a1bSJed Brown    FormJacobian - Evaluates nonlinear function's Jacobian
295c4762a1bSJed Brown 
296c4762a1bSJed Brown */
2979371c9d4SSatish Balay PetscErrorCode FormJacobian(TS ts, PetscReal ftime, Vec X, Mat A, Mat B, void *ptr) {
298c4762a1bSJed Brown   DM           da;
299c4762a1bSJed Brown   PetscInt     i, Mx, xs, xm;
300c4762a1bSJed Brown   MatStencil   row, cols[5];
301c4762a1bSJed Brown   PetscReal    hx, sx;
302c4762a1bSJed Brown   PetscScalar *x, vals[5];
303c4762a1bSJed Brown   Vec          localX;
304c4762a1bSJed Brown   UserCtx     *ctx = (UserCtx *)ptr;
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   PetscFunctionBegin;
3079566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
3089566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
3099566063dSJacob 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));
310c4762a1bSJed Brown 
3119371c9d4SSatish Balay   hx = 1.0 / (PetscReal)Mx;
3129371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /*
315c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
316c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
317c4762a1bSJed Brown      By placing code between these two statements, computations can be
318c4762a1bSJed Brown      done while messages are in transition.
319c4762a1bSJed Brown   */
3209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
3219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /*
324c4762a1bSJed Brown      Get pointers to vector data
325c4762a1bSJed Brown   */
3269566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localX, &x));
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   /*
329c4762a1bSJed Brown      Get local grid boundaries
330c4762a1bSJed Brown   */
3319566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   /*
334c4762a1bSJed Brown      Compute function over the locally owned part of the grid
335c4762a1bSJed Brown   */
336c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
337c4762a1bSJed Brown     row.i = i;
338c4762a1bSJed Brown     if (ctx->degenerate) {
339c4762a1bSJed Brown       /*PetscScalar c,r,l;
340c4762a1bSJed Brown       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
341c4762a1bSJed Brown       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
342c4762a1bSJed Brown       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */
343c4762a1bSJed Brown     } else {
3449371c9d4SSatish Balay       cols[0].i = i - 2;
3459371c9d4SSatish Balay       vals[0]   = -ctx->kappa * sx * sx;
3469371c9d4SSatish Balay       cols[1].i = i - 1;
3479371c9d4SSatish Balay       vals[1]   = 4.0 * ctx->kappa * sx * sx;
3489371c9d4SSatish Balay       cols[2].i = i;
3499371c9d4SSatish Balay       vals[2]   = -6.0 * ctx->kappa * sx * sx;
3509371c9d4SSatish Balay       cols[3].i = i + 1;
3519371c9d4SSatish Balay       vals[3]   = 4.0 * ctx->kappa * sx * sx;
3529371c9d4SSatish Balay       cols[4].i = i + 2;
3539371c9d4SSatish Balay       vals[4]   = -ctx->kappa * sx * sx;
354c4762a1bSJed Brown     }
3559566063dSJacob Faibussowitsch     PetscCall(MatSetValuesStencil(B, 1, &row, 5, cols, vals, INSERT_VALUES));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     if (ctx->cahnhillard) {
358c4762a1bSJed Brown       switch (ctx->energy) {
359c4762a1bSJed Brown       case 1: /* double well */
360c4762a1bSJed 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; */
361c4762a1bSJed Brown         break;
362c4762a1bSJed Brown       case 2: /* double obstacle */
363c4762a1bSJed Brown         /*        f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */
364c4762a1bSJed Brown         break;
3659371c9d4SSatish Balay       case 3: /* logarithmic + double well */ break;
3669371c9d4SSatish Balay       case 4: /* logarithmic + double obstacle */ break;
367c4762a1bSJed Brown       }
368c4762a1bSJed Brown     }
369c4762a1bSJed Brown   }
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /*
372c4762a1bSJed Brown      Restore vectors
373c4762a1bSJed Brown   */
3749566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
3759566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
3769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
378c4762a1bSJed Brown   if (A != B) {
3799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3809566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
381c4762a1bSJed Brown   }
382c4762a1bSJed Brown   PetscFunctionReturn(0);
383c4762a1bSJed Brown }
384c4762a1bSJed Brown /* ------------------------------------------------------------------- */
3859371c9d4SSatish Balay PetscErrorCode FormInitialSolution(DM da, Vec U) {
386c4762a1bSJed Brown   PetscInt           i, xs, xm, Mx, N, scale;
387c4762a1bSJed Brown   PetscScalar       *u;
388c4762a1bSJed Brown   PetscReal          r, hx, x;
389c4762a1bSJed Brown   const PetscScalar *f;
390c4762a1bSJed Brown   Vec                finesolution;
391c4762a1bSJed Brown   PetscViewer        viewer;
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   PetscFunctionBegin;
3949566063dSJacob 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));
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   hx = 1.0 / (PetscReal)Mx;
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /*
399c4762a1bSJed Brown      Get pointers to vector data
400c4762a1bSJed Brown   */
4019566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   /*
404c4762a1bSJed Brown      Get local grid boundaries
405c4762a1bSJed Brown   */
4069566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   /*
409c4762a1bSJed Brown       Seee heat.c for how to generate InitialSolution.heat
410c4762a1bSJed Brown   */
4119566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "InitialSolution.heat", FILE_MODE_READ, &viewer));
4129566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &finesolution));
4139566063dSJacob Faibussowitsch   PetscCall(VecLoad(finesolution, viewer));
4149566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
4159566063dSJacob Faibussowitsch   PetscCall(VecGetSize(finesolution, &N));
416c4762a1bSJed Brown   scale = N / Mx;
4179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(finesolution, &f));
418c4762a1bSJed Brown 
419c4762a1bSJed Brown   /*
420c4762a1bSJed Brown      Compute function over the locally owned part of the grid
421c4762a1bSJed Brown   */
422c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
423c4762a1bSJed Brown     x = i * hx;
424c4762a1bSJed Brown     r = PetscSqrtReal((x - .5) * (x - .5));
425c4762a1bSJed Brown     if (r < .125) u[i] = 1.0;
426c4762a1bSJed Brown     else u[i] = -.5;
427c4762a1bSJed Brown 
428c4762a1bSJed Brown     /* With the initial condition above the method is first order in space */
429c4762a1bSJed Brown     /* this is a smooth initial condition so the method becomes second order in space */
430c4762a1bSJed Brown     /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
431c4762a1bSJed Brown     u[i] = f[scale * i];
432c4762a1bSJed Brown   }
4339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(finesolution, &f));
4349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&finesolution));
435c4762a1bSJed Brown 
436c4762a1bSJed Brown   /*
437c4762a1bSJed Brown      Restore vectors
438c4762a1bSJed Brown   */
4399566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
440c4762a1bSJed Brown   PetscFunctionReturn(0);
441c4762a1bSJed Brown }
442c4762a1bSJed Brown 
443c4762a1bSJed Brown /*
444c4762a1bSJed Brown     This routine is not parallel
445c4762a1bSJed Brown */
4469371c9d4SSatish Balay PetscErrorCode MyMonitor(TS ts, PetscInt step, PetscReal time, Vec U, void *ptr) {
447c4762a1bSJed Brown   UserCtx            *ctx = (UserCtx *)ptr;
448c4762a1bSJed Brown   PetscDrawLG         lg;
449c4762a1bSJed Brown   PetscScalar        *u, l, r, c;
450c4762a1bSJed Brown   PetscInt            Mx, i, xs, xm, cnt;
451c4762a1bSJed Brown   PetscReal           x, y, hx, pause, sx, len, max, xx[4], yy[4], xx_netforce, yy_netforce, yup, ydown, y2, len2;
452c4762a1bSJed Brown   PetscDraw           draw;
453c4762a1bSJed Brown   Vec                 localU;
454c4762a1bSJed Brown   DM                  da;
455c4762a1bSJed Brown   int                 colors[] = {PETSC_DRAW_YELLOW, PETSC_DRAW_RED, PETSC_DRAW_BLUE, PETSC_DRAW_PLUM, PETSC_DRAW_BLACK};
456c4762a1bSJed Brown   /*
457c4762a1bSJed 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"}};
458c4762a1bSJed Brown    */
459c4762a1bSJed Brown   PetscDrawAxis       axis;
460c4762a1bSJed Brown   PetscDrawViewPorts *ports;
461c4762a1bSJed 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 */
462c4762a1bSJed Brown   PetscReal           vbounds[] = {-1.1, 1.1};
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   PetscFunctionBegin;
4659566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, vbounds));
4669566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 800, 600));
4679566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
4689566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
4699371c9d4SSatish 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));
4709566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
4719371c9d4SSatish Balay   hx = 1.0 / (PetscReal)Mx;
4729371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
4739566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
4749566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
4759566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &u));
476c4762a1bSJed Brown 
4779566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, &lg));
4789566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGGetDraw(lg, &draw));
4799566063dSJacob Faibussowitsch   PetscCall(PetscDrawCheckResizedWindow(draw));
480*48a46eb9SPierre Jolivet   if (!ctx->ports) PetscCall(PetscDrawViewPortsCreateRect(draw, 1, 3, &ctx->ports));
481c4762a1bSJed Brown   ports = ctx->ports;
4829566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGGetAxis(lg, &axis));
4839566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGReset(lg));
484c4762a1bSJed Brown 
4859371c9d4SSatish Balay   xx[0] = 0.0;
4869371c9d4SSatish Balay   xx[1] = 1.0;
4879371c9d4SSatish Balay   cnt   = 2;
4889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-zoom", xx, &cnt, NULL));
4899371c9d4SSatish Balay   xs = xx[0] / hx;
4909371c9d4SSatish Balay   xm = (xx[1] - xx[0]) / hx;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown   /*
493c4762a1bSJed Brown       Plot the  energies
494c4762a1bSJed Brown   */
4959566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetDimension(lg, 1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3)));
4969566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetColors(lg, colors + 1));
4979566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsSet(ports, 2));
498c4762a1bSJed Brown   x = hx * xs;
499c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
500c4762a1bSJed Brown     xx[0] = xx[1] = xx[2] = x;
501c4762a1bSJed 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);
502c4762a1bSJed Brown     else yy[0] = PetscRealPart(.25 * ctx->kappa * (u[i - 1] - u[i + 1]) * (u[i - 1] - u[i + 1]) * sx);
503c4762a1bSJed Brown 
504c4762a1bSJed Brown     if (ctx->cahnhillard) {
505c4762a1bSJed Brown       switch (ctx->energy) {
5069371c9d4SSatish Balay       case 1: /* double well */ yy[1] = .25 * PetscRealPart((1. - u[i] * u[i]) * (1. - u[i] * u[i])); break;
5079371c9d4SSatish Balay       case 2: /* double obstacle */ yy[1] = .5 * PetscRealPart(1. - u[i] * u[i]); break;
508c4762a1bSJed Brown       case 3: /* logarithm + double well */
509c4762a1bSJed Brown         yy[1] = .25 * PetscRealPart((1. - u[i] * u[i]) * (1. - u[i] * u[i]));
510c4762a1bSJed 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));
511c4762a1bSJed 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));
512c4762a1bSJed 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));
513c4762a1bSJed Brown         break;
514c4762a1bSJed Brown       case 4: /* logarithm + double obstacle */
515c4762a1bSJed Brown         yy[1] = .5 * theta_c * PetscRealPart(1.0 - u[i] * u[i]);
516c4762a1bSJed 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));
517c4762a1bSJed 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));
518c4762a1bSJed 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));
519c4762a1bSJed Brown         break;
5209371c9d4SSatish Balay       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "It will always be one of the values");
521c4762a1bSJed Brown       }
522c4762a1bSJed Brown     }
5239566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
524c4762a1bSJed Brown     x += hx;
525c4762a1bSJed Brown   }
5269566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetPause(draw, &pause));
5279566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, 0.0));
5289566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisSetLabels(axis, "Energy", "", ""));
5299566063dSJacob Faibussowitsch   /*  PetscCall(PetscDrawLGSetLegend(lg,legend[ctx->energy-1])); */
5309566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDraw(lg));
531c4762a1bSJed Brown 
532c4762a1bSJed Brown   /*
533c4762a1bSJed Brown       Plot the  forces
534c4762a1bSJed Brown   */
5359566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetDimension(lg, 0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3)));
5369566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetColors(lg, colors + 1));
5379566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsSet(ports, 1));
5389566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGReset(lg));
539c4762a1bSJed Brown   x   = xs * hx;
540c4762a1bSJed Brown   max = 0.;
541c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
542c4762a1bSJed Brown     xx[0] = xx[1] = xx[2] = xx[3] = x;
543c4762a1bSJed Brown     xx_netforce                   = x;
544c4762a1bSJed Brown     if (ctx->degenerate) {
545c4762a1bSJed Brown       c = (1. - u[i] * u[i]) * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx;
546c4762a1bSJed Brown       r = (1. - u[i + 1] * u[i + 1]) * (u[i] + u[i + 2] - 2.0 * u[i + 1]) * sx;
547c4762a1bSJed Brown       l = (1. - u[i - 1] * u[i - 1]) * (u[i - 2] + u[i] - 2.0 * u[i - 1]) * sx;
548c4762a1bSJed Brown     } else {
549c4762a1bSJed Brown       c = (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx;
550c4762a1bSJed Brown       r = (u[i] + u[i + 2] - 2.0 * u[i + 1]) * sx;
551c4762a1bSJed Brown       l = (u[i - 2] + u[i] - 2.0 * u[i - 1]) * sx;
552c4762a1bSJed Brown     }
553c4762a1bSJed Brown     yy[0]       = PetscRealPart(-ctx->kappa * (l + r - 2.0 * c) * sx);
554c4762a1bSJed Brown     yy_netforce = yy[0];
555c4762a1bSJed Brown     max         = PetscMax(max, PetscAbs(yy[0]));
556c4762a1bSJed Brown     if (ctx->cahnhillard) {
557c4762a1bSJed Brown       switch (ctx->energy) {
5589371c9d4SSatish Balay       case 1: /* double well */ 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); break;
5599371c9d4SSatish Balay       case 2: /* double obstacle */ yy[1] = -PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx; break;
560c4762a1bSJed Brown       case 3: /* logarithmic + double well */
561c4762a1bSJed 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);
562c4762a1bSJed Brown         if (ctx->truncation == 2) { /* quadratic */
563c4762a1bSJed 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;
564c4762a1bSJed 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;
565c4762a1bSJed 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);
566c4762a1bSJed Brown         } else { /* cubic */
567c4762a1bSJed Brown           a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol));
568c4762a1bSJed Brown           b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol);
569c4762a1bSJed 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);
570c4762a1bSJed 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);
571c4762a1bSJed 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);
572c4762a1bSJed Brown         }
573c4762a1bSJed Brown         break;
574c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
575c4762a1bSJed Brown         yy[1] = theta_c * PetscRealPart(-(u[i - 1] + u[i + 1] - 2.0 * u[i])) * sx;
576c4762a1bSJed Brown         if (ctx->truncation == 2) {
577c4762a1bSJed 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;
578c4762a1bSJed 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;
579c4762a1bSJed 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);
5809371c9d4SSatish Balay         } else {
581c4762a1bSJed Brown           a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol));
582c4762a1bSJed Brown           b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol);
583c4762a1bSJed 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);
584c4762a1bSJed 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);
585c4762a1bSJed 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);
586c4762a1bSJed Brown         }
587c4762a1bSJed Brown         break;
5889371c9d4SSatish Balay       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "It will always be one of the values");
589c4762a1bSJed Brown       }
590c4762a1bSJed Brown       if (ctx->energy < 3) {
591c4762a1bSJed Brown         max         = PetscMax(max, PetscAbs(yy[1]));
592c4762a1bSJed Brown         yy[2]       = yy[0] + yy[1];
593c4762a1bSJed Brown         yy_netforce = yy[2];
594c4762a1bSJed Brown       } else {
595c4762a1bSJed Brown         max         = PetscMax(max, PetscAbs(yy[1] + yy[2]));
596c4762a1bSJed Brown         yy[3]       = yy[0] + yy[1] + yy[2];
597c4762a1bSJed Brown         yy_netforce = yy[3];
598c4762a1bSJed Brown       }
599c4762a1bSJed Brown     }
600c4762a1bSJed Brown     if (ctx->netforce) {
6019566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGAddPoint(lg, &xx_netforce, &yy_netforce));
602c4762a1bSJed Brown     } else {
6039566063dSJacob Faibussowitsch       PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
604c4762a1bSJed Brown     }
605c4762a1bSJed Brown     x += hx;
606c4762a1bSJed Brown     /*if (max > 7200150000.0) */
607c4762a1bSJed Brown     /* printf("max very big when i = %d\n",i); */
608c4762a1bSJed Brown   }
6099566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisSetLabels(axis, "Right hand side", "", ""));
6109566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetLegend(lg, NULL));
6119566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDraw(lg));
612c4762a1bSJed Brown 
613c4762a1bSJed Brown   /*
614c4762a1bSJed Brown         Plot the solution
615c4762a1bSJed Brown   */
6169566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetDimension(lg, 1));
6179566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsSet(ports, 0));
6189566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGReset(lg));
619c4762a1bSJed Brown   x = hx * xs;
6209566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetLimits(lg, x, x + (xm - 1) * hx, -1.1, 1.1));
6219566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGSetColors(lg, colors));
622c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
623c4762a1bSJed Brown     xx[0] = x;
624c4762a1bSJed Brown     yy[0] = PetscRealPart(u[i]);
6259566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
626c4762a1bSJed Brown     x += hx;
627c4762a1bSJed Brown   }
6289566063dSJacob Faibussowitsch   PetscCall(PetscDrawAxisSetLabels(axis, "Solution", "", ""));
6299566063dSJacob Faibussowitsch   PetscCall(PetscDrawLGDraw(lg));
630c4762a1bSJed Brown 
631c4762a1bSJed Brown   /*
632c4762a1bSJed Brown       Print the  forces as arrows on the solution
633c4762a1bSJed Brown   */
634c4762a1bSJed Brown   x   = hx * xs;
635c4762a1bSJed Brown   cnt = xm / 60;
636c4762a1bSJed Brown   cnt = (!cnt) ? 1 : cnt;
637c4762a1bSJed Brown 
638c4762a1bSJed Brown   for (i = xs; i < xs + xm; i += cnt) {
639c4762a1bSJed Brown     y = yup = ydown = PetscRealPart(u[i]);
640c4762a1bSJed Brown     c               = (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx;
641c4762a1bSJed Brown     r               = (u[i] + u[i + 2] - 2.0 * u[i + 1]) * sx;
642c4762a1bSJed Brown     l               = (u[i - 2] + u[i] - 2.0 * u[i - 1]) * sx;
643c4762a1bSJed Brown     len             = -.5 * PetscRealPart(ctx->kappa * (l + r - 2.0 * c) * sx) / max;
6449566063dSJacob Faibussowitsch     PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_RED));
645c4762a1bSJed Brown     if (ctx->cahnhillard) {
646c4762a1bSJed Brown       if (len < 0.) ydown += len;
647c4762a1bSJed Brown       else yup += len;
648c4762a1bSJed Brown 
649c4762a1bSJed Brown       switch (ctx->energy) {
6509371c9d4SSatish Balay       case 1: /* double well */ 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; break;
6519371c9d4SSatish Balay       case 2: /* double obstacle */ len = -.5 * PetscRealPart(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max; break;
652c4762a1bSJed Brown       case 3: /* logarithmic + double well */
653c4762a1bSJed 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;
654c4762a1bSJed Brown         if (len < 0.) ydown += len;
655c4762a1bSJed Brown         else yup += len;
656c4762a1bSJed Brown 
657c4762a1bSJed Brown         if (ctx->truncation == 2) { /* quadratic */
658c4762a1bSJed 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;
659c4762a1bSJed 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;
660c4762a1bSJed 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);
661c4762a1bSJed Brown         } else { /* cubic */
662c4762a1bSJed Brown           a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol));
663c4762a1bSJed Brown           b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol);
664c4762a1bSJed 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);
665c4762a1bSJed 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);
666c4762a1bSJed 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);
667c4762a1bSJed Brown         }
668c4762a1bSJed Brown         y2 = len < 0 ? ydown : yup;
6699566063dSJacob Faibussowitsch         PetscCall(PetscDrawArrow(draw, x, y2, x, y2 + len2, PETSC_DRAW_PLUM));
670c4762a1bSJed Brown         break;
671c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
672c4762a1bSJed Brown         len = -.5 * theta_c * PetscRealPart(-(u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx / max);
673c4762a1bSJed Brown         if (len < 0.) ydown += len;
674c4762a1bSJed Brown         else yup += len;
675c4762a1bSJed Brown 
676c4762a1bSJed Brown         if (ctx->truncation == 2) { /* quadratic */
677c4762a1bSJed 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;
678c4762a1bSJed 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;
679c4762a1bSJed 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);
680c4762a1bSJed Brown         } else { /* cubic */
681c4762a1bSJed Brown           a = 2.0 * theta * (1.0 - 2.0 * tol) / (16.0 * tol * tol * (1.0 - tol) * (1.0 - tol));
682c4762a1bSJed Brown           b = theta / (4.0 * tol * (1.0 - tol)) - a * (1.0 - 2.0 * tol);
683c4762a1bSJed 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;
684c4762a1bSJed 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;
685c4762a1bSJed 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;
686c4762a1bSJed Brown         }
687c4762a1bSJed Brown         y2 = len < 0 ? ydown : yup;
6889566063dSJacob Faibussowitsch         PetscCall(PetscDrawArrow(draw, x, y2, x, y2 + len2, PETSC_DRAW_PLUM));
689c4762a1bSJed Brown         break;
690c4762a1bSJed Brown       }
6919566063dSJacob Faibussowitsch       PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_BLUE));
692c4762a1bSJed Brown     }
693c4762a1bSJed Brown     x += cnt * hx;
694c4762a1bSJed Brown   }
6959566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &x));
6969566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
6979566063dSJacob Faibussowitsch   PetscCall(PetscDrawStringSetSize(draw, .2, .2));
6989566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
6999566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetPause(draw, pause));
7009566063dSJacob Faibussowitsch   PetscCall(PetscDrawPause(draw));
701c4762a1bSJed Brown   PetscFunctionReturn(0);
702c4762a1bSJed Brown }
703c4762a1bSJed Brown 
7049371c9d4SSatish Balay PetscErrorCode MyDestroy(void **ptr) {
705c4762a1bSJed Brown   UserCtx *ctx = *(UserCtx **)ptr;
706c4762a1bSJed Brown 
707c4762a1bSJed Brown   PetscFunctionBegin;
7089566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsDestroy(ctx->ports));
709c4762a1bSJed Brown   PetscFunctionReturn(0);
710c4762a1bSJed Brown }
711c4762a1bSJed Brown 
712c4762a1bSJed Brown /*TEST
713c4762a1bSJed Brown 
714c4762a1bSJed Brown    test:
715c4762a1bSJed Brown      TODO: currently requires initial condition file generated by heat
716c4762a1bSJed Brown 
717c4762a1bSJed Brown TEST*/
718