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