Logo ROOT   6.18/05
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"
32
34
35
36/** \class TGraph2D
37 \ingroup Hist
38Graphics object made of three arrays X, Y and Z with the same number of points each.
39
40This 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
68Note that in any of these three cases, `SetPoint()` can be used to change a data
69point or add a new one. If the data point index (`i`) is greater than the
70current size of the internal arrays, they are automatically extended.
71
72Like TGraph the TGraph2D constructors do not have the TGraph2D title and name as parameters.
73A TGraph2D has the default title and name "Graph2D". To change the default title
74and name `SetTitle` and `SetName` should be called on the TGraph2D after its creation.
75
76Specific 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
90A TGraph2D can be also drawn with any options valid to draw a 2D histogram
91(like `COL`, `SURF`, `LEGO`, `CONT` etc..).
92
93When a TGraph2D is drawn with one of the 2D histogram drawing option,
94an intermediate 2D histogram is filled using the Delaunay triangles
95to interpolate the data set. The 2D histogram has equidistant bins along the X
96and Y directions. The number of bins along each direction can be change using
97`SetNpx()` and `SetNpy()`. Each bin is filled with the Z
98value found via a linear interpolation on the plane defined by the triangle above
99the (X,Y) coordinates of the bin center.
100
101The existing (X,Y,Z) points can be randomly scattered.
102The Delaunay triangles are build in the (X,Y) plane. These 2D triangles are then
103used to define flat planes in (X,Y,Z) over which the interpolation is done to fill
104the 2D histogram. The 3D triangles int takes build a 3D surface in
105the form of tessellating triangles at various angles. The triangles found can be
106drawn in 3D with one of the TGraph2D specific drawing options.
107
108The histogram generated by the Delaunay interpolation can be accessed using the
109`GetHistogram()` method.
110
111The axis settings (title, ranges etc ...) can be changed accessing the axis via
112the GetXaxis GetYaxis and GetZaxis methods. They access the histogram axis created
113at 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
125Example:
126
127Begin_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}
145End_Macro
146
1472D graphs can be fitted as shown by the following example:
148
149Begin_Macro(source)
150../../../tutorials/fit/graph2dfit.C
151End_Macro
152
153Example showing the PCOL option.
154
155Begin_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}
185End_Macro
186
187### Definition of Delaunay triangulation (After B. Delaunay)
188For a set S of points in the Euclidean plane, the unique triangulation DT(S)
189of S such that no point in S is inside the circumcircle of any triangle in
190DT(S). DT(S) is the dual of the Voronoi diagram of S.
191If n is the number of points in S, the Voronoi diagram of S is the partitioning
192of the plane containing S points into n convex polygons such that each polygon
193contains exactly one point and every point in a given polygon is closer to its
194central point than to any other. A Voronoi diagram is sometimes also known as
195a Dirichlet tessellation.
196
197\image html tgraph2d_delaunay.png
198
199[This applet](http://www.cs.cornell.edu/Info/People/chew/Delaunay.html)
200gives a nice practical view of Delaunay triangulation and Voronoi diagram.
201*/
202
203
204////////////////////////////////////////////////////////////////////////////////
205/// Graph2D default constructor
206
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;
228}
229
230
231////////////////////////////////////////////////////////////////////////////////
232/// Graph2D constructor with three vectors of ints as input.
233
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
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
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
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
337TGraph2D::TGraph2D(const char *name, const char *title,
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
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
379TGraph2D::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
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
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
526{
527 Clear();
528}
529
530
531////////////////////////////////////////////////////////////////////////////////
532/// Graph2D operator "="
533
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
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;
603
606 if (fDirectory) {
607 fDirectory->Append(this, kTRUE);
608 }
609 }
610}
611
612
613////////////////////////////////////////////////////////////////////////////////
614/// Browse
615
617{
618 Draw("p0");
619 gPad->Update();
620}
621
622
623////////////////////////////////////////////////////////////////////////////////
624/// Free all memory allocated by this object.
625
626void 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) {
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
661{
662 Bool_t addStatus = TH1::AddDirectoryStatus();
663 if (addStatus) {
664 SetDirectory(dir);
665 if (dir) {
667 }
668 }
669}
670
671
672////////////////////////////////////////////////////////////////////////////////
673/// Computes distance from point px,py to a graph
674
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
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
728{
729 if (fHistogram) fHistogram->ExecuteEvent(event, px, py);
730}
731
732
733////////////////////////////////////////////////////////////////////////////////
734/// search object named name in the list of functions
735
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
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
760TFitResultPtr 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/// option = "W" Set all weights to 1; ignore error bars
789/// = "U" Use a User specified fitting algorithm (via SetFCN)
790/// = "Q" Quiet mode (minimum printing)
791/// = "V" Verbose mode (default is between Q and V)
792/// = "R" Use the Range specified in the function range
793/// = "N" Do not store the graphics function, do not draw
794/// = "0" Do not plot the result of the fit. By default the fitted function
795/// is drawn unless the option "N" above is specified.
796/// = "+" Add this new fitted function to the list of fitted functions
797/// (by default, any previous function is deleted)
798/// = "C" In case of linear fitting, not calculate the chisquare
799/// (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
802/// coefficients (robust (resistant) regression), using
803/// the default fraction of good points
804/// "ROB=0.x" - compute the LTS regression coefficients, using
805/// 0.x as a fraction of good points
806/// = "S" The result of the fit is returned in the TFitResultPtr
807/// (see below Access to the Fit Result)
808///
809/// In order to use the Range option, one must first create a function
810/// with the expression to be fitted. For example, if your graph2d
811/// has a defined range between -4 and 4 and you want to fit a gaussian
812/// only in the interval 1 to 3, you can do:
813/// TF2 *f2 = new TF2("f2","gaus",1,3);
814/// graph2d->Fit("f2","R");
815///
816///
817/// Setting initial conditions
818/// ==========================
819/// Parameters must be initialized before invoking the Fit function.
820/// The setting of the parameter initial values is automatic for the
821/// predefined functions : poln, expo, gaus. One can however disable
822/// this automatic computation by specifying the option "B".
823/// You can specify boundary limits for some or all parameters via
824/// f2->SetParLimits(p_number, parmin, parmax);
825/// if parmin>=parmax, the parameter is fixed
826/// Note that you are not forced to fix the limits for all parameters.
827/// For example, if you fit a function with 6 parameters, you can do:
828/// func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
829/// func->SetParLimits(4,-10,-4);
830/// func->SetParLimits(5, 1,1);
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/// TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
852/// where MyFittingFunction is of type:
853/// extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
854///
855/// Associated functions
856/// ====================
857/// One or more object (typically a TF2*) can be added to the list
858/// of functions (fFunctions) associated to each graph.
859/// When TGraph::Fit is invoked, the fitted function is added to this list.
860/// Given a graph gr, one can retrieve an associated function
861/// with: TF2 *myfunc = gr->GetFunction("myfunc");
862///
863/// Access to the fit results
864/// =========================
865/// The function returns a TFitResultPtr which can hold a pointer to a TFitResult object.
866/// By default the TFitResultPtr contains only the status of the fit and it converts automatically to an
867/// integer. If the option "S" is instead used, TFitResultPtr contains the TFitResult and behaves as a smart
868/// pointer to it. For example one can do:
869/// TFitResultPtr r = graph->Fit("myFunc","S");
870/// TMatrixDSym cov = r->GetCovarianceMatrix(); // to access the covariance matrix
871/// Double_t par0 = r->Value(0); // retrieve the value for the parameter 0
872/// Double_t err0 = r->Error(0); // retrieve the error for the parameter 0
873/// r->Print("V"); // print full information of fit including covariance matrix
874/// r->Write(); // store the result in a file
875///
876/// The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also
877/// from the fitted function.
878/// If the graph is made persistent, the list of
879/// associated functions is also persistent. Given a pointer (see above)
880/// to an associated function myfunc, one can retrieve the function/fit
881/// parameters with calls such as:
882/// Double_t chi2 = myfunc->GetChisquare();
883/// Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
884/// Double_t err0 = myfunc->GetParError(0); //error on first parameter
885///
886/// Fit Statistics
887/// ==============
888/// You can change the statistics box to display the fit parameters with
889/// the TStyle::SetOptFit(mode) method. This mode has four digits.
890/// mode = pcev (default = 0111)
891/// v = 1; print name/values of parameters
892/// e = 1; print errors (if e=1, v must be 1)
893/// c = 1; print Chisquare/Number of degrees of freedom
894/// p = 1; print Probability
895///
896/// For example: gStyle->SetOptFit(1011);
897/// prints the fit probability, parameter names/values, and errors.
898/// You can change the position of the statistics box with these lines
899/// (where g is a pointer to the TGraph):
900///
901/// Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
902/// Root > st->SetX1NDC(newx1); //new x start position
903/// Root > st->SetX2NDC(newx2); //new x end position
904
906{
907 // internal graph2D fitting methods
908 Foption_t fitOption;
909 Option_t *goption = "";
911
912 // create range and minimizer options with default values
913 ROOT::Fit::DataRange range(2);
915 return ROOT::Fit::FitObject(this, f2 , fitOption , minOption, goption, range);
916}
917
918
919////////////////////////////////////////////////////////////////////////////////
920/// Display a GUI panel with all graph fit options.
921///
922/// See class TFitEditor for example
923
925{
926 if (!gPad)
927 gROOT->MakeDefCanvas();
928
929 if (!gPad) {
930 Error("FitPanel", "Unable to create a default canvas");
931 return;
932 }
933
934 // use plugin manager to create instance of TFitEditor
935 TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
936 if (handler && handler->LoadPlugin() != -1) {
937 if (handler->ExecPlugin(2, gPad, this) == 0)
938 Error("FitPanel", "Unable to crate the FitPanel");
939 } else
940 Error("FitPanel", "Unable to find the FitPanel plug-in");
941
942}
943
944
945////////////////////////////////////////////////////////////////////////////////
946/// Get x axis of the graph.
947
949{
950 TH1 *h = ((TGraph2D*)this)->GetHistogram("empty");
951 if (!h) return 0;
952 return h->GetXaxis();
953}
954
955
956////////////////////////////////////////////////////////////////////////////////
957/// Get y axis of the graph.
958
960{
961 TH1 *h = ((TGraph2D*)this)->GetHistogram("empty");
962 if (!h) return 0;
963 return h->GetYaxis();
964}
965
966
967////////////////////////////////////////////////////////////////////////////////
968/// Get z axis of the graph.
969
971{
972 TH1 *h = ((TGraph2D*)this)->GetHistogram("empty");
973 if (!h) return 0;
974 return h->GetZaxis();
975}
976
977
978////////////////////////////////////////////////////////////////////////////////
979/// Returns the X and Y graphs building a contour. A contour level may
980/// consist in several parts not connected to each other. This function
981/// returns them in a graphs' list.
982
984{
985 if (fNpoints <= 0) {
986 Error("GetContourList", "Empty TGraph2D");
987 return 0;
988 }
989
990 if (!fHistogram) GetHistogram("empty");
991
993
994 return fPainter->GetContourList(contour);
995}
996
997
998////////////////////////////////////////////////////////////////////////////////
999/// This function is called by Graph2DFitChisquare.
1000/// It always returns a negative value. Real implementation in TGraph2DErrors
1001
1003{
1004 return -1;
1005}
1006
1007
1008////////////////////////////////////////////////////////////////////////////////
1009/// This function is called by Graph2DFitChisquare.
1010/// It always returns a negative value. Real implementation in TGraph2DErrors
1011
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
1023{
1024 return -1;
1025}
1026
1027
1028////////////////////////////////////////////////////////////////////////////////
1029/// Add a TGraphDelaunay in the list of the fHistogram's functions
1030
1032{
1033
1035
1036 if (oldInterp) {
1037 TGraphDelaunay *dt = new TGraphDelaunay(this);
1038 dt->SetMaxIter(fMaxIter);
1040 fDelaunay = dt;
1042 if (!hl->FindObject("TGraphDelaunay")) hl->Add(fDelaunay);
1043 } else {
1044 TGraphDelaunay2D *dt = new TGraphDelaunay2D(this);
1046 fDelaunay = dt;
1048 if (!hl->FindObject("TGraphDelaunay2D")) hl->Add(fDelaunay);
1049 }
1050}
1051
1052////////////////////////////////////////////////////////////////////////////////
1053/// By default returns a pointer to the Delaunay histogram. If fHistogram
1054/// doesn't exist, books the 2D histogram fHistogram with a margin around
1055/// the hull. Calls TGraphDelaunay::Interpolate at each bin centre to build up
1056/// an interpolated 2D histogram.
1057/// If the "empty" option is selected, returns an empty histogram booked with
1058/// the limits of fX, fY and fZ. This option is used when the data set is
1059/// drawn with markers only. In that particular case there is no need to
1060/// find the Delaunay triangles.
1061/// By default use the new interpolation routine based on Triangles
1062/// If the option "old" the old interpolation is used
1063
1065{
1066 // for an empty graph create histogram in [0,1][0,1]
1067 if (fNpoints <= 0) {
1068 if (!fHistogram) {
1071 fHistogram = new TH2D(GetName(), GetTitle(), fNpx , 0., 1., fNpy, 0., 1.);
1072 TH1::AddDirectory(add);
1074 }
1075 return fHistogram;
1076 }
1077
1078 TString opt = option;
1079 opt.ToLower();
1080 Bool_t empty = opt.Contains("empty");
1081 Bool_t oldInterp = opt.Contains("old");
1082
1083 if (fHistogram) {
1084 if (!empty && fHistogram->GetEntries() == 0) {
1085 if (!fUserHisto) {
1086 delete fHistogram;
1087 fHistogram = nullptr;
1088 fDelaunay = nullptr;
1089 }
1090 } else if (fHistogram->GetEntries() == 0)
1091 {; }
1092 // check case if interpolation type has changed
1093 else if ( (TestBit(kOldInterpolation) && !oldInterp) || ( !TestBit(kOldInterpolation) && oldInterp ) ) {
1094 delete fHistogram;
1095 fHistogram = nullptr;
1096 fDelaunay = nullptr;
1097 }
1098 // normal case return existing histogram
1099 else {
1100 return fHistogram;
1101 }
1102 }
1103
1104 Double_t hxmax, hymax, hxmin, hymin;
1105
1106 // Book fHistogram if needed. It is not added in the current directory
1107 if (!fUserHisto) {
1114 hxmin = xmin - fMargin * (xmax - xmin);
1115 hymin = ymin - fMargin * (ymax - ymin);
1116 hxmax = xmax + fMargin * (xmax - xmin);
1117 hymax = ymax + fMargin * (ymax - ymin);
1118 if (TMath::Abs(hxmax - hxmin) < 0.0001) {
1119 if (TMath::Abs(hxmin) < 0.0001) {
1120 hxmin = -0.01;
1121 hxmax = 0.01;
1122 } else {
1123 hxmin = hxmin-TMath::Abs(hxmin)*0.01;
1124 hxmax = hxmax+TMath::Abs(hxmax)*0.01;
1125 }
1126 }
1127 if (TMath::Abs(hymax - hymin) < 0.0001) {
1128 if (TMath::Abs(hymin) < 0.0001) {
1129 hymin = -0.01;
1130 hymax = 0.01;
1131 } else {
1132 hymin = hymin-TMath::Abs(hymin)*0.01;
1133 hymax = hymax+TMath::Abs(hymax)*0.01;
1134 }
1135 }
1136 if (fHistogram) {
1137 fHistogram->GetXaxis()->SetLimits(hxmin, hxmax);
1138 fHistogram->GetYaxis()->SetLimits(hymin, hymax);
1139 } else {
1140 fHistogram = new TH2D(GetName(), GetTitle(),
1141 fNpx , hxmin, hxmax,
1142 fNpy, hymin, hymax);
1143 CreateInterpolator(oldInterp);
1144 }
1145 TH1::AddDirectory(add);
1147 } else {
1148 hxmin = fHistogram->GetXaxis()->GetXmin();
1149 hymin = fHistogram->GetYaxis()->GetXmin();
1150 hxmax = fHistogram->GetXaxis()->GetXmax();
1151 hymax = fHistogram->GetYaxis()->GetXmax();
1152 }
1153
1154 // Option "empty" is selected. An empty histogram is returned.
1155 if (empty) {
1156 Double_t hzmax, hzmin;
1157 if (fMinimum != -1111) {
1158 hzmin = fMinimum;
1159 } else {
1160 hzmin = GetZmin();
1161 }
1162 if (fMaximum != -1111) {
1163 hzmax = fMaximum;
1164 } else {
1165 hzmax = GetZmax();
1166 }
1167 if (hzmin == hzmax) {
1168 Double_t hz = hzmin;
1169 if (hz==0) hz = 1.;
1170 hzmin = hz - 0.01 * hz;
1171 hzmax = hz + 0.01 * hz;
1172 }
1173 fHistogram->SetMinimum(hzmin);
1174 fHistogram->SetMaximum(hzmax);
1175 return fHistogram;
1176 }
1177
1178 Double_t dx = (hxmax - hxmin) / fNpx;
1179 Double_t dy = (hymax - hymin) / fNpy;
1180
1181 Double_t x, y, z;
1182
1183 for (Int_t ix = 1; ix <= fNpx; ix++) {
1184 x = hxmin + (ix - 0.5) * dx;
1185 for (Int_t iy = 1; iy <= fNpy; iy++) {
1186 y = hymin + (iy - 0.5) * dy;
1187 // do interpolation
1188 if (oldInterp)
1189 z = ((TGraphDelaunay*)fDelaunay)->ComputeZ(x, y);
1190 else
1191 z = ((TGraphDelaunay2D*)fDelaunay)->ComputeZ(x, y);
1192
1193 fHistogram->Fill(x, y, z);
1194 }
1195 }
1196
1197
1198 if (fMinimum != -1111) fHistogram->SetMinimum(fMinimum);
1199 if (fMaximum != -1111) fHistogram->SetMaximum(fMaximum);
1200
1201 return fHistogram;
1202}
1203
1204
1205////////////////////////////////////////////////////////////////////////////////
1206/// Returns the X maximum
1207
1209{
1210 Double_t v = fX[0];
1211 for (Int_t i = 1; i < fNpoints; i++) if (fX[i] > v) v = fX[i];
1212 return v;
1213}
1214
1215
1216////////////////////////////////////////////////////////////////////////////////
1217/// Returns the X minimum
1218
1220{
1221 Double_t v = fX[0];
1222 for (Int_t i = 1; i < fNpoints; i++) if (fX[i] < v) v = fX[i];
1223 return v;
1224}
1225
1226
1227////////////////////////////////////////////////////////////////////////////////
1228/// Returns the Y maximum
1229
1231{
1232 Double_t v = fY[0];
1233 for (Int_t i = 1; i < fNpoints; i++) if (fY[i] > v) v = fY[i];
1234 return v;
1235}
1236
1237
1238////////////////////////////////////////////////////////////////////////////////
1239/// Returns the Y minimum
1240
1242{
1243 Double_t v = fY[0];
1244 for (Int_t i = 1; i < fNpoints; i++) if (fY[i] < v) v = fY[i];
1245 return v;
1246}
1247
1248
1249////////////////////////////////////////////////////////////////////////////////
1250/// Returns the Z maximum
1251
1253{
1254 Double_t v = fZ[0];
1255 for (Int_t i = 1; i < fNpoints; i++) if (fZ[i] > v) v = fZ[i];
1256 return v;
1257}
1258
1259
1260////////////////////////////////////////////////////////////////////////////////
1261/// Returns the Z minimum
1262
1264{
1265 Double_t v = fZ[0];
1266 for (Int_t i = 1; i < fNpoints; i++) if (fZ[i] < v) v = fZ[i];
1267 return v;
1268}
1269
1270////////////////////////////////////////////////////////////////////////////////
1271/// Get x, y and z values for point number i.
1272/// The function returns -1 in case of an invalid request or the point number otherwise
1273
1275{
1276 if (i < 0 || i >= fNpoints) return -1;
1277 if (!fX || !fY || !fZ) return -1;
1278 x = fX[i];
1279 y = fY[i];
1280 z = fZ[i];
1281 return i;
1282}
1283
1284////////////////////////////////////////////////////////////////////////////////
1285/// Finds the z value at the position (x,y) thanks to
1286/// the Delaunay interpolation.
1287
1289{
1290 if (fNpoints <= 0) {
1291 Error("Interpolate", "Empty TGraph2D");
1292 return 0;
1293 }
1294
1295 if (!fHistogram) GetHistogram("empty");
1296 if (!fDelaunay) {
1298 if (!TestBit(kOldInterpolation) ) {
1299 fDelaunay = hl->FindObject("TGraphDelaunay2D");
1300 if (!fDelaunay) fDelaunay = hl->FindObject("TGraphDelaunay");
1301 }
1302 else {
1303 // if using old implementation
1304 fDelaunay = hl->FindObject("TGraphDelaunay");
1305 if (!fDelaunay) fDelaunay = hl->FindObject("TGraphDelaunay2D");
1306 }
1307 }
1308
1309 if (!fDelaunay) return TMath::QuietNaN();
1310
1311 if (fDelaunay->IsA() == TGraphDelaunay2D::Class() )
1312 return ((TGraphDelaunay2D*)fDelaunay)->ComputeZ(x, y);
1313 else if (fDelaunay->IsA() == TGraphDelaunay::Class() )
1314 return ((TGraphDelaunay*)fDelaunay)->ComputeZ(x, y);
1315
1316 // cannot be here
1317 assert(false);
1318 return TMath::QuietNaN();
1319}
1320
1321
1322////////////////////////////////////////////////////////////////////////////////
1323/// Paints this 2D graph with its current attributes
1324
1326{
1327 if (fNpoints <= 0) {
1328 Error("Paint", "Empty TGraph2D");
1329 return;
1330 }
1331
1332 TString opt = option;
1333 opt.ToLower();
1334 if (opt.Contains("p") && !opt.Contains("tri")) {
1335 if (!opt.Contains("pol") &&
1336 !opt.Contains("sph") &&
1337 !opt.Contains("psr")) opt.Append("tri0");
1338 }
1339
1340 if (opt.Contains("line") && !opt.Contains("tri")) opt.Append("tri0");
1341
1342 if (opt.Contains("err") && !opt.Contains("tri")) opt.Append("tri0");
1343
1344 if (opt.Contains("tri0")) {
1345 GetHistogram("empty");
1346 } else if (opt.Contains("old")) {
1347 GetHistogram("old");
1348 } else {
1349 GetHistogram();
1350 }
1351
1360 fHistogram->Paint(opt.Data());
1361}
1362
1363
1364////////////////////////////////////////////////////////////////////////////////
1365/// Print 2D graph values.
1366
1368{
1369 for (Int_t i = 0; i < fNpoints; i++) {
1370 printf("x[%d]=%g, y[%d]=%g, z[%d]=%g\n", i, fX[i], i, fY[i], i, fZ[i]);
1371 }
1372}
1373
1374
1375////////////////////////////////////////////////////////////////////////////////
1376/// Projects a 2-d graph into 1 or 2-d histograms depending on the
1377/// option parameter
1378/// option may contain a combination of the characters x,y,z
1379/// option = "x" return the x projection into a TH1D histogram
1380/// option = "y" return the y projection into a TH1D histogram
1381/// option = "xy" return the x versus y projection into a TH2D histogram
1382/// option = "yx" return the y versus x projection into a TH2D histogram
1383
1385{
1386 if (fNpoints <= 0) {
1387 Error("Project", "Empty TGraph2D");
1388 return 0;
1389 }
1390
1391 TString opt = option;
1392 opt.ToLower();
1393
1394 Int_t pcase = 0;
1395 if (opt.Contains("x")) pcase = 1;
1396 if (opt.Contains("y")) pcase = 2;
1397 if (opt.Contains("xy")) pcase = 3;
1398 if (opt.Contains("yx")) pcase = 4;
1399
1400 // Create the projection histogram
1401 TH1D *h1 = 0;
1402 TH2D *h2 = 0;
1403 Int_t nch = strlen(GetName()) + opt.Length() + 2;
1404 char *name = new char[nch];
1405 snprintf(name, nch, "%s_%s", GetName(), option);
1406 nch = strlen(GetTitle()) + opt.Length() + 2;
1407 char *title = new char[nch];
1408 snprintf(title, nch, "%s_%s", GetTitle(), option);
1409
1410 Double_t hxmin = GetXmin();
1411 Double_t hxmax = GetXmax();
1412 Double_t hymin = GetYmin();
1413 Double_t hymax = GetYmax();
1414
1415 switch (pcase) {
1416 case 1:
1417 // "x"
1418 h1 = new TH1D(name, title, fNpx, hxmin, hxmax);
1419 break;
1420 case 2:
1421 // "y"
1422 h1 = new TH1D(name, title, fNpy, hymin, hymax);
1423 break;
1424 case 3:
1425 // "xy"
1426 h2 = new TH2D(name, title, fNpx, hxmin, hxmax, fNpy, hymin, hymax);
1427 break;
1428 case 4:
1429 // "yx"
1430 h2 = new TH2D(name, title, fNpy, hymin, hymax, fNpx, hxmin, hxmax);
1431 break;
1432 }
1433
1434 delete [] name;
1435 delete [] title;
1436 TH1 *h = h1;
1437 if (h2) h = h2;
1438 if (h == 0) return 0;
1439
1440 // Fill the projected histogram
1441 Double_t entries = 0;
1442 for (Int_t n = 0; n < fNpoints; n++) {
1443 switch (pcase) {
1444 case 1:
1445 // "x"
1446 h1->Fill(fX[n], fZ[n]);
1447 break;
1448 case 2:
1449 // "y"
1450 h1->Fill(fY[n], fZ[n]);
1451 break;
1452 case 3:
1453 // "xy"
1454 h2->Fill(fX[n], fY[n], fZ[n]);
1455 break;
1456 case 4:
1457 // "yx"
1458 h2->Fill(fY[n], fX[n], fZ[n]);
1459 break;
1460 }
1461 entries += fZ[n];
1462 }
1463 h->SetEntries(entries);
1464 return h;
1465}
1466
1467
1468////////////////////////////////////////////////////////////////////////////////
1469/// Deletes point number ipoint
1470
1472{
1473 if (ipoint < 0) return -1;
1474 if (ipoint >= fNpoints) return -1;
1475
1476 fNpoints--;
1477 Double_t *newX = new Double_t[fNpoints];
1478 Double_t *newY = new Double_t[fNpoints];
1479 Double_t *newZ = new Double_t[fNpoints];
1480 Int_t j = -1;
1481 for (Int_t i = 0; i < fNpoints + 1; i++) {
1482 if (i == ipoint) continue;
1483 j++;
1484 newX[j] = fX[i];
1485 newY[j] = fY[i];
1486 newZ[j] = fZ[i];
1487 }
1488 delete [] fX;
1489 delete [] fY;
1490 delete [] fZ;
1491 fX = newX;
1492 fY = newY;
1493 fZ = newZ;
1494 fSize = fNpoints;
1495 if (fHistogram) {
1496 delete fHistogram;
1497 fHistogram = nullptr;
1498 fDelaunay = nullptr;
1499 }
1500 return ipoint;
1501}
1502
1503
1504////////////////////////////////////////////////////////////////////////////////
1505/// Saves primitive as a C++ statement(s) on output stream out
1506
1507void TGraph2D::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1508{
1509 char quote = '"';
1510 out << " " << std::endl;
1511 if (gROOT->ClassSaved(TGraph2D::Class())) {
1512 out << " ";
1513 } else {
1514 out << " TGraph2D *";
1515 }
1516
1517 out << "graph2d = new TGraph2D(" << fNpoints << ");" << std::endl;
1518 out << " graph2d->SetName(" << quote << GetName() << quote << ");" << std::endl;
1519 out << " graph2d->SetTitle(" << quote << GetTitle() << quote << ");" << std::endl;
1520
1521 if (fDirectory == 0) {
1522 out << " graph2d->SetDirectory(0);" << std::endl;
1523 }
1524
1525 SaveFillAttributes(out, "graph2d", 0, 1001);
1526 SaveLineAttributes(out, "graph2d", 1, 1, 1);
1527 SaveMarkerAttributes(out, "graph2d", 1, 1, 1);
1528
1529 for (Int_t i = 0; i < fNpoints; i++) {
1530 out << " graph2d->SetPoint(" << i << "," << fX[i] << "," << fY[i] << "," << fZ[i] << ");" << std::endl;
1531 }
1532
1533 // save list of functions
1534 TIter next(fFunctions);
1535 TObject *obj;
1536 while ((obj = next())) {
1537 obj->SavePrimitive(out, "nodraw");
1538 out << " graph2d->GetListOfFunctions()->Add(" << obj->GetName() << ");" << std::endl;
1539 if (obj->InheritsFrom("TPaveStats")) {
1540 out << " ptstats->SetParent(graph2d->GetListOfFunctions());" << std::endl;
1541 } else if (obj->InheritsFrom("TF1")) {
1542 out << " " << obj->GetName() << "->SetParent(graph);\n";
1543 }
1544
1545 }
1546
1547 out << " graph2d->Draw(" << quote << option << quote << ");" << std::endl;
1548}
1549
1550
1551////////////////////////////////////////////////////////////////////////////////
1552/// Set number of points in the 2D graph.
1553/// Existing coordinates are preserved.
1554/// New coordinates above fNpoints are preset to 0.
1555
1557{
1558 if (n < 0) n = 0;
1559 if (n == fNpoints) return;
1560 if (n > fNpoints) SetPoint(n, 0, 0, 0);
1561 fNpoints = n;
1562}
1563
1564
1565////////////////////////////////////////////////////////////////////////////////
1566/// By default when an 2D graph is created, it is added to the list
1567/// of 2D graph objects in the current directory in memory.
1568/// This method removes reference to this 2D graph from current directory and add
1569/// reference to new directory dir. dir can be 0 in which case the
1570/// 2D graph does not belong to any directory.
1571
1573{
1574 if (fDirectory == dir) return;
1575 if (fDirectory) fDirectory->Remove(this);
1576 fDirectory = dir;
1577 if (fDirectory) fDirectory->Append(this);
1578}
1579
1580
1581////////////////////////////////////////////////////////////////////////////////
1582/// Sets the histogram to be filled.
1583/// If the 2D graph needs to be save in a TFile the following set should be
1584/// followed to read it back:
1585/// 1) Create TGraph2D
1586/// 2) Call g->SetHistogram(h), and do whatever you need to do
1587/// 3) Save g and h to the TFile, exit
1588/// 4) Open the TFile, retrieve g and h
1589/// 5) Call h->SetDirectory(0)
1590/// 6) Call g->SetHistogram(h) again
1591/// 7) Carry on as normal
1592
1594{
1595 fUserHisto = kTRUE;
1596 fHistogram = (TH2D*)h;
1597 fNpx = h->GetNbinsX();
1598 fNpy = h->GetNbinsY();
1600}
1601
1602
1603////////////////////////////////////////////////////////////////////////////////
1604/// Sets the extra space (in %) around interpolated area for the 2D histogram
1605
1607{
1608 if (m < 0 || m > 1) {
1609 Warning("SetMargin", "The margin must be >= 0 && <= 1, fMargin set to 0.1");
1610 fMargin = 0.1;
1611 } else {
1612 fMargin = m;
1613 }
1614 if (fHistogram) {
1615 delete fHistogram;
1616 fHistogram = nullptr;
1617 fDelaunay = nullptr;
1618 }
1619}
1620
1621
1622////////////////////////////////////////////////////////////////////////////////
1623/// Sets the histogram bin height for points lying outside the TGraphDelaunay
1624/// convex hull ie: the bins in the margin.
1625
1627{
1628 fZout = z;
1629 if (fHistogram) {
1630 delete fHistogram;
1631 fHistogram = nullptr;
1632 fDelaunay = nullptr;
1633 }
1634}
1635
1636
1637////////////////////////////////////////////////////////////////////////////////
1638/// Set maximum.
1639
1641{
1642 fMaximum = maximum;
1643 TH1 * h = GetHistogram();
1644 if (h) h->SetMaximum(maximum);
1645}
1646
1647
1648////////////////////////////////////////////////////////////////////////////////
1649/// Set minimum.
1650
1652{
1653 fMinimum = minimum;
1654 TH1 * h = GetHistogram();
1655 if (h) h->SetMinimum(minimum);
1656}
1657
1658
1659////////////////////////////////////////////////////////////////////////////////
1660/// Changes the name of this 2D graph
1661
1662void TGraph2D::SetName(const char *name)
1663{
1664 // 2D graphs are named objects in a THashList.
1665 // We must update the hashlist if we change the name
1666 if (fDirectory) fDirectory->Remove(this);
1667 fName = name;
1668 if (fDirectory) fDirectory->Append(this);
1669}
1670
1671
1672////////////////////////////////////////////////////////////////////////////////
1673/// Change the name and title of this 2D graph
1674///
1675
1676void TGraph2D::SetNameTitle(const char *name, const char *title)
1677{
1678 // 2D graphs are named objects in a THashList.
1679 // We must update the hashlist if we change the name
1680 if (fDirectory) fDirectory->Remove(this);
1681 fName = name;
1682 SetTitle(title);
1683 if (fDirectory) fDirectory->Append(this);
1684}
1685
1686
1687////////////////////////////////////////////////////////////////////////////////
1688/// Sets the number of bins along X used to draw the function
1689
1691{
1692 if (npx < 4) {
1693 Warning("SetNpx", "Number of points must be >4 && < 500, fNpx set to 4");
1694 fNpx = 4;
1695 } else if (npx > 500) {
1696 Warning("SetNpx", "Number of points must be >4 && < 500, fNpx set to 500");
1697 fNpx = 500;
1698 } else {
1699 fNpx = npx;
1700 }
1701 if (fHistogram) {
1702 delete fHistogram;
1703 fHistogram = nullptr;
1704 fDelaunay = nullptr;
1705 }
1706}
1707
1708
1709////////////////////////////////////////////////////////////////////////////////
1710/// Sets the number of bins along Y used to draw the function
1711
1713{
1714 if (npy < 4) {
1715 Warning("SetNpy", "Number of points must be >4 && < 500, fNpy set to 4");
1716 fNpy = 4;
1717 } else if (npy > 500) {
1718 Warning("SetNpy", "Number of points must be >4 && < 500, fNpy set to 500");
1719 fNpy = 500;
1720 } else {
1721 fNpy = npy;
1722 }
1723 if (fHistogram) {
1724 delete fHistogram;
1725 fHistogram = nullptr;
1726 fDelaunay = nullptr;
1727 }
1728}
1729
1730
1731////////////////////////////////////////////////////////////////////////////////
1732/// Sets point number n.
1733/// If n is greater than the current size, the arrays are automatically
1734/// extended.
1735
1737{
1738 if (n < 0) return;
1739
1740 if (!fX || !fY || !fZ || n >= fSize) {
1741 // re-allocate the object
1742 Int_t newN = TMath::Max(2 * fSize, n + 1);
1743 Double_t *savex = new Double_t [newN];
1744 Double_t *savey = new Double_t [newN];
1745 Double_t *savez = new Double_t [newN];
1746 if (fX && fSize) {
1747 memcpy(savex, fX, fSize * sizeof(Double_t));
1748 memset(&savex[fSize], 0, (newN - fSize)*sizeof(Double_t));
1749 delete [] fX;
1750 }
1751 if (fY && fSize) {
1752 memcpy(savey, fY, fSize * sizeof(Double_t));
1753 memset(&savey[fSize], 0, (newN - fSize)*sizeof(Double_t));
1754 delete [] fY;
1755 }
1756 if (fZ && fSize) {
1757 memcpy(savez, fZ, fSize * sizeof(Double_t));
1758 memset(&savez[fSize], 0, (newN - fSize)*sizeof(Double_t));
1759 delete [] fZ;
1760 }
1761 fX = savex;
1762 fY = savey;
1763 fZ = savez;
1764 fSize = newN;
1765 }
1766 fX[n] = x;
1767 fY[n] = y;
1768 fZ[n] = z;
1769 fNpoints = TMath::Max(fNpoints, n + 1);
1770}
1771
1772
1773////////////////////////////////////////////////////////////////////////////////
1774/// Sets the 2D graph title.
1775///
1776/// This method allows to change the global title and the axis' titles of a 2D
1777/// graph. If `g` is the 2D graph one can do:
1778///
1779/// ~~~ {.cpp}
1780/// g->SetTitle("Graph title; X axis title; Y axis title; Z axis title");
1781/// ~~~
1782
1783void TGraph2D::SetTitle(const char* title)
1784{
1785 fTitle = title;
1786 if (fHistogram) fHistogram->SetTitle(title);
1787}
1788
1789
1790////////////////////////////////////////////////////////////////////////////////
1791/// Stream a class object
1792
1793void TGraph2D::Streamer(TBuffer &b)
1794{
1795 if (b.IsReading()) {
1796 UInt_t R__s, R__c;
1797 Version_t R__v = b.ReadVersion(&R__s, &R__c);
1798 b.ReadClassBuffer(TGraph2D::Class(), this, R__v, R__s, R__c);
1799
1801 } else {
1802 b.WriteClassBuffer(TGraph2D::Class(), this);
1803 }
1804}
void Class()
Definition: Class.C:29
SVector< double, 2 > v
Definition: Dict.h:5
#define b(i)
Definition: RSha256.hxx:100
#define g(i)
Definition: RSha256.hxx:105
#define h(i)
Definition: RSha256.hxx:106
int Int_t
Definition: RtypesCore.h:41
short Version_t
Definition: RtypesCore.h:61
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
char * R__STRTOK_R(char *str, const char *delim, char **saveptr)
Definition: Rtypes.h:486
#define gDirectory
Definition: TDirectory.h:218
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float ymin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
float ymax
Definition: THbookFile.cxx:93
#define gROOT
Definition: TROOT.h:414
void Printf(const char *fmt,...)
R__EXTERN TSystem * gSystem
Definition: TSystem.h:560
#define gPad
Definition: TVirtualPad.h:286
#define snprintf
Definition: civetweb.c:1540
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition: DataRange.h:34
Fill Area Attributes class.
Definition: TAttFill.h:19
virtual Color_t GetFillColor() const
Return the fill area color.
Definition: TAttFill.h:30
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition: TAttFill.h:31
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
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:234
Line Attributes class.
Definition: TAttLine.h:18
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition: TAttLine.h:42
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
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:262
Marker Attributes class.
Definition: TAttMarker.h:19
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:245
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition: TAttMarker.h:32
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition: TAttMarker.h:31
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition: TAttMarker.h:33
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:41
Class to manage histogram axis.
Definition: TAxis.h:30
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
Double_t GetXmax() const
Definition: TAxis.h:134
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition: TAxis.h:154
Double_t GetXmin() const
Definition: TAxis.h:133
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
Describe directory structure in memory.
Definition: TDirectory.h:34
virtual void Append(TObject *obj, Bool_t replace=kFALSE)
Append object to this directory.
Definition: TDirectory.cxx:190
virtual TObject * Remove(TObject *)
Remove an object from the in-memory list.
A 2-Dim function with parameters.
Definition: TF2.h:29
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Definition: TFitResultPtr.h:31
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition: TGraph2D.h:40
Int_t fMaxIter
Maximum number of iterations to find Delaunay triangles.
Definition: TGraph2D.h:47
virtual Double_t GetYminE() const
Definition: TGraph2D.h:134
TGraph2D()
Graph2D default constructor.
Definition: TGraph2D.cxx:207
void Build(Int_t n)
Creates the 2D graph basic data structure.
Definition: TGraph2D.cxx:579
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:1288
Int_t fNpoints
Number of points in the data set.
Definition: TGraph2D.h:44
virtual void Print(Option_t *chopt="") const
Print 2D graph values.
Definition: TGraph2D.cxx:1367
virtual Double_t GetErrorZ(Int_t bin) const
This function is called by Graph2DFitChisquare.
Definition: TGraph2D.cxx:1022
virtual void FitPanel()
Display a GUI panel with all graph fit options.
Definition: TGraph2D.cxx:924
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:1626
Int_t fNpx
Number of bins along X in fHistogram.
Definition: TGraph2D.h:45
virtual Double_t GetYmaxE() const
Definition: TGraph2D.h:133
Double_t GetYmin() const
Returns the Y minimum.
Definition: TGraph2D.cxx:1241
Bool_t fUserHisto
Definition: TGraph2D.h:66
Double_t GetZmin() const
Returns the Z minimum.
Definition: TGraph2D.cxx:1263
Double_t * fZ
[fNpoints]
Definition: TGraph2D.h:51
virtual Double_t GetErrorX(Int_t bin) const
This function is called by Graph2DFitChisquare.
Definition: TGraph2D.cxx:1002
Double_t fMargin
Extra space (in %) around interpolated area for fHistogram.
Definition: TGraph2D.h:54
Double_t fMinimum
Minimum value for plotting along z.
Definition: TGraph2D.h:52
Double_t GetXmin() const
Returns the X minimum.
Definition: TGraph2D.cxx:1219
Int_t DistancetoPrimitive(Int_t px, Int_t py)
Computes distance from point px,py to a graph.
Definition: TGraph2D.cxx:675
virtual void Browse(TBrowser *)
Browse.
Definition: TGraph2D.cxx:616
TH2D * GetHistogram(Option_t *option="")
By default returns a pointer to the Delaunay histogram.
Definition: TGraph2D.cxx:1064
TVirtualHistPainter * fPainter
!Pointer to histogram painter
Definition: TGraph2D.h:60
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:760
virtual TObject * FindObject(const char *name) const
search object named name in the list of functions
Definition: TGraph2D.cxx:736
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Saves primitive as a C++ statement(s) on output stream out.
Definition: TGraph2D.cxx:1507
Double_t GetZmax() const
Returns the Z maximum.
Definition: TGraph2D.cxx:1252
Int_t RemovePoint(Int_t ipoint)
Deletes point number ipoint.
Definition: TGraph2D.cxx:1471
Double_t fMaximum
Maximum value for plotting along z.
Definition: TGraph2D.h:53
TAxis * GetZaxis() const
Get z axis of the graph.
Definition: TGraph2D.cxx:970
void SetMargin(Double_t m=0.1)
Sets the extra space (in %) around interpolated area for the 2D histogram.
Definition: TGraph2D.cxx:1606
void SetNpy(Int_t npx=40)
Sets the number of bins along Y used to draw the function.
Definition: TGraph2D.cxx:1712
Double_t fZout
fHistogram bin height for points lying outside the interpolated area
Definition: TGraph2D.h:55
TH2D * fHistogram
!2D histogram of z values linearly interpolated on the triangles
Definition: TGraph2D.h:57
virtual Double_t GetXminE() const
Definition: TGraph2D.h:132
TList * GetContourList(Double_t contour)
Returns the X and Y graphs building a contour.
Definition: TGraph2D.cxx:983
Double_t GetXmax() const
Returns the X maximum.
Definition: TGraph2D.cxx:1208
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y, Double_t &z) const
Get x, y and z values for point number i.
Definition: TGraph2D.cxx:1274
virtual void SetTitle(const char *title="")
Sets the 2D graph title.
Definition: TGraph2D.cxx:1783
TAxis * GetYaxis() const
Get y axis of the graph.
Definition: TGraph2D.cxx:959
virtual ~TGraph2D()
TGraph2D destructor.
Definition: TGraph2D.cxx:525
virtual Double_t GetErrorY(Int_t bin) const
This function is called by Graph2DFitChisquare.
Definition: TGraph2D.cxx:1012
TObject * fDelaunay
! Pointer to Delaunay interpolator object
Definition: TGraph2D.h:58
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:1384
virtual void SetHistogram(TH2 *h)
Sets the histogram to be filled.
Definition: TGraph2D.cxx:1593
virtual Double_t GetXmaxE() const
Definition: TGraph2D.h:131
virtual void Set(Int_t n)
Set number of points in the 2D graph.
Definition: TGraph2D.cxx:1556
void SetMinimum(Double_t minimum=-1111)
Set minimum.
Definition: TGraph2D.cxx:1651
TGraph2D & operator=(const TGraph2D &)
Graph2D operator "=".
Definition: TGraph2D.cxx:534
Double_t * fX
[fNpoints]
Definition: TGraph2D.h:49
Double_t * fY
[fNpoints] Data set to be plotted
Definition: TGraph2D.h:50
void SetMaximum(Double_t maximum=-1111)
Set maximum.
Definition: TGraph2D.cxx:1640
virtual void SetNameTitle(const char *name, const char *title)
Change the name and title of this 2D graph.
Definition: TGraph2D.cxx:1676
Double_t GetYmax() const
Returns the Y maximum.
Definition: TGraph2D.cxx:1230
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:1572
TDirectory * fDirectory
!Pointer to directory holding this 2D graph
Definition: TGraph2D.h:59
virtual void DirectoryAutoAdd(TDirectory *)
Perform the automatic addition of the graph to the given directory.
Definition: TGraph2D.cxx:660
void CreateInterpolator(Bool_t oldInterp)
Add a TGraphDelaunay in the list of the fHistogram's functions.
Definition: TGraph2D.cxx:1031
Int_t fNpy
Number of bins along Y in fHistogram.
Definition: TGraph2D.h:46
virtual void SetName(const char *name)
Changes the name of this 2D graph.
Definition: TGraph2D.cxx:1662
TList * fFunctions
Pointer to list of functions (fits and user)
Definition: TGraph2D.h:56
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Sets point number n.
Definition: TGraph2D.cxx:1736
void Paint(Option_t *option="")
Paints this 2D graph with its current attributes.
Definition: TGraph2D.cxx:1325
virtual void Draw(Option_t *option="P0")
Specific drawing options can be used to paint a TGraph2D:
Definition: TGraph2D.cxx:707
void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Executes action corresponding to one event.
Definition: TGraph2D.cxx:727
Int_t fSize
!Real size of fX, fY and fZ
Definition: TGraph2D.h:48
virtual void Clear(Option_t *option="")
Free all memory allocated by this object.
Definition: TGraph2D.cxx:626
@ kOldInterpolation
Definition: TGraph2D.h:69
TAxis * GetXaxis() const
Get x axis of the graph.
Definition: TGraph2D.cxx:948
void SetNpx(Int_t npx=40)
Sets the number of bins along X used to draw the function.
Definition: TGraph2D.cxx:1690
TGraphDelaunay2D generates a Delaunay triangulation of a TGraph2D.
void SetMarginBinsContent(Double_t z=0.)
TGraphDelaunay generates a Delaunay triangulation of a TGraph2D.
void SetMarginBinsContent(Double_t z=0.)
Sets the histogram bin height for points lying outside the convex hull ie: the bins in the margin.
void SetMaxIter(Int_t n=100000)
Defines the number of triangles tested for a Delaunay triangle (number of iterations) before abandoni...
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
The TH1 histogram class.
Definition: TH1.h:56
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6309
virtual Int_t GetNbinsY() const
Definition: TH1.h:293
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition: TH1.cxx:8476
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition: TH1.cxx:1225
@ kNoStats
don't draw stats box
Definition: TH1.h:160
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:316
TVirtualHistPainter * GetPainter(Option_t *option="")
Return pointer to painter.
Definition: TH1.cxx:4343
virtual Int_t GetNbinsX() const
Definition: TH1.h:292
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:394
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3258
TAxis * GetYaxis()
Definition: TH1.h:317
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:395
virtual Double_t GetEntries() const
Return the current number of entries.
Definition: TH1.cxx:4277
TList * GetListOfFunctions() const
Definition: TH1.h:239
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TH1.cxx:3154
virtual void Paint(Option_t *option="")
Control routine to paint any kind of histograms.
Definition: TH1.cxx:5809
static Bool_t AddDirectoryStatus()
Static function: cannot be inlined on Windows/NT.
Definition: TH1.cxx:705
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a line.
Definition: TH1.cxx:2719
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:289
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:292
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH2.h:82
A doubly linked list.
Definition: TList.h:44
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:575
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
TString fTitle
Definition: TNamed.h:33
TString fName
Definition: TNamed.h:32
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TObject.cxx:664
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
void MakeZombie()
Definition: TObject.h:49
void ResetBit(UInt_t f)
Definition: TObject.h:171
@ kCanDelete
if object in a list can be deleted
Definition: TObject.h:58
@ kInvalidObject
if object ctor succeeded but object should not be used
Definition: TObject.h:68
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition: TObject.h:60
Long_t ExecPlugin(int nargs, const T &... params)
Int_t LoadPlugin()
Load the plugin library for this handler.
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
Int_t Atoi() const
Return integer value of string.
Definition: TString.cxx:1921
Double_t Atof() const
Return floating-point value contained in string.
Definition: TString.cxx:1987
Bool_t IsFloat() const
Returns kTRUE if string contains a floating point or integer number.
Definition: TString.cxx:1791
const char * Data() const
Definition: TString.h:364
Bool_t IsDigit() const
Returns true if all characters in string are digits (0-9) or white spaces, i.e.
Definition: TString.cxx:1763
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
TString & Append(const char *cs)
Definition: TString.h:559
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:2311
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition: TSystem.cxx:1264
virtual TList * GetContourList(Double_t contour) const =0
TLine * line
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
TH1F * h1
Definition: legend1.C:5
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:964
void FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption)
Decode list of options into fitOption.
Definition: HFitImpl.cxx:681
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754
Definition: TMath.h:889
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * m
Definition: textangle.C:8