ROOT logo
// @(#)root/hist:$Id: TGraph2D.cxx,v 1.00
// Author: Olivier Couet

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

#include "Riostream.h"
#include "TROOT.h"
#include "TMath.h"
#include "TH2.h"
#include "TF2.h"
#include "TList.h"
#include "TGraph2D.h"
#include "TGraphDelaunay.h"
#include "TVirtualPad.h"
#include "TVirtualFitter.h"
#include "TPluginManager.h"
#include "TClass.h"
#include "TSystem.h"
#include <stdlib.h>

#include "HFitInterface.h"
#include "Fit/DataRange.h"
#include "Math/MinimizerOptions.h"

ClassImp(TGraph2D)


//______________________________________________________________________________
/* Begin_Html
<center><h2>Graph 2D class</h2></center>
A Graph2D is a graphics object made of three arrays X, Y and Z with the same
number of points each.
<p>
This class has different constructors:
<ol>
<p><li> With an array's dimension and three arrays x, y, and z:
<pre>
   TGraph2D *g = new TGraph2D(n, x, y, z);
</pre>
    x, y, z arrays can be doubles, floats, or ints.

<p><li> With an array's dimension only:
<pre>
   TGraph2D *g = new TGraph2D(n);
</pre>
   The internal arrays are then filled with SetPoint. The following line
   fills the the internal arrays at the position "i" with the values x,y,z.
<pre>
   g->SetPoint(i, x, y, z);
</pre>
<p><li> Without parameters:
<pre>
   TGraph2D *g = new TGraph2D();
</pre>
   again SetPoint must be used to fill the internal arrays.
<p><li> From a file:
<pre>
   TGraph2D *g = new TGraph2D("graph.dat");
</pre>
   Arrays are read from the ASCII file "graph.dat" according to a specifies
   format. The format's default value is "%lg %lg %lg"
</ol>

Note that in any of these three cases, SetPoint can be used to change a data
point or add a new one. If the data point index (i) is greater than the
current size of the internal arrays, they are automatically extended.
<p>
Specific drawing options can be used to paint a TGraph2D:

<table border=0>

<tr><th valign=top>"TRI"</th><td>
The Delaunay triangles are drawn using filled area.
An hidden surface drawing technique is used. The surface is
painted with the current fill area color. The edges of each
triangles are painted with the current line color.
</td></tr>

<tr><th valign=top>"TRIW</th><td>
The Delaunay triangles are drawn as wire frame
</td></tr>

<tr><th valign=top>"TRI1</th><td>
The Delaunay triangles are painted with color levels. The edges
of each triangles are painted with the current line color.
</td></tr>

<tr><th valign=top>"TRI2</th><td>
the Delaunay triangles are painted with color levels.
</td></tr>

<tr><th valign=top>"P"  </th><td>
Draw a marker at each vertex
</td></tr>

<tr><th valign=top>"P0" </th><td>
Draw a circle at each vertex. Each circle background is white.
</td></tr>

<tr><th valign=top>"PCOL" </th><td>
Draw a marker at each vertex. The color of each marker is
defined according to its Z position.
</td></tr>

<tr><th valign=top>"CONT" </th><td>
Draw contours.
</td></tr>

<tr><th valign=top>"LINE" </th><td>
Draw a 3D polyline.
</td></tr>

</table>

A TGraph2D can be also drawn with ANY options valid to draw a 2D histogram.
<p>
When a TGraph2D is drawn with one of the 2D histogram drawing option,
a intermediate 2D histogram is filled using the Delaunay triangles
technique to interpolate the data set.
<p>
TGraph2D linearly interpolate a Z value for any (X,Y) point given some
existing (X,Y,Z) points. The existing (X,Y,Z) points can be randomly
scattered. The algorithm works by joining the existing points to make
Delaunay triangles in (X,Y). These are then used to define flat planes
in (X,Y,Z) over which to interpolate. The interpolated surface thus takes
the form of tessellating triangles at various angles. Output can take the
form of a 2D histogram or a vector. The triangles found can be drawn in 3D.
<p>
This software cannot be guaranteed to work under all circumstances. They
were originally written to work with a few hundred points in an XY space
with similar X and Y ranges.
<p>
Example:

End_Html
Begin_Macro(source)
{
   TCanvas *c = new TCanvas("c","Graph2D example",0,0,600,400);
   Double_t x, y, z, P = 6.;
   Int_t np = 200;
   TGraph2D *dt = new TGraph2D();
   TRandom *r = new TRandom();
   for (Int_t N=0; N<np; N++) {
      x = 2*P*(r->Rndm(N))-P;
      y = 2*P*(r->Rndm(N))-P;
      z = (sin(x)/x)*(sin(y)/y)+0.2;
      dt->SetPoint(N,x,y,z);
   }
   gStyle->SetPalette(1);
   dt->Draw("surf1");
   return c;
}
End_Macro
Begin_Html

2D graphs can be fitted as shown by the following example:

End_Html
Begin_Macro(source)
../../../tutorials/fit/graph2dfit.C
End_Macro
Begin_Html

Example showing the PCOL option.

End_Html
Begin_Macro(source)
{
   TCanvas *c1 = new TCanvas("c1","Graph2D example",0,0,600,400);
   Double_t P = 5.;
   Int_t npx  = 20 ;
   Int_t npy  = 20 ;
   Double_t x = -P;
   Double_t y = -P;
   Double_t z;
   Int_t k = 0;
   Double_t dx = (2*P)/npx;
   Double_t dy = (2*P)/npy;
   TGraph2D *dt = new TGraph2D(npx*npy);
   dt->SetNpy(41);
   dt->SetNpx(40);
   for (Int_t i=0; i<npx; i++) {
      for (Int_t j=0; j<npy; j++) {
         z = sin(sqrt(x*x+y*y))+1;
         dt->SetPoint(k,x,y,z);
         k++;
         y = y+dx;
      }
      x = x+dx;
      y = -P;
   }
   gStyle->SetPalette(1);
   dt->SetMarkerStyle(20);
   dt->Draw("pcol");
   return c1;
}
End_Macro
Begin_Html

<h3>Definition of Delaunay triangulation (After B. Delaunay)</h3>
For a set S of points in the Euclidean plane, the unique triangulation DT(S)
of S such that no point in S is inside the circumcircle of any triangle in
DT(S). DT(S) is the dual of the Voronoi diagram of S.
If n is the number of points in S, the Voronoi diagram of S is the partitioning
of the plane containing S points into n convex polygons such that each polygon
contains exactly one point and every point in a given polygon is closer to its
central point than to any other. A Voronoi diagram is sometimes also known as
a Dirichlet tessellation.

<img src="gif/dtvd.gif">

<br>
<a href="http://www.cs.cornell.edu/Info/People/chew/Delaunay.html">This applet</a>
gives a nice practical view of Delaunay triangulation and Voronoi diagram.
End_Html */


//______________________________________________________________________________
TGraph2D::TGraph2D()
         : TNamed("Graph2D","Graph2D"), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(0)
{
   // Graph2D default constructor

   fSize      = 0;
   fMargin    = 0.;
   fNpx       = 40;
   fNpy       = 40;
   fDirectory = 0;
   fHistogram = 0;
   fMaximum   = -1111;
   fMinimum   = -1111;
   fX         = 0;
   fY         = 0;
   fZ         = 0;
   fZout      = 0;
   fMaxIter   = 100000;
   fPainter   = 0;
   fFunctions = new TList;
   fUserHisto = kFALSE;
}


//______________________________________________________________________________
TGraph2D::TGraph2D(Int_t n, Int_t *x, Int_t *y, Int_t *z)
         : TNamed("Graph2D","Graph2D"), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(n)
{
   // Graph2D constructor with three vectors of ints as input.

   Build(n);

   // Copy the input vectors into local arrays
   for (Int_t i=0; i<fNpoints; ++i) {
      fX[i] = (Double_t)x[i];
      fY[i] = (Double_t)y[i];
      fZ[i] = (Double_t)z[i];
   }
}


//______________________________________________________________________________
TGraph2D::TGraph2D(Int_t n, Float_t *x, Float_t *y, Float_t *z)
         : TNamed("Graph2D","Graph2D"), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(n)
{
   // Graph2D constructor with three vectors of floats as input.

   Build(n);

   // Copy the input vectors into local arrays
   for (Int_t i=0; i<fNpoints; ++i) {
      fX[i] = x[i];
      fY[i] = y[i];
      fZ[i] = z[i];
   }
}


//______________________________________________________________________________
TGraph2D::TGraph2D(Int_t n, Double_t *x, Double_t *y, Double_t *z)
         : TNamed("Graph2D","Graph2D"), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(n)
{
   // Graph2D constructor with three vectors of doubles as input.

   Build(n);

   // Copy the input vectors into local arrays
   for (Int_t i=0; i<fNpoints; ++i) {
      fX[i] = x[i];
      fY[i] = y[i];
      fZ[i] = z[i];
   }
}


//______________________________________________________________________________
TGraph2D::TGraph2D(TH2 *h2)
         : TNamed("Graph2D","Graph2D"), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(0)
{
   // Graph2D constructor with a TH2 (h2) as input.
   // Only the h2's bins within the X and Y axis ranges are used.
   // Empty bins, recognized when both content and errors are zero, are excluded.

   Build(10);

   SetName(h2->GetName());
   SetTitle(h2->GetTitle());

   TAxis *xaxis = h2->GetXaxis();
   TAxis *yaxis = h2->GetYaxis();
   Int_t xfirst = xaxis->GetFirst();
   Int_t xlast  = xaxis->GetLast();
   Int_t yfirst = yaxis->GetFirst();
   Int_t ylast  = yaxis->GetLast();


   Double_t x, y, z;
   Int_t k=0;

   for (Int_t i=xfirst; i<= xlast; i++) {
      for (Int_t j=yfirst; j<= ylast; j++) {
         x = xaxis->GetBinCenter(i);
         y = yaxis->GetBinCenter(j);
         z = h2->GetBinContent(i,j);
         Double_t ez = h2->GetBinError(i,j);
         if (z != 0. || ez != 0) {
            SetPoint(k, x, y, z);
            k++;
         }
      }
   }
}


//______________________________________________________________________________
TGraph2D::TGraph2D(const char *name,const char *title,
                   Int_t n, Double_t *x, Double_t *y, Double_t *z)
         : TNamed(name,title), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(n)
{
   // Graph2D constructor with name, title and three vectors of doubles as input.
   // name   : name of 2D graph (avoid blanks)
   // title  : 2D graph title
   //          if title is of the form "stringt;stringx;stringy;stringz"
   //          the 2D graph title is set to stringt, the x axis title to stringx,
   //          the y axis title to stringy,etc

   Build(n);

   // Copy the input vectors into local arrays
   for (Int_t i=0; i<fNpoints; ++i) {
      fX[i] = x[i];
      fY[i] = y[i];
      fZ[i] = z[i];
   }
}


//______________________________________________________________________________
TGraph2D::TGraph2D(Int_t n)
         : TNamed("Graph2D","Graph2D"), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(0)
{
   // Graph2D constructor. The arrays fX, fY and fZ should be filled via
   // calls to SetPoint

   Build(n);
}


//______________________________________________________________________________
TGraph2D::TGraph2D(const char *filename, const char *format, Option_t *)
         : TNamed("Graph2D",filename), TAttLine(1,1,1), TAttFill(0,1001),
           TAttMarker(), fNpoints(0)
{
   // Graph2D constructor reading input from filename
   // filename is assumed to contain at least three columns of numbers

   Build(100);

   Double_t x,y,z;
   FILE *fp = fopen(filename,"r");
   if (!fp) {
      MakeZombie();
      Error("TGraph2D", "Cannot open file: %s, TGraph2D is Zombie",filename);
      return;
   }
   char line[80];
   Int_t np = 0;
   while (fgets(line,80,fp)) {
      sscanf(&line[0],format,&x, &y, &z);
      SetPoint(np,x,y,z);
      np++;
   }

   fclose(fp);
}


//______________________________________________________________________________
TGraph2D::TGraph2D(const TGraph2D &g)
         : TNamed(g), TAttLine(g), TAttFill(g), TAttMarker(g)
{
   // Graph2D copy constructor.

   fNpoints = g.fNpoints;
   Build(fNpoints);

   for (Int_t n=0; n<fNpoints; n++) {
      fX[n] = g.fX[n];
      fY[n] = g.fY[n];
      fZ[n] = g.fZ[n];
   }
}


