xref: /petsc/src/ts/tutorials/phasefield/biharmonic.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*
5*c4762a1bSJed Brown   Solves the equation
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown     u_t = - kappa  \Delta \Delta u
8*c4762a1bSJed Brown     Periodic boundary conditions
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown Evolve the biharmonic heat equation:
11*c4762a1bSJed Brown ---------------
12*c4762a1bSJed 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
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown Evolve with the restriction that -1 <= u <= 1; i.e. as a variational inequality
15*c4762a1bSJed Brown ---------------
16*c4762a1bSJed 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
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown    u_t =  kappa \Delta \Delta u +   6.*u*(u_x)^2 + (3*u^2 - 12) \Delta u
19*c4762a1bSJed Brown     -1 <= u <= 1
20*c4762a1bSJed Brown     Periodic boundary conditions
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double well Initial hump shrinks then grows
23*c4762a1bSJed Brown ---------------
24*c4762a1bSJed 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
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown Initial hump neither shrinks nor grows when degenerate (otherwise similar solution)
27*c4762a1bSJed Brown 
28*c4762a1bSJed 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
29*c4762a1bSJed Brown 
30*c4762a1bSJed 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
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double obstacle
33*c4762a1bSJed Brown ---------------
34*c4762a1bSJed 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
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic + double well (never shrinks and then grows)
37*c4762a1bSJed Brown ---------------
38*c4762a1bSJed 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
39*c4762a1bSJed Brown 
40*c4762a1bSJed 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
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic +  double obstacle (never shrinks, never grows)
44*c4762a1bSJed Brown ---------------
45*c4762a1bSJed 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
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown */
51*c4762a1bSJed Brown #include <petscdm.h>
52*c4762a1bSJed Brown #include <petscdmda.h>
53*c4762a1bSJed Brown #include <petscts.h>
54*c4762a1bSJed Brown #include <petscdraw.h>
55*c4762a1bSJed Brown 
56*c4762a1bSJed 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*);
57*c4762a1bSJed Brown typedef struct {PetscBool cahnhillard;PetscBool degenerate;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta,theta_c;PetscInt truncation;PetscBool netforce; PetscDrawViewPorts *ports;} UserCtx;
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown int main(int argc,char **argv)
60*c4762a1bSJed Brown {
61*c4762a1bSJed Brown   TS             ts;                 /* nonlinear solver */
62*c4762a1bSJed Brown   Vec            x,r;                  /* solution, residual vectors */
63*c4762a1bSJed Brown   Mat            J;                    /* Jacobian matrix */
64*c4762a1bSJed Brown   PetscInt       steps,Mx;
65*c4762a1bSJed Brown   PetscErrorCode ierr;
66*c4762a1bSJed Brown   DM             da;
67*c4762a1bSJed Brown   PetscReal      dt;
68*c4762a1bSJed Brown   PetscBool      mymonitor;
69*c4762a1bSJed Brown   UserCtx        ctx;
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72*c4762a1bSJed Brown      Initialize program
73*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
75*c4762a1bSJed Brown   ctx.kappa       = 1.0;
76*c4762a1bSJed Brown   ierr            = PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ctx.degenerate  = PETSC_FALSE;
78*c4762a1bSJed Brown   ierr            = PetscOptionsGetBool(NULL,NULL,"-degenerate",&ctx.degenerate,NULL);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ctx.cahnhillard = PETSC_FALSE;
80*c4762a1bSJed Brown   ierr            = PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ctx.netforce    = PETSC_FALSE;
82*c4762a1bSJed Brown   ierr            = PetscOptionsGetBool(NULL,NULL,"-netforce",&ctx.netforce,NULL);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ctx.energy      = 1;
84*c4762a1bSJed Brown   ierr            = PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ctx.tol         = 1.0e-8;
86*c4762a1bSJed Brown   ierr            = PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ctx.theta       = .001;
88*c4762a1bSJed Brown   ctx.theta_c     = 1.0;
89*c4762a1bSJed Brown   ierr            = PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr            = PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ctx.truncation  = 1;
92*c4762a1bSJed Brown   ierr            = PetscOptionsGetInt(NULL,NULL,"-truncation",&ctx.truncation,NULL);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr            = PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
97*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,1,2,NULL,&da);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"Biharmonic heat equation: u");CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
103*c4762a1bSJed Brown   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106*c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
107*c4762a1bSJed Brown      vectors that are the same types
108*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113*c4762a1bSJed Brown      Create timestepping solver context
114*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,FormFunction,&ctx);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,J,J,FormJacobian,&ctx);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,.02);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126*c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
129*c4762a1bSJed Brown      routine. User can override with:
130*c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
131*c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
132*c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
133*c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
134*c4762a1bSJed Brown                          products within Newton-Krylov method
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138*c4762a1bSJed Brown      Customize nonlinear solver
139*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140*c4762a1bSJed Brown   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143*c4762a1bSJed Brown      Set initial conditions
144*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145*c4762a1bSJed Brown   ierr = FormInitialSolution(da,x);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   if (mymonitor) {
150*c4762a1bSJed Brown     ctx.ports = NULL;
151*c4762a1bSJed Brown     ierr      = TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy);CHKERRQ(ierr);
152*c4762a1bSJed Brown   }
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155*c4762a1bSJed Brown      Set runtime options
156*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
158*c4762a1bSJed Brown 
159*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160*c4762a1bSJed Brown      Solve nonlinear system
161*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162*c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr);
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
168*c4762a1bSJed Brown      are no longer needed.
169*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown   ierr = PetscFinalize();
177*c4762a1bSJed Brown   return ierr;
178*c4762a1bSJed Brown }
179*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
180*c4762a1bSJed Brown /*
181*c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown    Input Parameters:
184*c4762a1bSJed Brown .  ts - the TS context
185*c4762a1bSJed Brown .  X - input vector
186*c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown    Output Parameter:
189*c4762a1bSJed Brown .  F - function vector
190*c4762a1bSJed Brown  */
191*c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
192*c4762a1bSJed Brown {
193*c4762a1bSJed Brown   DM             da;
194*c4762a1bSJed Brown   PetscErrorCode ierr;
195*c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
196*c4762a1bSJed Brown   PetscReal      hx,sx;
197*c4762a1bSJed Brown   PetscScalar    *x,*f,c,r,l;
198*c4762a1bSJed Brown   Vec            localX;
199*c4762a1bSJed Brown   UserCtx        *ctx = (UserCtx*)ptr;
200*c4762a1bSJed 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 */
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown   PetscFunctionBegin;
203*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
205*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
208*c4762a1bSJed Brown 
209*c4762a1bSJed Brown   /*
210*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
211*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
212*c4762a1bSJed Brown      By placing code between these two statements, computations can be
213*c4762a1bSJed Brown      done while messages are in transition.
214*c4762a1bSJed Brown   */
215*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
217*c4762a1bSJed Brown 
218*c4762a1bSJed Brown   /*
219*c4762a1bSJed Brown      Get pointers to vector data
220*c4762a1bSJed Brown   */
221*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   /*
225*c4762a1bSJed Brown      Get local grid boundaries
226*c4762a1bSJed Brown   */
227*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   /*
230*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
231*c4762a1bSJed Brown   */
232*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
233*c4762a1bSJed Brown     if (ctx->degenerate) {
234*c4762a1bSJed Brown       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
235*c4762a1bSJed Brown       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
236*c4762a1bSJed Brown       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx;
237*c4762a1bSJed Brown     } else {
238*c4762a1bSJed Brown       c = (x[i-1] + x[i+1] - 2.0*x[i])*sx;
239*c4762a1bSJed Brown       r = (x[i] + x[i+2] - 2.0*x[i+1])*sx;
240*c4762a1bSJed Brown       l = (x[i-2] + x[i] - 2.0*x[i-1])*sx;
241*c4762a1bSJed Brown     }
242*c4762a1bSJed Brown     f[i] = -ctx->kappa*(l + r - 2.0*c)*sx;
243*c4762a1bSJed Brown     if (ctx->cahnhillard) {
244*c4762a1bSJed Brown       switch (ctx->energy) {
245*c4762a1bSJed Brown       case 1: /*  double well */
246*c4762a1bSJed 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;
247*c4762a1bSJed Brown         break;
248*c4762a1bSJed Brown       case 2: /* double obstacle */
249*c4762a1bSJed Brown         f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx;
250*c4762a1bSJed Brown         break;
251*c4762a1bSJed Brown       case 3: /* logarithmic + double well */
252*c4762a1bSJed 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;
253*c4762a1bSJed Brown         if (ctx->truncation==2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */
254*c4762a1bSJed 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;
255*c4762a1bSJed 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;
256*c4762a1bSJed 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;
257*c4762a1bSJed Brown         } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */
258*c4762a1bSJed Brown           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
259*c4762a1bSJed Brown           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
260*c4762a1bSJed 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;
261*c4762a1bSJed 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;
262*c4762a1bSJed 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;
263*c4762a1bSJed Brown         }
264*c4762a1bSJed Brown         break;
265*c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
266*c4762a1bSJed Brown         f[i] += -theta_c*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
267*c4762a1bSJed Brown         if (ctx->truncation==2) { /* quadratic */
268*c4762a1bSJed 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;
269*c4762a1bSJed 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;
270*c4762a1bSJed 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;
271*c4762a1bSJed Brown         } else { /* cubic */
272*c4762a1bSJed Brown           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
273*c4762a1bSJed Brown           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
274*c4762a1bSJed 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;
275*c4762a1bSJed 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;
276*c4762a1bSJed 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;
277*c4762a1bSJed Brown         }
278*c4762a1bSJed Brown         break;
279*c4762a1bSJed Brown       }
280*c4762a1bSJed Brown     }
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown   }
283*c4762a1bSJed Brown 
284*c4762a1bSJed Brown   /*
285*c4762a1bSJed Brown      Restore vectors
286*c4762a1bSJed Brown   */
287*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
289*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
290*c4762a1bSJed Brown   PetscFunctionReturn(0);
291*c4762a1bSJed Brown }
292*c4762a1bSJed Brown 
293*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
294*c4762a1bSJed Brown /*
295*c4762a1bSJed Brown    FormJacobian - Evaluates nonlinear function's Jacobian
296*c4762a1bSJed Brown 
297*c4762a1bSJed Brown */
298*c4762a1bSJed Brown PetscErrorCode FormJacobian(TS ts,PetscReal ftime,Vec X,Mat A,Mat B,void *ptr)
299*c4762a1bSJed Brown {
300*c4762a1bSJed Brown   DM             da;
301*c4762a1bSJed Brown   PetscErrorCode ierr;
302*c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
303*c4762a1bSJed Brown   MatStencil     row,cols[5];
304*c4762a1bSJed Brown   PetscReal      hx,sx;
305*c4762a1bSJed Brown   PetscScalar    *x,vals[5];
306*c4762a1bSJed Brown   Vec            localX;
307*c4762a1bSJed Brown   UserCtx        *ctx = (UserCtx*)ptr;
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown   PetscFunctionBegin;
310*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
311*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
312*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
313*c4762a1bSJed Brown 
314*c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
315*c4762a1bSJed Brown 
316*c4762a1bSJed Brown   /*
317*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
318*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
319*c4762a1bSJed Brown      By placing code between these two statements, computations can be
320*c4762a1bSJed Brown      done while messages are in transition.
321*c4762a1bSJed Brown   */
322*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
323*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown   /*
326*c4762a1bSJed Brown      Get pointers to vector data
327*c4762a1bSJed Brown   */
328*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
329*c4762a1bSJed Brown 
330*c4762a1bSJed Brown   /*
331*c4762a1bSJed Brown      Get local grid boundaries
332*c4762a1bSJed Brown   */
333*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown   /*
336*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
337*c4762a1bSJed Brown   */
338*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
339*c4762a1bSJed Brown     row.i = i;
340*c4762a1bSJed Brown     if (ctx->degenerate) {
341*c4762a1bSJed Brown       /*PetscScalar c,r,l;
342*c4762a1bSJed Brown       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
343*c4762a1bSJed Brown       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
344*c4762a1bSJed Brown       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */
345*c4762a1bSJed Brown     } else {
346*c4762a1bSJed Brown       cols[0].i = i - 2; vals[0] = -ctx->kappa*sx*sx;
347*c4762a1bSJed Brown       cols[1].i = i - 1; vals[1] =  4.0*ctx->kappa*sx*sx;
348*c4762a1bSJed Brown       cols[2].i = i    ; vals[2] = -6.0*ctx->kappa*sx*sx;
349*c4762a1bSJed Brown       cols[3].i = i + 1; vals[3] =  4.0*ctx->kappa*sx*sx;
350*c4762a1bSJed Brown       cols[4].i = i + 2; vals[4] = -ctx->kappa*sx*sx;
351*c4762a1bSJed Brown     }
352*c4762a1bSJed Brown     ierr = MatSetValuesStencil(B,1,&row,5,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
353*c4762a1bSJed Brown 
354*c4762a1bSJed Brown     if (ctx->cahnhillard) {
355*c4762a1bSJed Brown       switch (ctx->energy) {
356*c4762a1bSJed Brown       case 1: /* double well */
357*c4762a1bSJed 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; */
358*c4762a1bSJed Brown         break;
359*c4762a1bSJed Brown       case 2: /* double obstacle */
360*c4762a1bSJed Brown         /*        f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */
361*c4762a1bSJed Brown         break;
362*c4762a1bSJed Brown       case 3: /* logarithmic + double well */
363*c4762a1bSJed Brown         break;
364*c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
365*c4762a1bSJed Brown         break;
366*c4762a1bSJed Brown       }
367*c4762a1bSJed Brown     }
368*c4762a1bSJed Brown 
369*c4762a1bSJed Brown   }
370*c4762a1bSJed Brown 
371*c4762a1bSJed Brown   /*
372*c4762a1bSJed Brown      Restore vectors
373*c4762a1bSJed Brown   */
374*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
375*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
376*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
377*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
378*c4762a1bSJed Brown   if (A != B) {
379*c4762a1bSJed Brown     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
380*c4762a1bSJed Brown     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
381*c4762a1bSJed Brown   }
382*c4762a1bSJed Brown   PetscFunctionReturn(0);
383*c4762a1bSJed Brown }
384*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
385*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U)
386*c4762a1bSJed Brown {
387*c4762a1bSJed Brown   PetscErrorCode    ierr;
388*c4762a1bSJed Brown   PetscInt          i,xs,xm,Mx,N,scale;
389*c4762a1bSJed Brown   PetscScalar       *u;
390*c4762a1bSJed Brown   PetscReal         r,hx,x;
391*c4762a1bSJed Brown   const PetscScalar *f;
392*c4762a1bSJed Brown   Vec               finesolution;
393*c4762a1bSJed Brown   PetscViewer       viewer;
394*c4762a1bSJed Brown 
395*c4762a1bSJed Brown   PetscFunctionBegin;
396*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
397*c4762a1bSJed Brown 
398*c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx;
399*c4762a1bSJed Brown 
400*c4762a1bSJed Brown   /*
401*c4762a1bSJed Brown      Get pointers to vector data
402*c4762a1bSJed Brown   */
403*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
404*c4762a1bSJed Brown 
405*c4762a1bSJed Brown   /*
406*c4762a1bSJed Brown      Get local grid boundaries
407*c4762a1bSJed Brown   */
408*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
409*c4762a1bSJed Brown 
410*c4762a1bSJed Brown   /*
411*c4762a1bSJed Brown       Seee heat.c for how to generate InitialSolution.heat
412*c4762a1bSJed Brown   */
413*c4762a1bSJed Brown   ierr  = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer);CHKERRQ(ierr);
414*c4762a1bSJed Brown   ierr  = VecCreate(PETSC_COMM_WORLD,&finesolution);CHKERRQ(ierr);
415*c4762a1bSJed Brown   ierr  = VecLoad(finesolution,viewer);CHKERRQ(ierr);
416*c4762a1bSJed Brown   ierr  = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
417*c4762a1bSJed Brown   ierr  = VecGetSize(finesolution,&N);CHKERRQ(ierr);
418*c4762a1bSJed Brown   scale = N/Mx;
419*c4762a1bSJed Brown   ierr  = VecGetArrayRead(finesolution,&f);CHKERRQ(ierr);
420*c4762a1bSJed Brown 
421*c4762a1bSJed Brown   /*
422*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
423*c4762a1bSJed Brown   */
424*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
425*c4762a1bSJed Brown     x = i*hx;
426*c4762a1bSJed Brown     r = PetscSqrtReal((x-.5)*(x-.5));
427*c4762a1bSJed Brown     if (r < .125) u[i] = 1.0;
428*c4762a1bSJed Brown     else u[i] = -.5;
429*c4762a1bSJed Brown 
430*c4762a1bSJed Brown     /* With the initial condition above the method is first order in space */
431*c4762a1bSJed Brown     /* this is a smooth initial condition so the method becomes second order in space */
432*c4762a1bSJed Brown     /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
433*c4762a1bSJed Brown     u[i] = f[scale*i];
434*c4762a1bSJed Brown   }
435*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(finesolution,&f);CHKERRQ(ierr);
436*c4762a1bSJed Brown   ierr = VecDestroy(&finesolution);CHKERRQ(ierr);
437*c4762a1bSJed Brown 
438*c4762a1bSJed Brown   /*
439*c4762a1bSJed Brown      Restore vectors
440*c4762a1bSJed Brown   */
441*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
442*c4762a1bSJed Brown   PetscFunctionReturn(0);
443*c4762a1bSJed Brown }
444*c4762a1bSJed Brown 
445*c4762a1bSJed Brown /*
446*c4762a1bSJed Brown     This routine is not parallel
447*c4762a1bSJed Brown */
448*c4762a1bSJed Brown PetscErrorCode  MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
449*c4762a1bSJed Brown {
450*c4762a1bSJed Brown   UserCtx        *ctx = (UserCtx*)ptr;
451*c4762a1bSJed Brown   PetscDrawLG    lg;
452*c4762a1bSJed Brown   PetscErrorCode ierr;
453*c4762a1bSJed Brown   PetscScalar    *u,l,r,c;
454*c4762a1bSJed Brown   PetscInt       Mx,i,xs,xm,cnt;
455*c4762a1bSJed Brown   PetscReal      x,y,hx,pause,sx,len,max,xx[4],yy[4],xx_netforce,yy_netforce,yup,ydown,y2,len2;
456*c4762a1bSJed Brown   PetscDraw      draw;
457*c4762a1bSJed Brown   Vec            localU;
458*c4762a1bSJed Brown   DM             da;
459*c4762a1bSJed Brown   int            colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE,PETSC_DRAW_PLUM,PETSC_DRAW_BLACK};
460*c4762a1bSJed Brown   /*
461*c4762a1bSJed 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"}};
462*c4762a1bSJed Brown    */
463*c4762a1bSJed Brown   PetscDrawAxis      axis;
464*c4762a1bSJed Brown   PetscDrawViewPorts *ports;
465*c4762a1bSJed 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 */
466*c4762a1bSJed Brown   PetscReal          vbounds[] = {-1.1,1.1};
467*c4762a1bSJed Brown 
468*c4762a1bSJed Brown 
469*c4762a1bSJed Brown   PetscFunctionBegin;
470*c4762a1bSJed Brown   ierr = PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);CHKERRQ(ierr);
471*c4762a1bSJed Brown   ierr = PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),800,600);CHKERRQ(ierr);
472*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
473*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
474*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
475*c4762a1bSJed Brown                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
476*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
477*c4762a1bSJed Brown   hx   = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
478*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
479*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
480*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
481*c4762a1bSJed Brown 
482*c4762a1bSJed Brown   ierr = PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);CHKERRQ(ierr);
483*c4762a1bSJed Brown   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
484*c4762a1bSJed Brown   ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
485*c4762a1bSJed Brown   if (!ctx->ports) {
486*c4762a1bSJed Brown     ierr = PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);CHKERRQ(ierr);
487*c4762a1bSJed Brown   }
488*c4762a1bSJed Brown   ports = ctx->ports;
489*c4762a1bSJed Brown   ierr  = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr);
490*c4762a1bSJed Brown   ierr  = PetscDrawLGReset(lg);CHKERRQ(ierr);
491*c4762a1bSJed Brown 
492*c4762a1bSJed Brown   xx[0] = 0.0; xx[1] = 1.0; cnt = 2;
493*c4762a1bSJed Brown   ierr  = PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL);CHKERRQ(ierr);
494*c4762a1bSJed Brown   xs    = xx[0]/hx; xm = (xx[1] - xx[0])/hx;
495*c4762a1bSJed Brown 
496*c4762a1bSJed Brown   /*
497*c4762a1bSJed Brown       Plot the  energies
498*c4762a1bSJed Brown   */
499*c4762a1bSJed Brown   ierr = PetscDrawLGSetDimension(lg,1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3));CHKERRQ(ierr);
500*c4762a1bSJed Brown   ierr = PetscDrawLGSetColors(lg,colors+1);CHKERRQ(ierr);
501*c4762a1bSJed Brown   ierr = PetscDrawViewPortsSet(ports,2);CHKERRQ(ierr);
502*c4762a1bSJed Brown   x    = hx*xs;
503*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
504*c4762a1bSJed Brown     xx[0] = xx[1]  = xx[2] = x;
505*c4762a1bSJed 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);
506*c4762a1bSJed Brown     else                 yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);
507*c4762a1bSJed Brown 
508*c4762a1bSJed Brown     if (ctx->cahnhillard) {
509*c4762a1bSJed Brown       switch (ctx->energy) {
510*c4762a1bSJed Brown       case 1: /* double well */
511*c4762a1bSJed Brown         yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
512*c4762a1bSJed Brown         break;
513*c4762a1bSJed Brown       case 2: /* double obstacle */
514*c4762a1bSJed Brown         yy[1] = .5*PetscRealPart(1. - u[i]*u[i]);
515*c4762a1bSJed Brown         break;
516*c4762a1bSJed Brown       case 3: /* logarithm + double well */
517*c4762a1bSJed Brown         yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
518*c4762a1bSJed 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));
519*c4762a1bSJed 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));
520*c4762a1bSJed 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));
521*c4762a1bSJed Brown         break;
522*c4762a1bSJed Brown       case 4: /* logarithm + double obstacle */
523*c4762a1bSJed Brown         yy[1] = .5*theta_c*PetscRealPart(1.0-u[i]*u[i]);
524*c4762a1bSJed 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));
525*c4762a1bSJed 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));
526*c4762a1bSJed 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));
527*c4762a1bSJed Brown         break;
528*c4762a1bSJed Brown       default:
529*c4762a1bSJed Brown         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values");
530*c4762a1bSJed Brown         break;
531*c4762a1bSJed Brown       }
532*c4762a1bSJed Brown     }
533*c4762a1bSJed Brown     ierr = PetscDrawLGAddPoint(lg,xx,yy);CHKERRQ(ierr);
534*c4762a1bSJed Brown     x   += hx;
535*c4762a1bSJed Brown   }
536*c4762a1bSJed Brown   ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr);
537*c4762a1bSJed Brown   ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr);
538*c4762a1bSJed Brown   ierr = PetscDrawAxisSetLabels(axis,"Energy","","");CHKERRQ(ierr);
539*c4762a1bSJed Brown   /*  ierr = PetscDrawLGSetLegend(lg,legend[ctx->energy-1]);CHKERRQ(ierr); */
540*c4762a1bSJed Brown   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
541*c4762a1bSJed Brown 
542*c4762a1bSJed Brown   /*
543*c4762a1bSJed Brown       Plot the  forces
544*c4762a1bSJed Brown   */
545*c4762a1bSJed Brown   ierr = PetscDrawLGSetDimension(lg,0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3));CHKERRQ(ierr);
546*c4762a1bSJed Brown   ierr = PetscDrawLGSetColors(lg,colors+1);CHKERRQ(ierr);
547*c4762a1bSJed Brown   ierr = PetscDrawViewPortsSet(ports,1);CHKERRQ(ierr);
548*c4762a1bSJed Brown   ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
549*c4762a1bSJed Brown   x    = xs*hx;
550*c4762a1bSJed Brown   max  = 0.;
551*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
552*c4762a1bSJed Brown     xx[0] = xx[1] = xx[2] = xx[3] = x;
553*c4762a1bSJed Brown     xx_netforce = x;
554*c4762a1bSJed Brown     if (ctx->degenerate) {
555*c4762a1bSJed Brown       c = (1. - u[i]*u[i])*(u[i-1] + u[i+1] - 2.0*u[i])*sx;
556*c4762a1bSJed Brown       r = (1. - u[i+1]*u[i+1])*(u[i] + u[i+2] - 2.0*u[i+1])*sx;
557*c4762a1bSJed Brown       l = (1. - u[i-1]*u[i-1])*(u[i-2] + u[i] - 2.0*u[i-1])*sx;
558*c4762a1bSJed Brown     } else {
559*c4762a1bSJed Brown       c = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
560*c4762a1bSJed Brown       r = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
561*c4762a1bSJed Brown       l = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
562*c4762a1bSJed Brown     }
563*c4762a1bSJed Brown     yy[0]       = PetscRealPart(-ctx->kappa*(l + r - 2.0*c)*sx);
564*c4762a1bSJed Brown     yy_netforce = yy[0];
565*c4762a1bSJed Brown     max         = PetscMax(max,PetscAbs(yy[0]));
566*c4762a1bSJed Brown     if (ctx->cahnhillard) {
567*c4762a1bSJed Brown       switch (ctx->energy) {
568*c4762a1bSJed Brown       case 1: /* double well */
569*c4762a1bSJed 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);
570*c4762a1bSJed Brown         break;
571*c4762a1bSJed Brown       case 2: /* double obstacle */
572*c4762a1bSJed Brown         yy[1] = -PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx;
573*c4762a1bSJed Brown         break;
574*c4762a1bSJed Brown       case 3: /* logarithmic + double well */
575*c4762a1bSJed 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);
576*c4762a1bSJed Brown         if (ctx->truncation==2) { /* quadratic */
577*c4762a1bSJed 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;
578*c4762a1bSJed 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;
579*c4762a1bSJed 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);
580*c4762a1bSJed Brown         } else { /* cubic */
581*c4762a1bSJed Brown           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
582*c4762a1bSJed Brown           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
583*c4762a1bSJed 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);
584*c4762a1bSJed 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);
585*c4762a1bSJed 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);
586*c4762a1bSJed Brown         }
587*c4762a1bSJed Brown         break;
588*c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
589*c4762a1bSJed Brown         yy[1] = theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i]))*sx;
590*c4762a1bSJed Brown         if (ctx->truncation==2) {
591*c4762a1bSJed 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;
592*c4762a1bSJed 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;
593*c4762a1bSJed 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);
594*c4762a1bSJed Brown         }
595*c4762a1bSJed Brown         else {
596*c4762a1bSJed Brown           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
597*c4762a1bSJed Brown           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
598*c4762a1bSJed 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);
599*c4762a1bSJed 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);
600*c4762a1bSJed 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);
601*c4762a1bSJed Brown         }
602*c4762a1bSJed Brown         break;
603*c4762a1bSJed Brown       default:
604*c4762a1bSJed Brown         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values");
605*c4762a1bSJed Brown       break;
606*c4762a1bSJed Brown       }
607*c4762a1bSJed Brown       if (ctx->energy < 3) {
608*c4762a1bSJed Brown         max         = PetscMax(max,PetscAbs(yy[1]));
609*c4762a1bSJed Brown         yy[2]       = yy[0]+yy[1];
610*c4762a1bSJed Brown         yy_netforce = yy[2];
611*c4762a1bSJed Brown       } else {
612*c4762a1bSJed Brown         max         = PetscMax(max,PetscAbs(yy[1]+yy[2]));
613*c4762a1bSJed Brown         yy[3]       = yy[0]+yy[1]+yy[2];
614*c4762a1bSJed Brown         yy_netforce = yy[3];
615*c4762a1bSJed Brown       }
616*c4762a1bSJed Brown     }
617*c4762a1bSJed Brown     if (ctx->netforce) {
618*c4762a1bSJed Brown       ierr = PetscDrawLGAddPoint(lg,&xx_netforce,&yy_netforce);CHKERRQ(ierr);
619*c4762a1bSJed Brown     } else {
620*c4762a1bSJed Brown       ierr = PetscDrawLGAddPoint(lg,xx,yy);CHKERRQ(ierr);
621*c4762a1bSJed Brown     }
622*c4762a1bSJed Brown     x += hx;
623*c4762a1bSJed Brown     /*if (max > 7200150000.0) */
624*c4762a1bSJed Brown     /* printf("max very big when i = %d\n",i); */
625*c4762a1bSJed Brown   }
626*c4762a1bSJed Brown   ierr = PetscDrawAxisSetLabels(axis,"Right hand side","","");CHKERRQ(ierr);
627*c4762a1bSJed Brown   ierr = PetscDrawLGSetLegend(lg,NULL);CHKERRQ(ierr);
628*c4762a1bSJed Brown   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
629*c4762a1bSJed Brown 
630*c4762a1bSJed Brown   /*
631*c4762a1bSJed Brown         Plot the solution
632*c4762a1bSJed Brown   */
633*c4762a1bSJed Brown   ierr = PetscDrawLGSetDimension(lg,1);CHKERRQ(ierr);
634*c4762a1bSJed Brown   ierr = PetscDrawViewPortsSet(ports,0);CHKERRQ(ierr);
635*c4762a1bSJed Brown   ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
636*c4762a1bSJed Brown   x    = hx*xs;
637*c4762a1bSJed Brown   ierr = PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);CHKERRQ(ierr);
638*c4762a1bSJed Brown   ierr = PetscDrawLGSetColors(lg,colors);CHKERRQ(ierr);
639*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
640*c4762a1bSJed Brown     xx[0] = x;
641*c4762a1bSJed Brown     yy[0] = PetscRealPart(u[i]);
642*c4762a1bSJed Brown     ierr  = PetscDrawLGAddPoint(lg,xx,yy);CHKERRQ(ierr);
643*c4762a1bSJed Brown     x    += hx;
644*c4762a1bSJed Brown   }
645*c4762a1bSJed Brown   ierr = PetscDrawAxisSetLabels(axis,"Solution","","");CHKERRQ(ierr);
646*c4762a1bSJed Brown   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
647*c4762a1bSJed Brown 
648*c4762a1bSJed Brown   /*
649*c4762a1bSJed Brown       Print the  forces as arrows on the solution
650*c4762a1bSJed Brown   */
651*c4762a1bSJed Brown   x   = hx*xs;
652*c4762a1bSJed Brown   cnt = xm/60;
653*c4762a1bSJed Brown   cnt = (!cnt) ? 1 : cnt;
654*c4762a1bSJed Brown 
655*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i += cnt) {
656*c4762a1bSJed Brown     y    = yup = ydown = PetscRealPart(u[i]);
657*c4762a1bSJed Brown     c    = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
658*c4762a1bSJed Brown     r    = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
659*c4762a1bSJed Brown     l    = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
660*c4762a1bSJed Brown     len  = -.5*PetscRealPart(ctx->kappa*(l + r - 2.0*c)*sx)/max;
661*c4762a1bSJed Brown     ierr = PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);CHKERRQ(ierr);
662*c4762a1bSJed Brown     if (ctx->cahnhillard) {
663*c4762a1bSJed Brown       if (len < 0.) ydown += len;
664*c4762a1bSJed Brown       else yup += len;
665*c4762a1bSJed Brown 
666*c4762a1bSJed Brown       switch (ctx->energy) {
667*c4762a1bSJed Brown       case 1: /* double well */
668*c4762a1bSJed 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;
669*c4762a1bSJed Brown         break;
670*c4762a1bSJed Brown       case 2: /* double obstacle */
671*c4762a1bSJed Brown         len = -.5*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max;
672*c4762a1bSJed Brown         break;
673*c4762a1bSJed Brown       case 3: /* logarithmic + double well */
674*c4762a1bSJed 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;
675*c4762a1bSJed Brown         if (len < 0.) ydown += len;
676*c4762a1bSJed Brown         else yup += len;
677*c4762a1bSJed Brown 
678*c4762a1bSJed Brown         if (ctx->truncation==2) { /* quadratic */
679*c4762a1bSJed 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;
680*c4762a1bSJed 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;
681*c4762a1bSJed 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);
682*c4762a1bSJed Brown         } else { /* cubic */
683*c4762a1bSJed Brown           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
684*c4762a1bSJed Brown           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
685*c4762a1bSJed 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);
686*c4762a1bSJed 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);
687*c4762a1bSJed 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);
688*c4762a1bSJed Brown         }
689*c4762a1bSJed Brown         y2   = len < 0 ? ydown : yup;
690*c4762a1bSJed Brown         ierr = PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);CHKERRQ(ierr);
691*c4762a1bSJed Brown         break;
692*c4762a1bSJed Brown       case 4: /* logarithmic + double obstacle */
693*c4762a1bSJed Brown         len = -.5*theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i])*sx/max);
694*c4762a1bSJed Brown         if (len < 0.) ydown += len;
695*c4762a1bSJed Brown         else yup += len;
696*c4762a1bSJed Brown 
697*c4762a1bSJed Brown         if (ctx->truncation==2) { /* quadratic */
698*c4762a1bSJed 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;
699*c4762a1bSJed 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;
700*c4762a1bSJed 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);
701*c4762a1bSJed Brown         } else { /* cubic */
702*c4762a1bSJed Brown           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
703*c4762a1bSJed Brown           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
704*c4762a1bSJed 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;
705*c4762a1bSJed 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;
706*c4762a1bSJed 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;
707*c4762a1bSJed Brown         }
708*c4762a1bSJed Brown         y2   = len < 0 ? ydown : yup;
709*c4762a1bSJed Brown         ierr = PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);CHKERRQ(ierr);
710*c4762a1bSJed Brown         break;
711*c4762a1bSJed Brown       }
712*c4762a1bSJed Brown       ierr = PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);CHKERRQ(ierr);
713*c4762a1bSJed Brown     }
714*c4762a1bSJed Brown     x += cnt*hx;
715*c4762a1bSJed Brown   }
716*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&x);CHKERRQ(ierr);
717*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
718*c4762a1bSJed Brown   ierr = PetscDrawStringSetSize(draw,.2,.2);CHKERRQ(ierr);
719*c4762a1bSJed Brown   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
720*c4762a1bSJed Brown   ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr);
721*c4762a1bSJed Brown   ierr = PetscDrawPause(draw);CHKERRQ(ierr);
722*c4762a1bSJed Brown   PetscFunctionReturn(0);
723*c4762a1bSJed Brown }
724*c4762a1bSJed Brown 
725*c4762a1bSJed Brown PetscErrorCode  MyDestroy(void **ptr)
726*c4762a1bSJed Brown {
727*c4762a1bSJed Brown   UserCtx        *ctx = *(UserCtx**)ptr;
728*c4762a1bSJed Brown   PetscErrorCode ierr;
729*c4762a1bSJed Brown 
730*c4762a1bSJed Brown   PetscFunctionBegin;
731*c4762a1bSJed Brown   ierr = PetscDrawViewPortsDestroy(ctx->ports);CHKERRQ(ierr);
732*c4762a1bSJed Brown   PetscFunctionReturn(0);
733*c4762a1bSJed Brown }
734*c4762a1bSJed Brown 
735*c4762a1bSJed Brown 
736*c4762a1bSJed Brown 
737*c4762a1bSJed Brown /*TEST
738*c4762a1bSJed Brown 
739*c4762a1bSJed Brown    test:
740*c4762a1bSJed Brown      TODO: currently requires initial condition file generated by heat
741*c4762a1bSJed Brown 
742*c4762a1bSJed Brown TEST*/
743