xref: /petsc/src/ts/tutorials/phasefield/potentials.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode            ierr;
12   PetscInt                  Mx = 100,i;
13   PetscReal                 x,hx = .1/Mx,pause,xx[3],yy[3];
14   PetscDraw                 draw;
15   const char *const         legend[] = {"(1 - u^2)^2","1 - u^2","-(1 - u)log(1 - u)"};
16   PetscDrawAxis             axis;
17   PetscDrawViewPorts        *ports;
18 
19   PetscFunctionBegin;
20   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
21   CHKERRQ(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1200,800));
22   CHKERRQ(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),0,&lg));
23   CHKERRQ(PetscDrawLGGetDraw(lg,&draw));
24   CHKERRQ(PetscDrawCheckResizedWindow(draw));
25   CHKERRQ(PetscDrawViewPortsCreateRect(draw,1,2,&ports));
26   CHKERRQ(PetscDrawLGGetAxis(lg,&axis));
27   CHKERRQ(PetscDrawLGReset(lg));
28 
29   /*
30       Plot the  energies
31   */
32   CHKERRQ(PetscDrawLGSetDimension(lg,3));
33   CHKERRQ(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     CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy));
41     x    += hx;
42   }
43   CHKERRQ(PetscDrawGetPause(draw,&pause));
44   CHKERRQ(PetscDrawSetPause(draw,0.0));
45   CHKERRQ(PetscDrawAxisSetLabels(axis,"Energy","",""));
46   CHKERRQ(PetscDrawLGSetLegend(lg,legend));
47   CHKERRQ(PetscDrawLGDraw(lg));
48 
49   /*
50       Plot the  forces
51   */
52   CHKERRQ(PetscDrawViewPortsSet(ports,0));
53   CHKERRQ(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     CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy));
61     x    += hx;
62   }
63   CHKERRQ(PetscDrawAxisSetLabels(axis,"Derivative","",""));
64   CHKERRQ(PetscDrawLGSetLegend(lg,NULL));
65   CHKERRQ(PetscDrawLGDraw(lg));
66 
67   CHKERRQ(PetscDrawSetPause(draw,pause));
68   CHKERRQ(PetscDrawPause(draw));
69   CHKERRQ(PetscDrawViewPortsDestroy(ports));
70   ierr = PetscFinalize();
71   return ierr;
72 }
73 
74 /*TEST
75 
76    test:
77      requires: x
78 
79 TEST*/
80