Actual source code: biharmonic.c


  2: static char help[] = "Solves biharmonic equation in 1d.\n";

  4: /*
  5:   Solves the equation

  7:     u_t = - kappa  \Delta \Delta u
  8:     Periodic boundary conditions

 10: Evolve the biharmonic heat equation:
 11: ---------------
 12: ./biharmonic -ts_monitor -snes_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason  -draw_pause -2   -ts_type cn  -da_refine 5 -mymonitor

 14: Evolve with the restriction that -1 <= u <= 1; i.e. as a variational inequality
 15: ---------------
 16: ./biharmonic -ts_monitor -snes_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason  -draw_pause -2   -ts_type cn   -da_refine 5  -mymonitor

 18:    u_t =  kappa \Delta \Delta u +   6.*u*(u_x)^2 + (3*u^2 - 12) \Delta u
 19:     -1 <= u <= 1
 20:     Periodic boundary conditions

 22: Evolve the Cahn-Hillard equations: double well Initial hump shrinks then grows
 23: ---------------
 24: ./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

 26: Initial hump neither shrinks nor grows when degenerate (otherwise similar solution)

 28: ./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

 30: ./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

 32: Evolve the Cahn-Hillard equations: double obstacle
 33: ---------------
 34: ./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

 36: Evolve the Cahn-Hillard equations: logarithmic + double well (never shrinks and then grows)
 37: ---------------
 38: ./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

 40: ./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

 42: Evolve the Cahn-Hillard equations: logarithmic +  double obstacle (never shrinks, never grows)
 43: ---------------
 44: ./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: */
 47: #include <petscdm.h>
 48: #include <petscdmda.h>
 49: #include <petscts.h>
 50: #include <petscdraw.h>

 52: 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*);
 53: typedef struct {PetscBool cahnhillard;PetscBool degenerate;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta,theta_c;PetscInt truncation;PetscBool netforce; PetscDrawViewPorts *ports;} UserCtx;

 55: int main(int argc,char **argv)
 56: {
 57:   TS             ts;                 /* nonlinear solver */
 58:   Vec            x,r;                  /* solution, residual vectors */
 59:   Mat            J;                    /* Jacobian matrix */
 60:   PetscInt       steps,Mx;
 61:   DM             da;
 62:   PetscReal      dt;
 63:   PetscBool      mymonitor;
 64:   UserCtx        ctx;

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:      Initialize program
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 69:   PetscInitialize(&argc,&argv,(char*)0,help);
 70:   ctx.kappa       = 1.0;
 71:   PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
 72:   ctx.degenerate  = PETSC_FALSE;
 73:   PetscOptionsGetBool(NULL,NULL,"-degenerate",&ctx.degenerate,NULL);
 74:   ctx.cahnhillard = PETSC_FALSE;
 75:   PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
 76:   ctx.netforce    = PETSC_FALSE;
 77:   PetscOptionsGetBool(NULL,NULL,"-netforce",&ctx.netforce,NULL);
 78:   ctx.energy      = 1;
 79:   PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
 80:   ctx.tol         = 1.0e-8;
 81:   PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
 82:   ctx.theta       = .001;
 83:   ctx.theta_c     = 1.0;
 84:   PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
 85:   PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);
 86:   ctx.truncation  = 1;
 87:   PetscOptionsGetInt(NULL,NULL,"-truncation",&ctx.truncation,NULL);
 88:   PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91:      Create distributed array (DMDA) to manage parallel grid and vectors
 92:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 93:   DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,1,2,NULL,&da);
 94:   DMSetFromOptions(da);
 95:   DMSetUp(da);
 96:   DMDASetFieldName(da,0,"Biharmonic heat equation: u");
 97:   DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
 98:   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);

100:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Extract global vectors from DMDA; then duplicate for remaining
102:      vectors that are the same types
103:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104:   DMCreateGlobalVector(da,&x);
105:   VecDuplicate(x,&r);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:      Create timestepping solver context
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   TSCreate(PETSC_COMM_WORLD,&ts);
111:   TSSetDM(ts,da);
112:   TSSetProblemType(ts,TS_NONLINEAR);
113:   TSSetRHSFunction(ts,NULL,FormFunction,&ctx);
114:   DMSetMatType(da,MATAIJ);
115:   DMCreateMatrix(da,&J);
116:   TSSetRHSJacobian(ts,J,J,FormJacobian,&ctx);
117:   TSSetMaxTime(ts,.02);
118:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:      Create matrix data structure; set Jacobian evaluation routine

