Actual source code: biharmonic.c

petsc-3.8.3 2017-12-09
Report Typos and Errors

  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_vi_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason  -wait   -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_vi_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason  -wait   -ts_type cn   -da_refine 5 -vi -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_vi_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason   -wait   -ts_type cn    -da_refine 6 -vi  -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_vi_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason   -wait   -ts_type cn    -da_refine 6 -vi  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -degenerate -ts_monitor_draw_solution --mymonitor

 30: ./biharmonic -ts_monitor -snes_vi_monitor   -pc_type lu  -draw_pause .1 -snes_converged_reason   -wait   -ts_type cn    -da_refine 6 -vi  -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_vi_monitor  -pc_type lu  -draw_pause .1 -snes_converged_reason   -wait   -ts_type cn    -da_refine 5 -vi  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -energy 2 -snes_linesearch_monitor   -vi -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_vi_monitor  -pc_type lu  --snes_converged_reason  -wait   -ts_type cn    -da_refine 5 -vi  -kappa .0001 -ts_dt 5.96046e-06 -cahn-hillard -energy 3 -snes_linesearch_monitor -theta .00000001  -vi  -ts_monitor_draw_solution --ts_final_time 1. -mymonitor

 40: ./biharmonic -ts_monitor -snes_vi_monitor  -pc_type lu  --snes_converged_reason  -wait   -ts_type cn    -da_refine 5 -vi  -kappa .0001 -ts_dt 5.96046e-06 -cahn-hillard -energy 3 -snes_linesearch_monitor -theta .00000001  -vi  -ts_monitor_draw_solution --ts_final_time 1. -degenerate -mymonitor


 43: Evolve the Cahn-Hillard equations: logarithmic +  double obstacle (never shrinks, never grows)
 44: ---------------
 45: ./biharmonic -ts_monitor -snes_vi_monitor  -pc_type lu  --snes_converged_reason  -wait   -ts_type cn    -da_refine 5 -vi  -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -energy 4 -snes_linesearch_monitor -theta .00000001  -vi -ts_monitor_draw_solution --mymonitor




 50: */
 51:  #include <petscdm.h>
 52:  #include <petscdmda.h>
 53:  #include <petscts.h>
 54:  #include <petscdraw.h>

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

 59: int main(int argc,char **argv)
 60: {
 61:   TS             ts;                 /* nonlinear solver */
 62:   Vec            x,r;                  /* solution, residual vectors */
 63:   Mat            J;                    /* Jacobian matrix */
 64:   PetscInt       steps,Mx;
 66:   DM             da;
 67:   PetscReal      dt;
 68:   PetscReal      vbounds[] = {-1.1,1.1};
 69:   PetscBool      wait,vi = PETSC_FALSE,mymonitor;
 70:   Vec            ul,uh;
 71:   UserCtx        ctx;

 73:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 74:      Initialize program
 75:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 76:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 77:   ctx.kappa       = 1.0;
 78:   PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
 79:   ctx.degenerate  = PETSC_FALSE;
 80:   PetscOptionsGetBool(NULL,NULL,"-degenerate",&ctx.degenerate,NULL);
 81:   ctx.cahnhillard = PETSC_FALSE;
 82:   PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
 83:   PetscOptionsGetBool(NULL,NULL,"-vi",&vi,NULL);
 84:   ctx.netforce    = PETSC_FALSE;
 85:   PetscOptionsGetBool(NULL,NULL,"-netforce",&ctx.netforce,NULL);
 86:   ctx.energy      = 1;
 87:   PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
 88:   ctx.tol         = 1.0e-8;
 89:   PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
 90:   ctx.theta       = .001;
 91:   ctx.theta_c     = 1.0;
 92:   PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
 93:   PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);
 94:   ctx.truncation  = 1;
 95:   PetscOptionsGetInt(NULL,NULL,"-truncation",&ctx.truncation,NULL);
 96:   PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);
 97:   PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);
 98:   PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),800,600);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Create distributed array (DMDA) to manage parallel grid and vectors
102:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,1,2,NULL,&da);
104:   DMSetFromOptions(da);
105:   DMSetUp(da);
106:   DMDASetFieldName(da,0,"Biharmonic heat equation: u");
107:   DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);
108:   dt   = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx);

110:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:      Extract global vectors from DMDA; then duplicate for remaining
112:      vectors that are the same types
113:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114:   DMCreateGlobalVector(da,&x);
115:   VecDuplicate(x,&r);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:      Create timestepping solver context
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120:   TSCreate(PETSC_COMM_WORLD,&ts);
121:   TSSetDM(ts,da);
122:   TSSetProblemType(ts,TS_NONLINEAR);
123:   TSSetRHSFunction(ts,NULL,FormFunction,&ctx);
124:   DMSetMatType(da,MATAIJ);
125:   DMCreateMatrix(da,&J);
126:   TSSetRHSJacobian(ts,J,J,FormJacobian,&ctx);
127:   TSSetMaxTime(ts,.02);
128:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:      Create matrix data structure; set Jacobian evaluation routine

133:      Set Jacobian matrix data structure and default Jacobian evaluation
134:      routine. User can override with:
135:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
136:                 (unless user explicitly sets preconditioner)
137:      -snes_mf_operator : form preconditioning matrix as set by the user,
138:                          but use matrix-free approx for Jacobian-vector
139:                          products within Newton-Krylov method

141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: #if defined(f00)
143:   {
144:     SNES snes;
145:     DMCreateColoring(da,IS_COLORING_GLOBAL,&iscoloring);
146:     MatFDColoringCreate(J,iscoloring,&matfdcoloring);
147:     ISColoringDestroy(&iscoloring);
148:     MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
149:     MatFDColoringSetFromOptions(matfdcoloring);
150:      MatFDColoringSetUp(J,iscoloring,matfdcoloring);
151:     TSGetSNES(ts,&snes);
152:     SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
153:   }
154: #endif

156:   if (vi) {
157:     VecDuplicate(x,&ul);
158:     VecDuplicate(x,&uh);
159:     VecSet(ul,-1.0);
160:     VecSet(uh,1.0);
161:     TSVISetVariableBounds(ts,ul,uh);
162:   }

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Customize nonlinear solver
166:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167:   TSSetType(ts,TSCN);

169:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170:      Set initial conditions
171:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172:   FormInitialSolution(da,x);
173:   TSSetTimeStep(ts,dt);
174:   TSSetSolution(ts,x);

176:   if (mymonitor) {
177:     ctx.ports = NULL;
178:     TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy);
179:   }

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:      Set runtime options
183:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184:   TSSetFromOptions(ts);

186:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187:      Solve nonlinear system
188:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189:   TSSolve(ts,x);
190:   wait = PETSC_FALSE;
191:   PetscOptionsGetBool(NULL,NULL,"-wait",&wait,NULL);
192:   if (wait) {
193:     PetscSleep(-1);
194:   }
195:   TSGetStepNumber(ts,&steps);
196:   VecView(x,PETSC_VIEWER_BINARY_WORLD);

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:      Free work space.  All PETSc objects should be destroyed when they
200:      are no longer needed.
201:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202:   if (vi) {
203:     VecDestroy(&ul);
204:     VecDestroy(&uh);
205:   }
206:   MatDestroy(&J);
207: #if defined(f00)
208:   MatFDColoringDestroy(&matfdcoloring);
209: #endif
210:   VecDestroy(&x);
211:   VecDestroy(&r);
212:   TSDestroy(&ts);
213:   DMDestroy(&da);

215:   PetscFinalize();
216:   return ierr;
217: }
218: /* ------------------------------------------------------------------- */
219: /*
220:    FormFunction - Evaluates nonlinear function, F(x).

222:    Input Parameters:
223: .  ts - the TS context
224: .  X - input vector
225: .  ptr - optional user-defined context, as set by SNESSetFunction()

227:    Output Parameter:
228: .  F - function vector
229:  */
230: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
231: {
232:   DM             da;
234:   PetscInt       i,Mx,xs,xm;
235:   PetscReal      hx,sx;
236:   PetscScalar    *x,*f,c,r,l;
237:   Vec            localX;
238:   UserCtx        *ctx = (UserCtx*)ptr;
239:   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 */

242:   TSGetDM(ts,&da);
243:   DMGetLocalVector(da,&localX);
244:   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);

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