//______________________________________________________________________________
TGraph2D::~TGraph2D()
{
   // TGraph2D destructor.

   Clear();
}


//______________________________________________________________________________
TGraph2D& TGraph2D::operator=(const TGraph2D &g)
{
   // Graph2D operator "="

   if (this == &g) return *this;

   Clear();

   fNpoints = g.fNpoints;
   Build(fNpoints);

   for (Int_t n=0; n<fNpoints; n++) {
      fX[n] = g.fX[n];
      fY[n] = g.fY[n];
      fZ[n] = g.fZ[n];
   }
   return *this;
}

//______________________________________________________________________________
void TGraph2D::Build(Int_t n)
{
   // Creates the 2D graph basic data structure

   if (n <= 0) {
      Error("TGraph2D", "Invalid number of points (%d)", n);
      return;
   }

   fSize      = n;
   fMargin    = 0.;
   fNpx       = 40;
   fNpy       = 40;
   fDirectory = 0;
   fHistogram = 0;
   fMaximum   = -1111;
   fMinimum   = -1111;
   fX         = new Double_t[fSize];
   fY         = new Double_t[fSize];
   fZ         = new Double_t[fSize];
   fZout      = 0;
   fMaxIter   = 100000;
   fFunctions = new TList;
   fPainter   = 0;
   fUserHisto = kFALSE;

   if (TH1::AddDirectoryStatus()) {
      fDirectory = gDirectory;
      if (fDirectory) {
         fDirectory->Append(this,kTRUE);
      }
   }
}


//______________________________________________________________________________
void TGraph2D::Clear(Option_t * /*option = "" */)
{
   // Free all memory allocated by this object.

   delete [] fX; fX = 0;
   delete [] fY; fY = 0;
   delete [] fZ; fZ = 0;
   delete fHistogram; fHistogram = 0;
   if (fFunctions) {
      fFunctions->SetBit(kInvalidObject);
      fFunctions->Delete();
      delete fFunctions;
      fFunctions = 0;
   }
   if (fDirectory) {
      fDirectory->Remove(this);
      fDirectory = 0;
   }
   delete fPainter;
   fPainter   = 0;
}


//______________________________________________________________________________
void TGraph2D::DirectoryAutoAdd(TDirectory *dir)
{
   // Perform the automatic addition of the graph to the given directory
   //
   // Note this function is called in place when the semantic requires
   // this object to be added to a directory (I.e. when being read from
   // a TKey or being Cloned)

   Bool_t addStatus = TH1::AddDirectoryStatus();
   if (addStatus) {
      SetDirectory(dir);
   }
}


//______________________________________________________________________________
Int_t TGraph2D::DistancetoPrimitive(Int_t px, Int_t py)
{
   // Computes distance from point px,py to a graph

   Int_t distance = 9999;
   if (fHistogram) distance = fHistogram->DistancetoPrimitive(px,py);
   return distance;
}


//______________________________________________________________________________
void TGraph2D::Draw(Option_t *option)
{
   // Specific drawing options can be used to paint a TGraph2D:
   //
   //   "TRI"  : The Delaunay triangles are drawn using filled area.
   //            An hidden surface drawing technique is used. The surface is
   //            painted with the current fill area color. The edges of each
   //            triangles are painted with the current line color.
   //   "TRIW" : The Delaunay triangles are drawn as wire frame
   //   "TRI1" : The Delaunay triangles are painted with color levels. The edges
   //            of each triangles are painted with the current line color.
   //   "TRI2" : the Delaunay triangles are painted with color levels.
   //   "P"    : Draw a marker at each vertex
   //   "P0"   : Draw a circle at each vertex. Each circle background is white.
   //   "PCOL" : Draw a marker at each vertex. The color of each marker is
   //            defined according to its Z position.
   //   "CONT" : Draw contours
   //   "LINE" : Draw a 3D polyline
   //
   // A TGraph2D can be also drawn with ANY options valid to draw a 2D histogram.
   //
   // When a TGraph2D is drawn with one of the 2D histogram drawing option,
   // a intermediate 2D histogram is filled using the Delaunay triangles
   // technique to interpolate the data set.

   TString opt = option;
   opt.ToLower();
   if (gPad) {
      if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
      if (!opt.Contains("same")) {
         //the following statement is necessary in case one attempts to draw
         //a temporary histogram already in the current pad
         if (TestBit(kCanDelete)) gPad->GetListOfPrimitives()->Remove(this);
         gPad->Clear();
      }
   }
   AppendPad(opt.Data());
}


//______________________________________________________________________________
void TGraph2D::ExecuteEvent(Int_t event, Int_t px, Int_t py)
{
   // Executes action corresponding to one event

   if (fHistogram) fHistogram->ExecuteEvent(event, px, py);
}


//______________________________________________________________________________
TObject *TGraph2D::FindObject(const char *name) const
{
   // search object named name in the list of functions

   if (fFunctions) return fFunctions->FindObject(name);
   return 0;
}


//______________________________________________________________________________
TObject *TGraph2D::FindObject(const TObject *obj) const
{
   // search object obj in the list of functions

   if (fFunctions) return fFunctions->FindObject(obj);
   return 0;
}


//______________________________________________________________________________
TFitResultPtr TGraph2D::Fit(const char *fname, Option_t *option, Option_t *)
{
   // Fits this graph with function with name fname
   // Predefined functions such as gaus, expo and poln are automatically
   // created by ROOT.
   // fname can also be a formula, accepted by the linear fitter (linear parts divided
   // by "++" sign), for example "x++sin(y)" for fitting "[0]*x+[1]*sin(y)"


   char *linear;
   linear= (char*)strstr(fname, "++");
   TF2 *f2=0;
   if (linear)
      f2=new TF2(fname, fname);
   else{
      f2 = (TF2*)gROOT->GetFunction(fname);
      if (!f2) { Printf("Unknown function: %s",fname); return -1; }
   }
   return Fit(f2,option,"");

}


//______________________________________________________________________________
TFitResultPtr TGraph2D::Fit(TF2 *f2, Option_t *option, Option_t *)
{
   // Fits this 2D graph with function f2
   //
   //  f2 is an already predefined function created by TF2.
   //  Predefined functions such as gaus, expo and poln are automatically
   //  created by ROOT.
   //
   //  The list of fit options is given in parameter option.
   //     option = "W" Set all weights to 1; ignore error bars
   //            = "U" Use a User specified fitting algorithm (via SetFCN)
   //            = "Q" Quiet mode (minimum printing)
   //            = "V" Verbose mode (default is between Q and V)
   //            = "R" Use the Range specified in the function range
   //            = "N" Do not store the graphics function, do not draw
   //            = "0" Do not plot the result of the fit. By default the fitted function
   //                  is drawn unless the option "N" above is specified.
   //            = "+" Add this new fitted function to the list of fitted functions
   //                  (by default, any previous function is deleted)
   //            = "C" In case of linear fitting, not calculate the chisquare
   //                  (saves time)
   //            = "EX0" When fitting a TGraphErrors do not consider errors in the coordinate
   //            = "ROB" In case of linear fitting, compute the LTS regression
   //                     coefficients (robust (resistant) regression), using
   //                     the default fraction of good points
   //              "ROB=0.x" - compute the LTS regression coefficients, using
   //                           0.x as a fraction of good points
   //            = "S"  The result of the fit is returned in the TFitResultPtr 
   //                     (see below Access to the Fit Result) 
   //
   //  In order to use the Range option, one must first create a function
   //  with the expression to be fitted. For example, if your graph2d
   //  has a defined range between -4 and 4 and you want to fit a gaussian
   //  only in the interval 1 to 3, you can do:
   //       TF2 *f2 = new TF2("f2","gaus",1,3);
   //       graph2d->Fit("f2","R");
   //
   //
   //  Setting initial conditions
   //  ==========================
   //  Parameters must be initialized before invoking the Fit function.
   //  The setting of the parameter initial values is automatic for the
   //  predefined functions : poln, expo, gaus. One can however disable
   //  this automatic computation by specifying the option "B".
   //  You can specify boundary limits for some or all parameters via
   //       f2->SetParLimits(p_number, parmin, parmax);
   //  if parmin>=parmax, the parameter is fixed
   //  Note that you are not forced to fix the limits for all parameters.
   //  For example, if you fit a function with 6 parameters, you can do:
   //    func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
   //    func->SetParLimits(4,-10,-4);
   //    func->SetParLimits(5, 1,1);
   //  With this setup, parameters 0->3 can vary freely
   //  Parameter 4 has boundaries [-10,-4] with initial value -8
   //  Parameter 5 is fixed to 100.
   //
   //  Fit range
   //  =========
   //  The fit range can be specified in two ways:
   //    - specify rxmax > rxmin (default is rxmin=rxmax=0)
   //    - specify the option "R". In this case, the function will be taken
   //      instead of the full graph range.
   //
   //  Changing the fitting function
   //  =============================
   //  By default the fitting function Graph2DFitChisquare is used.
   //  To specify a User defined fitting function, specify option "U" and
   //  call the following functions:
   //    TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
   //  where MyFittingFunction is of type:
   //  extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
   //
   //  Associated functions
   //  ====================
   //  One or more object (typically a TF2*) can be added to the list
   //  of functions (fFunctions) associated to each graph.
   //  When TGraph::Fit is invoked, the fitted function is added to this list.
   //  Given a graph gr, one can retrieve an associated function
   //  with:  TF2 *myfunc = gr->GetFunction("myfunc");
   //
   //  Access to the fit results
   //  =========================
   //  The function returns a TFitResultPtr which can hold a  pointer to a TFitResult object.
   //  By default the TFitResultPtr contains only the status of the fit and it converts automatically to an
   //  integer. If the option "S" is instead used, TFitResultPtr contains the TFitResult and behaves as a smart 
   //  pointer to it. For example one can do: 
   //     TFitResult r    = graph->Fit("myFunc","S");
   //     TMatrixDSym cov = r->GetCovarianceMatrix();  //  to access the covariance matrix
   //     Double_t par0   = r->Value(0); // retrieve the value for the parameter 0 
   //     Double_t err0   = r->Error(0); // retrieve the error for the parameter 0 
   //     r->Print("V");     // print full information of fit including covariance matrix
   //     r->Write();        // store the result in a file
   // 
   //  The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also 
   //  from the fitted function. 
   //  If the graph is made persistent, the list of
   //  associated functions is also persistent. Given a pointer (see above)
   //  to an associated function myfunc, one can retrieve the function/fit
   //  parameters with calls such as:
   //    Double_t chi2 = myfunc->GetChisquare();
   //    Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
   //    Double_t err0 = myfunc->GetParError(0);  //error on first parameter
   //
   //  Fit Statistics
   //  ==============
   //  You can change the statistics box to display the fit parameters with
   //  the TStyle::SetOptFit(mode) method. This mode has four digits.
   //  mode = pcev  (default = 0111)
   //    v = 1;  print name/values of parameters
   //    e = 1;  print errors (if e=1, v must be 1)
   //    c = 1;  print Chisquare/Number of degress of freedom
   //    p = 1;  print Probability
   //
   //  For example: gStyle->SetOptFit(1011);
   //  prints the fit probability, parameter names/values, and errors.
   //  You can change the position of the statistics box with these lines
   //  (where g is a pointer to the TGraph):
   //
   //  Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
   //  Root > st->SetX1NDC(newx1); //new x start position
   //  Root > st->SetX2NDC(newx2); //new x end position

   // internal graph2D fitting methods
   Foption_t fitOption;
   Option_t *goption = "";
   ROOT::Fit::FitOptionsMake(option,fitOption);

   // create range and minimizer options with default values 
   ROOT::Fit::DataRange range(2); 
   ROOT::Math::MinimizerOptions minOption; 
   return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range); 
}


//______________________________________________________________________________
void TGraph2D::FitPanel()
{
   // Display a GUI panel with all graph fit options.
   //
   //   See class TFitEditor for example
   if (!gPad)
      gROOT->MakeDefCanvas();

   if (!gPad) {
      Error("FitPanel", "Unable to create a default canvas");
      return;
   }

   // use plugin manager to create instance of TFitEditor
   TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
   if (handler && handler->LoadPlugin() != -1) {
      if (handler->ExecPlugin(2, gPad, this) == 0)
         Error("FitPanel", "Unable to crate the FitPanel");
   }
   else
      Error("FitPanel", "Unable to find the FitPanel plug-in");

}


