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