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