Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TGraph2D.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id: TGraph2D.cxx,v 1.00
2 // Author: Olivier Couet
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "Riostream.h"
13 #include "TROOT.h"
14 #include "TMath.h"
15 #include "TH2.h"
16 #include "TF2.h"
17 #include "TList.h"
18 #include "TGraph2D.h"
19 #include "TGraphDelaunay.h"
20 #include "TGraphDelaunay2D.h"
21 #include "TVirtualPad.h"
22 #include "TVirtualFitter.h"
23 #include "TPluginManager.h"
24 #include "TClass.h"
25 #include "TSystem.h"
26 #include <stdlib.h>
27 #include <cassert>
28 
29 #include "HFitInterface.h"
30 #include "Fit/DataRange.h"
31 #include "Math/MinimizerOptions.h"
32 
33 ClassImp(TGraph2D);
34 
35 
36 /** \class TGraph2D
37  \ingroup Hist
38 Graphics object made of three arrays X, Y and Z with the same number of points each.
39 
40 This class has different constructors:
41 - With an array's dimension and three arrays x, y, and z:
42 ~~~ {.cpp}
43  TGraph2D *g = new TGraph2D(n, x, y, z);
44 ~~~
45  x, y, z arrays can be doubles, floats, or ints.
46 - With an array's dimension only:
47 ~~~ {.cpp}
48  TGraph2D *g = new TGraph2D(n);
49 ~~~
50  The internal arrays are then filled with `SetPoint()`. The following line
51  fills the internal arrays at the position `i` with the values
52  `x`, `y`, `z`.
53 ~~~ {.cpp}
54  g->SetPoint(i, x, y, z);
55 ~~~
56 - Without parameters:
57 ~~~ {.cpp}
58  TGraph2D *g = new TGraph2D();
59 ~~~
60  again `SetPoint()` must be used to fill the internal arrays.
61 - From a file:
62 ~~~ {.cpp}
63  TGraph2D *g = new TGraph2D("graph.dat");
64 ~~~
65  Arrays are read from the ASCII file "graph.dat" according to a specifies
66  format. The default format is `%%lg %%lg %%lg`
67 
68 Note that in any of these three cases, `SetPoint()` can be used to change a data
69 point or add a new one. If the data point index (`i`) is greater than the
70 current size of the internal arrays, they are automatically extended.
71 
72 Like TGraph the TGraph2D constructors do not have the TGraph2D title and name as parameters.
73 A TGraph2D has the default title and name "Graph2D". To change the default title
74 and name `SetTitle` and `SetName` should be called on the TGraph2D after its creation.
75 
76 Specific drawing options can be used to paint a TGraph2D:
77 
78 
79 | Option | Description |
80 |----------|-------------------------------------------------------------------|
81 | "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. |
82 | "TRIW" | The Delaunay triangles are drawn as wire frame. |
83 | "TRI1" | The Delaunay triangles are painted with color levels. The edges of each triangles are painted with the current line color. |
84 | "TRI2" | The Delaunay triangles are painted with color levels. |
85 | "P" | Draw a marker at each vertex. |
86 | "P0" | Draw a circle at each vertex. Each circle background is white. |
87 | "PCOL" | Draw a marker at each vertex. The color of each marker is defined according to its Z position. |
88 | "LINE" | Draw a 3D polyline. |
89 
90 A TGraph2D can be also drawn with any options valid to draw a 2D histogram
91 (like `COL`, `SURF`, `LEGO`, `CONT` etc..).
92 
93 When a TGraph2D is drawn with one of the 2D histogram drawing option,
94 an intermediate 2D histogram is filled using the Delaunay triangles
95 to interpolate the data set. The 2D histogram has equidistant bins along the X
96 and Y directions. The number of bins along each direction can be change using
97 `SetNpx()` and `SetNpy()`. Each bin is filled with the Z
98 value found via a linear interpolation on the plane defined by the triangle above
99 the (X,Y) coordinates of the bin center.
100 
101 The existing (X,Y,Z) points can be randomly scattered.
102 The Delaunay triangles are build in the (X,Y) plane. These 2D triangles are then
103 used to define flat planes in (X,Y,Z) over which the interpolation is done to fill
104 the 2D histogram. The 3D triangles int takes build a 3D surface in
105 the form of tessellating triangles at various angles. The triangles found can be
106 drawn in 3D with one of the TGraph2D specific drawing options.
107 
108 The histogram generated by the Delaunay interpolation can be accessed using the
109 `GetHistogram()` method.
110 
111 The axis settings (title, ranges etc ...) can be changed accessing the axis via
112 the GetXaxis GetYaxis and GetZaxis methods. They access the histogram axis created
113 at drawing time only. Therefore they should called after the TGraph2D is drawn:
114 
115 ~~~ {.cpp}
116  TGraph2D *g = new TGraph2D();
117 
118  [...]
119 
120  g->Draw("tri1");
121  gPad->Update();
122  g->GetXaxis()->SetTitle("X axis title");
123 ~~~
124 
125 Example:
126 
127 Begin_Macro(source)
128 {
129  TCanvas *c = new TCanvas("c","Graph2D example",0,0,600,400);
130  Double_t x, y, z, P = 6.;
131  Int_t np = 200;
132  TGraph2D *dt = new TGraph2D();
133  dt->SetTitle("Graph title; X axis title; Y axis title; Z axis title");
134  TRandom *r = new TRandom();
135  for (Int_t N=0; N<np; N++) {
136  x = 2*P*(r->Rndm(N))-P;
137  y = 2*P*(r->Rndm(N))-P;
138  z = (sin(x)/x)*(sin(y)/y)+0.2;
139  dt->SetPoint(N,x,y,z);
140  }
141  gStyle->SetPalette(1);
142  dt->Draw("surf1");
143  return c;
144 }
145 End_Macro
146 
147 2D graphs can be fitted as shown by the following example:
148 
149 Begin_Macro(source)
150 ../../../tutorials/fit/graph2dfit.C
151 End_Macro
152 
153 Example showing the PCOL option.
154 
155 Begin_Macro(source)
156 {
157  TCanvas *c1 = new TCanvas("c1","Graph2D example",0,0,600,400);
158  Double_t P = 5.;
159  Int_t npx = 20 ;
160  Int_t npy = 20 ;
161  Double_t x = -P;
162  Double_t y = -P;
163  Double_t z;
164  Int_t k = 0;
165  Double_t dx = (2*P)/npx;
166  Double_t dy = (2*P)/npy;
167  TGraph2D *dt = new TGraph2D(npx*npy);
168  dt->SetNpy(41);
169  dt->SetNpx(40);
170  for (Int_t i=0; i<npx; i++) {
171  for (Int_t j=0; j<npy; j++) {
172  z = sin(sqrt(x*x+y*y))+1;
173  dt->SetPoint(k,x,y,z);
174  k++;
175  y = y+dy;
176  }
177  x = x+dx;
178  y = -P;
179  }
180  gStyle->SetPalette(1);
181  dt->SetMarkerStyle(20);
182  dt->Draw("pcol");
183  return c1;
184 }
185 End_Macro
186 
187 ### Definition of Delaunay triangulation (After B. Delaunay)
188 For a set S of points in the Euclidean plane, the unique triangulation DT(S)
189 of S such that no point in S is inside the circumcircle of any triangle in
190 DT(S). DT(S) is the dual of the Voronoi diagram of S.
191 If n is the number of points in S, the Voronoi diagram of S is the partitioning
192 of the plane containing S points into n convex polygons such that each polygon
193 contains exactly one point and every point in a given polygon is closer to its
194 central point than to any other. A Voronoi diagram is sometimes also known as
195 a Dirichlet tessellation.
196 
197 \image html tgraph2d_delaunay.png
198 
199 [This applet](http://www.cs.cornell.edu/Info/People/chew/Delaunay.html)
200 gives a nice practical view of Delaunay triangulation and Voronoi diagram.
201 */
202 
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 /// Graph2D default constructor
206 
207 TGraph2D::TGraph2D()
208  : TNamed("Graph2D", "Graph2D"), TAttLine(1, 1, 1), TAttFill(0, 1001),
209  TAttMarker(), fNpoints(0)
210 {
211  fSize = 0;
212  fMargin = 0.;
213  fNpx = 40;
214  fNpy = 40;
215  fDirectory = 0;
216  fHistogram = 0;
217  fDelaunay = nullptr;
218  fMaximum = -1111;
219  fMinimum = -1111;
220  fX = 0;
221  fY = 0;
222  fZ = 0;
223  fZout = 0;
224  fMaxIter = 100000;
225  fPainter = 0;
226  fFunctions = new TList;
227  fUserHisto = kFALSE;
228 }
229 
230 
231 ////////////////////////////////////////////////////////////////////////////////
232 /// Graph2D constructor with three vectors of ints as input.
233 
234 TGraph2D::TGraph2D(Int_t n, Int_t *x, Int_t *y, Int_t *z)
235  : TNamed("Graph2D", "Graph2D"), TAttLine(1, 1, 1), TAttFill(0, 1001),
236  TAttMarker(), fNpoints(n)
237 {
238  Build(n);
239 
240  // Copy the input vectors into local arrays
241  for (Int_t i = 0; i < fNpoints; ++i) {
242  fX[i] = (Double_t)x[i];
243  fY[i] = (Double_t)y[i];
244  fZ[i] = (Double_t)z[i];
245  }
246 }
247 
248 
249 ////////////////////////////////////////////////////////////////////////////////
250 /// Graph2D constructor with three vectors of floats as input.
251 
252 TGraph2D::TGraph2D(Int_t n, Float_t *x, Float_t *y, Float_t *z)
253  : TNamed("Graph2D", "Graph2D"), TAttLine(1, 1, 1), TAttFill(0, 1001),
254  TAttMarker(), fNpoints(n)
255 {
256  Build(n);
257 
258  // Copy the input vectors into local arrays
259  for (Int_t i = 0; i < fNpoints; ++i) {
260  fX[i] = x[i];
261  fY[i] = y[i];
262  fZ[i] = z[i];
263  }
264 }
265 
266 
267 ////////////////////////////////////////////////////////////////////////////////
268 /// Graph2D constructor with three vectors of doubles as input.
269 
270 TGraph2D::TGraph2D(Int_t n, Double_t *x, Double_t *y, Double_t *z)
271  : TNamed("Graph2D", "Graph2D"), TAttLine(1, 1, 1), TAttFill(0, 1001),
272  TAttMarker(), fNpoints(n)
273 {
274  Build(n);
275 
276  // Copy the input vectors into local arrays
277  for (Int_t i = 0; i < fNpoints; ++i) {
278  fX[i] = x[i];
279  fY[i] = y[i];
280  fZ[i] = z[i];
281  }
282 }
283 
284 
285 ////////////////////////////////////////////////////////////////////////////////
286 /// Graph2D constructor with a TH2 (h2) as input.
287 /// Only the h2's bins within the X and Y axis ranges are used.
288 /// Empty bins, recognized when both content and errors are zero, are excluded.
289 
290 TGraph2D::TGraph2D(TH2 *h2)
291  : TNamed("Graph2D", "Graph2D"), TAttLine(1, 1, 1), TAttFill(0, 1001),
292  TAttMarker(), fNpoints(0)
293 {
294  Build(h2->GetNbinsX()*h2->GetNbinsY());
295 
296  TString gname = "Graph2D_from_" + TString(h2->GetName());
297  SetName(gname);
298  // need to call later because sets title in ref histogram
299  SetTitle(h2->GetTitle());
300 
301 
302 
303  TAxis *xaxis = h2->GetXaxis();
304  TAxis *yaxis = h2->GetYaxis();
305  Int_t xfirst = xaxis->GetFirst();
306  Int_t xlast = xaxis->GetLast();
307  Int_t yfirst = yaxis->GetFirst();
308  Int_t ylast = yaxis->GetLast();
309 
310 
311  Double_t x, y, z;
312  Int_t k = 0;
313 
314  for (Int_t i = xfirst; i <= xlast; i++) {
315  for (Int_t j = yfirst; j <= ylast; j++) {
316  x = xaxis->GetBinCenter(i);
317  y = yaxis->GetBinCenter(j);
318  z = h2->GetBinContent(i, j);
319  Double_t ez = h2->GetBinError(i, j);
320  if (z != 0. || ez != 0) {
321  SetPoint(k, x, y, z);
322  k++;
323  }
324  }
325  }
326 }
327 
328 
329 ////////////////////////////////////////////////////////////////////////////////
330 /// Graph2D constructor with name, title and three vectors of doubles as input.
331 /// name : name of 2D graph (avoid blanks)
332 /// title : 2D graph title
333 /// if title is of the form "stringt;stringx;stringy;stringz"
334 /// the 2D graph title is set to stringt, the x axis title to stringx,
335 /// the y axis title to stringy,etc
336 
337 TGraph2D::TGraph2D(const char *name, const char *title,
338  Int_t n, Double_t *x, Double_t *y, Double_t *z)
339  : TNamed(name, title), TAttLine(1, 1, 1), TAttFill(0, 1001),
340  TAttMarker(), fNpoints(n)
341 {
342  Build(n);
343 
344  // Copy the input vectors into local arrays
345  for (Int_t i = 0; i < fNpoints; ++i) {
346  fX[i] = x[i];
347  fY[i] = y[i];
348  fZ[i] = z[i];
349  }
350 }
351 
352 
353 ////////////////////////////////////////////////////////////////////////////////
354 /// Graph2D constructor. The arrays fX, fY and fZ should be filled via
355 /// calls to SetPoint
356 
357 TGraph2D::TGraph2D(Int_t n)
358  : TNamed("Graph2D", "Graph2D"), TAttLine(1, 1, 1), TAttFill(0, 1001),
359  TAttMarker(), fNpoints(n)
360 {
361  Build(n);
362  for (Int_t i = 0; i < fNpoints; i++) {
363  fX[i] = 0.;
364  fY[i] = 0.;
365  fZ[i] = 0.;
366  }
367 }
368 
369 
370 ////////////////////////////////////////////////////////////////////////////////
371 /// Graph2D constructor reading input from filename
372 /// filename is assumed to contain at least three columns of numbers.
373 /// For files separated by a specific delimiter different from ' ' and '\t' (e.g. ';' in csv files)
374 /// you can avoid using %*s to bypass this delimiter by explicitly specify the "option" argument,
375 /// e.g. option=" \t,;" for columns of figures separated by any of these characters (' ', '\t', ',', ';')
376 /// used once (e.g. "1;1") or in a combined way (" 1;,;; 1").
377 /// Note in that case, the instantiation is about 2 times slower.
378 
379 TGraph2D::TGraph2D(const char *filename, const char *format, Option_t *option)
380  : TNamed("Graph2D", filename), TAttLine(1, 1, 1), TAttFill(0, 1001),
381  TAttMarker(), fNpoints(0)
382 {
383  Double_t x, y, z;
384  TString fname = filename;
385  gSystem->ExpandPathName(fname);
386 
387  std::ifstream infile(fname.Data());
388  if (!infile.good()) {
389  MakeZombie();
390  Error("TGraph2D", "Cannot open file: %s, TGraph2D is Zombie", filename);
391  return;
392  } else {
393  Build(100);
394  }
395  std::string line;
396  Int_t np = 0;
397 
398  if (strcmp(option, "") == 0) { // No delimiters specified (standard constructor).
399 
400  while (std::getline(infile, line, '\n')) {
401  if (3 != sscanf(line.c_str(), format, &x, &y, &z)) {
402  continue; // skip empty and ill-formed lines
403  }
404  SetPoint(np, x, y, z);
405  np++;
406  }
407 
408  } else { // A delimiter has been specified in "option"
409 
410  // Checking format and creating its boolean equivalent
411  TString format_ = TString(format) ;
412  format_.ReplaceAll(" ", "") ;
413  format_.ReplaceAll("\t", "") ;
414  format_.ReplaceAll("lg", "") ;
415  format_.ReplaceAll("s", "") ;
416  format_.ReplaceAll("%*", "0") ;
417  format_.ReplaceAll("%", "1") ;
418  if (!format_.IsDigit()) {
419  Error("TGraph2D", "Incorrect input format! Allowed format tags are {\"%%lg\",\"%%*lg\" or \"%%*s\"}");
420  return;
421  }
422  Int_t ntokens = format_.Length() ;
423  if (ntokens < 3) {
424  Error("TGraph2D", "Incorrect input format! Only %d tag(s) in format whereas 3 \"%%lg\" tags are expected!", ntokens);
425  return;
426  }
427  Int_t ntokensToBeSaved = 0 ;
428  Bool_t * isTokenToBeSaved = new Bool_t [ntokens] ;
429  for (Int_t idx = 0; idx < ntokens; idx++) {
430  isTokenToBeSaved[idx] = TString::Format("%c", format_[idx]).Atoi() ; //atoi(&format_[idx]) does not work for some reason...
431  if (isTokenToBeSaved[idx] == 1) {
432  ntokensToBeSaved++ ;
433  }
434  }
435  if (ntokens >= 3 && ntokensToBeSaved != 3) { //first condition not to repeat the previous error message
436  Error("TGraph2D", "Incorrect input format! There are %d \"%%lg\" tag(s) in format whereas 3 and only 3 are expected!", ntokensToBeSaved);
437  delete [] isTokenToBeSaved ;
438  return;
439  }
440 
441  // Initializing loop variables
442  Bool_t isLineToBeSkipped = kFALSE ; //empty and ill-formed lines
443  char * token = NULL ;
444  TString token_str = "" ;
445  Int_t token_idx = 0 ;
446  Double_t * value = new Double_t [3] ; //x,y,z buffers
447  Int_t value_idx = 0 ;
448 
449  // Looping
450  char *rest;
451  while (std::getline(infile, line, '\n')) {
452  if (line != "") {
453  if (line[line.size() - 1] == char(13)) { // removing DOS CR character
454  line.erase(line.end() - 1, line.end()) ;
455  }
456  token = R__STRTOK_R(const_cast<char*>(line.c_str()), option, &rest);
457  while (token != NULL && value_idx < 3) {
458  if (isTokenToBeSaved[token_idx]) {
459  token_str = TString(token) ;
460  token_str.ReplaceAll("\t", "") ;
461  if (!token_str.IsFloat()) {
462  isLineToBeSkipped = kTRUE ;
463  break ;
464  } else {
465  value[value_idx] = token_str.Atof() ;
466  value_idx++ ;
467  }
468  }
469  token = R__STRTOK_R(NULL, option, &rest); // next token
470  token_idx++ ;
471  }
472  if (!isLineToBeSkipped && value_idx == 3) {
473  x = value[0] ;
474  y = value[1] ;
475  z = value[2] ;
476  SetPoint(np, x, y, z) ;
477  np++ ;
478  }
479  }
480  isLineToBeSkipped = kFALSE ;
481  token = NULL ;
482  token_idx = 0 ;
483  value_idx = 0 ;
484  }
485 
486  // Cleaning
487  delete [] isTokenToBeSaved ;
488  delete [] value ;
489  delete token ;
490  }
491  infile.close();
492 }
493 
494 
495 ////////////////////////////////////////////////////////////////////////////////
496 /// Graph2D copy constructor.
497 /// copy everything apart from the list of contained functions
498 
499 TGraph2D::TGraph2D(const TGraph2D &g)
500 : TNamed(g), TAttLine(g), TAttFill(g), TAttMarker(g),
501  fX(0), fY(0), fZ(0),
502  fHistogram(0), fDirectory(0), fPainter(0)
503 {
504  fFunctions = new TList(); // do not copy the functions
505 
506  // use operator=
507  (*this) = g;
508 
509  // append TGraph2D to gdirectory
510  if (TH1::AddDirectoryStatus()) {
511  fDirectory = gDirectory;
512  if (fDirectory) {
513  // append without replacing existing objects
514  fDirectory->Append(this);
515  }
516  }
517 
518 
519 }
520 
521 
522 ////////////////////////////////////////////////////////////////////////////////
523 /// TGraph2D destructor.
524 
525 TGraph2D::~TGraph2D()
526 {
527  Clear();
528 }
529 
530 
531 ////////////////////////////////////////////////////////////////////////////////
532 /// Graph2D operator "="
533 
534 TGraph2D& TGraph2D::operator=(const TGraph2D &g)
535 {
536  if (this == &g) return *this;
537 
538  // delete before existing contained objects
539  if (fX) delete [] fX;
540  if (fY) delete [] fY;
541  if (fZ) delete [] fZ;
542  if (fHistogram && !fUserHisto) {
543  delete fHistogram;
544  fHistogram = nullptr;
545  fDelaunay = nullptr;
546  }
547  // copy everything except the function list
548  fNpoints = g.fNpoints;
549  fNpx = g.fNpx;
550  fNpy = g.fNpy;
551  fMaxIter = g.fMaxIter;
552  fSize = fNpoints; // force size to be the same of npoints
553  fX = (fSize > 0) ? new Double_t[fSize] : 0;
554  fY = (fSize > 0) ? new Double_t[fSize] : 0;
555  fZ = (fSize > 0) ? new Double_t[fSize] : 0;
556  fMinimum = g.fMinimum;
557  fMaximum = g.fMaximum;
558  fMargin = g.fMargin;
559  fZout = g.fZout;
560  fUserHisto = g.fUserHisto;
561  if (g.fHistogram)
562  fHistogram = (fUserHisto ) ? g.fHistogram : new TH2D(*g.fHistogram);
563 
564 
565 
566  // copy the points
567  for (Int_t n = 0; n < fSize; n++) {
568  fX[n] = g.fX[n];
569  fY[n] = g.fY[n];
570  fZ[n] = g.fZ[n];
571  }
572 
573  return *this;
574 }
575 
576 ////////////////////////////////////////////////////////////////////////////////
577 /// Creates the 2D graph basic data structure
578 
579 void TGraph2D::Build(Int_t n)
580 {
581  if (n <= 0) {
582  Error("TGraph2D", "Invalid number of points (%d)", n);
583  return;
584  }
585 
586  fSize = n;
587  fMargin = 0.;
588  fNpx = 40;
589  fNpy = 40;
590  fDirectory = 0;
591  fHistogram = 0;
592  fDelaunay = nullptr;
593  fMaximum = -1111;
594  fMinimum = -1111;
595  fX = new Double_t[fSize];
596  fY = new Double_t[fSize];
597  fZ = new Double_t[fSize];
598  fZout = 0;
599  fMaxIter = 100000;
600  fFunctions = new TList;
601  fPainter = 0;
602  fUserHisto = kFALSE;
603 
604  if (TH1::AddDirectoryStatus()) {
605  fDirectory = gDirectory;
606  if (fDirectory) {
607  fDirectory->Append(this, kTRUE);
608  }
609  }
610 }
611 
612 
613 ////////////////////////////////////////////////////////////////////////////////
614 /// Browse
615 
616 void TGraph2D::Browse(TBrowser *)
617 {
618  Draw("p0");
619  gPad->Update();
620 }
621 
622 
623 ////////////////////////////////////////////////////////////////////////////////
624 /// Free all memory allocated by this object.
625 
626 void TGraph2D::Clear(Option_t * /*option = "" */)
627 {
628  if (fX) delete [] fX;
629  fX = 0;
630  if (fY) delete [] fY;
631  fY = 0;
632  if (fZ) delete [] fZ;
633  fZ = 0;
634  fSize = fNpoints = 0;
635  if (fHistogram && !fUserHisto) {
636  delete fHistogram;
637  fHistogram = nullptr;
638  fDelaunay = nullptr;
639  }
640  if (fFunctions) {
641  fFunctions->SetBit(kInvalidObject);
642  fFunctions->Delete();
643  delete fFunctions;
644  fFunctions = 0;
645  }
646  if (fDirectory) {
647  fDirectory->Remove(this);
648  fDirectory = 0;
649  }
650 }
651 
652 
653 ////////////////////////////////////////////////////////////////////////////////
654 /// Perform the automatic addition of the graph to the given directory
655 ///
656 /// Note this function is called in place when the semantic requires
657 /// this object to be added to a directory (I.e. when being read from
658 /// a TKey or being Cloned)
659 
660 void TGraph2D::DirectoryAutoAdd(TDirectory *dir)
661 {
662  Bool_t addStatus = TH1::AddDirectoryStatus();
663  if (addStatus) {
664  SetDirectory(dir);
665  if (dir) {
666  ResetBit(kCanDelete);
667  }
668  }
669 }
670 
671 
672 ////////////////////////////////////////////////////////////////////////////////
673 /// Computes distance from point px,py to a graph
674 
675 Int_t TGraph2D::DistancetoPrimitive(Int_t px, Int_t py)
676 {
677  Int_t distance = 9999;
678  if (fHistogram) distance = fHistogram->DistancetoPrimitive(px, py);
679  return distance;
680 }
681 
682 
683 ////////////////////////////////////////////////////////////////////////////////
684 /// Specific drawing options can be used to paint a TGraph2D:
685 ///
686 /// - "TRI" : The Delaunay triangles are drawn using filled area.
687 /// An hidden surface drawing technique is used. The surface is
688 /// painted with the current fill area color. The edges of each
689 /// triangles are painted with the current line color.
690 /// - "TRIW" : The Delaunay triangles are drawn as wire frame
691 /// - "TRI1" : The Delaunay triangles are painted with color levels. The edges
692 /// of each triangles are painted with the current line color.
693 /// - "TRI2" : the Delaunay triangles are painted with color levels.
694 /// - "P" : Draw a marker at each vertex
695 /// - "P0" : Draw a circle at each vertex. Each circle background is white.
696 /// - "PCOL" : Draw a marker at each vertex. The color of each marker is
697 /// defined according to its Z position.
698 /// - "CONT" : Draw contours
699 /// - "LINE" : Draw a 3D polyline
700 ///
701 /// A TGraph2D can be also drawn with ANY options valid to draw a 2D histogram.
702 ///
703 /// When a TGraph2D is drawn with one of the 2D histogram drawing option,
704 /// a intermediate 2D histogram is filled using the Delaunay triangles
705 /// technique to interpolate the data set.
706 
707 void TGraph2D::Draw(Option_t *option)
708 {
709  TString opt = option;
710  opt.ToLower();
711  if (gPad) {
712  if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
713  if (!opt.Contains("same")) {
714  //the following statement is necessary in case one attempts to draw
715  //a temporary histogram already in the current pad
716  if (TestBit(kCanDelete)) gPad->GetListOfPrimitives()->Remove(this);
717  gPad->Clear();
718  }
719  }
720  AppendPad(opt.Data());
721 }
722 
723 
724 ////////////////////////////////////////////////////////////////////////////////
725 /// Executes action corresponding to one event
726 
727 void TGraph2D::ExecuteEvent(Int_t event, Int_t px, Int_t py)
728 {
729  if (fHistogram) fHistogram->ExecuteEvent(event, px, py);
730 }
731 
732 
733 ////////////////////////////////////////////////////////////////////////////////
734 /// search object named name in the list of functions
735 
736 TObject *TGraph2D::FindObject(const char *name) const
737 {
738  if (fFunctions) return fFunctions->FindObject(name);
739  return 0;
740 }
741 
742 
743 ////////////////////////////////////////////////////////////////////////////////
744 /// search object obj in the list of functions
745 
746 TObject *TGraph2D::FindObject(const TObject *obj) const
747 {
748  if (fFunctions) return fFunctions->FindObject(obj);
749  return 0;
750 }
751 
752 
753 ////////////////////////////////////////////////////////////////////////////////
754 /// Fits this graph with function with name fname
755 /// Predefined functions such as gaus, expo and poln are automatically
756 /// created by ROOT.
757 /// fname can also be a formula, accepted by the linear fitter (linear parts divided
758 /// by "++" sign), for example "x++sin(y)" for fitting "[0]*x+[1]*sin(y)"
759 
760 TFitResultPtr TGraph2D::Fit(const char *fname, Option_t *option, Option_t *)
761 {
762 
763  char *linear;
764  linear = (char*)strstr(fname, "++");
765 
766  if (linear) {
767  TF2 f2(fname, fname);
768  return Fit(&f2, option, "");
769  }
770  TF2 * f2 = (TF2*)gROOT->GetFunction(fname);
771  if (!f2) {
772  Printf("Unknown function: %s", fname);
773  return -1;
774  }
775  return Fit(f2, option, "");
776 
777 }
778 
779 
780 ////////////////////////////////////////////////////////////////////////////////
781 /// Fits this 2D graph with function f2
782 ///
783 /// f2 is an already predefined function created by TF2.
784 /// Predefined functions such as gaus, expo and poln are automatically
785 /// created by ROOT.
786 ///
787 /// The list of fit options is given in parameter option:
788 ///
789 /// | Option | Description |
790 /// |----------|-------------------------------------------------------------------|
791 /// | "W" | Set all weights to 1; ignore error bars |
792 /// | "U" | Use a User specified fitting algorithm (via SetFCN) |
793 /// | "Q" | Quiet mode (minimum printing) |
794 /// | "V" | Verbose mode (default is between Q and V) |
795 /// | "R" | Use the Range specified in the function range |
796 /// | "N" | Do not store the graphics function, do not draw |
797 /// | "0" | Do not plot the result of the fit. By default the fitted function is drawn unless the option "N" above is specified. |
798 /// | "+" | Add this new fitted function to the list of fitted functions (by default, any previous function is deleted) |
799 /// | "C" | In case of linear fitting, not calculate the chisquare (saves time) |
800 /// | "EX0" | When fitting a TGraphErrors do not consider errors in the coordinate |
801 /// | "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 |
802 /// | "S" | The result of the fit is returned in the TFitResultPtr (see below Access to the Fit Result) |
803 ///
804 /// In order to use the Range option, one must first create a function
805 /// with the expression to be fitted. For example, if your graph2d
806 /// has a defined range between -4 and 4 and you want to fit a gaussian
807 /// only in the interval 1 to 3, you can do:
808 /// ~~~ {.cpp}
809 /// TF2 *f2 = new TF2("f2","gaus",1,3);
810 /// graph2d->Fit("f2","R");
811 /// ~~~
812 ///
813 /// ### Setting initial conditions
814 ///
815 /// Parameters must be initialized before invoking the Fit function.
816 /// The setting of the parameter initial values is automatic for the
817 /// predefined functions : poln, expo, gaus. One can however disable
818 /// this automatic computation by specifying the option "B".
819 /// You can specify boundary limits for some or all parameters via
820 /// ~~~ {.cpp}
821 /// f2->SetParLimits(p_number, parmin, parmax);
822 /// ~~~
823 /// if parmin>=parmax, the parameter is fixed
824 /// Note that you are not forced to fix the limits for all parameters.
825 /// For example, if you fit a function with 6 parameters, you can do:
826 /// ~~~ {.cpp}
827 /// func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
828 /// func->SetParLimits(4,-10,-4);
829 /// func->SetParLimits(5, 1,1);
830 /// ~~~
831 /// With this setup, parameters 0->3 can vary freely
832 /// Parameter 4 has boundaries [-10,-4] with initial value -8
833 /// Parameter 5 is fixed to 100.
834 ///
835 /// ### Fit range
836 ///
837 /// The fit range can be specified in two ways:
838 /// - specify rxmax > rxmin (default is rxmin=rxmax=0)
839 /// - specify the option "R". In this case, the function will be taken
840 /// instead of the full graph range.
841 ///
842 /// ### Changing the fitting function
843 ///
844 /// By default a chi2 fitting function is used for fitting a TGraph.
845 /// The function is implemented in FitUtil::EvaluateChi2.
846 /// In case of TGraph2DErrors an effective chi2 is used
847 /// (see TGraphErrors fit in TGraph::Fit) and is implemented in
848 /// FitUtil::EvaluateChi2Effective
849 /// To specify a User defined fitting function, specify option "U" and
850 /// call the following functions:
851 /// ~~~ {.cpp}
852 /// TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
853 /// ~~~
854 /// where MyFittingFunction is of type:
855 /// ~~~ {.cpp}
856 /// extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
857 /// ~~~
858 ///
859 /// ### Associated functions
860 ///
861 /// One or more object (typically a TF2*) can be added to the list
862 /// of functions (fFunctions) associated to each graph.
863 /// When TGraph::Fit is invoked, the fitted function is added to this list.
864 /// Given a graph gr, one can retrieve an associated function
865 /// with: TF2 *myfunc = gr->GetFunction("myfunc");
866 ///
867 /// ### Access to the fit results
868 ///
869 /// The function returns a TFitResultPtr which can hold a pointer to a TFitResult object.
870 /// By default the TFitResultPtr contains only the status of the fit and it converts automatically to an
871 /// integer. If the option "S" is instead used, TFitResultPtr contains the TFitResult and behaves as a smart
872 /// pointer to it. For example one can do:
873 /// ~~~ {.cpp}
874 /// TFitResultPtr r = graph->Fit("myFunc","S");
875 /// TMatrixDSym cov = r->GetCovarianceMatrix(); // to access the covariance matrix
876 /// Double_t par0 = r->Value(0); // retrieve the value for the parameter 0
877 /// Double_t err0 = r->Error(0); // retrieve the error for the parameter 0
878 /// r->Print("V"); // print full information of fit including covariance matrix
879 /// r->Write(); // store the result in a file
880 /// ~~~
881 ///
882 /// The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also
883 /// from the fitted function.
884 /// If the graph is made persistent, the list of
885 /// associated functions is also persistent. Given a pointer (see above)
886 /// to an associated function myfunc, one can retrieve the function/fit
887 /// parameters with calls such as:
888 /// ~~~ {.cpp}
889 /// Double_t chi2 = myfunc->GetChisquare();
890 /// Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
891 /// Double_t err0 = myfunc->GetParError(0); //error on first parameter
892 /// ~~~
893 ///
894 /// ### Fit Statistics
895 ///
896 /// You can change the statistics box to display the fit parameters with
897 /// the TStyle::SetOptFit(mode) method. This mode has four digits.
898 /// mode = pcev (default = 0111)
899 /// - v = 1; print name/values of parameters
900 /// - e = 1; print errors (if e=1, v must be 1)
901 /// - c = 1; print Chisquare/Number of degrees of freedom
902 /// - p = 1; print Probability
903 ///
904 /// For example: gStyle->SetOptFit(1011);
905 /// prints the fit probability, parameter names/values, and errors.
906 /// You can change the position of the statistics box with these lines
907 /// (where g is a pointer to the TGraph):
908 ///
909 /// ~~~ {.cpp}
910 /// Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
911 /// Root > st->SetX1NDC(newx1); //new x start position
912 /// Root > st->SetX2NDC(newx2); //new x end position
913 /// ~~~
914 
915 TFitResultPtr TGraph2D::Fit(TF2 *f2, Option_t *option, Option_t *)
916 {
917  // internal graph2D fitting methods
918  Foption_t fitOption;
919  Option_t *goption = "";
920  ROOT::Fit::FitOptionsMake(ROOT::Fit::kGraph, option, fitOption);
921 
922  // create range and minimizer options with default values
923  ROOT::Fit::DataRange range(2);
924  ROOT::Math::MinimizerOptions minOption;
925  return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range);
926 }
927 
928 
929 ////////////////////////////////////////////////////////////////////////////////
930 /// Display a GUI panel with all graph fit options.
931 ///
932 /// See class TFitEditor for example
933 
934 void TGraph2D::FitPanel()
935 {
936  if (!gPad)
937  gROOT->MakeDefCanvas();
938 
939  if (!gPad) {
940  Error("FitPanel", "Unable to create a default canvas");
941  return;
942  }
943 
944  // use plugin manager to create instance of TFitEditor
945  TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
946  if (handler && handler->LoadPlugin() != -1) {
947  if (handler->ExecPlugin(2, gPad, this) == 0)
948  Error("FitPanel", "Unable to crate the FitPanel");
949  } else
950  Error("FitPanel", "Unable to find the FitPanel plug-in");
951 
952 }
953 
954 
955 ////////////////////////////////////////////////////////////////////////////////
956 /// Get x axis of the graph.
957 
958 TAxis *TGraph2D::GetXaxis() const
959 {
960  TH1 *h = ((TGraph2D*)this)->GetHistogram("empty");
961  if (!h) return 0;
962  return h->GetXaxis();
963 }
964 
965 
966 ////////////////////////////////////////////////////////////////////////////////
967 /// Get y axis of the graph.
968 
969 TAxis *TGraph2D::GetYaxis() const
970 {
971  TH1 *h = ((TGraph2D*)this)->GetHistogram("empty");
972  if (!h) return 0;
973  return h->GetYaxis();
974 }
975 
976 
977 ////////////////////////////////////////////////////////////////////////////////
978 /// Get z axis of the graph.
979 
980 TAxis *TGraph2D::GetZaxis() const
981 {
982  TH1 *h = ((TGraph2D*)this)->GetHistogram("empty");
983  if (!h) return 0;
984  return h->GetZaxis();
985 }
986 
987 
988 ////////////////////////////////////////////////////////////////////////////////
989 /// Returns the X and Y graphs building a contour. A contour level may
990 /// consist in several parts not connected to each other. This function
991 /// returns them in a graphs' list.
992 
993 TList *TGraph2D::GetContourList(Double_t contour)
994 {
995  if (fNpoints <= 0) {
996  Error("GetContourList", "Empty TGraph2D");
997  return 0;
998  }
999 
1000  if (!fHistogram) GetHistogram("empty");
1001 
1002  if (!fPainter) fPainter = fHistogram->GetPainter();
1003 
1004  return fPainter->GetContourList(contour);
1005 }
1006 
1007 
1008 ////////////////////////////////////////////////////////////////////////////////
1009 /// This function is called by Graph2DFitChisquare.
1010 /// It always returns a negative value. Real implementation in TGraph2DErrors
1011 
1012 Double_t TGraph2D::GetErrorX(Int_t) const
1013 {
1014  return -1;
1015 }
1016 
1017 
1018 ////////////////////////////////////////////////////////////////////////////////
1019 /// This function is called by Graph2DFitChisquare.
1020 /// It always returns a negative value. Real implementation in TGraph2DErrors
1021 
1022 Double_t TGraph2D::GetErrorY(Int_t) const
1023 {
1024  return -1;
1025 }
1026 
1027 
1028 ////////////////////////////////////////////////////////////////////////////////
1029 /// This function is called by Graph2DFitChisquare.
1030 /// It always returns a negative value. Real implementation in TGraph2DErrors
1031 
1032 Double_t TGraph2D::GetErrorZ(Int_t) const
1033 {
1034  return -1;
1035 }
1036 
1037 
1038 ////////////////////////////////////////////////////////////////////////////////
1039 /// Add a TGraphDelaunay in the list of the fHistogram's functions
1040 
1041 void TGraph2D::CreateInterpolator(Bool_t oldInterp)
1042 {
1043 
1044  TList *hl = fHistogram->GetListOfFunctions();
1045 
1046  if (oldInterp) {
1047  TGraphDelaunay *dt = new TGraphDelaunay(this);
1048  dt->SetMaxIter(fMaxIter);
1049  dt->SetMarginBinsContent(fZout);
1050  fDelaunay = dt;
1051  SetBit(kOldInterpolation);
1052  if (!hl->FindObject("TGraphDelaunay")) hl->Add(fDelaunay);
1053  } else {
1054  TGraphDelaunay2D *dt = new TGraphDelaunay2D(this);
1055  dt->SetMarginBinsContent(fZout);
1056  fDelaunay = dt;
1057  ResetBit(kOldInterpolation);
1058  if (!hl->FindObject("TGraphDelaunay2D")) hl->Add(fDelaunay);
1059  }
1060 }
1061 
1062 ////////////////////////////////////////////////////////////////////////////////
1063 /// By default returns a pointer to the Delaunay histogram. If fHistogram
1064 /// doesn't exist, books the 2D histogram fHistogram with a margin around
1065 /// the hull. Calls TGraphDelaunay::Interpolate at each bin centre to build up
1066 /// an interpolated 2D histogram.
1067 ///
1068 /// If the "empty" option is selected, returns an empty histogram booked with
1069 /// the limits of fX, fY and fZ. This option is used when the data set is
1070 /// drawn with markers only. In that particular case there is no need to
1071 /// find the Delaunay triangles.
1072 ///
1073 /// By default use the new interpolation routine based on Triangles
1074 /// If the option "old" the old interpolation is used
1075 
1076 TH2D *TGraph2D::GetHistogram(Option_t *option)
1077 {
1078  // for an empty graph create histogram in [0,1][0,1]
1079  if (fNpoints <= 0) {
1080  if (!fHistogram) {
1081  Bool_t add = TH1::AddDirectoryStatus();
1082  TH1::AddDirectory(kFALSE);
1083  fHistogram = new TH2D(GetName(), GetTitle(), fNpx , 0., 1., fNpy, 0., 1.);
1084  TH1::AddDirectory(add);
1085  fHistogram->SetBit(TH1::kNoStats);
1086  }
1087  return fHistogram;
1088  }
1089 
1090  TString opt = option;
1091  opt.ToLower();
1092  Bool_t empty = opt.Contains("empty");
1093  Bool_t oldInterp = opt.Contains("old");
1094 
1095  if (fHistogram) {
1096  if (!empty && fHistogram->GetEntries() == 0) {
1097  if (!fUserHisto) {
1098  delete fHistogram;
1099  fHistogram = nullptr;
1100  fDelaunay = nullptr;
1101  }
1102  } else if (fHistogram->GetEntries() == 0)
1103  {; }
1104  // check case if interpolation type has changed
1105  else if ( (TestBit(kOldInterpolation) && !oldInterp) || ( !TestBit(kOldInterpolation) && oldInterp ) ) {
1106  delete fHistogram;
1107  fHistogram = nullptr;
1108  fDelaunay = nullptr;
1109  }
1110  // normal case return existing histogram
1111  else {
1112  return fHistogram;
1113  }
1114  }
1115 
1116  Double_t hxmax, hymax, hxmin, hymin;
1117 
1118  // Book fHistogram if needed. It is not added in the current directory
1119  if (!fUserHisto) {
1120  Bool_t add = TH1::AddDirectoryStatus();
1121  TH1::AddDirectory(kFALSE);
1122  Double_t xmax = GetXmaxE();
1123  Double_t ymax = GetYmaxE();
1124  Double_t xmin = GetXminE();
1125  Double_t ymin = GetYminE();
1126  hxmin = xmin - fMargin * (xmax - xmin);
1127  hymin = ymin - fMargin * (ymax - ymin);
1128  hxmax = xmax + fMargin * (xmax - xmin);
1129  hymax = ymax + fMargin * (ymax - ymin);
1130  if (TMath::Abs(hxmax - hxmin) < 0.0001) {
1131  if (TMath::Abs(hxmin) < 0.0001) {
1132  hxmin = -0.01;
1133  hxmax = 0.01;
1134  } else {
1135  hxmin = hxmin-TMath::Abs(hxmin)*0.01;
1136  hxmax = hxmax+TMath::Abs(hxmax)*0.01;
1137  }
1138  }
1139  if (TMath::Abs(hymax - hymin) < 0.0001) {
1140  if (TMath::Abs(hymin) < 0.0001) {
1141  hymin = -0.01;
1142  hymax = 0.01;
1143  } else {
1144  hymin = hymin-TMath::Abs(hymin)*0.01;
1145  hymax = hymax+TMath::Abs(hymax)*0.01;
1146  }
1147  }
1148  if (fHistogram) {
1149  fHistogram->GetXaxis()->SetLimits(hxmin, hxmax);
1150  fHistogram->GetYaxis()->SetLimits(hymin, hymax);
1151  } else {
1152  fHistogram = new TH2D(GetName(), GetTitle(),
1153  fNpx , hxmin, hxmax,
1154  fNpy, hymin, hymax);
1155  CreateInterpolator(oldInterp);
1156  }
1157  TH1::AddDirectory(add);
1158  fHistogram->SetBit(TH1::kNoStats);
1159  } else {
1160  hxmin = fHistogram->GetXaxis()->GetXmin();
1161  hymin = fHistogram->GetYaxis()->GetXmin();
1162  hxmax = fHistogram->GetXaxis()->GetXmax();
1163  hymax = fHistogram->GetYaxis()->GetXmax();
1164  }
1165 
1166  // Option "empty" is selected. An empty histogram is returned.
1167  if (empty) {
1168  Double_t hzmax, hzmin;
1169  if (fMinimum != -1111) {
1170  hzmin = fMinimum;
1171  } else {
1172  hzmin = GetZmin();
1173  }
1174  if (fMaximum != -1111) {
1175  hzmax = fMaximum;
1176  } else {
1177  hzmax = GetZmax();
1178  }
1179  if (hzmin == hzmax) {
1180  Double_t hz = hzmin;
1181  if (hz==0) hz = 1.;
1182  hzmin = hz - 0.01 * hz;
1183  hzmax = hz + 0.01 * hz;
1184  }
1185  fHistogram->SetMinimum(hzmin);
1186  fHistogram->SetMaximum(hzmax);
1187  return fHistogram;
1188  }
1189 
1190  Double_t dx = (hxmax - hxmin) / fNpx;
1191  Double_t dy = (hymax - hymin) / fNpy;
1192 
1193  Double_t x, y, z;
1194 
1195  for (Int_t ix = 1; ix <= fNpx; ix++) {
1196  x = hxmin + (ix - 0.5) * dx;
1197  for (Int_t iy = 1; iy <= fNpy; iy++) {
1198  y = hymin + (iy - 0.5) * dy;
1199  // do interpolation
1200  if (oldInterp)
1201  z = ((TGraphDelaunay*)fDelaunay)->ComputeZ(x, y);
1202  else
1203  z = ((TGraphDelaunay2D*)fDelaunay)->ComputeZ(x, y);
1204 
1205  fHistogram->Fill(x, y, z);
1206  }
1207  }
1208 
1209 
1210  if (fMinimum != -1111) fHistogram->SetMinimum(fMinimum);
1211  if (fMaximum != -1111) fHistogram->SetMaximum(fMaximum);
1212 
1213  return fHistogram;
1214 }
1215 
1216 
1217 ////////////////////////////////////////////////////////////////////////////////
1218 /// Returns the X maximum
1219 
1220 Double_t TGraph2D::GetXmax() const
1221 {
1222  Double_t v = fX[0];
1223  for (Int_t i = 1; i < fNpoints; i++) if (fX[i] > v) v = fX[i];
1224  return v;
1225 }
1226 
1227 
1228 ////////////////////////////////////////////////////////////////////////////////
1229 /// Returns the X minimum
1230 
1231 Double_t TGraph2D::GetXmin() const
1232 {
1233  Double_t v = fX[0];
1234  for (Int_t i = 1; i < fNpoints; i++) if (fX[i] < v) v = fX[i];
1235  return v;
1236 }
1237 
1238 
1239 ////////////////////////////////////////////////////////////////////////////////
1240 /// Returns the Y maximum
1241 
1242 Double_t TGraph2D::GetYmax() const
1243 {
1244  Double_t v = fY[0];
1245  for (Int_t i = 1; i < fNpoints; i++) if (fY[i] > v) v = fY[i];
1246  return v;
1247 }
1248 
1249 
1250 ////////////////////////////////////////////////////////////////////////////////
1251 /// Returns the Y minimum
1252 
1253 Double_t TGraph2D::GetYmin() const
1254 {
1255  Double_t v = fY[0];
1256  for (Int_t i = 1; i < fNpoints; i++) if (fY[i] < v) v = fY[i];
1257  return v;
1258 }
1259 
1260 
1261 ////////////////////////////////////////////////////////////////////////////////
1262 /// Returns the Z maximum
1263 
1264 Double_t TGraph2D::GetZmax() const
1265 {
1266  Double_t v = fZ[0];
1267  for (Int_t i = 1; i < fNpoints; i++) if (fZ[i] > v) v = fZ[i];
1268  return v;
1269 }
1270 
1271 
1272 ////////////////////////////////////////////////////////////////////////////////
1273 /// Returns the Z minimum
1274 
1275 Double_t TGraph2D::GetZmin() const
1276 {
1277  Double_t v = fZ[0];
1278  for (Int_t i = 1; i < fNpoints; i++) if (fZ[i] < v) v = fZ[i];
1279  return v;
1280 }
1281 
1282 ////////////////////////////////////////////////////////////////////////////////
1283 /// Get x, y and z values for point number i.
1284 /// The function returns -1 in case of an invalid request or the point number otherwise
1285 
1286 Int_t TGraph2D::GetPoint(Int_t i, Double_t &x, Double_t &y, Double_t &z) const
1287 {
1288  if (i < 0 || i >= fNpoints) return -1;
1289  if (!fX || !fY || !fZ) return -1;
1290  x = fX[i];
1291  y = fY[i];
1292  z = fZ[i];
1293  return i;
1294 }
1295 
1296 ////////////////////////////////////////////////////////////////////////////////
1297 /// Finds the z value at the position (x,y) thanks to
1298 /// the Delaunay interpolation.
1299 
1300 Double_t TGraph2D::Interpolate(Double_t x, Double_t y)
1301 {
1302  if (fNpoints <= 0) {
1303  Error("Interpolate", "Empty TGraph2D");
1304  return 0;
1305  }
1306 
1307  if (!fHistogram) GetHistogram("empty");
1308  if (!fDelaunay) {
1309  TList *hl = fHistogram->GetListOfFunctions();
1310  if (!TestBit(kOldInterpolation) ) {
1311  fDelaunay = hl->FindObject("TGraphDelaunay2D");
1312  if (!fDelaunay) fDelaunay = hl->FindObject("TGraphDelaunay");
1313  }
1314  else {
1315  // if using old implementation
1316  fDelaunay = hl->FindObject("TGraphDelaunay");
1317  if (!fDelaunay) fDelaunay = hl->FindObject("TGraphDelaunay2D");
1318  }
1319  }
1320 
1321  if (!fDelaunay) return TMath::QuietNaN();
1322 
1323  if (fDelaunay->IsA() == TGraphDelaunay2D::Class() )
1324  return ((TGraphDelaunay2D*)fDelaunay)->ComputeZ(x, y);
1325  else if (fDelaunay->IsA() == TGraphDelaunay::Class() )
1326  return ((TGraphDelaunay*)fDelaunay)->ComputeZ(x, y);
1327 
1328  // cannot be here
1329  assert(false);
1330  return TMath::QuietNaN();
1331 }
1332 
1333 
1334 ////////////////////////////////////////////////////////////////////////////////
1335 /// Paints this 2D graph with its current attributes
1336 
1337 void TGraph2D::Paint(Option_t *option)
1338 {
1339  if (fNpoints <= 0) {
1340  Error("Paint", "Empty TGraph2D");
1341  return;
1342  }
1343 
1344  TString opt = option;
1345  opt.ToLower();
1346  if (opt.Contains("p") && !opt.Contains("tri")) {
1347  if (!opt.Contains("pol") &&
1348  !opt.Contains("sph") &&
1349  !opt.Contains("psr")) opt.Append("tri0");
1350  }
1351 
1352  if (opt.Contains("line") && !opt.Contains("tri")) opt.Append("tri0");
1353 
1354  if (opt.Contains("err") && !opt.Contains("tri")) opt.Append("tri0");
1355 
1356  if (opt.Contains("tri0")) {
1357  GetHistogram("empty");
1358  } else if (opt.Contains("old")) {
1359  GetHistogram("old");
1360  } else {
1361  GetHistogram();
1362  }
1363 
1364  fHistogram->SetLineColor(GetLineColor());
1365  fHistogram->SetLineStyle(GetLineStyle());
1366  fHistogram->SetLineWidth(GetLineWidth());
1367  fHistogram->SetFillColor(GetFillColor());
1368  fHistogram->SetFillStyle(GetFillStyle());
1369  fHistogram->SetMarkerColor(GetMarkerColor());
1370  fHistogram->SetMarkerStyle(GetMarkerStyle());
1371  fHistogram->SetMarkerSize(GetMarkerSize());
1372  fHistogram->Paint(opt.Data());
1373 }
1374 
1375 
1376 ////////////////////////////////////////////////////////////////////////////////
1377 /// Print 2D graph values.
1378 
1379 void TGraph2D::Print(Option_t *) const
1380 {
1381  for (Int_t i = 0; i < fNpoints; i++) {
1382  printf("x[%d]=%g, y[%d]=%g, z[%d]=%g\n", i, fX[i], i, fY[i], i, fZ[i]);
1383  }
1384 }
1385 
1386 
1387 ////////////////////////////////////////////////////////////////////////////////
1388 /// Projects a 2-d graph into 1 or 2-d histograms depending on the option parameter.
1389 /// option may contain a combination of the characters x,y,z:
1390 ///
1391 /// - option = "x" return the x projection into a TH1D histogram
1392 /// - option = "y" return the y projection into a TH1D histogram
1393 /// - option = "xy" return the x versus y projection into a TH2D histogram
1394 /// - option = "yx" return the y versus x projection into a TH2D histogram
1395 
1396 TH1 *TGraph2D::Project(Option_t *option) const
1397 {
1398  if (fNpoints <= 0) {
1399  Error("Project", "Empty TGraph2D");
1400  return 0;
1401  }
1402 
1403  TString opt = option;
1404  opt.ToLower();
1405 
1406  Int_t pcase = 0;
1407  if (opt.Contains("x")) pcase = 1;
1408  if (opt.Contains("y")) pcase = 2;
1409  if (opt.Contains("xy")) pcase = 3;
1410  if (opt.Contains("yx")) pcase = 4;
1411 
1412  // Create the projection histogram
1413  TH1D *h1 = 0;
1414  TH2D *h2 = 0;
1415  Int_t nch = strlen(GetName()) + opt.Length() + 2;
1416  char *name = new char[nch];
1417  snprintf(name, nch, "%s_%s", GetName(), option);
1418  nch = strlen(GetTitle()) + opt.Length() + 2;
1419  char *title = new char[nch];
1420  snprintf(title, nch, "%s_%s", GetTitle(), option);
1421 
1422  Double_t hxmin = GetXmin();
1423  Double_t hxmax = GetXmax();
1424  Double_t hymin = GetYmin();
1425  Double_t hymax = GetYmax();
1426 
1427  switch (pcase) {
1428  case 1:
1429  // "x"
1430  h1 = new TH1D(name, title, fNpx, hxmin, hxmax);
1431  break;
1432  case 2:
1433  // "y"
1434  h1 = new TH1D(name, title, fNpy, hymin, hymax);
1435  break;
1436  case 3:
1437  // "xy"
1438  h2 = new TH2D(name, title, fNpx, hxmin, hxmax, fNpy, hymin, hymax);
1439  break;
1440  case 4:
1441  // "yx"
1442  h2 = new TH2D(name, title, fNpy, hymin, hymax, fNpx, hxmin, hxmax);
1443  break;
1444  }
1445 
1446  delete [] name;
1447  delete [] title;
1448  TH1 *h = h1;
1449  if (h2) h = h2;
1450  if (h == 0) return 0;
1451 
1452  // Fill the projected histogram
1453  Double_t entries = 0;
1454  for (Int_t n = 0; n < fNpoints; n++) {
1455  switch (pcase) {
1456  case 1:
1457  // "x"
1458  h1->Fill(fX[n], fZ[n]);
1459  break;
1460  case 2:
1461  // "y"
1462  h1->Fill(fY[n], fZ[n]);
1463  break;
1464  case 3:
1465  // "xy"
1466  h2->Fill(fX[n], fY[n], fZ[n]);
1467  break;
1468  case 4:
1469  // "yx"
1470  h2->Fill(fY[n], fX[n], fZ[n]);
1471  break;
1472  }
1473  entries += fZ[n];
1474  }
1475  h->SetEntries(entries);
1476  return h;
1477 }
1478 
1479 
1480 ////////////////////////////////////////////////////////////////////////////////
1481 /// Deletes point number ipoint
1482 
1483 Int_t TGraph2D::RemovePoint(Int_t ipoint)
1484 {
1485  if (ipoint < 0) return -1;
1486  if (ipoint >= fNpoints) return -1;
1487 
1488  fNpoints--;
1489  Double_t *newX = new Double_t[fNpoints];
1490  Double_t *newY = new Double_t[fNpoints];
1491  Double_t *newZ = new Double_t[fNpoints];
1492  Int_t j = -1;
1493  for (Int_t i = 0; i < fNpoints + 1; i++) {
1494  if (i == ipoint) continue;
1495  j++;
1496  newX[j] = fX[i];
1497  newY[j] = fY[i];
1498  newZ[j] = fZ[i];
1499  }
1500  delete [] fX;
1501  delete [] fY;
1502  delete [] fZ;
1503  fX = newX;
1504  fY = newY;
1505  fZ = newZ;
1506  fSize = fNpoints;
1507  if (fHistogram) {
1508  delete fHistogram;
1509  fHistogram = nullptr;
1510  fDelaunay = nullptr;
1511  }
1512  return ipoint;
1513 }
1514 
1515 
1516 ////////////////////////////////////////////////////////////////////////////////
1517 /// Saves primitive as a C++ statement(s) on output stream out
1518 
1519 void TGraph2D::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1520 {
1521  char quote = '"';
1522  out << " " << std::endl;
1523  if (gROOT->ClassSaved(TGraph2D::Class())) {
1524  out << " ";
1525  } else {
1526  out << " TGraph2D *";
1527  }
1528 
1529  out << "graph2d = new TGraph2D(" << fNpoints << ");" << std::endl;
1530  out << " graph2d->SetName(" << quote << GetName() << quote << ");" << std::endl;
1531  out << " graph2d->SetTitle(" << quote << GetTitle() << ";"
1532  << GetXaxis()->GetTitle() << ";"
1533  << GetYaxis()->GetTitle() << ";"
1534  << GetZaxis()->GetTitle() << quote << ");" << std::endl;
1535 
1536  if (fDirectory == 0) {
1537  out << " graph2d->SetDirectory(0);" << std::endl;
1538  }
1539 
1540  SaveFillAttributes(out, "graph2d", 0, 1001);
1541  SaveLineAttributes(out, "graph2d", 1, 1, 1);
1542  SaveMarkerAttributes(out, "graph2d", 1, 1, 1);
1543 
1544  for (Int_t i = 0; i < fNpoints; i++) {
1545  out << " graph2d->SetPoint(" << i << "," << fX[i] << "," << fY[i] << "," << fZ[i] << ");" << std::endl;
1546  }
1547 
1548  // save list of functions
1549  TIter next(fFunctions);
1550  TObject *obj;
1551  while ((obj = next())) {
1552  obj->SavePrimitive(out, "nodraw");
1553  out << " graph2d->GetListOfFunctions()->Add(" << obj->GetName() << ");" << std::endl;
1554  if (obj->InheritsFrom("TPaveStats")) {
1555  out << " ptstats->SetParent(graph2d->GetListOfFunctions());" << std::endl;
1556  } else if (obj->InheritsFrom("TF1")) {
1557  out << " " << obj->GetName() << "->SetParent(graph);\n";
1558  }
1559 
1560  }
1561 
1562  out << " graph2d->Draw(" << quote << option << quote << ");" << std::endl;
1563 }
1564 
1565 
1566 ////////////////////////////////////////////////////////////////////////////////
1567 /// Set number of points in the 2D graph.
1568 /// Existing coordinates are preserved.
1569 /// New coordinates above fNpoints are preset to 0.
1570 
1571 void TGraph2D::Set(Int_t n)
1572 {
1573  if (n < 0) n = 0;
1574  if (n == fNpoints) return;
1575  if (n > fNpoints) SetPoint(n, 0, 0, 0);
1576  fNpoints = n;
1577 }
1578 
1579 
1580 ////////////////////////////////////////////////////////////////////////////////
1581 /// By default when an 2D graph is created, it is added to the list
1582 /// of 2D graph objects in the current directory in memory.
1583 /// This method removes reference to this 2D graph from current directory and add
1584 /// reference to new directory dir. dir can be 0 in which case the
1585 /// 2D graph does not belong to any directory.
1586 
1587 void TGraph2D::SetDirectory(TDirectory *dir)
1588 {
1589  if (fDirectory == dir) return;
1590  if (fDirectory) fDirectory->Remove(this);
1591  fDirectory = dir;
1592  if (fDirectory) fDirectory->Append(this);
1593 }
1594 
1595 
1596 ////////////////////////////////////////////////////////////////////////////////
1597 /// Sets the histogram to be filled.
1598 /// If the 2D graph needs to be save in a TFile the following set should be
1599 /// followed to read it back:
1600 /// 1. Create TGraph2D
1601 /// 2. Call g->SetHistogram(h), and do whatever you need to do
1602 /// 3. Save g and h to the TFile, exit
1603 /// 4. Open the TFile, retrieve g and h
1604 /// 5. Call h->SetDirectory(0)
1605 /// 6. Call g->SetHistogram(h) again
1606 /// 7. Carry on as normal
1607 
1608 void TGraph2D::SetHistogram(TH2 *h)
1609 {
1610  fUserHisto = kTRUE;
1611  fHistogram = (TH2D*)h;
1612  fNpx = h->GetNbinsX();
1613  fNpy = h->GetNbinsY();
1614  CreateInterpolator(kTRUE);
1615 }
1616 
1617 
1618 ////////////////////////////////////////////////////////////////////////////////
1619 /// Sets the extra space (in %) around interpolated area for the 2D histogram
1620 
1621 void TGraph2D::SetMargin(Double_t m)
1622 {
1623  if (m < 0 || m > 1) {
1624  Warning("SetMargin", "The margin must be >= 0 && <= 1, fMargin set to 0.1");
1625  fMargin = 0.1;
1626  } else {
1627  fMargin = m;
1628  }
1629  if (fHistogram) {
1630  delete fHistogram;
1631  fHistogram = nullptr;
1632  fDelaunay = nullptr;
1633  }
1634 }
1635 
1636 
1637 ////////////////////////////////////////////////////////////////////////////////
1638 /// Sets the histogram bin height for points lying outside the TGraphDelaunay
1639 /// convex hull ie: the bins in the margin.
1640 
1641 void TGraph2D::SetMarginBinsContent(Double_t z)
1642 {
1643  fZout = z;
1644  if (fHistogram) {
1645  delete fHistogram;
1646  fHistogram = nullptr;
1647  fDelaunay = nullptr;
1648  }
1649 }
1650 
1651 
1652 ////////////////////////////////////////////////////////////////////////////////
1653 /// Set maximum.
1654 
1655 void TGraph2D::SetMaximum(Double_t maximum)
1656 {
1657  fMaximum = maximum;
1658  TH1 * h = GetHistogram();
1659  if (h) h->SetMaximum(maximum);
1660 }
1661 
1662 
1663 ////////////////////////////////////////////////////////////////////////////////
1664 /// Set minimum.
1665 
1666 void TGraph2D::SetMinimum(Double_t minimum)
1667 {
1668  fMinimum = minimum;
1669  TH1 * h = GetHistogram();
1670  if (h) h->SetMinimum(minimum);
1671 }
1672 
1673 
1674 ////////////////////////////////////////////////////////////////////////////////
1675 /// Changes the name of this 2D graph
1676 
1677 void TGraph2D::SetName(const char *name)
1678 {
1679  // 2D graphs are named objects in a THashList.
1680  // We must update the hashlist if we change the name
1681  if (fDirectory) fDirectory->Remove(this);
1682  fName = name;
1683  if (fDirectory) fDirectory->Append(this);
1684 }
1685 
1686 
1687 ////////////////////////////////////////////////////////////////////////////////
1688 /// Change the name and title of this 2D graph
1689 ///
1690 
1691 void TGraph2D::SetNameTitle(const char *name, const char *title)
1692 {
1693  // 2D graphs are named objects in a THashList.
1694  // We must update the hashlist if we change the name
1695  if (fDirectory) fDirectory->Remove(this);
1696  fName = name;
1697  SetTitle(title);
1698  if (fDirectory) fDirectory->Append(this);
1699 }
1700 
1701 
1702 ////////////////////////////////////////////////////////////////////////////////
1703 /// Sets the number of bins along X used to draw the function
1704 
1705 void TGraph2D::SetNpx(Int_t npx)
1706 {
1707  if (npx < 4) {
1708  Warning("SetNpx", "Number of points must be >4 && < 500, fNpx set to 4");
1709  fNpx = 4;
1710  } else if (npx > 500) {
1711  Warning("SetNpx", "Number of points must be >4 && < 500, fNpx set to 500");
1712  fNpx = 500;
1713  } else {
1714  fNpx = npx;
1715  }
1716  if (fHistogram) {
1717  delete fHistogram;
1718  fHistogram = nullptr;
1719  fDelaunay = nullptr;
1720  }
1721 }
1722 
1723 
1724 ////////////////////////////////////////////////////////////////////////////////
1725 /// Sets the number of bins along Y used to draw the function
1726 
1727 void TGraph2D::SetNpy(Int_t npy)
1728 {
1729  if (npy < 4) {
1730  Warning("SetNpy", "Number of points must be >4 && < 500, fNpy set to 4");
1731  fNpy = 4;
1732  } else if (npy > 500) {
1733  Warning("SetNpy", "Number of points must be >4 && < 500, fNpy set to 500");
1734  fNpy = 500;
1735  } else {
1736  fNpy = npy;
1737  }
1738  if (fHistogram) {
1739  delete fHistogram;
1740  fHistogram = nullptr;
1741  fDelaunay = nullptr;
1742  }
1743 }
1744 
1745 
1746 ////////////////////////////////////////////////////////////////////////////////
1747 /// Sets point number n.
1748 /// If n is greater than the current size, the arrays are automatically
1749 /// extended.
1750 
1751 void TGraph2D::SetPoint(Int_t n, Double_t x, Double_t y, Double_t z)
1752 {
1753  if (n < 0) return;
1754 
1755  if (!fX || !fY || !fZ || n >= fSize) {
1756  // re-allocate the object
1757  Int_t newN = TMath::Max(2 * fSize, n + 1);
1758  Double_t *savex = new Double_t [newN];
1759  Double_t *savey = new Double_t [newN];
1760  Double_t *savez = new Double_t [newN];
1761  if (fX && fSize) {
1762  memcpy(savex, fX, fSize * sizeof(Double_t));
1763  memset(&savex[fSize], 0, (newN - fSize)*sizeof(Double_t));
1764  delete [] fX;
1765  }
1766  if (fY && fSize) {
1767  memcpy(savey, fY, fSize * sizeof(Double_t));
1768  memset(&savey[fSize], 0, (newN - fSize)*sizeof(Double_t));
1769  delete [] fY;
1770  }
1771  if (fZ && fSize) {
1772  memcpy(savez, fZ, fSize * sizeof(Double_t));
1773  memset(&savez[fSize], 0, (newN - fSize)*sizeof(Double_t));
1774  delete [] fZ;
1775  }
1776  fX = savex;
1777  fY = savey;
1778  fZ = savez;
1779  fSize = newN;
1780  }
1781  fX[n] = x;
1782  fY[n] = y;
1783  fZ[n] = z;
1784  fNpoints = TMath::Max(fNpoints, n + 1);
1785 }
1786 
1787 
1788 ////////////////////////////////////////////////////////////////////////////////
1789 /// Sets the 2D graph title.
1790 ///
1791 /// This method allows to change the global title and the axis' titles of a 2D
1792 /// graph. If `g` is the 2D graph one can do:
1793 ///
1794 /// ~~~ {.cpp}
1795 /// g->SetTitle("Graph title; X axis title; Y axis title; Z axis title");
1796 /// ~~~
1797 
1798 void TGraph2D::SetTitle(const char* title)
1799 {
1800  fTitle = title;
1801  if (fHistogram) fHistogram->SetTitle(title);
1802 }
1803 
1804 
1805 ////////////////////////////////////////////////////////////////////////////////
1806 /// Stream a class object
1807 
1808 void TGraph2D::Streamer(TBuffer &b)
1809 {
1810  if (b.IsReading()) {
1811  UInt_t R__s, R__c;
1812  Version_t R__v = b.ReadVersion(&R__s, &R__c);
1813  b.ReadClassBuffer(TGraph2D::Class(), this, R__v, R__s, R__c);
1814 
1815  ResetBit(kMustCleanup);
1816  } else {
1817  b.WriteClassBuffer(TGraph2D::Class(), this);
1818  }
1819 }