248:   /*
249:      Scatter ghost points to local vector,using the 2-step process
250:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
251:      By placing code between these two statements, computations can be
252:      done while messages are in transition.
253:   */
254:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
255:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

257:   /*
258:      Get pointers to vector data
259:   */
260:   DMDAVecGetArrayRead(da,localX,&x);
261:   DMDAVecGetArray(da,F,&f);

263:   /*
264:      Get local grid boundaries
265:   */
266:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

268:   /*
269:      Compute function over the locally owned part of the grid
270:   */
271:   for (i=xs; i<xs+xm; i++) {
272:     if (ctx->degenerate) {
273:       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
274:       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
275:       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx;
276:     } else {
277:       c = (x[i-1] + x[i+1] - 2.0*x[i])*sx;
278:       r = (x[i] + x[i+2] - 2.0*x[i+1])*sx;
279:       l = (x[i-2] + x[i] - 2.0*x[i-1])*sx;
280:     }
281:     f[i] = -ctx->kappa*(l + r - 2.0*c)*sx;
282:     if (ctx->cahnhillard) {
283:       switch (ctx->energy) {
284:       case 1: /*  double well */
285:         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;
286:         break;
287:       case 2: /* double obstacle */
288:         f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx;
289:         break;
290:       case 3: /* logarithmic + double well */
291:         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;
292:         if (ctx->truncation==2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */
293:           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;
294:           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;
295:           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;
296:         } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */
297:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
298:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
299:           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;
300:           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;
301:           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;
302:         }
303:         break;
304:       case 4: /* logarithmic + double obstacle */
305:         f[i] += -theta_c*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
306:         if (ctx->truncation==2) { /* quadratic */
307:           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;
308:           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;
309:           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;
310:         } else { /* cubic */
311:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
312:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
313:           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;
314:           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;
315:           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;
316:         }
317:         break;
318:       }
319:     }

321:   }

323:   /*
324:      Restore vectors
325:   */
326:   DMDAVecRestoreArrayRead(da,localX,&x);
327:   DMDAVecRestoreArray(da,F,&f);
328:   DMRestoreLocalVector(da,&localX);
329:   return(0);
330: }

332: /* ------------------------------------------------------------------- */
333: /*
334:    FormJacobian - Evaluates nonlinear function's Jacobian

336: */
337: PetscErrorCode FormJacobian(TS ts,PetscReal ftime,Vec X,Mat A,Mat B,void *ptr)
338: {
339:   DM             da;
341:   PetscInt       i,Mx,xs,xm;
342:   MatStencil     row,cols[5];
343:   PetscReal      hx,sx;
344:   PetscScalar    *x,vals[5];
345:   Vec            localX;
346:   UserCtx        *ctx = (UserCtx*)ptr;

349:   TSGetDM(ts,&da);
350:   DMGetLocalVector(da,&localX);
351:   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);

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

355:   /*
356:      Scatter ghost points to local vector,using the 2-step process
357:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
358:      By placing code between these two statements, computations can be
359:      done while messages are in transition.
360:   */
361:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
362:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

364:   /*
365:      Get pointers to vector data
366:   */
367:   DMDAVecGetArrayRead(da,localX,&x);

369:   /*
370:      Get local grid boundaries
371:   */
372:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

374:   /*
375:      Compute function over the locally owned part of the grid
376:   */
377:   for (i=xs; i<xs+xm; i++) {
378:     row.i = i;
379:     if (ctx->degenerate) {
380:       /*PetscScalar c,r,l;
381:       c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
382:       r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
383:       l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */
384:     } else {
385:       cols[0].i = i - 2; vals[0] = -ctx->kappa*sx*sx;
386:       cols[1].i = i - 1; vals[1] =  4.0*ctx->kappa*sx*sx;
387:       cols[2].i = i    ; vals[2] = -6.0*ctx->kappa*sx*sx;
388:       cols[3].i = i + 1; vals[3] =  4.0*ctx->kappa*sx*sx;
389:       cols[4].i = i + 2; vals[4] = -ctx->kappa*sx*sx;
390:     }
391:     MatSetValuesStencil(B,1,&row,5,cols,vals,INSERT_VALUES);

393:     if (ctx->cahnhillard) {
394:       switch (ctx->energy) {
395:       case 1: /* double well */
396:         /*  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; */
397:         break;
398:       case 2: /* double obstacle */
399:         /*        f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */
400:         break;
401:       case 3: /* logarithmic + double well */
402:         break;
403:       case 4: /* logarithmic + double obstacle */
404:         break;
405:       }
406:     }

408:   }

