xref: /petsc/src/ts/tutorials/phasefield/heat.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Solves heat equation in 1d.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*
5*c4762a1bSJed Brown   Solves the equation
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown     u_t = kappa  \Delta u
8*c4762a1bSJed Brown        Periodic boundary conditions
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown Evolve the  heat equation:
11*c4762a1bSJed Brown ---------------
12*c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor  -pc_type lu  -draw_pause .1 -snes_converged_reason  -ts_type cn  -da_refine 5 -mymonitor
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown Evolve the  Allen-Cahn equation:
15*c4762a1bSJed Brown ---------------
16*c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor  -pc_type lu  -draw_pause .1 -snes_converged_reason  -ts_type cn  -da_refine 5   -allen-cahn -kappa .001 -ts_max_time 5  -mymonitor
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown Evolve the  Allen-Cahn equation: zoom in on part of the domain
19*c4762a1bSJed Brown ---------------
20*c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor  -pc_type lu   -snes_converged_reason     -ts_type cn  -da_refine 5   -allen-cahn -kappa .001 -ts_max_time 5  -zoom .25,.45  -mymonitor
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown The option -square_initial indicates it should use a square wave initial condition otherwise it loads the file InitialSolution.heat as the initial solution. You should run with
24*c4762a1bSJed Brown ./heat -square_initial -ts_monitor -snes_monitor  -pc_type lu   -snes_converged_reason    -ts_type cn  -da_refine 9 -ts_max_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25  -ts_max_steps 15
25*c4762a1bSJed Brown to generate InitialSolution.heat
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown */
28*c4762a1bSJed Brown #include <petscdm.h>
29*c4762a1bSJed Brown #include <petscdmda.h>
30*c4762a1bSJed Brown #include <petscts.h>
31*c4762a1bSJed Brown #include <petscdraw.h>
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown /*
34*c4762a1bSJed Brown    User-defined routines
35*c4762a1bSJed Brown */
36*c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec),MyMonitor(TS,PetscInt,PetscReal,Vec,void*),MyDestroy(void**);
37*c4762a1bSJed Brown typedef struct {PetscReal kappa;PetscBool allencahn;PetscDrawViewPorts *ports;} UserCtx;
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown int main(int argc,char **argv)
40*c4762a1bSJed Brown {
41*c4762a1bSJed Brown   TS             ts;                           /* time integrator */
42*c4762a1bSJed Brown   Vec            x,r;                          /* solution, residual vectors */
43*c4762a1bSJed Brown   PetscInt       steps,Mx;
44*c4762a1bSJed Brown   PetscErrorCode ierr;
45*c4762a1bSJed Brown   DM             da;
46*c4762a1bSJed Brown   PetscReal      dt;
47*c4762a1bSJed Brown   UserCtx        ctx;
48*c4762a1bSJed Brown   PetscBool      mymonitor;
49*c4762a1bSJed Brown   PetscViewer    viewer;
50*c4762a1bSJed Brown   PetscBool      flg;
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53*c4762a1bSJed Brown      Initialize program
54*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
56*c4762a1bSJed Brown   ctx.kappa     = 1.0;
57*c4762a1bSJed Brown   ierr          = PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ctx.allencahn = PETSC_FALSE;
59*c4762a1bSJed Brown   ierr          = PetscOptionsHasName(NULL,NULL,"-allen-cahn",&ctx.allencahn);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr          = PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor);CHKERRQ(ierr);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
64*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,1,2,NULL,&da);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"Heat equation: u");CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
70*c4762a1bSJed Brown   dt   = 1.0/(ctx.kappa*Mx*Mx);
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73*c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
74*c4762a1bSJed Brown      vectors that are the same types
75*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80*c4762a1bSJed Brown      Create timestepping solver context
81*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,FormFunction,&ctx);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88*c4762a1bSJed Brown      Customize nonlinear solver
89*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90*c4762a1bSJed Brown   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93*c4762a1bSJed Brown      Set initial conditions
94*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95*c4762a1bSJed Brown   ierr = FormInitialSolution(da,x);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,.02);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   if (mymonitor) {
103*c4762a1bSJed Brown     ctx.ports = NULL;
104*c4762a1bSJed Brown     ierr      = TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy);CHKERRQ(ierr);
105*c4762a1bSJed Brown   }
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108*c4762a1bSJed Brown      Set runtime options
109*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113*c4762a1bSJed Brown      Solve nonlinear system
114*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115*c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-square_initial",&flg);CHKERRQ(ierr);
118*c4762a1bSJed Brown   if (flg) {
119*c4762a1bSJed Brown     ierr  = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
120*c4762a1bSJed Brown     ierr  = VecView(x,viewer);CHKERRQ(ierr);
121*c4762a1bSJed Brown     ierr  = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
122*c4762a1bSJed Brown   }
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
126*c4762a1bSJed Brown      are no longer needed.
127*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   ierr = PetscFinalize();
134*c4762a1bSJed Brown   return ierr;
135*c4762a1bSJed Brown }
136*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
137*c4762a1bSJed Brown /*
138*c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown    Input Parameters:
141*c4762a1bSJed Brown .  ts - the TS context
142*c4762a1bSJed Brown .  X - input vector
143*c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown    Output Parameter:
146*c4762a1bSJed Brown .  F - function vector
147*c4762a1bSJed Brown  */
148*c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
149*c4762a1bSJed Brown {
150*c4762a1bSJed Brown   DM             da;
151*c4762a1bSJed Brown   PetscErrorCode ierr;
152*c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
153*c4762a1bSJed Brown   PetscReal      hx,sx;
154*c4762a1bSJed Brown   PetscScalar    *x,*f;
155*c4762a1bSJed Brown   Vec            localX;
156*c4762a1bSJed Brown   UserCtx        *ctx = (UserCtx*)ptr;
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown   PetscFunctionBegin;
159*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown   /*
166*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
167*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
168*c4762a1bSJed Brown      By placing code between these two statements, computations can be
169*c4762a1bSJed Brown      done while messages are in transition.
170*c4762a1bSJed Brown   */
171*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown   /*
175*c4762a1bSJed Brown      Get pointers to vector data
176*c4762a1bSJed Brown   */
177*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /*
181*c4762a1bSJed Brown      Get local grid boundaries
182*c4762a1bSJed Brown   */
183*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown   /*
186*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
187*c4762a1bSJed Brown   */
188*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
189*c4762a1bSJed Brown     f[i] = ctx->kappa*(x[i-1] + x[i+1] - 2.0*x[i])*sx;
190*c4762a1bSJed Brown     if (ctx->allencahn) f[i] += (x[i] - x[i]*x[i]*x[i]);
191*c4762a1bSJed Brown   }
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown   /*
194*c4762a1bSJed Brown      Restore vectors
195*c4762a1bSJed Brown   */
196*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
197*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
198*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
199*c4762a1bSJed Brown   PetscFunctionReturn(0);
200*c4762a1bSJed Brown }
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
203*c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U)
204*c4762a1bSJed Brown {
205*c4762a1bSJed Brown   PetscErrorCode    ierr;
206*c4762a1bSJed Brown   PetscInt          i,xs,xm,Mx,scale=1,N;
207*c4762a1bSJed Brown   PetscScalar       *u;
208*c4762a1bSJed Brown   const PetscScalar *f;
209*c4762a1bSJed Brown   PetscReal         hx,x,r;
210*c4762a1bSJed Brown   Vec               finesolution;
211*c4762a1bSJed Brown   PetscViewer       viewer;
212*c4762a1bSJed Brown   PetscBool         flg;
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown   PetscFunctionBegin;
215*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown   hx = 1.0/(PetscReal)Mx;
218*c4762a1bSJed Brown 
219*c4762a1bSJed Brown   /*
220*c4762a1bSJed Brown      Get pointers to vector data
221*c4762a1bSJed Brown   */
222*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   /*
225*c4762a1bSJed Brown      Get local grid boundaries
226*c4762a1bSJed Brown   */
227*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   /*  InitialSolution is obtained with
230*c4762a1bSJed Brown       ./heat -ts_monitor -snes_monitor  -pc_type lu   -snes_converged_reason    -ts_type cn  -da_refine 9 -ts_max_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25  -ts_max_steps 15
231*c4762a1bSJed Brown   */
232*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-square_initial",&flg);CHKERRQ(ierr);
233*c4762a1bSJed Brown   if (!flg) {
234*c4762a1bSJed Brown     ierr  = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer);CHKERRQ(ierr);
235*c4762a1bSJed Brown     ierr  = VecCreate(PETSC_COMM_WORLD,&finesolution);CHKERRQ(ierr);
236*c4762a1bSJed Brown     ierr  = VecLoad(finesolution,viewer);CHKERRQ(ierr);
237*c4762a1bSJed Brown     ierr  = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
238*c4762a1bSJed Brown     ierr  = VecGetSize(finesolution,&N);CHKERRQ(ierr);
239*c4762a1bSJed Brown     scale = N/Mx;
240*c4762a1bSJed Brown     ierr  = VecGetArrayRead(finesolution,&f);CHKERRQ(ierr);
241*c4762a1bSJed Brown   }
242*c4762a1bSJed Brown 
243*c4762a1bSJed Brown   /*
244*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
245*c4762a1bSJed Brown   */
246*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
247*c4762a1bSJed Brown     x = i*hx;
248*c4762a1bSJed Brown     r = PetscSqrtScalar((x-.5)*(x-.5));
249*c4762a1bSJed Brown     if (r < .125) u[i] = 1.0;
250*c4762a1bSJed Brown     else u[i] = -.5;
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown     /* With the initial condition above the method is first order in space */
253*c4762a1bSJed Brown     /* this is a smooth initial condition so the method becomes second order in space */
254*c4762a1bSJed Brown     /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
255*c4762a1bSJed Brown     /*  u[i] = f[scale*i];*/
256*c4762a1bSJed Brown     if (!flg) u[i] = f[scale*i];
257*c4762a1bSJed Brown   }
258*c4762a1bSJed Brown   if (!flg) {
259*c4762a1bSJed Brown     ierr = VecRestoreArrayRead(finesolution,&f);CHKERRQ(ierr);
260*c4762a1bSJed Brown     ierr = VecDestroy(&finesolution);CHKERRQ(ierr);
261*c4762a1bSJed Brown   }
262*c4762a1bSJed Brown 
263*c4762a1bSJed Brown   /*
264*c4762a1bSJed Brown      Restore vectors
265*c4762a1bSJed Brown   */
266*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
267*c4762a1bSJed Brown   PetscFunctionReturn(0);
268*c4762a1bSJed Brown }
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown /*
271*c4762a1bSJed Brown     This routine is not parallel
272*c4762a1bSJed Brown */
273*c4762a1bSJed Brown PetscErrorCode  MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr)
274*c4762a1bSJed Brown {
275*c4762a1bSJed Brown   UserCtx            *ctx = (UserCtx*)ptr;
276*c4762a1bSJed Brown   PetscDrawLG        lg;
277*c4762a1bSJed Brown   PetscErrorCode     ierr;
278*c4762a1bSJed Brown   PetscScalar        *u;
279*c4762a1bSJed Brown   PetscInt           Mx,i,xs,xm,cnt;
280*c4762a1bSJed Brown   PetscReal          x,y,hx,pause,sx,len,max,xx[2],yy[2];
281*c4762a1bSJed Brown   PetscDraw          draw;
282*c4762a1bSJed Brown   Vec                localU;
283*c4762a1bSJed Brown   DM                 da;
284*c4762a1bSJed Brown   int                colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE};
285*c4762a1bSJed Brown   const char*const   legend[] = {"-kappa (\\grad u,\\grad u)","(1 - u^2)^2"};
286*c4762a1bSJed Brown   PetscDrawAxis      axis;
287*c4762a1bSJed Brown   PetscDrawViewPorts *ports;
288*c4762a1bSJed Brown   PetscReal          vbounds[] = {-1.1,1.1};
289*c4762a1bSJed Brown 
290*c4762a1bSJed Brown   PetscFunctionBegin;
291*c4762a1bSJed Brown   ierr = PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds);CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1200,800);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
294*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
295*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
296*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
297*c4762a1bSJed Brown   hx   = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx);
298*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
299*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
300*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown   ierr = PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg);CHKERRQ(ierr);
303*c4762a1bSJed Brown   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
304*c4762a1bSJed Brown   ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
305*c4762a1bSJed Brown   if (!ctx->ports) {
306*c4762a1bSJed Brown     ierr = PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports);CHKERRQ(ierr);
307*c4762a1bSJed Brown   }
308*c4762a1bSJed Brown   ports = ctx->ports;
309*c4762a1bSJed Brown   ierr  = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr);
310*c4762a1bSJed Brown   ierr  = PetscDrawLGReset(lg);CHKERRQ(ierr);
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown   xx[0] = 0.0; xx[1] = 1.0; cnt = 2;
313*c4762a1bSJed Brown   ierr  = PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL);CHKERRQ(ierr);
314*c4762a1bSJed Brown   xs    = xx[0]/hx; xm = (xx[1] - xx[0])/hx;
315*c4762a1bSJed Brown 
316*c4762a1bSJed Brown   /*
317*c4762a1bSJed Brown       Plot the  energies
318*c4762a1bSJed Brown   */
319*c4762a1bSJed Brown   ierr = PetscDrawLGSetDimension(lg,1 + (ctx->allencahn ? 1 : 0));CHKERRQ(ierr);
320*c4762a1bSJed Brown   ierr = PetscDrawLGSetColors(lg,colors+1);CHKERRQ(ierr);
321*c4762a1bSJed Brown   ierr = PetscDrawViewPortsSet(ports,2);CHKERRQ(ierr);
322*c4762a1bSJed Brown   x    = hx*xs;
323*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
324*c4762a1bSJed Brown     xx[0] = xx[1] = x;
325*c4762a1bSJed Brown     yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx);
326*c4762a1bSJed Brown     if (ctx->allencahn) yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i]));
327*c4762a1bSJed Brown     ierr = PetscDrawLGAddPoint(lg,xx,yy);CHKERRQ(ierr);
328*c4762a1bSJed Brown     x   += hx;
329*c4762a1bSJed Brown   }
330*c4762a1bSJed Brown   ierr = PetscDrawGetPause(draw,&pause);CHKERRQ(ierr);
331*c4762a1bSJed Brown   ierr = PetscDrawSetPause(draw,0.0);CHKERRQ(ierr);
332*c4762a1bSJed Brown   ierr = PetscDrawAxisSetLabels(axis,"Energy","","");CHKERRQ(ierr);
333*c4762a1bSJed Brown   ierr = PetscDrawLGSetLegend(lg,legend);CHKERRQ(ierr);
334*c4762a1bSJed Brown   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
335*c4762a1bSJed Brown 
336*c4762a1bSJed Brown   /*
337*c4762a1bSJed Brown       Plot the  forces
338*c4762a1bSJed Brown   */
339*c4762a1bSJed Brown   ierr = PetscDrawViewPortsSet(ports,1);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
341*c4762a1bSJed Brown   x    = xs*hx;
342*c4762a1bSJed Brown   max  = 0.;
343*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
344*c4762a1bSJed Brown     xx[0] = xx[1] = x;
345*c4762a1bSJed Brown     yy[0] = PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx);
346*c4762a1bSJed Brown     max   = PetscMax(max,PetscAbs(yy[0]));
347*c4762a1bSJed Brown     if (ctx->allencahn) {
348*c4762a1bSJed Brown       yy[1] = PetscRealPart(u[i] - u[i]*u[i]*u[i]);
349*c4762a1bSJed Brown       max   = PetscMax(max,PetscAbs(yy[1]));
350*c4762a1bSJed Brown     }
351*c4762a1bSJed Brown     ierr = PetscDrawLGAddPoint(lg,xx,yy);CHKERRQ(ierr);
352*c4762a1bSJed Brown     x   += hx;
353*c4762a1bSJed Brown   }
354*c4762a1bSJed Brown   ierr = PetscDrawAxisSetLabels(axis,"Right hand side","","");CHKERRQ(ierr);
355*c4762a1bSJed Brown   ierr = PetscDrawLGSetLegend(lg,NULL);CHKERRQ(ierr);
356*c4762a1bSJed Brown   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
357*c4762a1bSJed Brown 
358*c4762a1bSJed Brown   /*
359*c4762a1bSJed Brown         Plot the solution
360*c4762a1bSJed Brown   */
361*c4762a1bSJed Brown   ierr = PetscDrawLGSetDimension(lg,1);CHKERRQ(ierr);
362*c4762a1bSJed Brown   ierr = PetscDrawViewPortsSet(ports,0);CHKERRQ(ierr);
363*c4762a1bSJed Brown   ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
364*c4762a1bSJed Brown   x    = hx*xs;
365*c4762a1bSJed Brown   ierr = PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1);CHKERRQ(ierr);
366*c4762a1bSJed Brown   ierr = PetscDrawLGSetColors(lg,colors);CHKERRQ(ierr);
367*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
368*c4762a1bSJed Brown     xx[0] = x;
369*c4762a1bSJed Brown     yy[0] = PetscRealPart(u[i]);
370*c4762a1bSJed Brown     ierr  = PetscDrawLGAddPoint(lg,xx,yy);CHKERRQ(ierr);
371*c4762a1bSJed Brown     x    += hx;
372*c4762a1bSJed Brown   }
373*c4762a1bSJed Brown   ierr = PetscDrawAxisSetLabels(axis,"Solution","","");CHKERRQ(ierr);
374*c4762a1bSJed Brown   ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
375*c4762a1bSJed Brown 
376*c4762a1bSJed Brown   /*
377*c4762a1bSJed Brown       Print the  forces as arrows on the solution
378*c4762a1bSJed Brown   */
379*c4762a1bSJed Brown   x   = hx*xs;
380*c4762a1bSJed Brown   cnt = xm/60;
381*c4762a1bSJed Brown   cnt = (!cnt) ? 1 : cnt;
382*c4762a1bSJed Brown 
383*c4762a1bSJed Brown   for (i=xs; i<xs+xm; i += cnt) {
384*c4762a1bSJed Brown     y    = PetscRealPart(u[i]);
385*c4762a1bSJed Brown     len  = .5*PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max;
386*c4762a1bSJed Brown     ierr = PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED);CHKERRQ(ierr);
387*c4762a1bSJed Brown     if (ctx->allencahn) {
388*c4762a1bSJed Brown       len  = .5*PetscRealPart(u[i] - u[i]*u[i]*u[i])/max;
389*c4762a1bSJed Brown       ierr = PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE);CHKERRQ(ierr);
390*c4762a1bSJed Brown     }
391*c4762a1bSJed Brown     x += cnt*hx;
392*c4762a1bSJed Brown   }
393*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&x);CHKERRQ(ierr);
394*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
395*c4762a1bSJed Brown   ierr = PetscDrawStringSetSize(draw,.2,.2);CHKERRQ(ierr);
396*c4762a1bSJed Brown   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
397*c4762a1bSJed Brown   ierr = PetscDrawSetPause(draw,pause);CHKERRQ(ierr);
398*c4762a1bSJed Brown   ierr = PetscDrawPause(draw);CHKERRQ(ierr);
399*c4762a1bSJed Brown   PetscFunctionReturn(0);
400*c4762a1bSJed Brown }
401*c4762a1bSJed Brown 
402*c4762a1bSJed Brown PetscErrorCode  MyDestroy(void **ptr)
403*c4762a1bSJed Brown {
404*c4762a1bSJed Brown   UserCtx        *ctx = *(UserCtx**)ptr;
405*c4762a1bSJed Brown   PetscErrorCode ierr;
406*c4762a1bSJed Brown 
407*c4762a1bSJed Brown   PetscFunctionBegin;
408*c4762a1bSJed Brown   ierr = PetscDrawViewPortsDestroy(ctx->ports);CHKERRQ(ierr);
409*c4762a1bSJed Brown   PetscFunctionReturn(0);
410*c4762a1bSJed Brown }
411*c4762a1bSJed Brown 
412*c4762a1bSJed Brown /*TEST
413*c4762a1bSJed Brown 
414*c4762a1bSJed Brown    test:
415*c4762a1bSJed Brown      args: -ts_monitor -snes_monitor  -pc_type lu  -snes_converged_reason   -ts_type cn  -da_refine 2 -square_initial
416*c4762a1bSJed Brown 
417*c4762a1bSJed Brown    test:
418*c4762a1bSJed Brown      suffix: 2
419*c4762a1bSJed Brown      args: -ts_monitor -snes_monitor  -pc_type lu  -snes_converged_reason   -ts_type cn  -da_refine 5 -mymonitor -square_initial -allen-cahn -kappa .001
420*c4762a1bSJed Brown      requires: x
421*c4762a1bSJed Brown 
422*c4762a1bSJed Brown TEST*/
423