1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Plots the various potentials used in the examples.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdmda.h> 5c4762a1bSJed Brown #include <petscts.h> 6c4762a1bSJed Brown #include <petscdraw.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **argv) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscDrawLG lg; 11c4762a1bSJed Brown PetscInt Mx = 100,i; 12c4762a1bSJed Brown PetscReal x,hx = .1/Mx,pause,xx[3],yy[3]; 13c4762a1bSJed Brown PetscDraw draw; 14c4762a1bSJed Brown const char *const legend[] = {"(1 - u^2)^2","1 - u^2","-(1 - u)log(1 - u)"}; 15c4762a1bSJed Brown PetscDrawAxis axis; 16c4762a1bSJed Brown PetscDrawViewPorts *ports; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscFunctionBegin; 19*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,0,help)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1200,800)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),0,&lg)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGGetDraw(lg,&draw)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawCheckResizedWindow(draw)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsCreateRect(draw,1,2,&ports)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGGetAxis(lg,&axis)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGReset(lg)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown Plot the energies 30c4762a1bSJed Brown */ 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetDimension(lg,3)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsSet(ports,1)); 33c4762a1bSJed Brown x = .9; 34c4762a1bSJed Brown for (i=0; i<Mx; i++) { 35c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = x; 36c4762a1bSJed Brown yy[0] = (1.-x*x)*(1. - x*x); 37c4762a1bSJed Brown yy[1] = (1. - x*x); 38c4762a1bSJed Brown yy[2] = -(1.-x)*PetscLogReal(1.-x); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy)); 40c4762a1bSJed Brown x += hx; 41c4762a1bSJed Brown } 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawGetPause(draw,&pause)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetPause(draw,0.0)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Energy","","")); 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetLegend(lg,legend)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGDraw(lg)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* 49c4762a1bSJed Brown Plot the forces 50c4762a1bSJed Brown */ 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsSet(ports,0)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGReset(lg)); 53c4762a1bSJed Brown x = .9; 54c4762a1bSJed Brown for (i=0; i<Mx; i++) { 55c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = x; 56c4762a1bSJed Brown yy[0] = x*x*x - x; 57c4762a1bSJed Brown yy[1] = -x; 58c4762a1bSJed Brown yy[2] = 1.0 + PetscLogReal(1. - x); 595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy)); 60c4762a1bSJed Brown x += hx; 61c4762a1bSJed Brown } 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Derivative","","")); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetLegend(lg,NULL)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGDraw(lg)); 65c4762a1bSJed Brown 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetPause(draw,pause)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawPause(draw)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsDestroy(ports)); 69*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 70*b122ec5aSJacob Faibussowitsch return 0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /*TEST 74c4762a1bSJed Brown 75c4762a1bSJed Brown test: 76c4762a1bSJed Brown requires: x 77c4762a1bSJed Brown 78c4762a1bSJed Brown TEST*/ 79