//______________________________________________________________________________
TAxis *TGraph2D::GetXaxis() const
{
   // Get x axis of the graph.

   //if (!gPad) return 0;
   TH1 *h = ((TGraph2D*)this)->GetHistogram();
   if (!h) return 0;
   return h->GetXaxis();
}


//______________________________________________________________________________
TAxis *TGraph2D::GetYaxis() const
{
   // Get y axis of the graph.

   //if (!gPad) return 0;
   TH1 *h = ((TGraph2D*)this)->GetHistogram();
   if (!h) return 0;
   return h->GetYaxis();
}


//______________________________________________________________________________
TAxis *TGraph2D::GetZaxis() const
{
   // Get z axis of the graph.

   //if (!gPad) return 0;
   TH1 *h = ((TGraph2D*)this)->GetHistogram();
   if (!h) return 0;
   return h->GetZaxis();
}


//______________________________________________________________________________
TList *TGraph2D::GetContourList(Double_t contour)
{
   // Returns the X and Y graphs building a contour. A contour level may
   // consist in several parts not connected to each other. This function
   // returns them in a graphs' list.

   if (fNpoints <= 0) {
      Error("GetContourList", "Empty TGraph2D");
      return 0;
   }

   if(!fHistogram) GetHistogram("empty");

   if (!fPainter) fPainter = fHistogram->GetPainter();

   return fPainter->GetContourList(contour);
}


//______________________________________________________________________________
Double_t TGraph2D::GetErrorX(Int_t) const
{
   // This function is called by Graph2DFitChisquare.
   // It always returns a negative value. Real implementation in TGraph2DErrors

   return -1;
}


//______________________________________________________________________________
Double_t TGraph2D::GetErrorY(Int_t) const
{
   // This function is called by Graph2DFitChisquare.
   // It always returns a negative value. Real implementation in TGraph2DErrors

   return -1;
}


//______________________________________________________________________________
Double_t TGraph2D::GetErrorZ(Int_t) const
{
   // This function is called by Graph2DFitChisquare.
   // It always returns a negative value. Real implementation in TGraph2DErrors

   return -1;
}


//______________________________________________________________________________
TH2D *TGraph2D::GetHistogram(Option_t *option)
{
   // By default returns a pointer to the Delaunay histogram. If fHistogram
   // doesn't exist, books the 2D histogram fHistogram with a margin around
   // the hull. Calls TGraphDelaunay::Interpolate at each bin centre to build up
   // an interpolated 2D histogram.
   // If the "empty" option is selected, returns an empty histogram booked with
   // the limits of fX, fY and fZ. This option is used when the data set is
   // drawn with markers only. In that particular case there is no need to
   // find the Delaunay triangles.

   if (fNpoints <= 0) {
      Error("GetHistogram", "Empty TGraph2D");
      return 0;
   }

   TString opt = option;
   opt.ToLower();
   Bool_t empty = opt.Contains("empty");

   if (fHistogram) {
      if (!empty && fHistogram->GetEntries()==0) {
         if (!fUserHisto) {
            delete fHistogram;
            fHistogram = 0;
         }
      } else {
         return fHistogram;
      }
   }

   Double_t hxmax, hymax, hxmin, hymin;

   // Book fHistogram if needed. It is not added in the current directory
   if (!fUserHisto) {
      Bool_t add = TH1::AddDirectoryStatus();
      TH1::AddDirectory(kFALSE);
      Double_t xmax  = GetXmax();
      Double_t ymax  = GetYmax();
      Double_t xmin  = GetXmin();
      Double_t ymin  = GetYmin();
      hxmin = xmin-fMargin*(xmax-xmin);
      hymin = ymin-fMargin*(ymax-ymin);
      hxmax = xmax+fMargin*(xmax-xmin);
      hymax = ymax+fMargin*(ymax-ymin);
      fHistogram = new TH2D(GetName(),GetTitle(),
                            fNpx ,hxmin, hxmax,
                            fNpy, hymin, hymax);
      TH1::AddDirectory(add);
      fHistogram->SetBit(TH1::kNoStats);
   } else {
      hxmin = fHistogram->GetXaxis()->GetXmin();
      hymin = fHistogram->GetYaxis()->GetXmin();
      hxmax = fHistogram->GetXaxis()->GetXmax();
      hymax = fHistogram->GetYaxis()->GetXmax();
   }

   // Add a TGraphDelaunay in the list of the fHistogram's functions
   TGraphDelaunay *dt = new TGraphDelaunay(this);
   dt->SetMaxIter(fMaxIter);
   dt->SetMarginBinsContent(fZout);
   TList *hl = fHistogram->GetListOfFunctions();
   hl->Add(dt);

   // Option "empty" is selected. An empty histogram is returned.
   if (empty) {
      if (fMinimum != -1111) {
         fHistogram->SetMinimum(fMinimum);
      } else {
         fHistogram->SetMinimum(GetZmin());
      }
      if (fMaximum != -1111) {
         fHistogram->SetMaximum(fMaximum);
      } else {
         fHistogram->SetMaximum(GetZmax());
      }
      return fHistogram;
   }

   Double_t dx = (hxmax-hxmin)/fNpx;
   Double_t dy = (hymax-hymin)/fNpy;

   Double_t x, y, z;
   for (Int_t ix=1; ix<=fNpx; ix++) {
      x  = hxmin+(ix-0.5)*dx;
      for (Int_t iy=1; iy<=fNpy; iy++) {
         y  = hymin+(iy-0.5)*dy;
         z  = dt->ComputeZ(x, y);
         fHistogram->Fill(x, y, z);
      }
   }

   if (fMinimum != -1111) fHistogram->SetMinimum(fMinimum);
   if (fMaximum != -1111) fHistogram->SetMaximum(fMaximum);

   return fHistogram;
}


//______________________________________________________________________________
Double_t TGraph2D::GetXmax() const
{
   // Returns the X maximum

   Double_t v = fX[0];
   for (Int_t i=1; i<fNpoints; i++) if (fX[i]>v) v=fX[i];
   return v;
}


//______________________________________________________________________________
Double_t TGraph2D::GetXmin() const
{
   // Returns the X minimum

   Double_t v = fX[0];
   for (Int_t i=1; i<fNpoints; i++) if (fX[i]<v) v=fX[i];
   return v;
}


//______________________________________________________________________________
Double_t TGraph2D::GetYmax() const
{
   // Returns the Y maximum

   Double_t v = fY[0];
   for (Int_t i=1; i<fNpoints; i++) if (fY[i]>v) v=fY[i];
   return v;
}


//______________________________________________________________________________
Double_t TGraph2D::GetYmin() const
{
   // Returns the Y minimum

   Double_t v = fY[0];
   for (Int_t i=1; i<fNpoints; i++) if (fY[i]<v) v=fY[i];
   return v;
}


//______________________________________________________________________________
Double_t TGraph2D::GetZmax() const
{
   // Returns the Z maximum

   Double_t v = fZ[0];
   for (Int_t i=1; i<fNpoints; i++) if (fZ[i]>v) v=fZ[i];
   return v;
}


//______________________________________________________________________________
Double_t TGraph2D::GetZmin() const
{
   // Returns the Z minimum

   Double_t v = fZ[0];
   for (Int_t i=1; i<fNpoints; i++) if (fZ[i]<v) v=fZ[i];
   return v;
}


//______________________________________________________________________________
Double_t TGraph2D::Interpolate(Double_t x, Double_t y)
{
   // Finds the z value at the position (x,y) thanks to
   // the Delaunay interpolation.

   if (fNpoints <= 0) {
      Error("Interpolate", "Empty TGraph2D");
      return 0;
   }

   TGraphDelaunay *dt;

   if(!fHistogram) GetHistogram("empty");

   TList *hl = fHistogram->GetListOfFunctions();
   dt = (TGraphDelaunay*)hl->FindObject("TGraphDelaunay");

   return dt->ComputeZ(x, y);
}


//______________________________________________________________________________
void TGraph2D::Paint(Option_t *option)
{
   // Paints this 2D graph with its current attributes

   if (fNpoints <= 0) {
      Error("Paint", "Empty TGraph2D");
      return;
   }

   TString opt = option;
   opt.ToLower();
   if (opt.Contains("p") && !opt.Contains("tri")) {
      if (!opt.Contains("pol") &&
          !opt.Contains("sph") &&
          !opt.Contains("psr")) opt.Append("tri0");
   }

   if (opt.Contains("line") && !opt.Contains("tri")) opt.Append("tri0");

   if (opt.Contains("tri0")) {
      GetHistogram("empty");
   } else {
      GetHistogram();
   }

   fHistogram->SetLineColor(GetLineColor());
   fHistogram->SetLineStyle(GetLineStyle());
   fHistogram->SetLineWidth(GetLineWidth());
   fHistogram->SetFillColor(GetFillColor());
   fHistogram->SetFillStyle(GetFillStyle());
   fHistogram->SetMarkerColor(GetMarkerColor());
   fHistogram->SetMarkerStyle(GetMarkerStyle());
   fHistogram->SetMarkerSize(GetMarkerSize());
   fHistogram->Paint(opt.Data());
}


//______________________________________________________________________________
TH1 *TGraph2D::Project(Option_t *option) const
{
   // Projects a 2-d graph into 1 or 2-d histograms depending on the
   // option parameter
   // option may contain a combination of the characters x,y,z
   // option = "x" return the x projection into a TH1D histogram
   // option = "y" return the y projection into a TH1D histogram
   // option = "xy" return the x versus y projection into a TH2D histogram
   // option = "yx" return the y versus x projection into a TH2D histogram

   if (fNpoints <= 0) {
      Error("Project", "Empty TGraph2D");
      return 0;
   }

   TString opt = option; opt.ToLower();

   Int_t pcase = 0;
   if (opt.Contains("x"))  pcase = 1;
   if (opt.Contains("y"))  pcase = 2;
   if (opt.Contains("xy")) pcase = 3;
   if (opt.Contains("yx")) pcase = 4;

   // Create the projection histogram
   TH1D *h1 = 0;
   TH2D *h2 = 0;
   Int_t nch = strlen(GetName()) +opt.Length() +2;
   char *name = new char[nch];
   sprintf(name,"%s_%s",GetName(),option);
   nch = strlen(GetTitle()) +opt.Length() +2;
   char *title = new char[nch];
   sprintf(title,"%s_%s",GetTitle(),option);

   Double_t hxmin = GetXmin();
   Double_t hxmax = GetXmax();
   Double_t hymin = GetYmin();
   Double_t hymax = GetYmax();

   switch (pcase) {
      case 1:
         // "x"
         h1 = new TH1D(name,title,fNpx,hxmin,hxmax);
         break;
      case 2:
         // "y"
         h1 = new TH1D(name,title,fNpy,hymin,hymax);
         break;
      case 3:
         // "xy"
         h2 = new TH2D(name,title,fNpx,hxmin,hxmax,fNpy,hymin,hymax);
         break;
      case 4:
         // "yx"
         h2 = new TH2D(name,title,fNpy,hymin,hymax,fNpx,hxmin,hxmax);
         break;
   }

   delete [] name;
   delete [] title;
   TH1 *h = h1;
   if (h2) h = h2;
   if (h == 0) return 0;

   // Fill the projected histogram
   Double_t entries = 0;
   for (Int_t n=0; n<fNpoints; n++) {
      switch (pcase) {
         case 1:
            // "x"
            h1->Fill(fX[n], fZ[n]);
            break;
         case 2:
            // "y"
            h1->Fill(fY[n], fZ[n]);
            break;
         case 3:
            // "xy"
            h2->Fill(fX[n], fY[n], fZ[n]);
            break;
         case 4:
            // "yx"
            h2->Fill(fY[n], fX[n], fZ[n]);
            break;
      }
      entries += fZ[n];
   }
   h->SetEntries(entries);
   return h;
}