410:   /*
411:      Restore vectors
412:   */
413:   DMDAVecRestoreArrayRead(da,localX,&x);
414:   DMRestoreLocalVector(da,&localX);
415:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
416:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
417:   if (A != B) {
418:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
419:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
420:   }
421:   return(0);
422: }
423: /* ------------------------------------------------------------------- */
424: PetscErrorCode FormInitialSolution(DM da,Vec U)
425: {
426:   PetscErrorCode    ierr;
427:   PetscInt          i,xs,xm,Mx,N,scale;
428:   PetscScalar       *u;
429:   PetscReal         r,hx,x;
430:   const PetscScalar *f;
431:   Vec               finesolution;
432:   PetscViewer       viewer;

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

437:   hx = 1.0/(PetscReal)Mx;

439:   /*
440:      Get pointers to vector data
441:   */
442:   DMDAVecGetArray(da,U,&u);

444:   /*
445:      Get local grid boundaries
446:   */
447:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

449:   /*  InitialSolution.biharmonic is obtained by running
450:        ./heat -square_initial -ts_monitor -snes_monitor  -pc_type lu   -snes_converged_reason    -ts_type cn  -da_refine 9 -ts_final_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25  -ts_max_steps 30
451:        After changing the initial grid spacing to 10 and the stencil width to 2 in the DMDA create.
452:     */
453:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.biharmonic",FILE_MODE_READ,&viewer);
454:   VecCreate(PETSC_COMM_WORLD,&finesolution);
455:   VecLoad(finesolution,viewer);
456:   PetscViewerDestroy(&viewer);
457:   VecGetSize(finesolution,&N);
458:   scale = N/Mx;
459:   VecGetArrayRead(finesolution,&f);

461:   /*
462:      Compute function over the locally owned part of the grid
463:   */
464:   for (i=xs; i<xs+xm; i++) {
465:     x = i*hx;
466:     r = PetscSqrtReal((x-.5)*(x-.5));
467:     if (r < .125) u[i] = 1.0;
468:     else u[i] = -.5;

470:     /* With the initial condition above the method is first order in space */
471:     /* this is a smooth initial condition so the method becomes second order in space */
472:     /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
473:     u[i] = f[scale*i];
474:   }
475:   VecRestoreArrayRead(finesolution,&f);
476:   VecDestroy(&finesolution);

478:   /*
479:      Restore vectors
480:   */
481:   DMDAVecRestoreArray(da,U,&u);
482:   return(0);
483: }

485: /*
486:     This routine is not parallel
487: */
488: PetscErrorCode  MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
489: {
490:   UserCtx        *ctx = (UserCtx*)ptr;
491:   PetscDrawLG    lg;
493:   PetscScalar    *u,l,r,c;
494:   PetscInt       Mx,i,xs,xm,cnt;
495:   PetscReal      x,y,hx,pause,sx,len,max,xx[4],yy[4],xx_netforce,yy_netforce,yup,ydown,y2,len2;
496:   PetscDraw      draw;
497:   Vec            localU;
498:   DM             da;
499:   int            colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE,PETSC_DRAW_PLUM,PETSC_DRAW_BLACK};
500:   /*
501:   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"}};
502:    */
503:   PetscDrawAxis      axis;
504:   PetscDrawViewPorts *ports;
505:   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 */


509:   TSGetDM(ts,&da);
510:   DMGetLocalVector(da,&localU);
511:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
512:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
513:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
514:   hx   = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
515:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
516:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
517:   DMDAVecGetArrayRead(da,localU,&u);

