1c4762a1bSJed Brown static char help[] = "Solves heat equation in 1d.\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown Solves the equation
5c4762a1bSJed Brown
6c4762a1bSJed Brown u_t = kappa \Delta u
7c4762a1bSJed Brown Periodic boundary conditions
8c4762a1bSJed Brown
9c4762a1bSJed Brown Evolve the heat equation:
10c4762a1bSJed Brown ---------------
11c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type cn -da_refine 5 -mymonitor
12c4762a1bSJed Brown
13c4762a1bSJed Brown Evolve the Allen-Cahn equation:
14c4762a1bSJed Brown ---------------
15c4762a1bSJed 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
16c4762a1bSJed Brown
17c4762a1bSJed Brown Evolve the Allen-Cahn equation: zoom in on part of the domain
18c4762a1bSJed Brown ---------------
19c4762a1bSJed 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
20c4762a1bSJed Brown
21c4762a1bSJed 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
22188af4bfSBarry Smith ./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_time_step .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25 -ts_max_steps 15
23c4762a1bSJed Brown to generate InitialSolution.heat
24c4762a1bSJed Brown
25c4762a1bSJed Brown */
26c4762a1bSJed Brown #include <petscdm.h>
27c4762a1bSJed Brown #include <petscdmda.h>
28c4762a1bSJed Brown #include <petscts.h>
29c4762a1bSJed Brown #include <petscdraw.h>
30c4762a1bSJed Brown
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown User-defined routines
33c4762a1bSJed Brown */
34*2a8381b2SBarry Smith extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, PetscCtx), FormInitialSolution(DM, Vec), MyMonitor(TS, PetscInt, PetscReal, Vec, PetscCtx);
35*2a8381b2SBarry Smith extern PetscCtxDestroyFn MyDestroy;
369371c9d4SSatish Balay typedef struct {
379371c9d4SSatish Balay PetscReal kappa;
389371c9d4SSatish Balay PetscBool allencahn;
399371c9d4SSatish Balay PetscDrawViewPorts *ports;
40*2a8381b2SBarry Smith } AppCtx;
41c4762a1bSJed Brown
main(int argc,char ** argv)42d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown TS ts; /* time integrator */
45c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */
46c4762a1bSJed Brown PetscInt steps, Mx;
47c4762a1bSJed Brown DM da;
48c4762a1bSJed Brown PetscReal dt;
49*2a8381b2SBarry Smith AppCtx ctx;
50c4762a1bSJed Brown PetscBool mymonitor;
51c4762a1bSJed Brown PetscViewer viewer;
52c4762a1bSJed Brown PetscBool flg;
53c4762a1bSJed Brown
54c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55c4762a1bSJed Brown Initialize program
56c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57327415f7SBarry Smith PetscFunctionBeginUser;
58c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
59c4762a1bSJed Brown ctx.kappa = 1.0;
609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &ctx.kappa, NULL));
61c4762a1bSJed Brown ctx.allencahn = PETSC_FALSE;
629566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-allen-cahn", &ctx.allencahn));
639566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-mymonitor", &mymonitor));
64c4762a1bSJed Brown
65c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
67c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
689566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, 1, 2, NULL, &da));
699566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
709566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
719566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "Heat equation: u"));
729566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
73c4762a1bSJed Brown dt = 1.0 / (ctx.kappa * Mx * Mx);
74c4762a1bSJed Brown
75c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining
77c4762a1bSJed Brown vectors that are the same types
78c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
799566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r));
81c4762a1bSJed Brown
82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown Create timestepping solver context
84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
859566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
869566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
879566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
889566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &ctx));
89c4762a1bSJed Brown
90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91c4762a1bSJed Brown Customize nonlinear solver
92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
939566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
94c4762a1bSJed Brown
95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown Set initial conditions
97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
989566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x));
999566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
1009566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, .02));
1019566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_INTERPOLATE));
1029566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x));
103c4762a1bSJed Brown
104c4762a1bSJed Brown if (mymonitor) {
105c4762a1bSJed Brown ctx.ports = NULL;
1069566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyMonitor, &ctx, MyDestroy));
107c4762a1bSJed Brown }
108c4762a1bSJed Brown
109c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown Set runtime options
111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1129566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
113c4762a1bSJed Brown
114c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown Solve nonlinear system
116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1179566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
1189566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
1199566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-square_initial", &flg));
120c4762a1bSJed Brown if (flg) {
1219566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "InitialSolution.heat", FILE_MODE_WRITE, &viewer));
1229566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer));
1239566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
124c4762a1bSJed Brown }
125c4762a1bSJed Brown
126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
128c4762a1bSJed Brown are no longer needed.
129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
1329566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1339566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
134c4762a1bSJed Brown
1359566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
136b122ec5aSJacob Faibussowitsch return 0;
137c4762a1bSJed Brown }
138c4762a1bSJed Brown /* ------------------------------------------------------------------- */
139c4762a1bSJed Brown /*
140c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x).
141c4762a1bSJed Brown
142c4762a1bSJed Brown Input Parameters:
143c4762a1bSJed Brown . ts - the TS context
144c4762a1bSJed Brown . X - input vector
145c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction()
146c4762a1bSJed Brown
147c4762a1bSJed Brown Output Parameter:
148c4762a1bSJed Brown . F - function vector
149c4762a1bSJed Brown */
FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,PetscCtx Ctx)150*2a8381b2SBarry Smith PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, PetscCtx Ctx)
151d71ae5a4SJacob Faibussowitsch {
152c4762a1bSJed Brown DM da;
153c4762a1bSJed Brown PetscInt i, Mx, xs, xm;
154c4762a1bSJed Brown PetscReal hx, sx;
155c4762a1bSJed Brown PetscScalar *x, *f;
156c4762a1bSJed Brown Vec localX;
157*2a8381b2SBarry Smith AppCtx *ctx = (AppCtx *)Ctx;
158c4762a1bSJed Brown
159c4762a1bSJed Brown PetscFunctionBegin;
1609566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
1619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX));
1629566063dSJacob Faibussowitsch PetscCall(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));
163c4762a1bSJed Brown
1649371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx;
1659371c9d4SSatish Balay sx = 1.0 / (hx * hx);
166c4762a1bSJed Brown
167c4762a1bSJed Brown /*
168c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
169c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
170c4762a1bSJed Brown By placing code between these two statements, computations can be
171c4762a1bSJed Brown done while messages are in transition.
172c4762a1bSJed Brown */
1739566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1749566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
175c4762a1bSJed Brown
176c4762a1bSJed Brown /*
177c4762a1bSJed Brown Get pointers to vector data
178c4762a1bSJed Brown */
1799566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x));
1809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
181c4762a1bSJed Brown
182c4762a1bSJed Brown /*
183c4762a1bSJed Brown Get local grid boundaries
184c4762a1bSJed Brown */
1859566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
186c4762a1bSJed Brown
187c4762a1bSJed Brown /*
188c4762a1bSJed Brown Compute function over the locally owned part of the grid
189c4762a1bSJed Brown */
190c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
191c4762a1bSJed Brown f[i] = ctx->kappa * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx;
192c4762a1bSJed Brown if (ctx->allencahn) f[i] += (x[i] - x[i] * x[i] * x[i]);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown
195c4762a1bSJed Brown /*
196c4762a1bSJed Brown Restore vectors
197c4762a1bSJed Brown */
1989566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
1999566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
2009566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX));
2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
202c4762a1bSJed Brown }
203c4762a1bSJed Brown
204c4762a1bSJed Brown /* ------------------------------------------------------------------- */
FormInitialSolution(DM da,Vec U)205d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialSolution(DM da, Vec U)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown PetscInt i, xs, xm, Mx, scale = 1, N;
208c4762a1bSJed Brown PetscScalar *u;
209c4762a1bSJed Brown const PetscScalar *f;
210c4762a1bSJed Brown PetscReal hx, x, r;
211c4762a1bSJed Brown Vec finesolution;
212c4762a1bSJed Brown PetscViewer viewer;
213c4762a1bSJed Brown PetscBool flg;
214c4762a1bSJed Brown
215c4762a1bSJed Brown PetscFunctionBegin;
2169566063dSJacob Faibussowitsch PetscCall(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));
217c4762a1bSJed Brown
218c4762a1bSJed Brown hx = 1.0 / (PetscReal)Mx;
219c4762a1bSJed Brown
220c4762a1bSJed Brown /*
221c4762a1bSJed Brown Get pointers to vector data
222c4762a1bSJed Brown */
2239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
224c4762a1bSJed Brown
225c4762a1bSJed Brown /*
226c4762a1bSJed Brown Get local grid boundaries
227c4762a1bSJed Brown */
2289566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
229c4762a1bSJed Brown
230c4762a1bSJed Brown /* InitialSolution is obtained with
231188af4bfSBarry Smith ./heat -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 9 -ts_max_time 1.e-4 -ts_time_step .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25 -ts_max_steps 15
232c4762a1bSJed Brown */
2339566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-square_initial", &flg));
234c4762a1bSJed Brown if (!flg) {
2359566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "InitialSolution.heat", FILE_MODE_READ, &viewer));
2369566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &finesolution));
2379566063dSJacob Faibussowitsch PetscCall(VecLoad(finesolution, viewer));
2389566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
2399566063dSJacob Faibussowitsch PetscCall(VecGetSize(finesolution, &N));
240c4762a1bSJed Brown scale = N / Mx;
2419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(finesolution, &f));
242c4762a1bSJed Brown }
243c4762a1bSJed Brown
244c4762a1bSJed Brown /*
245c4762a1bSJed Brown Compute function over the locally owned part of the grid
246c4762a1bSJed Brown */
247c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
248c4762a1bSJed Brown x = i * hx;
2499fa27a79SStefano Zampini r = PetscSqrtReal((x - .5) * (x - .5));
250c4762a1bSJed Brown if (r < .125) u[i] = 1.0;
251c4762a1bSJed Brown else u[i] = -.5;
252c4762a1bSJed Brown
253c4762a1bSJed Brown /* With the initial condition above the method is first order in space */
254c4762a1bSJed Brown /* this is a smooth initial condition so the method becomes second order in space */
255c4762a1bSJed Brown /*u[i] = PetscSinScalar(2*PETSC_PI*x); */
256c4762a1bSJed Brown /* u[i] = f[scale*i];*/
257c4762a1bSJed Brown if (!flg) u[i] = f[scale * i];
258c4762a1bSJed Brown }
259c4762a1bSJed Brown if (!flg) {
2609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(finesolution, &f));
2619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&finesolution));
262c4762a1bSJed Brown }
263c4762a1bSJed Brown
264c4762a1bSJed Brown /*
265c4762a1bSJed Brown Restore vectors
266c4762a1bSJed Brown */
2679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown
271c4762a1bSJed Brown /*
272c4762a1bSJed Brown This routine is not parallel
273c4762a1bSJed Brown */
MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,PetscCtx Ctx)274*2a8381b2SBarry Smith PetscErrorCode MyMonitor(TS ts, PetscInt step, PetscReal time, Vec U, PetscCtx Ctx)
275d71ae5a4SJacob Faibussowitsch {
276*2a8381b2SBarry Smith AppCtx *ctx = (AppCtx *)Ctx;
277c4762a1bSJed Brown PetscDrawLG lg;
278c4762a1bSJed Brown PetscScalar *u;
279c4762a1bSJed Brown PetscInt Mx, i, xs, xm, cnt;
280c4762a1bSJed Brown PetscReal x, y, hx, pause, sx, len, max, xx[2], yy[2];
281c4762a1bSJed Brown PetscDraw draw;
282c4762a1bSJed Brown Vec localU;
283c4762a1bSJed Brown DM da;
284c4762a1bSJed Brown int colors[] = {PETSC_DRAW_YELLOW, PETSC_DRAW_RED, PETSC_DRAW_BLUE};
285c4762a1bSJed Brown const char *const legend[] = {"-kappa (\\grad u,\\grad u)", "(1 - u^2)^2"};
286c4762a1bSJed Brown PetscDrawAxis axis;
287c4762a1bSJed Brown PetscDrawViewPorts *ports;
288c4762a1bSJed Brown PetscReal vbounds[] = {-1.1, 1.1};
289c4762a1bSJed Brown
290c4762a1bSJed Brown PetscFunctionBegin;
2919566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, vbounds));
2929566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1200, 800));
2939566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
2949566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
2959566063dSJacob Faibussowitsch PetscCall(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));
2969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
2979371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx;
2989371c9d4SSatish Balay sx = 1.0 / (hx * hx);
2999566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
3009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
3019566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
302c4762a1bSJed Brown
3039566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, &lg));
3049566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetDraw(lg, &draw));
3059566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw));
30648a46eb9SPierre Jolivet if (!ctx->ports) PetscCall(PetscDrawViewPortsCreateRect(draw, 1, 3, &ctx->ports));
307c4762a1bSJed Brown ports = ctx->ports;
3089566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lg, &axis));
3099566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg));
310c4762a1bSJed Brown
3119371c9d4SSatish Balay xx[0] = 0.0;
3129371c9d4SSatish Balay xx[1] = 1.0;
3139371c9d4SSatish Balay cnt = 2;
3149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-zoom", xx, &cnt, NULL));
3159371c9d4SSatish Balay xs = xx[0] / hx;
3169371c9d4SSatish Balay xm = (xx[1] - xx[0]) / hx;
317c4762a1bSJed Brown
318c4762a1bSJed Brown /*
319c4762a1bSJed Brown Plot the energies
320c4762a1bSJed Brown */
3219566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 1 + (ctx->allencahn ? 1 : 0)));
3229566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors + 1));
3239566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 2));
324c4762a1bSJed Brown x = hx * xs;
325c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
326c4762a1bSJed Brown xx[0] = xx[1] = x;
327c4762a1bSJed Brown yy[0] = PetscRealPart(.25 * ctx->kappa * (u[i - 1] - u[i + 1]) * (u[i - 1] - u[i + 1]) * sx);
328c4762a1bSJed Brown if (ctx->allencahn) yy[1] = .25 * PetscRealPart((1. - u[i] * u[i]) * (1. - u[i] * u[i]));
3299566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
330c4762a1bSJed Brown x += hx;
331c4762a1bSJed Brown }
3329566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause));
3339566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0));
3349566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Energy", "", ""));
3359566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, legend));
3369566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg));
337c4762a1bSJed Brown
338c4762a1bSJed Brown /*
339c4762a1bSJed Brown Plot the forces
340c4762a1bSJed Brown */
3419566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 1));
3429566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg));
343c4762a1bSJed Brown x = xs * hx;
344c4762a1bSJed Brown max = 0.;
345c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
346c4762a1bSJed Brown xx[0] = xx[1] = x;
347c4762a1bSJed Brown yy[0] = PetscRealPart(ctx->kappa * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx);
348c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[0]));
349c4762a1bSJed Brown if (ctx->allencahn) {
350c4762a1bSJed Brown yy[1] = PetscRealPart(u[i] - u[i] * u[i] * u[i]);
351c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[1]));
352c4762a1bSJed Brown }
3539566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
354c4762a1bSJed Brown x += hx;
355c4762a1bSJed Brown }
3569566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Right hand side", "", ""));
3579566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, NULL));
3589566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg));
359c4762a1bSJed Brown
360c4762a1bSJed Brown /*
361c4762a1bSJed Brown Plot the solution
362c4762a1bSJed Brown */
3639566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 1));
3649566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 0));
3659566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg));
366c4762a1bSJed Brown x = hx * xs;
3679566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLimits(lg, x, x + (xm - 1) * hx, -1.1, 1.1));
3689566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors));
369c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
370c4762a1bSJed Brown xx[0] = x;
371c4762a1bSJed Brown yy[0] = PetscRealPart(u[i]);
3729566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
373c4762a1bSJed Brown x += hx;
374c4762a1bSJed Brown }
3759566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution", "", ""));
3769566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg));
377c4762a1bSJed Brown
378c4762a1bSJed Brown /*
379c4762a1bSJed Brown Print the forces as arrows on the solution
380c4762a1bSJed Brown */
381c4762a1bSJed Brown x = hx * xs;
382c4762a1bSJed Brown cnt = xm / 60;
383c4762a1bSJed Brown cnt = (!cnt) ? 1 : cnt;
384c4762a1bSJed Brown
385c4762a1bSJed Brown for (i = xs; i < xs + xm; i += cnt) {
386c4762a1bSJed Brown y = PetscRealPart(u[i]);
387c4762a1bSJed Brown len = .5 * PetscRealPart(ctx->kappa * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max;
3889566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_RED));
389c4762a1bSJed Brown if (ctx->allencahn) {
390c4762a1bSJed Brown len = .5 * PetscRealPart(u[i] - u[i] * u[i] * u[i]) / max;
3919566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_BLUE));
392c4762a1bSJed Brown }
393c4762a1bSJed Brown x += cnt * hx;
394c4762a1bSJed Brown }
3959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &x));
3969566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
3979566063dSJacob Faibussowitsch PetscCall(PetscDrawStringSetSize(draw, .2, .2));
3989566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw));
3999566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause));
4009566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw));
4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
402c4762a1bSJed Brown }
403c4762a1bSJed Brown
MyDestroy(PetscCtxRt Ctx)404*2a8381b2SBarry Smith PetscErrorCode MyDestroy(PetscCtxRt Ctx)
405d71ae5a4SJacob Faibussowitsch {
406*2a8381b2SBarry Smith AppCtx *ctx = *(AppCtx **)Ctx;
407c4762a1bSJed Brown
408c4762a1bSJed Brown PetscFunctionBegin;
4099566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ctx->ports));
4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
411c4762a1bSJed Brown }
412c4762a1bSJed Brown
413c4762a1bSJed Brown /*TEST
414c4762a1bSJed Brown
415c4762a1bSJed Brown test:
416c4762a1bSJed Brown args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 2 -square_initial
417c4762a1bSJed Brown
418c4762a1bSJed Brown test:
419c4762a1bSJed Brown suffix: 2
420c4762a1bSJed 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
421c4762a1bSJed Brown requires: x
422c4762a1bSJed Brown
423c4762a1bSJed Brown TEST*/
424