//______________________________________________________________________________
Int_t TGraph2D::RemovePoint(Int_t ipoint)
{
   // Deletes point number ipoint

   if (ipoint < 0) return -1;
   if (ipoint >= fNpoints) return -1;

   fNpoints--;
   Double_t *newX = new Double_t[fNpoints];
   Double_t *newY = new Double_t[fNpoints];
   Double_t *newZ = new Double_t[fNpoints];
   Int_t j = -1;
   for (Int_t i=0;i<fNpoints+1;i++) {
      if (i == ipoint) continue;
      j++;
      newX[j] = fX[i];
      newY[j] = fY[i];
      newZ[j] = fZ[i];
   }
   delete [] fX;
   delete [] fY;
   delete [] fZ;
   fX = newX;
   fY = newY;
   fZ = newZ;
   fSize = fNpoints;
   if (fHistogram) {delete fHistogram; fHistogram = 0;}
   return ipoint;
}


//______________________________________________________________________________
void TGraph2D::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
   // Saves primitive as a C++ statement(s) on output stream out

   char quote = '"';
   out<<"   "<<endl;
   if (gROOT->ClassSaved(TGraph2D::Class())) {
      out<<"   ";
   } else {
      out<<"   TGraph2D *";
   }

   out<<"graph2d = new TGraph2D("<<fNpoints<<");"<<endl;
   out<<"   graph2d->SetName("<<quote<<GetName()<<quote<<");"<<endl;
   out<<"   graph2d->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;

   if (fDirectory == 0) {
      out<<"   "<<GetName()<<"->SetDirectory(0);"<<endl;
   }

   SaveFillAttributes(out,"graph2d",0,1001);
   SaveLineAttributes(out,"graph2d",1,1,1);
   SaveMarkerAttributes(out,"graph2d",1,1,1);

   for (Int_t i=0;i<fNpoints;i++) {
      out<<"   graph2d->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<","<<fZ[i]<<");"<<endl;
   }

   // save list of functions
   TIter next(fFunctions);
   TObject *obj;
   while ((obj=next())) {
      obj->SavePrimitive(out,"nodraw");
      out<<"   graph2d->GetListOfFunctions()->Add("<<obj->GetName()<<");"<<endl;
      if (obj->InheritsFrom("TPaveStats")) {
         out<<"   ptstats->SetParent(graph2d->GetListOfFunctions());"<<endl;
      }
   }

   out<<"   graph2d->Draw("<<quote<<option<<quote<<");"<<endl;
}


//______________________________________________________________________________
void TGraph2D::Set(Int_t n)
{
   // Set number of points in the 2D graph.
   // Existing coordinates are preserved.
   // New coordinates above fNpoints are preset to 0.

   if (n < 0) n = 0;
   if (n == fNpoints) return;
   if (n >  fNpoints) SetPoint(n,0,0,0);
   fNpoints = n;
}


//______________________________________________________________________________
void TGraph2D::SetDirectory(TDirectory *dir)
{
   // By default when an 2D graph is created, it is added to the list
   // of 2D graph objects in the current directory in memory.
   // Remove reference to this 2D graph from current directory and add
   // reference to new directory dir. dir can be 0 in which case the
   // 2D graph does not belong to any directory.

   if (fDirectory == dir) return;
   if (fDirectory) fDirectory->Remove(this);
   fDirectory = dir;
   if (fDirectory) fDirectory->Append(this);
}


//______________________________________________________________________________
void TGraph2D::SetHistogram(TH2 *h)
{
   // Sets the histogram to be filled

   fUserHisto = kTRUE;
   fHistogram = (TH2D*)h;
   fNpx       = h->GetNbinsX();
   fNpy       = h->GetNbinsY();
}


//______________________________________________________________________________
void TGraph2D::SetMargin(Double_t m)
{
   // Sets the extra space (in %) around interpolated area for the 2D histogram

   if (m<0 || m>1) {
      Warning("SetMargin","The margin must be >= 0 && <= 1, fMargin set to 0.1");
      fMargin = 0.1;
   } else {
      fMargin = m;
   }
   if (fHistogram) {delete fHistogram; fHistogram = 0;}
}


//______________________________________________________________________________
void TGraph2D::SetMarginBinsContent(Double_t z)
{
   // Sets the histogram bin height for points lying outside the TGraphDelaunay
   // convex hull ie: the bins in the margin.

   fZout = z;
   if (fHistogram) {delete fHistogram; fHistogram = 0;}
}


//______________________________________________________________________________
void TGraph2D::SetMaximum(Double_t maximum)
{
   // Set maximum.

   fMaximum = maximum;
   GetHistogram()->SetMaximum(maximum);
}


//______________________________________________________________________________
void TGraph2D::SetMinimum(Double_t minimum)
{
   // Set minimum.

   fMinimum = minimum;
   GetHistogram()->SetMinimum(minimum);
}


//______________________________________________________________________________
void TGraph2D::SetName(const char *name)
{
   // Changes the name of this 2D graph

   //  2D graphs are named objects in a THashList.
   //  We must update the hashlist if we change the name
   if (fDirectory) fDirectory->Remove(this);
   fName = name;
   if (fDirectory) fDirectory->Append(this);
}


//______________________________________________________________________________
void TGraph2D::SetNameTitle(const char *name, const char *title)
{
   // Change the name and title of this 2D graph
   //

   //  2D graphs are named objects in a THashList.
   //  We must update the hashlist if we change the name
   if (fDirectory) fDirectory->Remove(this);
   fName  = name;
   SetTitle(title);
   if (fDirectory) fDirectory->Append(this);
}


//______________________________________________________________________________
void TGraph2D::SetNpx(Int_t npx)
{
   // Sets the number of bins along X used to draw the function

   if (npx < 4) {
      Warning("SetNpx","Number of points must be >4 && < 500, fNpx set to 4");
      fNpx = 4;
   } else if(npx > 500) {
      Warning("SetNpx","Number of points must be >4 && < 500, fNpx set to 500");
      fNpx = 500;
   } else {
      fNpx = npx;
   }
   if (fHistogram) {delete fHistogram; fHistogram = 0;}
}


//______________________________________________________________________________
void TGraph2D::SetNpy(Int_t npy)
{
   // Sets the number of bins along Y used to draw the function

   if (npy < 4) {
      Warning("SetNpy","Number of points must be >4 && < 500, fNpy set to 4");
      fNpy = 4;
   } else if(npy > 500) {
      Warning("SetNpy","Number of points must be >4 && < 500, fNpy set to 500");
      fNpy = 500;
   } else {
      fNpy = npy;
   }
   if (fHistogram) {delete fHistogram; fHistogram = 0;}
}


//______________________________________________________________________________
void TGraph2D::SetPoint(Int_t n, Double_t x, Double_t y, Double_t z)
{
   // Sets point number n.
   // If n is greater than the current size, the arrays are automatically
   // extended.

   if (n < 0) return;

   if (!fX || !fY || !fZ || n >= fSize) {
   // re-allocate the object
      Int_t newN = TMath::Max(2*fSize,n+1);
      Double_t *savex = new Double_t [newN];
      Double_t *savey = new Double_t [newN];
      Double_t *savez = new Double_t [newN];
      if (fX && fSize) {
         memcpy(savex,fX,fSize*sizeof(Double_t));
         memset(&savex[fSize],0,(newN-fSize)*sizeof(Double_t));
         delete [] fX;
      }
      if (fY && fSize) {
         memcpy(savey,fY,fSize*sizeof(Double_t));
         memset(&savey[fSize],0,(newN-fSize)*sizeof(Double_t));
         delete [] fY;
      }
      if (fZ && fSize) {
         memcpy(savez,fZ,fSize*sizeof(Double_t));
         memset(&savez[fSize],0,(newN-fSize)*sizeof(Double_t));
         delete [] fZ;
      }
      fX    = savex;
      fY    = savey;
      fZ    = savez;
      fSize = newN;
   }
   fX[n]    = x;
   fY[n]    = y;
   fZ[n]    = z;
   fNpoints = TMath::Max(fNpoints,n+1);
}


//______________________________________________________________________________
void TGraph2D::SetTitle(const char* title)
{
   // Sets graph title

   fTitle = title;
   if (fHistogram) fHistogram->SetTitle(title);
}