519:   PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);
520:   PetscDrawLGGetDraw(lg,&draw);
521:   PetscDrawCheckResizedWindow(draw);
522:   if (!ctx->ports) {
523:     PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);
524:   }
525:   ports = ctx->ports;
526:   PetscDrawLGGetAxis(lg,&axis);
527:   PetscDrawLGReset(lg);

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

533:   /*
534:       Plot the  energies
535:   */
536:   PetscDrawLGSetDimension(lg,1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3));
537:   PetscDrawLGSetColors(lg,colors+1);
538:   PetscDrawViewPortsSet(ports,2);
539:   x    = hx*xs;
540:   for (i=xs; i<xs+xm; i++) {
541:     xx[0] = xx[1]  = xx[2] = x;
542:     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);
543:     else                 yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);

545:     if (ctx->cahnhillard) {
546:       switch (ctx->energy) {
547:       case 1: /* double well */
548:         yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
549:         break;
550:       case 2: /* double obstacle */
551:         yy[1] = .5*PetscRealPart(1. - u[i]*u[i]);
552:         break;
553:       case 3: /* logarithm + double well */
554:         yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
555:         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));
556:         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));
557:         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));
558:         break;
559:       case 4: /* logarithm + double obstacle */
560:         yy[1] = .5*theta_c*PetscRealPart(1.0-u[i]*u[i]);
561:         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));
562:         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));
563:         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));
564:         break;
565:       default:
566:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values");
567:         break;
568:       }
569:     }
570:     PetscDrawLGAddPoint(lg,xx,yy);
571:     x   += hx;
572:   }
573:   PetscDrawGetPause(draw,&pause);
574:   PetscDrawSetPause(draw,0.0);
575:   PetscDrawAxisSetLabels(axis,"Energy","","");
576:   /*  PetscDrawLGSetLegend(lg,legend[ctx->energy-1]); */
577:   PetscDrawLGDraw(lg);

579:   /*
580:       Plot the  forces
581:   */
582:   PetscDrawLGSetDimension(lg,0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3));
583:   PetscDrawLGSetColors(lg,colors+1);
584:   PetscDrawViewPortsSet(ports,1);
585:   PetscDrawLGReset(lg);
586:   x    = xs*hx;
587:   max  = 0.;
588:   for (i=xs; i<xs+xm; i++) {
589:     xx[0] = xx[1] = xx[2] = xx[3] = x;
590:     xx_netforce = x;
591:     if (ctx->degenerate) {
592:       c = (1. - u[i]*u[i])*(u[i-1] + u[i+1] - 2.0*u[i])*sx;
593:       r = (1. - u[i+1]*u[i+1])*(u[i] + u[i+2] - 2.0*u[i+1])*sx;
594:       l = (1. - u[i-1]*u[i-1])*(u[i-2] + u[i] - 2.0*u[i-1])*sx;
595:     } else {
596:       c = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
597:       r = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
598:       l = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
599:     }
600:     yy[0]       = PetscRealPart(-ctx->kappa*(l + r - 2.0*c)*sx);
601:     yy_netforce = yy[0];
602:     max         = PetscMax(max,PetscAbs(yy[0]));
603:     if (ctx->cahnhillard) {
604:       switch (ctx->energy) {
605:       case 1: /* double well */
606:         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);
607:         break;
608:       case 2: /* double obstacle */
609:         yy[1] = -PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx;
610:         break;
611:       case 3: /* logarithmic + double well */
612:         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);
613:         if (ctx->truncation==2) { /* quadratic */
614:           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;
615:           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;
616:           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);
617:         } else { /* cubic */
618:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
619:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
620:           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);
621:           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);
622:           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);
623:         }
624:         break;
625:       case 4: /* logarithmic + double obstacle */
626:         yy[1] = theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i]))*sx;
627:         if (ctx->truncation==2) {
628:           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;
629:           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;
630:           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);
631:         }
632:         else {
633:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
634:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
635:           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);
636:           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);
637:           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);
638:         }
639:         break;
640:       default:
641:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values");
642:       break;
643:       }
644:       if (ctx->energy < 3) {
645:         max         = PetscMax(max,PetscAbs(yy[1]));
646:         yy[2]       = yy[0]+yy[1];
647:         yy_netforce = yy[2];
648:       } else {
649:         max         = PetscMax(max,PetscAbs(yy[1]+yy[2]));
650:         yy[3]       = yy[0]+yy[1]+yy[2];
651:         yy_netforce = yy[3];
652:       }
653:     }
654:     if (ctx->netforce) {
655:       PetscDrawLGAddPoint(lg,&xx_netforce,&yy_netforce);
656:     } else {
657:       PetscDrawLGAddPoint(lg,xx,yy);
658:     }
659:     x += hx;
660:     /*if (max > 7200150000.0) */
661:     /* printf("max very big when i = %d\n",i); */
662:   }
663:   PetscDrawAxisSetLabels(axis,"Right hand side","","");
664:   PetscDrawLGSetLegend(lg,NULL);
665:   PetscDrawLGDraw(lg);

