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 8*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 9*d71ae5a4SJacob Faibussowitsch { 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; 19327415f7SBarry Smith PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 219566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1200, 800)); 229566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 0, &lg)); 239566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetDraw(lg, &draw)); 249566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 259566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreateRect(draw, 1, 2, &ports)); 269566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lg, &axis)); 279566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* 30c4762a1bSJed Brown Plot the energies 31c4762a1bSJed Brown */ 329566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 3)); 339566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 1)); 34c4762a1bSJed Brown x = .9; 35c4762a1bSJed Brown for (i = 0; i < Mx; i++) { 36c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = x; 37c4762a1bSJed Brown yy[0] = (1. - x * x) * (1. - x * x); 38c4762a1bSJed Brown yy[1] = (1. - x * x); 39c4762a1bSJed Brown yy[2] = -(1. - x) * PetscLogReal(1. - x); 409566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 41c4762a1bSJed Brown x += hx; 42c4762a1bSJed Brown } 439566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause)); 449566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 459566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Energy", "", "")); 469566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, legend)); 479566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Plot the forces 51c4762a1bSJed Brown */ 529566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 0)); 539566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 54c4762a1bSJed Brown x = .9; 55c4762a1bSJed Brown for (i = 0; i < Mx; i++) { 56c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = x; 57c4762a1bSJed Brown yy[0] = x * x * x - x; 58c4762a1bSJed Brown yy[1] = -x; 59c4762a1bSJed Brown yy[2] = 1.0 + PetscLogReal(1. - x); 609566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 61c4762a1bSJed Brown x += hx; 62c4762a1bSJed Brown } 639566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Derivative", "", "")); 649566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause)); 689566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 699566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ports)); 709566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 71b122ec5aSJacob Faibussowitsch return 0; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /*TEST 75c4762a1bSJed Brown 76c4762a1bSJed Brown test: 77c4762a1bSJed Brown requires: x 78c4762a1bSJed Brown 79c4762a1bSJed Brown TEST*/ 80