//_______________________________________________________________________
void TGraph2D::Streamer(TBuffer &b)
{
   // Stream a class object

   if (b.IsReading()) {
      UInt_t R__s, R__c;
      Version_t R__v = b.ReadVersion(&R__s, &R__c);
      b.ReadClassBuffer(TGraph2D::Class(), this, R__v, R__s, R__c);

      ResetBit(kCanDelete);
      ResetBit(kMustCleanup);
   } else {
      b.WriteClassBuffer(TGraph2D::Class(),this);
   }
}
 TGraph2D.cxx:1
 TGraph2D.cxx:2
 TGraph2D.cxx:3
 TGraph2D.cxx:4
 TGraph2D.cxx:5
 TGraph2D.cxx:6
 TGraph2D.cxx:7
 TGraph2D.cxx:8
 TGraph2D.cxx:9
 TGraph2D.cxx:10
 TGraph2D.cxx:11
 TGraph2D.cxx:12
 TGraph2D.cxx:13
 TGraph2D.cxx:14
 TGraph2D.cxx:15
 TGraph2D.cxx:16
 TGraph2D.cxx:17
 TGraph2D.cxx:18
 TGraph2D.cxx:19
 TGraph2D.cxx:20
 TGraph2D.cxx:21
 TGraph2D.cxx:22
 TGraph2D.cxx:23
 TGraph2D.cxx:24
 TGraph2D.cxx:25
 TGraph2D.cxx:26
 TGraph2D.cxx:27
 TGraph2D.cxx:28
 TGraph2D.cxx:29
 TGraph2D.cxx:30
 TGraph2D.cxx:31
 TGraph2D.cxx:32
 TGraph2D.cxx:33
 TGraph2D.cxx:34
 TGraph2D.cxx:35
 TGraph2D.cxx:36
 TGraph2D.cxx:37
 TGraph2D.cxx:38
 TGraph2D.cxx:39
 TGraph2D.cxx:40
 TGraph2D.cxx:41
 TGraph2D.cxx:42
 TGraph2D.cxx:43
 TGraph2D.cxx:44
 TGraph2D.cxx:45
 TGraph2D.cxx:46
 TGraph2D.cxx:47
 TGraph2D.cxx:48
 TGraph2D.cxx:49
 TGraph2D.cxx:50
 TGraph2D.cxx:51
 TGraph2D.cxx:52
 TGraph2D.cxx:53
 TGraph2D.cxx:54
 TGraph2D.cxx:55
 TGraph2D.cxx:56
 TGraph2D.cxx:57
 TGraph2D.cxx:58
 TGraph2D.cxx:59
 TGraph2D.cxx:60
 TGraph2D.cxx:61
 TGraph2D.cxx:62
 TGraph2D.cxx:63
 TGraph2D.cxx:64
 TGraph2D.cxx:65
 TGraph2D.cxx:66
 TGraph2D.cxx:67
 TGraph2D.cxx:68
 TGraph2D.cxx:69
 TGraph2D.cxx:70
 TGraph2D.cxx:71
 TGraph2D.cxx:72
 TGraph2D.cxx:73
 TGraph2D.cxx:74
 TGraph2D.cxx:75
 TGraph2D.cxx:76
 TGraph2D.cxx:77
 TGraph2D.cxx:78
 TGraph2D.cxx:79
 TGraph2D.cxx:80
 TGraph2D.cxx:81
 TGraph2D.cxx:82
 TGraph2D.cxx:83
 TGraph2D.cxx:84
 TGraph2D.cxx:85
 TGraph2D.cxx:86
 TGraph2D.cxx:87
 TGraph2D.cxx:88
 TGraph2D.cxx:89
 TGraph2D.cxx:90
 TGraph2D.cxx:91
 TGraph2D.cxx:92
 TGraph2D.cxx:93
 TGraph2D.cxx:94
 TGraph2D.cxx:95
 TGraph2D.cxx:96
 TGraph2D.cxx:97
 TGraph2D.cxx:98
 TGraph2D.cxx:99
 TGraph2D.cxx:100
 TGraph2D.cxx:101
 TGraph2D.cxx:102
 TGraph2D.cxx:103
 TGraph2D.cxx:104
 TGraph2D.cxx:105
 TGraph2D.cxx:106
 TGraph2D.cxx:107
 TGraph2D.cxx:108
 TGraph2D.cxx:109
 TGraph2D.cxx:110
 TGraph2D.cxx:111
 TGraph2D.cxx:112
 TGraph2D.cxx:113
 TGraph2D.cxx:114
 TGraph2D.cxx:115
 TGraph2D.cxx:116
 TGraph2D.cxx:117
 TGraph2D.cxx:118
 TGraph2D.cxx:119
 TGraph2D.cxx:120
 TGraph2D.cxx:121
 TGraph2D.cxx:122
 TGraph2D.cxx:123
 TGraph2D.cxx:124
 TGraph2D.cxx:125
 TGraph2D.cxx:126
 TGraph2D.cxx:127
 TGraph2D.cxx:128
 TGraph2D.cxx:129
 TGraph2D.cxx:130
 TGraph2D.cxx:131
 TGraph2D.cxx:132
 TGraph2D.cxx:133
 TGraph2D.cxx:134
 TGraph2D.cxx:135
 TGraph2D.cxx:136
 TGraph2D.cxx:137
 TGraph2D.cxx:138
 TGraph2D.cxx:139
 TGraph2D.cxx:140
 TGraph2D.cxx:141
 TGraph2D.cxx:142
 TGraph2D.cxx:143
 TGraph2D.cxx:144
 TGraph2D.cxx:145
 TGraph2D.cxx:146
 TGraph2D.cxx:147
 TGraph2D.cxx:148
 TGraph2D.cxx:149
 TGraph2D.cxx:150
 TGraph2D.cxx:151
 TGraph2D.cxx:152
 TGraph2D.cxx:153
 TGraph2D.cxx:154
 TGraph2D.cxx:155
 TGraph2D.cxx:156
 TGraph2D.cxx:157
 TGraph2D.cxx:158
 TGraph2D.cxx:159
 TGraph2D.cxx:160
 TGraph2D.cxx:161
 TGraph2D.cxx:162
 TGraph2D.cxx:163
 TGraph2D.cxx:164
 TGraph2D.cxx:165
 TGraph2D.cxx:166
 TGraph2D.cxx:167
 TGraph2D.cxx:168
 TGraph2D.cxx:169
 TGraph2D.cxx:170
 TGraph2D.cxx:171
 TGraph2D.cxx:172
 TGraph2D.cxx:173
 TGraph2D.cxx:174
 TGraph2D.cxx:175
 TGraph2D.cxx:176
 TGraph2D.cxx:177
 TGraph2D.cxx:178
 TGraph2D.cxx:179
 TGraph2D.cxx:180
 TGraph2D.cxx:181
 TGraph2D.cxx:182
 TGraph2D.cxx:183
 TGraph2D.cxx:184
 TGraph2D.cxx:185
 TGraph2D.cxx:186
 TGraph2D.cxx:187
 TGraph2D.cxx:188
 TGraph2D.cxx:189
 TGraph2D.cxx:190
 TGraph2D.cxx:191
 TGraph2D.cxx:192
 TGraph2D.cxx:193
 TGraph2D.cxx:194
 TGraph2D.cxx:195
 TGraph2D.cxx:196
 TGraph2D.cxx:197
 TGraph2D.cxx:198
 TGraph2D.cxx:199
 TGraph2D.cxx:200
 TGraph2D.cxx:201
 TGraph2D.cxx:202
 TGraph2D.cxx:203
 TGraph2D.cxx:204
 TGraph2D.cxx:205
 TGraph2D.cxx:206
 TGraph2D.cxx:207
 TGraph2D.cxx:208
 TGraph2D.cxx:209
 TGraph2D.cxx:210
 TGraph2D.cxx:211
 TGraph2D.cxx:212
 TGraph2D.cxx:213
 TGraph2D.cxx:214
 TGraph2D.cxx:215
 TGraph2D.cxx:216
 TGraph2D.cxx:217
 TGraph2D.cxx:218
 TGraph2D.cxx:219
 TGraph2D.cxx:220
 TGraph2D.cxx:221
 TGraph2D.cxx:222
 TGraph2D.cxx:223
 TGraph2D.cxx:224
 TGraph2D.cxx:225
 TGraph2D.cxx:226
 TGraph2D.cxx:227
 TGraph2D.cxx:228
 TGraph2D.cxx:229
 TGraph2D.cxx:230
 TGraph2D.cxx:231
 TGraph2D.cxx:232
 TGraph2D.cxx:233
 TGraph2D.cxx:234
 TGraph2D.cxx:235
 TGraph2D.cxx:236
 TGraph2D.cxx:237
 TGraph2D.cxx:238
 TGraph2D.cxx:239
 TGraph2D.cxx:240
 TGraph2D.cxx:241
 TGraph2D.cxx:242
 TGraph2D.cxx:243
 TGraph2D.cxx:244
 TGraph2D.cxx:245
 TGraph2D.cxx:246
 TGraph2D.cxx:247
 TGraph2D.cxx:248
 TGraph2D.cxx:249
 TGraph2D.cxx:250
 TGraph2D.cxx:251
 TGraph2D.cxx:252
 TGraph2D.cxx:253
 TGraph2D.cxx:254
 TGraph2D.cxx:255
 TGraph2D.cxx:256
 TGraph2D.cxx:257
 TGraph2D.cxx:258
 TGraph2D.cxx:259
 TGraph2D.cxx:260
 TGraph2D.cxx:261
 TGraph2D.cxx:262
 TGraph2D.cxx:263
 TGraph2D.cxx:264
 TGraph2D.cxx:265
 TGraph2D.cxx:266
 TGraph2D.cxx:267
 TGraph2D.cxx:268
 TGraph2D.cxx:269
 TGraph2D.cxx:270
 TGraph2D.cxx:271
 TGraph2D.cxx:272
 TGraph2D.cxx:273
 TGraph2D.cxx:274
 TGraph2D.cxx:275
 TGraph2D.cxx:276
 TGraph2D.cxx:277
 TGraph2D.cxx:278
 TGraph2D.cxx:279
 TGraph2D.cxx:280
 TGraph2D.cxx:281
 TGraph2D.cxx:282
 TGraph2D.cxx:283
 TGraph2D.cxx:284
 TGraph2D.cxx:285
 TGraph2D.cxx:286
 TGraph2D.cxx:287
 TGraph2D.cxx:288
 TGraph2D.cxx:289
 TGraph2D.cxx:290
 TGraph2D.cxx:291
 TGraph2D.cxx:292
 TGraph2D.cxx:293
 TGraph2D.cxx:294
 TGraph2D.cxx:295
 TGraph2D.cxx:296
 TGraph2D.cxx:297
 TGraph2D.cxx:298
 TGraph2D.cxx:299
 TGraph2D.cxx:300
 TGraph2D.cxx:301
 TGraph2D.cxx:302
 TGraph2D.cxx:303
 TGraph2D.cxx:304
 TGraph2D.cxx:305
 TGraph2D.cxx:306
 TGraph2D.cxx:307
 TGraph2D.cxx:308
 TGraph2D.cxx:309
 TGraph2D.cxx:310
 TGraph2D.cxx:311
 TGraph2D.cxx:312
 TGraph2D.cxx:313
 TGraph2D.cxx:314
 TGraph2D.cxx:315
 TGraph2D.cxx:316
 TGraph2D.cxx:317
 TGraph2D.cxx:318
 TGraph2D.cxx:319
 TGraph2D.cxx:320
 TGraph2D.cxx:321
 TGraph2D.cxx:322
 TGraph2D.cxx:323
 TGraph2D.cxx:324
 TGraph2D.cxx:325
 TGraph2D.cxx:326
 TGraph2D.cxx:327
 TGraph2D.cxx:328
 TGraph2D.cxx:329
 TGraph2D.cxx:330
 TGraph2D.cxx:331
 TGraph2D.cxx:332
 TGraph2D.cxx:333
 TGraph2D.cxx:334
 TGraph2D.cxx:335
 TGraph2D.cxx:336
 TGraph2D.cxx:337
 TGraph2D.cxx:338
 TGraph2D.cxx:339
 TGraph2D.cxx:340
 TGraph2D.cxx:341
 TGraph2D.cxx:342
 TGraph2D.cxx:343
 TGraph2D.cxx:344
 TGraph2D.cxx:345
 TGraph2D.cxx:346
 TGraph2D.cxx:347
 TGraph2D.cxx:348
 TGraph2D.cxx:349
 TGraph2D.cxx:350
 TGraph2D.cxx:351
 TGraph2D.cxx:352
 TGraph2D.cxx:353
 TGraph2D.cxx:354
 TGraph2D.cxx:355
 TGraph2D.cxx:356
 TGraph2D.cxx:357
 TGraph2D.cxx:358
 TGraph2D.cxx:359
 TGraph2D.cxx:360
 TGraph2D.cxx:361
 TGraph2D.cxx:362
 TGraph2D.cxx:363
 TGraph2D.cxx:364
 TGraph2D.cxx:365
 TGraph2D.cxx:366
 TGraph2D.cxx:367
 TGraph2D.cxx:368
 TGraph2D.cxx:369
 TGraph2D.cxx:370
 TGraph2D.cxx:371
 TGraph2D.cxx:372
 TGraph2D.cxx:373
 TGraph2D.cxx:374
 TGraph2D.cxx:375
 TGraph2D.cxx:376
 TGraph2D.cxx:377
 TGraph2D.cxx:378
 TGraph2D.cxx:379
 TGraph2D.cxx:380
 TGraph2D.cxx:381
 TGraph2D.cxx:382
 TGraph2D.cxx:383
 TGraph2D.cxx:384
 TGraph2D.cxx:385
 TGraph2D.cxx:386
 TGraph2D.cxx:387
 TGraph2D.cxx:388
 TGraph2D.cxx:389
 TGraph2D.cxx:390
 TGraph2D.cxx:391
 TGraph2D.cxx:392
 TGraph2D.cxx:393
 TGraph2D.cxx:394
 TGraph2D.cxx:395
 TGraph2D.cxx:396
 TGraph2D.cxx:397
 TGraph2D.cxx:398
 TGraph2D.cxx:399
 TGraph2D.cxx:400
 TGraph2D.cxx:401
 TGraph2D.cxx:402
 TGraph2D.cxx:403
 TGraph2D.cxx:404
 TGraph2D.cxx:405
 TGraph2D.cxx:406
 TGraph2D.cxx:407
 TGraph2D.cxx:408
 TGraph2D.cxx:409
 TGraph2D.cxx:410
 TGraph2D.cxx:411
 TGraph2D.cxx:412
 TGraph2D.cxx:413
 TGraph2D.cxx:414
 TGraph2D.cxx:415
 TGraph2D.cxx:416
 TGraph2D.cxx:417
 TGraph2D.cxx:418
 TGraph2D.cxx:419
 TGraph2D.cxx:420
 TGraph2D.cxx:421
 TGraph2D.cxx:422
 TGraph2D.cxx:423
 TGraph2D.cxx:424
 TGraph2D.cxx:425
 TGraph2D.cxx:426
 TGraph2D.cxx:427
 TGraph2D.cxx:428
 TGraph2D.cxx:429
 TGraph2D.cxx:430
 TGraph2D.cxx:431
 TGraph2D.cxx:432
 TGraph2D.cxx:433
 TGraph2D.cxx:434
 TGraph2D.cxx:435
 TGraph2D.cxx:436
 TGraph2D.cxx:437
 TGraph2D.cxx:438
 TGraph2D.cxx:439
 TGraph2D.cxx:440
 TGraph2D.cxx:441
 TGraph2D.cxx:442
 TGraph2D.cxx:443
 TGraph2D.cxx:444
 TGraph2D.cxx:445
 TGraph2D.cxx:446
 TGraph2D.cxx:447
 TGraph2D.cxx:448
 TGraph2D.cxx:449
 TGraph2D.cxx:450
 TGraph2D.cxx:451
 TGraph2D.cxx:452
 TGraph2D.cxx:453
 TGraph2D.cxx:454
 TGraph2D.cxx:455
 TGraph2D.cxx:456
 TGraph2D.cxx:457
 TGraph2D.cxx:458
 TGraph2D.cxx:459
 TGraph2D.cxx:460
 TGraph2D.cxx:461
 TGraph2D.cxx:462
 TGraph2D.cxx:463
 TGraph2D.cxx:464
 TGraph2D.cxx:465
 TGraph2D.cxx:466
 TGraph2D.cxx:467
 TGraph2D.cxx:468
 TGraph2D.cxx:469
 TGraph2D.cxx:470
 TGraph2D.cxx:471
 TGraph2D.cxx:472
 TGraph2D.cxx:473
 TGraph2D.cxx:474
 TGraph2D.cxx:475
 TGraph2D.cxx:476
 TGraph2D.cxx:477
 TGraph2D.cxx:478
 TGraph2D.cxx:479
 TGraph2D.cxx:480
 TGraph2D.cxx:481
 TGraph2D.cxx:482
 TGraph2D.cxx:483
 TGraph2D.cxx:484
 TGraph2D.cxx:485
 TGraph2D.cxx:486
 TGraph2D.cxx:487
 TGraph2D.cxx:488
 TGraph2D.cxx:489
 TGraph2D.cxx:490
 TGraph2D.cxx:491
 TGraph2D.cxx:492
 TGraph2D.cxx:493
 TGraph2D.cxx:494
 TGraph2D.cxx:495
 TGraph2D.cxx:496
 TGraph2D.cxx:497
 TGraph2D.cxx:498
 TGraph2D.cxx:499
 TGraph2D.cxx:500
 TGraph2D.cxx:501
 TGraph2D.cxx:502
 TGraph2D.cxx:503
 TGraph2D.cxx:504
 TGraph2D.cxx:505
 TGraph2D.cxx:506
 TGraph2D.cxx:507
 TGraph2D.cxx:508
 TGraph2D.cxx:509
 TGraph2D.cxx:510
 TGraph2D.cxx:511
 TGraph2D.cxx:512
 TGraph2D.cxx:513
 TGraph2D.cxx:514
 TGraph2D.cxx:515
 TGraph2D.cxx:516
 TGraph2D.cxx:517
 TGraph2D.cxx:518
 TGraph2D.cxx:519
 TGraph2D.cxx:520
 TGraph2D.cxx:521
 TGraph2D.cxx:522
 TGraph2D.cxx:523
 TGraph2D.cxx:524
 TGraph2D.cxx:525
 TGraph2D.cxx:526
 TGraph2D.cxx:527
 TGraph2D.cxx:528
 TGraph2D.cxx:529
 TGraph2D.cxx:530
 TGraph2D.cxx:531
 TGraph2D.cxx:532
 TGraph2D.cxx:533
 TGraph2D.cxx:534
 TGraph2D.cxx:535
 TGraph2D.cxx:536
 TGraph2D.cxx:537
 TGraph2D.cxx:538
 TGraph2D.cxx:539
 TGraph2D.cxx:540
 TGraph2D.cxx:541
 TGraph2D.cxx:542
 TGraph2D.cxx:543
 TGraph2D.cxx:544
 TGraph2D.cxx:545
 TGraph2D.cxx:546
 TGraph2D.cxx:547
 TGraph2D.cxx:548
 TGraph2D.cxx:549
 TGraph2D.cxx:550
 TGraph2D.cxx:551
 TGraph2D.cxx:552
 TGraph2D.cxx:553
 TGraph2D.cxx:554
 TGraph2D.cxx:555
 TGraph2D.cxx:556
 TGraph2D.cxx:557
 TGraph2D.cxx:558
 TGraph2D.cxx:559
 TGraph2D.cxx:560
 TGraph2D.cxx:561
 TGraph2D.cxx:562
 TGraph2D.cxx:563
 TGraph2D.cxx:564
 TGraph2D.cxx:565
 TGraph2D.cxx:566
 TGraph2D.cxx:567
 TGraph2D.cxx:568
 TGraph2D.cxx:569
 TGraph2D.cxx:570
 TGraph2D.cxx:571
 TGraph2D.cxx:572
 TGraph2D.cxx:573
 TGraph2D.cxx:574
 TGraph2D.cxx:575
 TGraph2D.cxx:576
 TGraph2D.cxx:577
 TGraph2D.cxx:578
 TGraph2D.cxx:579
 TGraph2D.cxx:580
 TGraph2D.cxx:581
 TGraph2D.cxx:582
 TGraph2D.cxx:583
 TGraph2D.cxx:584
 TGraph2D.cxx:585
 TGraph2D.cxx:586
 TGraph2D.cxx:587
 TGraph2D.cxx:588
 TGraph2D.cxx:589
 TGraph2D.cxx:590
 TGraph2D.cxx:591
 TGraph2D.cxx:592
 TGraph2D.cxx:593
 TGraph2D.cxx:594
 TGraph2D.cxx:595
 TGraph2D.cxx:596
 TGraph2D.cxx:597
 TGraph2D.cxx:598
 TGraph2D.cxx:599
 TGraph2D.cxx:600
 TGraph2D.cxx:601
 TGraph2D.cxx:602
 TGraph2D.cxx:603
 TGraph2D.cxx:604
 TGraph2D.cxx:605
 TGraph2D.cxx:606
 TGraph2D.cxx:607
 TGraph2D.cxx:608
 TGraph2D.cxx:609
 TGraph2D.cxx:610
 TGraph2D.cxx:611
 TGraph2D.cxx:612
 TGraph2D.cxx:613
 TGraph2D.cxx:614
 TGraph2D.cxx:615
 TGraph2D.cxx:616
 TGraph2D.cxx:617
 TGraph2D.cxx:618
 TGraph2D.cxx:619
 TGraph2D.cxx:620
 TGraph2D.cxx:621
 TGraph2D.cxx:622
 TGraph2D.cxx:623
 TGraph2D.cxx:624
 TGraph2D.cxx:625
 TGraph2D.cxx:626
 TGraph2D.cxx:627
 TGraph2D.cxx:628
 TGraph2D.cxx:629
 TGraph2D.cxx:630
 TGraph2D.cxx:631
 TGraph2D.cxx:632
 TGraph2D.cxx:633
 TGraph2D.cxx:634
 TGraph2D.cxx:635
 TGraph2D.cxx:636
 TGraph2D.cxx:637
 TGraph2D.cxx:638
 TGraph2D.cxx:639
 TGraph2D.cxx:640
 TGraph2D.cxx:641
 TGraph2D.cxx:642
 TGraph2D.cxx:643
 TGraph2D.cxx:644
 TGraph2D.cxx:645
 TGraph2D.cxx:646
 TGraph2D.cxx:647
 TGraph2D.cxx:648
 TGraph2D.cxx:649
 TGraph2D.cxx:650
 TGraph2D.cxx:651
 TGraph2D.cxx:652
 TGraph2D.cxx:653
 TGraph2D.cxx:654
 TGraph2D.cxx:655
 TGraph2D.cxx:656
 TGraph2D.cxx:657
 TGraph2D.cxx:658
 TGraph2D.cxx:659
 TGraph2D.cxx:660
 TGraph2D.cxx:661
 TGraph2D.cxx:662
 TGraph2D.cxx:663
 TGraph2D.cxx:664
 TGraph2D.cxx:665
 TGraph2D.cxx:666
 TGraph2D.cxx:667
 TGraph2D.cxx:668
 TGraph2D.cxx:669
 TGraph2D.cxx:670
 TGraph2D.cxx:671
 TGraph2D.cxx:672
 TGraph2D.cxx:673
 TGraph2D.cxx:674
 TGraph2D.cxx:675
 TGraph2D.cxx:676
 TGraph2D.cxx:677
 TGraph2D.cxx:678
 TGraph2D.cxx:679
 TGraph2D.cxx:680
 TGraph2D.cxx:681
 TGraph2D.cxx:682
 TGraph2D.cxx:683
 TGraph2D.cxx:684
 TGraph2D.cxx:685
 TGraph2D.cxx:686
 TGraph2D.cxx:687
 TGraph2D.cxx:688
 TGraph2D.cxx:689
 TGraph2D.cxx:690
 TGraph2D.cxx:691
 TGraph2D.cxx:692
 TGraph2D.cxx:693
 TGraph2D.cxx:694
 TGraph2D.cxx:695
 TGraph2D.cxx:696
 TGraph2D.cxx:697
 TGraph2D.cxx:698
 TGraph2D.cxx:699
 TGraph2D.cxx:700
 TGraph2D.cxx:701
 TGraph2D.cxx:702
 TGraph2D.cxx:703
 TGraph2D.cxx:704
 TGraph2D.cxx:705
 TGraph2D.cxx:706
 TGraph2D.cxx:707
 TGraph2D.cxx:708
 TGraph2D.cxx:709
 TGraph2D.cxx:710
 TGraph2D.cxx:711
 TGraph2D.cxx:712
 TGraph2D.cxx:713
 TGraph2D.cxx:714
 TGraph2D.cxx:715
 TGraph2D.cxx:716
 TGraph2D.cxx:717
 TGraph2D.cxx:718
 TGraph2D.cxx:719
 TGraph2D.cxx:720
 TGraph2D.cxx:721
 TGraph2D.cxx:722
 TGraph2D.cxx:723
 TGraph2D.cxx:724
 TGraph2D.cxx:725
 TGraph2D.cxx:726
 TGraph2D.cxx:727
 TGraph2D.cxx:728
 TGraph2D.cxx:729
 TGraph2D.cxx:730
 TGraph2D.cxx:731
 TGraph2D.cxx:732
 TGraph2D.cxx:733
 TGraph2D.cxx:734
 TGraph2D.cxx:735
 TGraph2D.cxx:736
 TGraph2D.cxx:737
 TGraph2D.cxx:738
 TGraph2D.cxx:739
 TGraph2D.cxx:740
 TGraph2D.cxx:741
 TGraph2D.cxx:742
 TGraph2D.cxx:743
 TGraph2D.cxx:744
 TGraph2D.cxx:745
 TGraph2D.cxx:746
 TGraph2D.cxx:747
 TGraph2D.cxx:748
 TGraph2D.cxx:749
 TGraph2D.cxx:750
 TGraph2D.cxx:751
 TGraph2D.cxx:752
 TGraph2D.cxx:753
 TGraph2D.cxx:754
 TGraph2D.cxx:755
 TGraph2D.cxx:756
 TGraph2D.cxx:757
 TGraph2D.cxx:758
 TGraph2D.cxx:759
 TGraph2D.cxx:760
 TGraph2D.cxx:761
 TGraph2D.cxx:762
 TGraph2D.cxx:763
 TGraph2D.cxx:764
 TGraph2D.cxx:765
 TGraph2D.cxx:766
 TGraph2D.cxx:767
 TGraph2D.cxx:768
 TGraph2D.cxx:769
 TGraph2D.cxx:770
 TGraph2D.cxx:771
 TGraph2D.cxx:772
 TGraph2D.cxx:773
 TGraph2D.cxx:774
 TGraph2D.cxx:775
 TGraph2D.cxx:776
 TGraph2D.cxx:777
 TGraph2D.cxx:778
 TGraph2D.cxx:779
 TGraph2D.cxx:780
 TGraph2D.cxx:781
 TGraph2D.cxx:782
 TGraph2D.cxx:783
 TGraph2D.cxx:784
 TGraph2D.cxx:785
 TGraph2D.cxx:786
 TGraph2D.cxx:787
 TGraph2D.cxx:788
 TGraph2D.cxx:789
 TGraph2D.cxx:790
 TGraph2D.cxx:791
 TGraph2D.cxx:792
 TGraph2D.cxx:793
 TGraph2D.cxx:794
 TGraph2D.cxx:795
 TGraph2D.cxx:796
 TGraph2D.cxx:797
 TGraph2D.cxx:798
 TGraph2D.cxx:799
 TGraph2D.cxx:800
 TGraph2D.cxx:801
 TGraph2D.cxx:802
 TGraph2D.cxx:803
 TGraph2D.cxx:804
 TGraph2D.cxx:805
 TGraph2D.cxx:806
 TGraph2D.cxx:807
 TGraph2D.cxx:808
 TGraph2D.cxx:809
 TGraph2D.cxx:810
 TGraph2D.cxx:811
 TGraph2D.cxx:812
 TGraph2D.cxx:813
 TGraph2D.cxx:814
 TGraph2D.cxx:815
 TGraph2D.cxx:816
 TGraph2D.cxx:817
 TGraph2D.cxx:818
 TGraph2D.cxx:819
 TGraph2D.cxx:820
 TGraph2D.cxx:821
 TGraph2D.cxx:822
 TGraph2D.cxx:823
 TGraph2D.cxx:824
 TGraph2D.cxx:825
 TGraph2D.cxx:826
 TGraph2D.cxx:827
 TGraph2D.cxx:828
 TGraph2D.cxx:829
 TGraph2D.cxx:830
 TGraph2D.cxx:831
 TGraph2D.cxx:832
 TGraph2D.cxx:833
 TGraph2D.cxx:834
 TGraph2D.cxx:835
 TGraph2D.cxx:836
 TGraph2D.cxx:837
 TGraph2D.cxx:838
 TGraph2D.cxx:839
 TGraph2D.cxx:840
 TGraph2D.cxx:841
 TGraph2D.cxx:842
 TGraph2D.cxx:843
 TGraph2D.cxx:844
 TGraph2D.cxx:845
 TGraph2D.cxx:846
 TGraph2D.cxx:847
 TGraph2D.cxx:848
 TGraph2D.cxx:849
 TGraph2D.cxx:850
 TGraph2D.cxx:851
 TGraph2D.cxx:852
 TGraph2D.cxx:853
 TGraph2D.cxx:854
 TGraph2D.cxx:855
 TGraph2D.cxx:856
 TGraph2D.cxx:857
 TGraph2D.cxx:858
 TGraph2D.cxx:859
 TGraph2D.cxx:860
 TGraph2D.cxx:861
 TGraph2D.cxx:862
 TGraph2D.cxx:863
 TGraph2D.cxx:864
 TGraph2D.cxx:865
 TGraph2D.cxx:866
 TGraph2D.cxx:867
 TGraph2D.cxx:868
 TGraph2D.cxx:869
 TGraph2D.cxx:870
 TGraph2D.cxx:871
 TGraph2D.cxx:872
 TGraph2D.cxx:873
 TGraph2D.cxx:874
 TGraph2D.cxx:875
 TGraph2D.cxx:876
 TGraph2D.cxx:877
 TGraph2D.cxx:878
 TGraph2D.cxx:879
 TGraph2D.cxx:880
 TGraph2D.cxx:881
 TGraph2D.cxx:882
 TGraph2D.cxx:883
 TGraph2D.cxx:884
 TGraph2D.cxx:885
 TGraph2D.cxx:886
 TGraph2D.cxx:887
 TGraph2D.cxx:888
 TGraph2D.cxx:889
 TGraph2D.cxx:890
 TGraph2D.cxx:891
 TGraph2D.cxx:892
 TGraph2D.cxx:893
 TGraph2D.cxx:894
 TGraph2D.cxx:895
 TGraph2D.cxx:896
 TGraph2D.cxx:897
 TGraph2D.cxx:898
 TGraph2D.cxx:899
 TGraph2D.cxx:900
 TGraph2D.cxx:901
 TGraph2D.cxx:902
 TGraph2D.cxx:903
 TGraph2D.cxx:904
 TGraph2D.cxx:905
 TGraph2D.cxx:906
 TGraph2D.cxx:907
 TGraph2D.cxx:908
 TGraph2D.cxx:909
 TGraph2D.cxx:910
 TGraph2D.cxx:911
 TGraph2D.cxx:912
 TGraph2D.cxx:913
 TGraph2D.cxx:914
 TGraph2D.cxx:915
 TGraph2D.cxx:916
 TGraph2D.cxx:917
 TGraph2D.cxx:918
 TGraph2D.cxx:919
 TGraph2D.cxx:920
 TGraph2D.cxx:921
 TGraph2D.cxx:922
 TGraph2D.cxx:923
 TGraph2D.cxx:924
 TGraph2D.cxx:925
 TGraph2D.cxx:926
 TGraph2D.cxx:927
 TGraph2D.cxx:928
 TGraph2D.cxx:929
 TGraph2D.cxx:930
 TGraph2D.cxx:931
 TGraph2D.cxx:932
 TGraph2D.cxx:933
 TGraph2D.cxx:934
 TGraph2D.cxx:935
 TGraph2D.cxx:936
 TGraph2D.cxx:937
 TGraph2D.cxx:938
 TGraph2D.cxx:939
 TGraph2D.cxx:940
 TGraph2D.cxx:941
 TGraph2D.cxx:942
 TGraph2D.cxx:943
 TGraph2D.cxx:944
 TGraph2D.cxx:945
 TGraph2D.cxx:946
 TGraph2D.cxx:947
 TGraph2D.cxx:948
 TGraph2D.cxx:949
 TGraph2D.cxx:950
 TGraph2D.cxx:951
 TGraph2D.cxx:952
 TGraph2D.cxx:953
 TGraph2D.cxx:954
 TGraph2D.cxx:955
 TGraph2D.cxx:956
 TGraph2D.cxx:957
 TGraph2D.cxx:958
 TGraph2D.cxx:959
 TGraph2D.cxx:960
 TGraph2D.cxx:961
 TGraph2D.cxx:962
 TGraph2D.cxx:963
 TGraph2D.cxx:964
 TGraph2D.cxx:965
 TGraph2D.cxx:966
 TGraph2D.cxx:967
 TGraph2D.cxx:968
 TGraph2D.cxx:969
 TGraph2D.cxx:970
 TGraph2D.cxx:971
 TGraph2D.cxx:972
 TGraph2D.cxx:973
 TGraph2D.cxx:974
 TGraph2D.cxx:975
 TGraph2D.cxx:976
 TGraph2D.cxx:977
 TGraph2D.cxx:978
 TGraph2D.cxx:979
 TGraph2D.cxx:980
 TGraph2D.cxx:981
 TGraph2D.cxx:982
 TGraph2D.cxx:983
 TGraph2D.cxx:984
 TGraph2D.cxx:985
 TGraph2D.cxx:986
 TGraph2D.cxx:987
 TGraph2D.cxx:988
 TGraph2D.cxx:989
 TGraph2D.cxx:990
 TGraph2D.cxx:991
 TGraph2D.cxx:992
 TGraph2D.cxx:993
 TGraph2D.cxx:994
 TGraph2D.cxx:995
 TGraph2D.cxx:996
 TGraph2D.cxx:997
 TGraph2D.cxx:998
 TGraph2D.cxx:999
 TGraph2D.cxx:1000
 TGraph2D.cxx:1001
 TGraph2D.cxx:1002
 TGraph2D.cxx:1003
 TGraph2D.cxx:1004
 TGraph2D.cxx:1005
 TGraph2D.cxx:1006
 TGraph2D.cxx:1007
 TGraph2D.cxx:1008
 TGraph2D.cxx:1009
 TGraph2D.cxx:1010
 TGraph2D.cxx:1011
 TGraph2D.cxx:1012
 TGraph2D.cxx:1013
 TGraph2D.cxx:1014
 TGraph2D.cxx:1015
 TGraph2D.cxx:1016
 TGraph2D.cxx:1017
 TGraph2D.cxx:1018
 TGraph2D.cxx:1019
 TGraph2D.cxx:1020
 TGraph2D.cxx:1021
 TGraph2D.cxx:1022
 TGraph2D.cxx:1023
 TGraph2D.cxx:1024
 TGraph2D.cxx:1025
 TGraph2D.cxx:1026
 TGraph2D.cxx:1027
 TGraph2D.cxx:1028
 TGraph2D.cxx:1029
 TGraph2D.cxx:1030
 TGraph2D.cxx:1031
 TGraph2D.cxx:1032
 TGraph2D.cxx:1033
 TGraph2D.cxx:1034
 TGraph2D.cxx:1035
 TGraph2D.cxx:1036
 TGraph2D.cxx:1037
 TGraph2D.cxx:1038
 TGraph2D.cxx:1039
 TGraph2D.cxx:1040
 TGraph2D.cxx:1041
 TGraph2D.cxx:1042
 TGraph2D.cxx:1043
 TGraph2D.cxx:1044
 TGraph2D.cxx:1045
 TGraph2D.cxx:1046
 TGraph2D.cxx:1047
 TGraph2D.cxx:1048
 TGraph2D.cxx:1049
 TGraph2D.cxx:1050
 TGraph2D.cxx:1051
 TGraph2D.cxx:1052
 TGraph2D.cxx:1053
 TGraph2D.cxx:1054
 TGraph2D.cxx:1055
 TGraph2D.cxx:1056
 TGraph2D.cxx:1057
 TGraph2D.cxx:1058
 TGraph2D.cxx:1059
 TGraph2D.cxx:1060
 TGraph2D.cxx:1061
 TGraph2D.cxx:1062
 TGraph2D.cxx:1063
 TGraph2D.cxx:1064
 TGraph2D.cxx:1065
 TGraph2D.cxx:1066
 TGraph2D.cxx:1067
 TGraph2D.cxx:1068
 TGraph2D.cxx:1069
 TGraph2D.cxx:1070
 TGraph2D.cxx:1071
 TGraph2D.cxx:1072
 TGraph2D.cxx:1073
 TGraph2D.cxx:1074
 TGraph2D.cxx:1075
 TGraph2D.cxx:1076
 TGraph2D.cxx:1077
 TGraph2D.cxx:1078
 TGraph2D.cxx:1079
 TGraph2D.cxx:1080
 TGraph2D.cxx:1081
 TGraph2D.cxx:1082
 TGraph2D.cxx:1083
 TGraph2D.cxx:1084
 TGraph2D.cxx:1085
 TGraph2D.cxx:1086
 TGraph2D.cxx:1087
 TGraph2D.cxx:1088
 TGraph2D.cxx:1089
 TGraph2D.cxx:1090
 TGraph2D.cxx:1091
 TGraph2D.cxx:1092
 TGraph2D.cxx:1093
 TGraph2D.cxx:1094
 TGraph2D.cxx:1095
 TGraph2D.cxx:1096
 TGraph2D.cxx:1097
 TGraph2D.cxx:1098
 TGraph2D.cxx:1099
 TGraph2D.cxx:1100
 TGraph2D.cxx:1101
 TGraph2D.cxx:1102
 TGraph2D.cxx:1103
 TGraph2D.cxx:1104
 TGraph2D.cxx:1105
 TGraph2D.cxx:1106
 TGraph2D.cxx:1107
 TGraph2D.cxx:1108
 TGraph2D.cxx:1109
 TGraph2D.cxx:1110
 TGraph2D.cxx:1111
 TGraph2D.cxx:1112
 TGraph2D.cxx:1113
 TGraph2D.cxx:1114
 TGraph2D.cxx:1115
 TGraph2D.cxx:1116
 TGraph2D.cxx:1117
 TGraph2D.cxx:1118
 TGraph2D.cxx:1119
 TGraph2D.cxx:1120
 TGraph2D.cxx:1121
 TGraph2D.cxx:1122
 TGraph2D.cxx:1123
 TGraph2D.cxx:1124
 TGraph2D.cxx:1125
 TGraph2D.cxx:1126
 TGraph2D.cxx:1127
 TGraph2D.cxx:1128
 TGraph2D.cxx:1129
 TGraph2D.cxx:1130
 TGraph2D.cxx:1131
 TGraph2D.cxx:1132
 TGraph2D.cxx:1133
 TGraph2D.cxx:1134
 TGraph2D.cxx:1135
 TGraph2D.cxx:1136
 TGraph2D.cxx:1137
 TGraph2D.cxx:1138
 TGraph2D.cxx:1139
 TGraph2D.cxx:1140
 TGraph2D.cxx:1141
 TGraph2D.cxx:1142
 TGraph2D.cxx:1143
 TGraph2D.cxx:1144
 TGraph2D.cxx:1145
 TGraph2D.cxx:1146
 TGraph2D.cxx:1147
 TGraph2D.cxx:1148
 TGraph2D.cxx:1149
 TGraph2D.cxx:1150
 TGraph2D.cxx:1151
 TGraph2D.cxx:1152
 TGraph2D.cxx:1153
 TGraph2D.cxx:1154
 TGraph2D.cxx:1155
 TGraph2D.cxx:1156
 TGraph2D.cxx:1157
 TGraph2D.cxx:1158
 TGraph2D.cxx:1159
 TGraph2D.cxx:1160
 TGraph2D.cxx:1161
 TGraph2D.cxx:1162
 TGraph2D.cxx:1163
 TGraph2D.cxx:1164
 TGraph2D.cxx:1165
 TGraph2D.cxx:1166
 TGraph2D.cxx:1167
 TGraph2D.cxx:1168
 TGraph2D.cxx:1169
 TGraph2D.cxx:1170
 TGraph2D.cxx:1171
 TGraph2D.cxx:1172
 TGraph2D.cxx:1173
 TGraph2D.cxx:1174
 TGraph2D.cxx:1175
 TGraph2D.cxx:1176
 TGraph2D.cxx:1177
 TGraph2D.cxx:1178
 TGraph2D.cxx:1179
 TGraph2D.cxx:1180
 TGraph2D.cxx:1181
 TGraph2D.cxx:1182
 TGraph2D.cxx:1183
 TGraph2D.cxx:1184
 TGraph2D.cxx:1185
 TGraph2D.cxx:1186
 TGraph2D.cxx:1187
 TGraph2D.cxx:1188
 TGraph2D.cxx:1189
 TGraph2D.cxx:1190
 TGraph2D.cxx:1191
 TGraph2D.cxx:1192
 TGraph2D.cxx:1193
 TGraph2D.cxx:1194
 TGraph2D.cxx:1195
 TGraph2D.cxx:1196
 TGraph2D.cxx:1197
 TGraph2D.cxx:1198
 TGraph2D.cxx:1199
 TGraph2D.cxx:1200
 TGraph2D.cxx:1201
 TGraph2D.cxx:1202
 TGraph2D.cxx:1203
 TGraph2D.cxx:1204
 TGraph2D.cxx:1205
 TGraph2D.cxx:1206
 TGraph2D.cxx:1207
 TGraph2D.cxx:1208
 TGraph2D.cxx:1209
 TGraph2D.cxx:1210
 TGraph2D.cxx:1211
 TGraph2D.cxx:1212
 TGraph2D.cxx:1213
 TGraph2D.cxx:1214
 TGraph2D.cxx:1215
 TGraph2D.cxx:1216
 TGraph2D.cxx:1217
 TGraph2D.cxx:1218
 TGraph2D.cxx:1219
 TGraph2D.cxx:1220
 TGraph2D.cxx:1221
 TGraph2D.cxx:1222
 TGraph2D.cxx:1223
 TGraph2D.cxx:1224
 TGraph2D.cxx:1225
 TGraph2D.cxx:1226
 TGraph2D.cxx:1227
 TGraph2D.cxx:1228
 TGraph2D.cxx:1229
 TGraph2D.cxx:1230
 TGraph2D.cxx:1231
 TGraph2D.cxx:1232
 TGraph2D.cxx:1233
 TGraph2D.cxx:1234
 TGraph2D.cxx:1235
 TGraph2D.cxx:1236
 TGraph2D.cxx:1237
 TGraph2D.cxx:1238
 TGraph2D.cxx:1239
 TGraph2D.cxx:1240
 TGraph2D.cxx:1241
 TGraph2D.cxx:1242
 TGraph2D.cxx:1243
 TGraph2D.cxx:1244
 TGraph2D.cxx:1245
 TGraph2D.cxx:1246
 TGraph2D.cxx:1247
 TGraph2D.cxx:1248
 TGraph2D.cxx:1249
 TGraph2D.cxx:1250
 TGraph2D.cxx:1251
 TGraph2D.cxx:1252
 TGraph2D.cxx:1253
 TGraph2D.cxx:1254
 TGraph2D.cxx:1255
 TGraph2D.cxx:1256
 TGraph2D.cxx:1257
 TGraph2D.cxx:1258
 TGraph2D.cxx:1259
 TGraph2D.cxx:1260
 TGraph2D.cxx:1261
 TGraph2D.cxx:1262
 TGraph2D.cxx:1263
 TGraph2D.cxx:1264
 TGraph2D.cxx:1265
 TGraph2D.cxx:1266
 TGraph2D.cxx:1267
 TGraph2D.cxx:1268
 TGraph2D.cxx:1269
 TGraph2D.cxx:1270
 TGraph2D.cxx:1271
 TGraph2D.cxx:1272
 TGraph2D.cxx:1273
 TGraph2D.cxx:1274
 TGraph2D.cxx:1275
 TGraph2D.cxx:1276
 TGraph2D.cxx:1277
 TGraph2D.cxx:1278
 TGraph2D.cxx:1279
 TGraph2D.cxx:1280
 TGraph2D.cxx:1281
 TGraph2D.cxx:1282
 TGraph2D.cxx:1283
 TGraph2D.cxx:1284
 TGraph2D.cxx:1285
 TGraph2D.cxx:1286
 TGraph2D.cxx:1287
 TGraph2D.cxx:1288
 TGraph2D.cxx:1289
 TGraph2D.cxx:1290
 TGraph2D.cxx:1291
 TGraph2D.cxx:1292
 TGraph2D.cxx:1293
 TGraph2D.cxx:1294
 TGraph2D.cxx:1295
 TGraph2D.cxx:1296
 TGraph2D.cxx:1297
 TGraph2D.cxx:1298
 TGraph2D.cxx:1299
 TGraph2D.cxx:1300
 TGraph2D.cxx:1301
 TGraph2D.cxx:1302
 TGraph2D.cxx:1303
 TGraph2D.cxx:1304
 TGraph2D.cxx:1305
 TGraph2D.cxx:1306
 TGraph2D.cxx:1307
 TGraph2D.cxx:1308
 TGraph2D.cxx:1309
 TGraph2D.cxx:1310
 TGraph2D.cxx:1311
 TGraph2D.cxx:1312
 TGraph2D.cxx:1313
 TGraph2D.cxx:1314
 TGraph2D.cxx:1315
 TGraph2D.cxx:1316
 TGraph2D.cxx:1317
 TGraph2D.cxx:1318
 TGraph2D.cxx:1319
 TGraph2D.cxx:1320
 TGraph2D.cxx:1321
 TGraph2D.cxx:1322
 TGraph2D.cxx:1323
 TGraph2D.cxx:1324
 TGraph2D.cxx:1325
 TGraph2D.cxx:1326
 TGraph2D.cxx:1327
 TGraph2D.cxx:1328
 TGraph2D.cxx:1329
 TGraph2D.cxx:1330
 TGraph2D.cxx:1331
 TGraph2D.cxx:1332
 TGraph2D.cxx:1333
 TGraph2D.cxx:1334
 TGraph2D.cxx:1335
 TGraph2D.cxx:1336
 TGraph2D.cxx:1337
 TGraph2D.cxx:1338
 TGraph2D.cxx:1339
 TGraph2D.cxx:1340
 TGraph2D.cxx:1341
 TGraph2D.cxx:1342
 TGraph2D.cxx:1343
 TGraph2D.cxx:1344
 TGraph2D.cxx:1345
 TGraph2D.cxx:1346
 TGraph2D.cxx:1347
 TGraph2D.cxx:1348
 TGraph2D.cxx:1349
 TGraph2D.cxx:1350
 TGraph2D.cxx:1351
 TGraph2D.cxx:1352
 TGraph2D.cxx:1353
 TGraph2D.cxx:1354
 TGraph2D.cxx:1355
 TGraph2D.cxx:1356
 TGraph2D.cxx:1357
 TGraph2D.cxx:1358
 TGraph2D.cxx:1359
 TGraph2D.cxx:1360
 TGraph2D.cxx:1361
 TGraph2D.cxx:1362
 TGraph2D.cxx:1363
 TGraph2D.cxx:1364
 TGraph2D.cxx:1365
 TGraph2D.cxx:1366
 TGraph2D.cxx:1367
 TGraph2D.cxx:1368
 TGraph2D.cxx:1369
 TGraph2D.cxx:1370
 TGraph2D.cxx:1371
 TGraph2D.cxx:1372
 TGraph2D.cxx:1373
 TGraph2D.cxx:1374
 TGraph2D.cxx:1375
 TGraph2D.cxx:1376
 TGraph2D.cxx:1377
 TGraph2D.cxx:1378
 TGraph2D.cxx:1379
 TGraph2D.cxx:1380
 TGraph2D.cxx:1381
 TGraph2D.cxx:1382
 TGraph2D.cxx:1383
 TGraph2D.cxx:1384
 TGraph2D.cxx:1385
 TGraph2D.cxx:1386
 TGraph2D.cxx:1387
 TGraph2D.cxx:1388
 TGraph2D.cxx:1389
 TGraph2D.cxx:1390
 TGraph2D.cxx:1391
 TGraph2D.cxx:1392
 TGraph2D.cxx:1393
 TGraph2D.cxx:1394
 TGraph2D.cxx:1395
 TGraph2D.cxx:1396
 TGraph2D.cxx:1397
 TGraph2D.cxx:1398
 TGraph2D.cxx:1399
 TGraph2D.cxx:1400
 TGraph2D.cxx:1401
 TGraph2D.cxx:1402
 TGraph2D.cxx:1403
 TGraph2D.cxx:1404
 TGraph2D.cxx:1405
 TGraph2D.cxx:1406
 TGraph2D.cxx:1407
 TGraph2D.cxx:1408
 TGraph2D.cxx:1409
 TGraph2D.cxx:1410
 TGraph2D.cxx:1411
 TGraph2D.cxx:1412
 TGraph2D.cxx:1413
 TGraph2D.cxx:1414
 TGraph2D.cxx:1415
 TGraph2D.cxx:1416
 TGraph2D.cxx:1417
 TGraph2D.cxx:1418
 TGraph2D.cxx:1419
 TGraph2D.cxx:1420
 TGraph2D.cxx:1421
 TGraph2D.cxx:1422
 TGraph2D.cxx:1423
 TGraph2D.cxx:1424
 TGraph2D.cxx:1425
 TGraph2D.cxx:1426
 TGraph2D.cxx:1427
 TGraph2D.cxx:1428
 TGraph2D.cxx:1429
 TGraph2D.cxx:1430
 TGraph2D.cxx:1431
 TGraph2D.cxx:1432
 TGraph2D.cxx:1433
 TGraph2D.cxx:1434
 TGraph2D.cxx:1435
 TGraph2D.cxx:1436
 TGraph2D.cxx:1437
 TGraph2D.cxx:1438
 TGraph2D.cxx:1439
 TGraph2D.cxx:1440
 TGraph2D.cxx:1441
 TGraph2D.cxx:1442
 TGraph2D.cxx:1443
 TGraph2D.cxx:1444
 TGraph2D.cxx:1445
 TGraph2D.cxx:1446
 TGraph2D.cxx:1447
 TGraph2D.cxx:1448
 TGraph2D.cxx:1449
 TGraph2D.cxx:1450
 TGraph2D.cxx:1451
 TGraph2D.cxx:1452
 TGraph2D.cxx:1453
 TGraph2D.cxx:1454
 TGraph2D.cxx:1455
 TGraph2D.cxx:1456
 TGraph2D.cxx:1457
 TGraph2D.cxx:1458
 TGraph2D.cxx:1459
 TGraph2D.cxx:1460
 TGraph2D.cxx:1461
 TGraph2D.cxx:1462
 TGraph2D.cxx:1463
 TGraph2D.cxx:1464
 TGraph2D.cxx:1465
 TGraph2D.cxx:1466
 TGraph2D.cxx:1467
 TGraph2D.cxx:1468
 TGraph2D.cxx:1469
 TGraph2D.cxx:1470
 TGraph2D.cxx:1471
 TGraph2D.cxx:1472
 TGraph2D.cxx:1473
 TGraph2D.cxx:1474
 TGraph2D.cxx:1475
 TGraph2D.cxx:1476
 TGraph2D.cxx:1477
 TGraph2D.cxx:1478
 TGraph2D.cxx:1479
 TGraph2D.cxx:1480
 TGraph2D.cxx:1481
 TGraph2D.cxx:1482
 TGraph2D.cxx:1483
 TGraph2D.cxx:1484
 TGraph2D.cxx:1485
 TGraph2D.cxx:1486
 TGraph2D.cxx:1487
 TGraph2D.cxx:1488
 TGraph2D.cxx:1489
 TGraph2D.cxx:1490
 TGraph2D.cxx:1491
 TGraph2D.cxx:1492
 TGraph2D.cxx:1493
 TGraph2D.cxx:1494
 TGraph2D.cxx:1495
 TGraph2D.cxx:1496