123:      Set Jacobian matrix data structure and default Jacobian evaluation
124:      routine. User can override with:
125:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
126:                 (unless user explicitly sets preconditioner)
127:      -snes_mf_operator : form preconditioning matrix as set by the user,
128:                          but use matrix-free approx for Jacobian-vector
129:                          products within Newton-Krylov method

131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Customize nonlinear solver
134:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   TSSetType(ts,TSCN);

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Set initial conditions
139:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   FormInitialSolution(da,x);
141:   TSSetTimeStep(ts,dt);
142:   TSSetSolution(ts,x);

144:   if (mymonitor) {
145:     ctx.ports = NULL;
146:     TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy);
147:   }

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:      Set runtime options
151:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152:   TSSetFromOptions(ts);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:      Solve nonlinear system
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   TSSolve(ts,x);
158:   TSGetStepNumber(ts,&steps);
159:   VecView(x,PETSC_VIEWER_BINARY_WORLD);

161:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162:      Free work space.  All PETSc objects should be destroyed when they
163:      are no longer needed.
164:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   MatDestroy(&J);
166:   VecDestroy(&x);
167:   VecDestroy(&r);
168:   TSDestroy(&ts);
169:   DMDestroy(&da);

171:   PetscFinalize();
172:   return 0;
173: }
174: /* ------------------------------------------------------------------- */
175: /*
176:    FormFunction - Evaluates nonlinear function, F(x).

178:    Input Parameters:
179: .  ts - the TS context
180: .  X - input vector
181: .  ptr - optional user-defined context, as set by SNESSetFunction()

183:    Output Parameter:
184: .  F - function vector
185:  */
186: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
187: {
188:   DM             da;
189:   PetscInt       i,Mx,xs,xm;
190:   PetscReal      hx,sx;
191:   PetscScalar    *x,*f,c,r,l;
192:   Vec            localX;
193:   UserCtx        *ctx = (UserCtx*)ptr;
194:   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 */

196:   TSGetDM(ts,&da);
197:   DMGetLocalVector(da,&localX);
198:   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);

200:   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);

202:   /*
203:      Scatter ghost points to local vector,using the 2-step process
204:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
205:      By placing code between these two statements, computations can be
206:      done while messages are in transition.
207:   */
208:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
209:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

211:   /*
212:      Get pointers to vector data
213:   */
214:   DMDAVecGetArrayRead(da,localX,&x);
215:   DMDAVecGetArray(da,F,&f);

217:   /*
218:      Get local grid boundaries
219:   */
220:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

222:   /*
223:      Compute function over the locally owned part of the grid
224:   */
225:   for (i=xs; i<xs+xm; i++) {
226:     if (ctx->degenerate) {
227:       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
228:       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
229:       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx;
230:     } else {
231:       c = (x[i-1] + x[i+1] - 2.0*x[i])*sx;
232:       r = (x[i] + x[i+2] - 2.0*x[i+1])*sx;
233:       l = (x[i-2] + x[i] - 2.0*x[i-1])*sx;
234:     }
235:     f[i] = -ctx->kappa*(l + r - 2.0*c)*sx;
236:     if (ctx->cahnhillard) {
237:       switch (ctx->energy) {
238:       case 1: /*  double well */
239:         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;
240:         break;
241:       case 2: /* double obstacle */
242:         f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx;
243:         break;
244:       case 3: /* logarithmic + double well */
245:         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;
246:         if (ctx->truncation==2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */
247:           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;
248:           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;
249:           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;
250:         } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */
251:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
252:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
253:           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;
254:           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;
255:           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;
256:         }
257:         break;
258:       case 4: /* logarithmic + double obstacle */
259:         f[i] += -theta_c*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
260:         if (ctx->truncation==2) { /* quadratic */
261:           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;
262:           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;
263:           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;
264:         } else { /* cubic */
265:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
266:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
267:           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;
268:           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;
269:           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;
270:         }
271:         break;
272:       }
273:     }

275:   }

277:   /*
278:      Restore vectors
279:   */
280:   DMDAVecRestoreArrayRead(da,localX,&x);
281:   DMDAVecRestoreArray(da,F,&f);
282:   DMRestoreLocalVector(da,&localX);
283:   return 0;
284: }

286: /* ------------------------------------------------------------------- */
287: /*
288:    FormJacobian - Evaluates nonlinear function's Jacobian

290: */
291: PetscErrorCode FormJacobian(TS ts,PetscReal ftime,Vec X,Mat A,Mat B,void *ptr)
292: {
293:   DM             da;
294:   PetscInt       i,Mx,xs,xm;
295:   MatStencil     row,cols[5];
296:   PetscReal      hx,sx;
297:   PetscScalar    *x,vals[5];
298:   Vec            localX;
299:   UserCtx        *ctx = (UserCtx*)ptr;

301:   TSGetDM(ts,&da);
302:   DMGetLocalVector(da,&localX);
303:   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);