667:   /*
668:         Plot the solution
669:   */
670:   PetscDrawLGSetDimension(lg,1);
671:   PetscDrawViewPortsSet(ports,0);
672:   PetscDrawLGReset(lg);
673:   x    = hx*xs;
674:   PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);
675:   PetscDrawLGSetColors(lg,colors);
676:   for (i=xs; i<xs+xm; i++) {
677:     xx[0] = x;
678:     yy[0] = PetscRealPart(u[i]);
679:     PetscDrawLGAddPoint(lg,xx,yy);
680:     x    += hx;
681:   }
682:   PetscDrawAxisSetLabels(axis,"Solution","","");
683:   PetscDrawLGDraw(lg);

685:   /*
686:       Print the  forces as arrows on the solution
687:   */
688:   x   = hx*xs;
689:   cnt = xm/60;
690:   cnt = (!cnt) ? 1 : cnt;

692:   for (i=xs; i<xs+xm; i += cnt) {
693:     y    = yup = ydown = PetscRealPart(u[i]);
694:     c    = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
695:     r    = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
696:     l    = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
697:     len  = -.5*PetscRealPart(ctx->kappa*(l + r - 2.0*c)*sx)/max;
698:     PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);
699:     if (ctx->cahnhillard) {
700:       if (len < 0.) ydown += len;
701:       else yup += len;

703:       switch (ctx->energy) {
704:       case 1: /* double well */
705:         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;
706:         break;
707:       case 2: /* double obstacle */
708:         len = -.5*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max;
709:         break;
710:       case 3: /* logarithmic + double well */
711:         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;
712:         if (len < 0.) ydown += len;
713:         else yup += len;

715:         if (ctx->truncation==2) { /* quadratic */
716:           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;
717:           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;
718:           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);
719:         } else { /* cubic */
720:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
721:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
722:           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);
723:           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);
724:           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);
725:         }
726:         y2   = len < 0 ? ydown : yup;
727:         PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);
728:         break;
729:       case 4: /* logarithmic + double obstacle */
730:         len = -.5*theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i])*sx/max);
731:         if (len < 0.) ydown += len;
732:         else yup += len;

734:         if (ctx->truncation==2) { /* quadratic */
735:           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;
736:           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;
737:           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);
738:         } else { /* cubic */
739:           a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
740:           b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
741:           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;
742:           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;
743:           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;
744:         }
745:         y2   = len < 0 ? ydown : yup;
746:         PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);
747:         break;
748:       }
749:       PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);
750:     }
751:     x += cnt*hx;
752:   }
753:   DMDAVecRestoreArrayRead(da,localU,&x);
754:   DMRestoreLocalVector(da,&localU);
755:   PetscDrawStringSetSize(draw,.2,.2);
756:   PetscDrawFlush(draw);
757:   PetscDrawSetPause(draw,pause);
758:   PetscDrawPause(draw);
759:   return(0);
760: }

762: PetscErrorCode  MyDestroy(void **ptr)
763: {
764:   UserCtx        *ctx = *(UserCtx**)ptr;

768:   PetscDrawViewPortsDestroy(ctx->ports);
769:   return(0);
770: }