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