305:   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);

307:   /*
308:      Scatter ghost points to local vector,using the 2-step process
309:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
310:      By placing code between these two statements, computations can be
311:      done while messages are in transition.
312:   */
313:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
314:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

316:   /*
317:      Get pointers to vector data
318:   */
319:   DMDAVecGetArrayRead(da,localX,&x);

321:   /*
322:      Get local grid boundaries
323:   */
324:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

326:   /*
327:      Compute function over the locally owned part of the grid
328:   */
329:   for (i=xs; i<xs+xm; i++) {
330:     row.i = i;
331:     if (ctx->degenerate) {
332:       /*PetscScalar c,r,l;
333:       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
334:       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
335:       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */
336:     } else {
337:       cols[0].i = i - 2; vals[0] = -ctx->kappa*sx*sx;
338:       cols[1].i = i - 1; vals[1] =  4.0*ctx->kappa*sx*sx;
339:       cols[2].i = i    ; vals[2] = -6.0*ctx->kappa*sx*sx;
340:       cols[3].i = i + 1; vals[3] =  4.0*ctx->kappa*sx*sx;
341:       cols[4].i = i + 2; vals[4] = -ctx->kappa*sx*sx;
342:     }
343:     MatSetValuesStencil(B,1,&row,5,cols,vals,INSERT_VALUES);

345:     if (ctx->cahnhillard) {
346:       switch (ctx->energy) {
347:       case 1: /* double well */
348:         /*  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; */
349:         break;
350:       case 2: /* double obstacle */
351:         /*        f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */
352:         break;
353:       case 3: /* logarithmic + double well */
354:         break;
355:       case 4: /* logarithmic + double obstacle */
356:         break;
357:       }
358:     }

360:   }

362:   /*
363:      Restore vectors
364:   */
365:   DMDAVecRestoreArrayRead(da,localX,&x);
366:   DMRestoreLocalVector(da,&localX);
367:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
368:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
369:   if (A != B) {
370:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
371:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
372:   }
373:   return 0;
374: }
375: /* ------------------------------------------------------------------- */
376: PetscErrorCode FormInitialSolution(DM da,Vec U)
377: {
378:   PetscInt          i,xs,xm,Mx,N,scale;
379:   PetscScalar       *u;
380:   PetscReal         r,hx,x;
381:   const PetscScalar *f;
382:   Vec               finesolution;
383:   PetscViewer       viewer;

385:   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);

387:   hx = 1.0/(PetscReal)Mx;

389:   /*
390:      Get pointers to vector data
391:   */
392:   DMDAVecGetArray(da,U,&u);

394:   /*
395:      Get local grid boundaries
396:   */
397:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

399:   /*
400:       Seee heat.c for how to generate InitialSolution.heat
401:   */
402:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer);
403:   VecCreate(PETSC_COMM_WORLD,&finesolution);
404:   VecLoad(finesolution,viewer);
405:   PetscViewerDestroy(&viewer);
406:   VecGetSize(finesolution,&N);
407:   scale = N/Mx;
408:   VecGetArrayRead(finesolution,&f);

410:   /*
411:      Compute function over the locally owned part of the grid
412:   */
413:   for (i=xs; i<xs+xm; i++) {
414:     x = i*hx;
415:     r = PetscSqrtReal((x-.5)*(x-.5));
416:     if (r < .125) u[i] = 1.0;
417:     else u[i] = -.5;

419:     /* With the initial condition above the method is first order in space */
420:     /* this is a smooth initial condition so the method becomes second order in space */
421:     /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
422:     u[i] = f[scale*i];
423:   }
424:   VecRestoreArrayRead(finesolution,&f);
425:   VecDestroy(&finesolution);

427:   /*
428:      Restore vectors
429:   */
430:   DMDAVecRestoreArray(da,U,&u);
431:   return 0;
432: }

434: /*
435:     This routine is not parallel
436: */
437: PetscErrorCode  MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
438: {
439:   UserCtx        *ctx = (UserCtx*)ptr;
440:   PetscDrawLG    lg;
441:   PetscScalar    *u,l,r,c;
442:   PetscInt       Mx,i,xs,xm,cnt;
443:   PetscReal      x,y,hx,pause,sx,len,max,xx[4],yy[4],xx_netforce,yy_netforce,yup,ydown,y2,len2;
444:   PetscDraw      draw;
445:   Vec            localU;
446:   DM             da;
447:   int            colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE,PETSC_DRAW_PLUM,PETSC_DRAW_BLACK};
448:   /*
449:   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"}};
450:    */
451:   PetscDrawAxis      axis;
452:   PetscDrawViewPorts *ports;
453:   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 */
454:   PetscReal          vbounds[] = {-1.1,1.1};

