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