xref: /petsc/src/ts/tutorials/phasefield/potentials.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1 
2 static char help[] = "Plots the various potentials used in the examples.\n";
3 
4 #include <petscdmda.h>
5 #include <petscts.h>
6 #include <petscdraw.h>
7 
8 int main(int argc, char **argv)
9 {
10   PetscDrawLG         lg;
11   PetscInt            Mx = 100, i;
12   PetscReal           x, hx = .1 / Mx, pause, xx[3], yy[3];
13   PetscDraw           draw;
14   const char *const   legend[] = {"(1 - u^2)^2", "1 - u^2", "-(1 - u)log(1 - u)"};
15   PetscDrawAxis       axis;
16   PetscDrawViewPorts *ports;
17 
18   PetscFunctionBegin;
19   PetscFunctionBeginUser;
20   PetscCall(PetscInitialize(&argc, &argv, 0, help));
21   PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1200, 800));
22   PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 0, &lg));
23   PetscCall(PetscDrawLGGetDraw(lg, &draw));
24   PetscCall(PetscDrawCheckResizedWindow(draw));
25   PetscCall(PetscDrawViewPortsCreateRect(draw, 1, 2, &ports));
26   PetscCall(PetscDrawLGGetAxis(lg, &axis));
27   PetscCall(PetscDrawLGReset(lg));
28 
29   /*
30       Plot the  energies
31   */
32   PetscCall(PetscDrawLGSetDimension(lg, 3));
33   PetscCall(PetscDrawViewPortsSet(ports, 1));
34   x = .9;
35   for (i = 0; i < Mx; i++) {
36     xx[0] = xx[1] = xx[2] = x;
37     yy[0]                 = (1. - x * x) * (1. - x * x);
38     yy[1]                 = (1. - x * x);
39     yy[2]                 = -(1. - x) * PetscLogReal(1. - x);
40     PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
41     x += hx;
42   }
43   PetscCall(PetscDrawGetPause(draw, &pause));
44   PetscCall(PetscDrawSetPause(draw, 0.0));
45   PetscCall(PetscDrawAxisSetLabels(axis, "Energy", "", ""));
46   PetscCall(PetscDrawLGSetLegend(lg, legend));
47   PetscCall(PetscDrawLGDraw(lg));
48 
49   /*
50       Plot the  forces
51   */
52   PetscCall(PetscDrawViewPortsSet(ports, 0));
53   PetscCall(PetscDrawLGReset(lg));
54   x = .9;
55   for (i = 0; i < Mx; i++) {
56     xx[0] = xx[1] = xx[2] = x;
57     yy[0]                 = x * x * x - x;
58     yy[1]                 = -x;
59     yy[2]                 = 1.0 + PetscLogReal(1. - x);
60     PetscCall(PetscDrawLGAddPoint(lg, xx, yy));
61     x += hx;
62   }
63   PetscCall(PetscDrawAxisSetLabels(axis, "Derivative", "", ""));
64   PetscCall(PetscDrawLGSetLegend(lg, NULL));
65   PetscCall(PetscDrawLGDraw(lg));
66 
67   PetscCall(PetscDrawSetPause(draw, pause));
68   PetscCall(PetscDrawPause(draw));
69   PetscCall(PetscDrawViewPortsDestroy(ports));
70   PetscCall(PetscFinalize());
71   return 0;
72 }
73 
74 /*TEST
75 
76    test:
77      requires: x
78 
79 TEST*/
80