456:   PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);
457:   PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),800,600);
458:   TSGetDM(ts,&da);
459:   DMGetLocalVector(da,&localU);
460:   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
461:                       PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
462:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
463:   hx   = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
464:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
465:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
466:   DMDAVecGetArrayRead(da,localU,&u);

468:   PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);
469:   PetscDrawLGGetDraw(lg,&draw);
470:   PetscDrawCheckResizedWindow(draw);
471:   if (!ctx->ports) {
472:     PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);
473:   }
474:   ports = ctx->ports;
475:   PetscDrawLGGetAxis(lg,&axis);
476:   PetscDrawLGReset(lg);

478:   xx[0] = 0.0; xx[1] = 1.0; cnt = 2;
479:   PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL);
480:   xs    = xx[0]/hx; xm = (xx[1] - xx[0])/hx;

482:   /*
483:       Plot the  energies
484:   */
485:   PetscDrawLGSetDimension(lg,1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3));
486:   PetscDrawLGSetColors(lg,colors+1);
487:   PetscDrawViewPortsSet(ports,2);
488:   x    = hx*xs;
489:   for (i=xs; i<xs+xm; i++) {
490:     xx[0] = xx[1]  = xx[2] = x;
491:     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);
492:     else                 yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);

494:     if (ctx->cahnhillard) {
495:       switch (ctx->energy) {
496:       case 1: /* double well */
497:         yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
498:         break;
499:       case 2: /* double obstacle */
500:         yy[1] = .5*PetscRealPart(1. - u[i]*u[i]);
501:         break;
502:       case 3: /* logarithm + double well */
503:         yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
504:         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));
505:         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));
506:         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));
507:         break;
508:       case 4: /* logarithm + double obstacle */
509:         yy[1] = .5*theta_c*PetscRealPart(1.0-u[i]*u[i]);
510:         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));
511:         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));
512:         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));
513:         break;
514:       default:
515:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values");
516:       }
517:     }
518:     PetscDrawLGAddPoint(lg,xx,yy);
519:     x   += hx;
520:   }
521:   PetscDrawGetPause(draw,&pause);
522:   PetscDrawSetPause(draw,0.0);
523:   PetscDrawAxisSetLabels(axis,"Energy","","");
524:   /*  PetscDrawLGSetLegend(lg,legend[ctx->energy-1]); */
525:   PetscDrawLGDraw(lg);

527:   /*
528:       Plot the  forces
529:   */
530:   PetscDrawLGSetDimension(lg,0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3));
531:   PetscDrawLGSetColors(lg,colors+1);
532:   PetscDrawViewPortsSet(ports,1);
533:   PetscDrawLGReset(lg);
534:   x    = xs*hx;
535:   max  = 0.;
536:   for (i=xs; i<xs+xm; i++) {
537:     xx[0] = xx[1] = xx[2] = xx[3] = x;
538:     xx_netforce = x;
539:     if (ctx->degenerate) {
540:       c = (1. - u[i]*u[i])*(u[i-1] + u[i+1] - 2.0*u[i])*sx;
541:       r = (1. - u[i+1]*u[i+1])*(u[i] + u[i+2] - 2.0*u[i+1])*sx;
542:       l = (1. - u[i-1]*u[i-1])*(u[i-2] + u[i] - 2.0*u[i-1])*sx;
543:     } else {
544:       c = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
545:       r = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
546:       l = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
547:     }
548:     yy[0]       = PetscRealPart(-ctx->kappa*(l + r - 2.0*c)*sx);
549:     yy_netforce = yy[0];
550:     max         = PetscMax(max,PetscAbs(yy[0]));
551:     if (ctx->cahnhillard) {
552:       switch (ctx->energy) {
553:       case 1: /* double well */
554:         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);
555:         break;
556:       case 2: /* double obstacle */
557:         yy[1] = -PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx;
558:         break;
559:       case 3: /* logarithmic + double well */
560:         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);
561:         if (ctx->truncation==2) { /* quadratic */
562:           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;
563:           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;
564:           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);
565:         } else { /* cubic */
566:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
567:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
568:           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);
569:           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);
570:           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);
571:         }
572:         break;
573:       case 4: /* logarithmic + double obstacle */
574:         yy[1] = theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i]))*sx;
575:         if (ctx->truncation==2) {
576:           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;
577:           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;
578:           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);
579:         }
580:         else {
581:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
582:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
583:           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:           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:           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:         }
587:         break;
588:       default:
589:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values");
590:       }
591:       if (ctx->energy < 3) {
592:         max         = PetscMax(max,PetscAbs(yy[1]));
593:         yy[2]       = yy[0]+yy[1];
594:         yy_netforce = yy[2];
595:       } else {
596:         max         = PetscMax(max,PetscAbs(yy[1]+yy[2]));
597:         yy[3]       = yy[0]+yy[1]+yy[2];
598:         yy_netforce = yy[3];
599:       }
600:     }
601:     if (ctx->netforce) {
602:       PetscDrawLGAddPoint(lg,&xx_netforce,&yy_netforce);
603:     } else {
604:       PetscDrawLGAddPoint(lg,xx,yy);
605:     }
606:     x += hx;
607:     /*if (max > 7200150000.0) */
608:     /* printf("max very big when i = %d\n",i); */
609:   }
610:   PetscDrawAxisSetLabels(axis,"Right hand side","","");
611:   PetscDrawLGSetLegend(lg,NULL);
612:   PetscDrawLGDraw(lg);

