Actual source code: biharmonic.c
petsc-3.7.1 2016-05-15
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;
61: int main(int argc,char **argv)
62: {
63: TS ts; /* nonlinear solver */
64: Vec x,r; /* solution, residual vectors */
65: Mat J; /* Jacobian matrix */
66: PetscInt steps,Mx,maxsteps = 10000000;
68: DM da;
69: PetscReal dt;
70: PetscReal vbounds[] = {-1.1,1.1};
71: PetscBool wait,vi = PETSC_FALSE,mymonitor;
72: Vec ul,uh;
73: UserCtx ctx;
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Initialize program
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: PetscInitialize(&argc,&argv,(char*)0,help);
79: ctx.kappa = 1.0;
80: PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);
81: ctx.degenerate = PETSC_FALSE;
82: PetscOptionsGetBool(NULL,NULL,"-degenerate",&ctx.degenerate,NULL);
83: ctx.cahnhillard = PETSC_FALSE;
84: PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL);
85: PetscOptionsGetBool(NULL,NULL,"-vi",&vi,NULL);
86: ctx.netforce = PETSC_FALSE;
87: PetscOptionsGetBool(NULL,NULL,"-netforce",&ctx.netforce,NULL);
88: ctx.energy = 1;
89: PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL);
90: ctx.tol = 1.0e-8;
91: PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL);
92: ctx.theta = .001;
93: ctx.theta_c = 1.0;
94: PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL);
95: PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL);
96: ctx.truncation = 1;
97: PetscOptionsGetInt(NULL,NULL,"-truncation",&ctx.truncation,NULL);
98: PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);
99: PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);
100: PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),800,600);
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Create distributed array (DMDA) to manage parallel grid and vectors
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, -10,1,2,NULL,&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: TSSetDuration(ts,maxsteps,.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: TSSetInitialTimeStep(ts,0.0,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: TSGetTimeStepNumber(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(0);
217: }
218: /* ------------------------------------------------------------------- */
221: /*
222: FormFunction - Evaluates nonlinear function, F(x).
224: Input Parameters:
225: . ts - the TS context
226: . X - input vector
227: . ptr - optional user-defined context, as set by SNESSetFunction()
229: Output Parameter:
230: . F - function vector
231: */
232: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
233: {
234: DM da;
236: PetscInt i,Mx,xs,xm;
237: PetscReal hx,sx;
238: PetscScalar *x,*f,c,r,l;
239: Vec localX;
240: UserCtx *ctx = (UserCtx*)ptr;
241: 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 */
244: TSGetDM(ts,&da);
245: DMGetLocalVector(da,&localX);
246: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
247: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
249: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
251: /*
252: Scatter ghost points to local vector,using the 2-step process
253: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
254: By placing code between these two statements, computations can be
255: done while messages are in transition.
256: */
257: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
258: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
260: /*
261: Get pointers to vector data
262: */
263: DMDAVecGetArrayRead(da,localX,&x);
264: DMDAVecGetArray(da,F,&f);
266: /*
267: Get local grid boundaries
268: */
269: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
271: /*
272: Compute function over the locally owned part of the grid
273: */
274: for (i=xs; i<xs+xm; i++) {
275: if (ctx->degenerate) {
276: c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
277: r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
278: l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx;
279: } else {
280: c = (x[i-1] + x[i+1] - 2.0*x[i])*sx;
281: r = (x[i] + x[i+2] - 2.0*x[i+1])*sx;
282: l = (x[i-2] + x[i] - 2.0*x[i-1])*sx;
283: }
284: f[i] = -ctx->kappa*(l + r - 2.0*c)*sx;
285: if (ctx->cahnhillard) {
286: switch (ctx->energy) {
287: case 1: /* double well */
288: 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;
289: break;
290: case 2: /* double obstacle */
291: f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx;
292: break;
293: case 3: /* logarithmic + double well */
294: 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;
295: if (ctx->truncation==2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */
296: 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;
297: 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;
298: 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;
299: } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */
300: a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
301: b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
302: 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;
303: 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;
304: 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;
305: }
306: break;
307: case 4: /* logarithmic + double obstacle */
308: f[i] += -theta_c*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
309: if (ctx->truncation==2) { /* quadratic */
310: 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;
311: 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;
312: 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;
313: } else { /* cubic */
314: a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
315: b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
316: 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;
317: 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;
318: 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;
319: }
320: break;
321: }
322: }
324: }
326: /*
327: Restore vectors
328: */
329: DMDAVecRestoreArrayRead(da,localX,&x);
330: DMDAVecRestoreArray(da,F,&f);
331: DMRestoreLocalVector(da,&localX);
332: return(0);
333: }
335: /* ------------------------------------------------------------------- */
338: /*
339: FormJacobian - Evaluates nonlinear function's Jacobian
341: */
342: PetscErrorCode FormJacobian(TS ts,PetscReal ftime,Vec X,Mat A,Mat B,void *ptr)
343: {
344: DM da;
346: PetscInt i,Mx,xs,xm;
347: MatStencil row,cols[5];
348: PetscReal hx,sx;
349: PetscScalar *x,vals[5];
350: Vec localX;
351: UserCtx *ctx = (UserCtx*)ptr;
354: TSGetDM(ts,&da);
355: DMGetLocalVector(da,&localX);
356: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
357: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
359: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
361: /*
362: Scatter ghost points to local vector,using the 2-step process
363: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
364: By placing code between these two statements, computations can be
365: done while messages are in transition.
366: */
367: DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
368: DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
370: /*
371: Get pointers to vector data
372: */
373: DMDAVecGetArrayRead(da,localX,&x);
375: /*
376: Get local grid boundaries
377: */
378: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
380: /*
381: Compute function over the locally owned part of the grid
382: */
383: for (i=xs; i<xs+xm; i++) {
384: row.i = i;
385: if (ctx->degenerate) {
386: /*PetscScalar c,r,l;
387: c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
388: r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx;
389: l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */
390: } else {
391: cols[0].i = i - 2; vals[0] = -ctx->kappa*sx*sx;
392: cols[1].i = i - 1; vals[1] = 4.0*ctx->kappa*sx*sx;
393: cols[2].i = i ; vals[2] = -6.0*ctx->kappa*sx*sx;
394: cols[3].i = i + 1; vals[3] = 4.0*ctx->kappa*sx*sx;
395: cols[4].i = i + 2; vals[4] = -ctx->kappa*sx*sx;
396: }
397: MatSetValuesStencil(B,1,&row,5,cols,vals,INSERT_VALUES);
399: if (ctx->cahnhillard) {
400: switch (ctx->energy) {
401: case 1: /* double well */
402: /* 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; */
403: break;
404: case 2: /* double obstacle */
405: /* f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */
406: break;
407: case 3: /* logarithmic + double well */
408: break;
409: case 4: /* logarithmic + double obstacle */
410: break;
411: }
412: }
414: }
416: /*
417: Restore vectors
418: */
419: DMDAVecRestoreArrayRead(da,localX,&x);
420: DMRestoreLocalVector(da,&localX);
421: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
422: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
423: if (A != B) {
424: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
425: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
426: }
427: return(0);
428: }
429: /* ------------------------------------------------------------------- */
432: PetscErrorCode FormInitialSolution(DM da,Vec U)
433: {
434: PetscErrorCode ierr;
435: PetscInt i,xs,xm,Mx,N,scale;
436: PetscScalar *u;
437: PetscReal r,hx,x;
438: const PetscScalar *f;
439: Vec finesolution;
440: PetscViewer viewer;
443: 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);
445: hx = 1.0/(PetscReal)Mx;
447: /*
448: Get pointers to vector data
449: */
450: DMDAVecGetArray(da,U,&u);
452: /*
453: Get local grid boundaries
454: */
455: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
457: /* InitialSolution.biharmonic is obtained by running
458: ./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
459: After changing the initial grid spacing to 10 and the stencil width to 2 in the DMDA create.
460: */
461: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.biharmonic",FILE_MODE_READ,&viewer);
462: VecCreate(PETSC_COMM_WORLD,&finesolution);
463: VecLoad(finesolution,viewer);
464: PetscViewerDestroy(&viewer);
465: VecGetSize(finesolution,&N);
466: scale = N/Mx;
467: VecGetArrayRead(finesolution,&f);
469: /*
470: Compute function over the locally owned part of the grid
471: */
472: for (i=xs; i<xs+xm; i++) {
473: x = i*hx;
474: r = PetscSqrtReal((x-.5)*(x-.5));
475: if (r < .125) u[i] = 1.0;
476: else u[i] = -.5;
478: /* With the initial condition above the method is first order in space */
479: /* this is a smooth initial condition so the method becomes second order in space */
480: /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
481: u[i] = f[scale*i];
482: }
483: VecRestoreArrayRead(finesolution,&f);
484: VecDestroy(&finesolution);
486: /*
487: Restore vectors
488: */
489: DMDAVecRestoreArray(da,U,&u);
490: return(0);
491: }
495: /*
496: This routine is not parallel
497: */
498: PetscErrorCode MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
499: {
500: UserCtx *ctx = (UserCtx*)ptr;
501: PetscDrawLG lg;
503: PetscScalar *u,l,r,c;
504: PetscInt Mx,i,xs,xm,cnt;
505: PetscReal x,y,hx,pause,sx,len,max,xx[4],yy[4],xx_netforce,yy_netforce,yup,ydown,y2,len2;
506: PetscDraw draw;
507: Vec localU;
508: DM da;
509: int colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE,PETSC_DRAW_PLUM,PETSC_DRAW_BLACK};
510: /*
511: 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"}};
512: */
513: PetscDrawAxis axis;
514: PetscDrawViewPorts *ports;
515: 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 */
519: TSGetDM(ts,&da);
520: DMGetLocalVector(da,&localU);
521: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
522: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
523: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
524: hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
525: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
526: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
527: DMDAVecGetArrayRead(da,localU,&u);
529: PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);
530: PetscDrawLGGetDraw(lg,&draw);
531: PetscDrawCheckResizedWindow(draw);
532: if (!ctx->ports) {
533: PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);
534: }
535: ports = ctx->ports;
536: PetscDrawLGGetAxis(lg,&axis);
537: PetscDrawLGReset(lg);
539: xx[0] = 0.0; xx[1] = 1.0; cnt = 2;
540: PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL);
541: xs = xx[0]/hx; xm = (xx[1] - xx[0])/hx;
543: /*
544: Plot the energies
545: */
546: PetscDrawLGSetDimension(lg,1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3));
547: PetscDrawLGSetColors(lg,colors+1);
548: PetscDrawViewPortsSet(ports,2);
549: x = hx*xs;
550: for (i=xs; i<xs+xm; i++) {
551: xx[0] = xx[1] = xx[2] = x;
552: 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);
553: else yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);
555: if (ctx->cahnhillard) {
556: switch (ctx->energy) {
557: case 1: /* double well */
558: yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
559: break;
560: case 2: /* double obstacle */
561: yy[1] = .5*PetscRealPart(1. - u[i]*u[i]);
562: break;
563: case 3: /* logarithm + double well */
564: yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
565: 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));
566: 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));
567: 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));
568: break;
569: case 4: /* logarithm + double obstacle */
570: yy[1] = .5*theta_c*PetscRealPart(1.0-u[i]*u[i]);
571: 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));
572: 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));
573: 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));
574: break;
575: }
576: }
577: PetscDrawLGAddPoint(lg,xx,yy);
578: x += hx;
579: }
580: PetscDrawGetPause(draw,&pause);
581: PetscDrawSetPause(draw,0.0);
582: PetscDrawAxisSetLabels(axis,"Energy","","");
583: /* PetscDrawLGSetLegend(lg,legend[ctx->energy-1]); */
584: PetscDrawLGDraw(lg);
586: /*
587: Plot the forces
588: */
589: PetscDrawLGSetDimension(lg,0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3));
590: PetscDrawLGSetColors(lg,colors+1);
591: PetscDrawViewPortsSet(ports,1);
592: PetscDrawLGReset(lg);
593: x = xs*hx;
594: max = 0.;
595: for (i=xs; i<xs+xm; i++) {
596: xx[0] = xx[1] = xx[2] = xx[3] = x;
597: xx_netforce = x;
598: if (ctx->degenerate) {
599: c = (1. - u[i]*u[i])*(u[i-1] + u[i+1] - 2.0*u[i])*sx;
600: r = (1. - u[i+1]*u[i+1])*(u[i] + u[i+2] - 2.0*u[i+1])*sx;
601: l = (1. - u[i-1]*u[i-1])*(u[i-2] + u[i] - 2.0*u[i-1])*sx;
602: } else {
603: c = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
604: r = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
605: l = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
606: }
607: yy[0] = PetscRealPart(-ctx->kappa*(l + r - 2.0*c)*sx);
608: yy_netforce = yy[0];
609: max = PetscMax(max,PetscAbs(yy[0]));
610: if (ctx->cahnhillard) {
611: switch (ctx->energy) {
612: case 1: /* double well */
613: 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);
614: break;
615: case 2: /* double obstacle */
616: yy[1] = -PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx;
617: break;
618: case 3: /* logarithmic + double well */
619: 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);
620: if (ctx->truncation==2) { /* quadratic */
621: 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;
622: 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;
623: 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);
624: } else { /* cubic */
625: a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
626: b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
627: 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);
628: 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);
629: 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);
630: }
631: break;
632: case 4: /* logarithmic + double obstacle */
633: yy[1] = theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i]))*sx;
634: if (ctx->truncation==2) {
635: 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;
636: 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;
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: else {
640: a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
641: b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
642: 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);
643: 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);
644: 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);
645: }
646: break;
647: }
648: if (ctx->energy < 3) {
649: max = PetscMax(max,PetscAbs(yy[1]));
650: yy[2] = yy[0]+yy[1];
651: yy_netforce = yy[2];
652: } else {
653: max = PetscMax(max,PetscAbs(yy[1]+yy[2]));
654: yy[3] = yy[0]+yy[1]+yy[2];
655: yy_netforce = yy[3];
656: }
657: }
658: if (ctx->netforce) {
659: PetscDrawLGAddPoint(lg,&xx_netforce,&yy_netforce);
660: } else {
661: PetscDrawLGAddPoint(lg,xx,yy);
662: }
663: x += hx;
664: /*if (max > 7200150000.0) */
665: /* printf("max very big when i = %d\n",i); */
666: }
667: PetscDrawAxisSetLabels(axis,"Right hand side","","");
668: PetscDrawLGSetLegend(lg,NULL);
669: PetscDrawLGDraw(lg);
671: /*
672: Plot the solution
673: */
674: PetscDrawLGSetDimension(lg,1);
675: PetscDrawViewPortsSet(ports,0);
676: PetscDrawLGReset(lg);
677: x = hx*xs;
678: PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);
679: PetscDrawLGSetColors(lg,colors);
680: for (i=xs; i<xs+xm; i++) {
681: xx[0] = x;
682: yy[0] = PetscRealPart(u[i]);
683: PetscDrawLGAddPoint(lg,xx,yy);
684: x += hx;
685: }
686: PetscDrawAxisSetLabels(axis,"Solution","","");
687: PetscDrawLGDraw(lg);
689: /*
690: Print the forces as arrows on the solution
691: */
692: x = hx*xs;
693: cnt = xm/60;
694: cnt = (!cnt) ? 1 : cnt;
696: for (i=xs; i<xs+xm; i += cnt) {
697: y = yup = ydown = PetscRealPart(u[i]);
698: c = (u[i-1] + u[i+1] - 2.0*u[i])*sx;
699: r = (u[i] + u[i+2] - 2.0*u[i+1])*sx;
700: l = (u[i-2] + u[i] - 2.0*u[i-1])*sx;
701: len = -.5*PetscRealPart(ctx->kappa*(l + r - 2.0*c)*sx)/max;
702: PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);
703: if (ctx->cahnhillard) {
704: if (len < 0.) ydown += len;
705: else yup += len;
707: switch (ctx->energy) {
708: case 1: /* double well */
709: 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;
710: break;
711: case 2: /* double obstacle */
712: len = -.5*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max;
713: break;
714: case 3: /* logarithmic + double well */
715: 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;
716: if (len < 0.) ydown += len;
717: else yup += len;
719: if (ctx->truncation==2) { /* quadratic */
720: 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;
721: 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;
722: 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);
723: } else { /* cubic */
724: a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
725: b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
726: 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);
727: 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);
728: 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);
729: }
730: y2 = len < 0 ? ydown : yup;
731: PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);
732: break;
733: case 4: /* logarithmic + double obstacle */
734: len = -.5*theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i])*sx/max);
735: if (len < 0.) ydown += len;
736: else yup += len;
738: if (ctx->truncation==2) { /* quadratic */
739: 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;
740: 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;
741: 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);
742: } else { /* cubic */
743: a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol));
744: b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol);
745: 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;
746: 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;
747: 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;
748: }
749: y2 = len < 0 ? ydown : yup;
750: PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM);
751: break;
752: }
753: PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);
754: }
755: x += cnt*hx;
756: }
757: DMDAVecRestoreArrayRead(da,localU,&x);
758: DMRestoreLocalVector(da,&localU);
759: PetscDrawStringSetSize(draw,.2,.2);
760: PetscDrawFlush(draw);
761: PetscDrawSetPause(draw,pause);
762: PetscDrawPause(draw);
763: return(0);
764: }
768: PetscErrorCode MyDestroy(void **ptr)
769: {
770: UserCtx *ctx = *(UserCtx**)ptr;
774: PetscDrawViewPortsDestroy(ctx->ports);
775: return(0);
776: }