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