614:   /*
615:         Plot the solution
616:   */
617:   PetscDrawLGSetDimension(lg,1);
618:   PetscDrawViewPortsSet(ports,0);
619:   PetscDrawLGReset(lg);
620:   x    = hx*xs;
621:   PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);
622:   PetscDrawLGSetColors(lg,colors);
623:   for (i=xs; i<xs+xm; i++) {
624:     xx[0] = x;
625:     yy[0] = PetscRealPart(u[i]);
626:     PetscDrawLGAddPoint(lg,xx,yy);
627:     x    += hx;
628:   }
629:   PetscDrawAxisSetLabels(axis,"Solution","","");
630:   PetscDrawLGDraw(lg);

632:   /*
633:       Print the  forces as arrows on the solution
634:   */
635:   x   = hx*xs;
636:   cnt = xm/60;
637:   cnt = (!cnt) ? 1 : cnt;

639:   for (i=xs; i<xs+xm; i += cnt) {
640:     y    = yup = ydown = PetscRealPart(u[i]);
641:     c    = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
642:     r    = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
643:     l    = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
644:     len  = -.5*PetscRealPart(ctx->kappa*(l + r - 2.0*c)*sx)/max;
645:     PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);
646:     if (ctx->cahnhillard) {
647:       if (len < 0.) ydown += len;
648:       else yup += len;

650:       switch (ctx->energy) {
651:       case 1: /* double well */
652:         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;
653:         break;
654:       case 2: /* double obstacle */
655:         len = -.5*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max;
656:         break;
657:       case 3: /* logarithmic + double well */
658:         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;
659:         if (len < 0.) ydown += len;
660:         else yup += len;

662:         if (ctx->truncation==2) { /* quadratic */
663:           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;
664:           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;
665:           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);
666:         } else { /* cubic */
667:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
668:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
669:           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);
670:           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);
671:           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);
672:         }
673:         y2   = len < 0 ? ydown : yup;
674:         PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);
675:         break;
676:       case 4: /* logarithmic + double obstacle */
677:         len = -.5*theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i])*sx/max);
678:         if (len < 0.) ydown += len;
679:         else yup += len;

681:         if (ctx->truncation==2) { /* quadratic */
682:           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;
683:           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;
684:           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);
685:         } else { /* cubic */
686:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
687:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
688:           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;
689:           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;
690:           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;
691:         }
692:         y2   = len < 0 ? ydown : yup;
693:         PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);
694:         break;
695:       }
696:       PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);
697:     }
698:     x += cnt*hx;
699:   }
700:   DMDAVecRestoreArrayRead(da,localU,&x);
701:   DMRestoreLocalVector(da,&localU);
702:   PetscDrawStringSetSize(draw,.2,.2);
703:   PetscDrawFlush(draw);
704:   PetscDrawSetPause(draw,pause);
705:   PetscDrawPause(draw);
706:   return 0;
707: }

709: PetscErrorCode  MyDestroy(void **ptr)
710: {
711:   UserCtx        *ctx = *(UserCtx**)ptr;

713:   PetscDrawViewPortsDestroy(ctx->ports);
714:   return 0;
715: }

717: /*TEST

719:    test:
720:      TODO: currently requires initial condition file generated by heat

722: TEST*/