Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGraph.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Rene Brun, Olivier Couet 12/12/94
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
13#include "TROOT.h"
14#include "TBuffer.h"
15#include "TEnv.h"
16#include "TGraph.h"
17#include "TH1.h"
18#include "TF1.h"
19#include "TStyle.h"
20#include "TMath.h"
21#include "TVectorD.h"
22#include "Foption.h"
23#include "TRandom.h"
24#include "TSpline.h"
25#include "TVirtualFitter.h"
26#include "TVirtualPad.h"
28#include "TBrowser.h"
29#include "TSystem.h"
30#include "TPluginManager.h"
31#include "strtok.h"
32
33#include <cstdlib>
34#include <string>
35#include <cassert>
36#include <iostream>
37#include <fstream>
38#include <cstring>
39#include <numeric>
40#include <limits>
41#include <iomanip>
42
43#include "HFitInterface.h"
44#include "Fit/DataRange.h"
46
48
49
50////////////////////////////////////////////////////////////////////////////////
51
52/** \class TGraph
53 \ingroup Graphs
54A TGraph is an object made of two arrays X and Y with npoints each.
55The TGraph painting is performed thanks to the TGraphPainter
56class. All details about the various painting options are given in this class.
57
58#### Notes
59
60 - Unlike histogram or tree (or even TGraph2D), TGraph objects
61 are not automatically attached to the current TFile, in order to keep the
62 management and size of the TGraph as small as possible.
63 - The TGraph constructors do not have the TGraph title and name as parameters.
64 A TGraph has the default title and name "Graph". To change the default title
65 and name `SetTitle` and `SetName` should be called on the TGraph after its creation.
66 TGraph was a light weight object to start with, like TPolyline or TPolyMarker.
67 That’s why it did not have any title and name parameters in the constructors.
68
69#### Example
70
71The picture below gives an example:
72
73Begin_Macro(source)
74{
75 double x[100], y[100];
76 int n = 20;
77 for (int i=0;i<n;i++) {
78 x[i] = i*0.1;
79 y[i] = 10*sin(x[i]+0.2);
80 }
81 auto g = new TGraph(n,x,y);
82 g->SetTitle("Graph title;X title;Y title");
83 g->Draw("AC*");
84}
85End_Macro
86
87#### Default X-Points
88
89If one doesn't specify the points in the x-axis, they will get the default values 0, 1, 2, 3, (etc. depending
90on the length of the y-points):
91
92Begin_Macro(source)
93{
94 double y[6] = {3, 8, 1, 10, 5, 7};
95 auto g = new TGraph(6,y);
96 g->SetTitle("A Graph with default X points");
97 g->Draw();
98}
99End_Macro
100
101#### X-axis zooming
102
103The underlying x axis of a TGraph is based on a virtual fixed binwidth histogram, which means that one can not infinitely zoom on it, just down to one bin.
104The number of virtual bins is the maximum between 100 and the number of points in the TGraph. If you find a case where you would like to zoom deeper than allowed,
105you can either start canvas in web mode, or keep in classic mode but predefine the histogram before the first drawing of the graph:
106
107\code{.cpp}
108TGraph gr;
109for(auto i = 0; i < 100; ++i) gr.AddPoint(i,i);
110gr.AddPoint(10000, 0);
111// If one draws here, one can not zoom between x = 0 and x = 50, since first bin goes from x = 0 to 100.
112auto h1 = new TH1F("hist", "hist", 10000, -10., 11000.); // Define underlying hist with 10000 instead of 101 bins
113h1->SetMaximum(120);
114h1->SetStats(0);
115h1->SetDirectory(nullptr);
116gr.SetHistogram(h1);
117gr.Draw("A*")
118\endcode
119
120*/
121
122////////////////////////////////////////////////////////////////////////////////
123/// Graph default constructor.
124
126{
127 fNpoints = -1; //will be reset to 0 in CtorAllocate
128 if (!CtorAllocate()) return;
129}
130
131////////////////////////////////////////////////////////////////////////////////
132/// Constructor with only the number of points set
133/// the arrays x and y will be set later
134
136 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
137{
138 fNpoints = n;
139 if (!CtorAllocate()) return;
140 FillZero(0, fNpoints);
141}
142
143////////////////////////////////////////////////////////////////////////////////
144/// Graph normal constructor with ints.
145
147 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
148{
149 if (!x || !y) {
150 fNpoints = 0;
151 } else {
152 fNpoints = n;
153 }
154 if (!CtorAllocate()) return;
155 for (Int_t i = 0; i < n; i++) {
156 fX[i] = (Double_t)x[i];
157 fY[i] = (Double_t)y[i];
158 }
159}
160
161////////////////////////////////////////////////////////////////////////////////
162/// Graph normal constructor with floats.
163
165 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
166{
167 if (!x || !y) {
168 fNpoints = 0;
169 } else {
170 fNpoints = n;
171 }
172 if (!CtorAllocate()) return;
173 for (Int_t i = 0; i < n; i++) {
174 fX[i] = x[i];
175 fY[i] = y[i];
176 }
177}
178
179////////////////////////////////////////////////////////////////////////////////
180/// Default X-Points constructor. The points along the x-axis get the default
181/// values `start`, `start+step`, `start+2*step`, `start+3*step`, etc ...
182
184 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
185{
186 if (!y) {
187 fNpoints = 0;
188 } else {
189 fNpoints = n;
190 }
191 if (!CtorAllocate()) return;
192 for (Int_t i = 0; i < n; i++) {
193 fX[i] = start+i*step;
194 fY[i] = y[i];
195 }
196}
197
198////////////////////////////////////////////////////////////////////////////////
199/// Graph normal constructor with doubles.
200
202 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
203{
204 if (!x || !y) {
205 fNpoints = 0;
206 } else {
207 fNpoints = n;
208 }
209 if (!CtorAllocate()) return;
210 n = fNpoints * sizeof(Double_t);
211 memcpy(fX, x, n);
212 memcpy(fY, y, n);
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Copy constructor for this graph
217
220{
224 else fFunctions = new TList;
225 if (gr.fHistogram) {
227 fHistogram->SetDirectory(nullptr);
228 } else {
229 fHistogram = nullptr;
230 }
233 if (!fMaxSize) {
234 fX = fY = nullptr;
235 return;
236 } else {
237 fX = new Double_t[fMaxSize];
238 fY = new Double_t[fMaxSize];
239 }
240
241 Int_t n = gr.GetN() * sizeof(Double_t);
242 memcpy(fX, gr.fX, n);
243 memcpy(fY, gr.fY, n);
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Equal operator for this graph
248
250{
251 if (this != &gr) {
253 TAttLine::operator=(gr);
254 TAttFill::operator=(gr);
255 TAttMarker::operator=(gr);
256
259
260 // delete list of functions and their contents before copying it
261 if (fFunctions) {
262 // delete previous lists of functions
263 if (!fFunctions->IsEmpty()) {
265 // use TList::Remove to take into account the case the same object is
266 // added multiple times in the list
267 TObject *obj;
268 while ((obj = fFunctions->First())) {
269 while (fFunctions->Remove(obj)) { }
270 delete obj;
271 }
272 }
273 delete fFunctions;
274 }
275
277 else fFunctions = new TList;
278
279 if (fHistogram) delete fHistogram;
280 if (gr.fHistogram) {
281 fHistogram = new TH1F(*(gr.fHistogram));
282 fHistogram->SetDirectory(nullptr);
283 } else {
284 fHistogram = nullptr;
285 }
286
289 if (fX) delete [] fX;
290 if (fY) delete [] fY;
291 if (!fMaxSize) {
292 fX = fY = nullptr;
293 return *this;
294 } else {
295 fX = new Double_t[fMaxSize];
296 fY = new Double_t[fMaxSize];
297 }
298
299 Int_t n = gr.GetN() * sizeof(Double_t);
300 if (n > 0) {
301 memcpy(fX, gr.fX, n);
302 memcpy(fY, gr.fY, n);
303 }
304 }
305 return *this;
306}
307
308////////////////////////////////////////////////////////////////////////////////
309/// Graph constructor with two vectors of floats in input
310/// A graph is build with the X coordinates taken from vx and Y coord from vy
311/// The number of points in the graph is the minimum of number of points
312/// in vx and vy.
313
315 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
316{
317 fNpoints = TMath::Min(vx.GetNrows(), vy.GetNrows());
318 if (!CtorAllocate()) return;
319 Int_t ivxlow = vx.GetLwb();
320 Int_t ivylow = vy.GetLwb();
321 for (Int_t i = 0; i < fNpoints; i++) {
322 fX[i] = vx(i + ivxlow);
323 fY[i] = vy(i + ivylow);
324 }
325}
326
327////////////////////////////////////////////////////////////////////////////////
328/// Graph constructor with two vectors of doubles in input
329/// A graph is build with the X coordinates taken from vx and Y coord from vy
330/// The number of points in the graph is the minimum of number of points
331/// in vx and vy.
332
334 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
335{
336 fNpoints = TMath::Min(vx.GetNrows(), vy.GetNrows());
337 if (!CtorAllocate()) return;
338 Int_t ivxlow = vx.GetLwb();
339 Int_t ivylow = vy.GetLwb();
340 for (Int_t i = 0; i < fNpoints; i++) {
341 fX[i] = vx(i + ivxlow);
342 fY[i] = vy(i + ivylow);
343 }
344}
345
346////////////////////////////////////////////////////////////////////////////////
347/// Graph constructor importing its parameters from the TH1 object passed as argument
348
350 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
351{
352 if (!h) {
353 Error("TGraph", "Pointer to histogram is null");
354 fNpoints = 0;
355 return;
356 }
357 if (h->GetDimension() != 1) {
358 Error("TGraph", "Histogram must be 1-D; h %s is %d-D", h->GetName(), h->GetDimension());
359 fNpoints = 0;
360 } else {
361 fNpoints = ((TH1*)h)->GetXaxis()->GetNbins();
362 }
363
364 if (!CtorAllocate()) return;
365
366 TAxis *xaxis = ((TH1*)h)->GetXaxis();
367 for (Int_t i = 0; i < fNpoints; i++) {
368 fX[i] = xaxis->GetBinCenter(i + 1);
369 fY[i] = h->GetBinContent(i + 1);
370 }
371 h->TAttLine::Copy(*this);
372 h->TAttFill::Copy(*this);
373 h->TAttMarker::Copy(*this);
374
375 std::string gname = "Graph_from_" + std::string(h->GetName());
376 SetName(gname.c_str());
377 SetTitle(h->GetTitle());
378}
379
380////////////////////////////////////////////////////////////////////////////////
381/// Graph constructor importing its parameters from the TF1 object passed as argument
382/// - if option =="" (default), a TGraph is created with points computed
383/// at the fNpx points of f.
384/// - if option =="d", a TGraph is created with points computed with the derivatives
385/// at the fNpx points of f.
386/// - if option =="i", a TGraph is created with points computed with the integral
387/// at the fNpx points of f.
388/// - if option =="I", a TGraph is created with points computed with the integral
389/// at the fNpx+1 points of f and the integral is normalized to 1.
390
392 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
393{
394 char coption = ' ';
395 if (!f) {
396 Error("TGraph", "Pointer to function is null");
397 fNpoints = 0;
398 } else {
399 fNpoints = f->GetNpx();
400 if (option) coption = *option;
401 if (coption == 'i' || coption == 'I') fNpoints++;
402 }
403 if (!CtorAllocate()) return;
404
405 Double_t xmin = f->GetXmin();
406 Double_t xmax = f->GetXmax();
407 Double_t dx = (xmax - xmin) / fNpoints;
408 Double_t integ = 0;
409 Int_t i;
410 for (i = 0; i < fNpoints; i++) {
411 if (coption == 'i' || coption == 'I') {
412 fX[i] = xmin + i * dx;
413 if (i == 0) fY[i] = 0;
414 else fY[i] = integ + ((TF1*)f)->Integral(fX[i] - dx, fX[i]);
415 integ = fY[i];
416 } else if (coption == 'd' || coption == 'D') {
417 fX[i] = xmin + (i + 0.5) * dx;
418 fY[i] = ((TF1*)f)->Derivative(fX[i]);
419 } else {
420 fX[i] = xmin + (i + 0.5) * dx;
421 fY[i] = ((TF1*)f)->Eval(fX[i]);
422 }
423 }
424 if (integ != 0 && coption == 'I') {
425 for (i = 1; i < fNpoints; i++) fY[i] /= integ;
426 }
427
428 f->TAttLine::Copy(*this);
429 f->TAttFill::Copy(*this);
430 f->TAttMarker::Copy(*this);
431
432 SetName(f->GetName());
433 SetTitle(f->GetTitle());
434}
435
436////////////////////////////////////////////////////////////////////////////////
437/// Graph constructor reading input from filename.
438///
439/// `filename` is assumed to contain at least two columns of numbers.
440/// The string format is by default `"%lg %lg"`.
441/// This is a standard c formatting for `scanf()`.
442/// For example, set format to `"%lg,%lg"` for a comma-separated file.
443///
444/// If columns of numbers should be skipped, a `"%*lg"` or `"%*s"` for each column
445/// can be added, e.g. `"%lg %*lg %lg"` would read x-values from the first and
446/// y-values from the third column.
447///
448/// For files separated by a specific delimiter different from ' ' and '\\t' (e.g.
449/// ';' in csv files) you can avoid using `%*s` to bypass this delimiter by explicitly
450/// specify the `option` argument,
451/// e.g. option=`" \\t,;"` for columns of figures separated by any of these characters
452/// (' ', '\\t', ',', ';')
453/// used once (e.g. `"1;1"`) or in a combined way (`" 1;,;; 1"`).
454/// Note in that case, the instantiation is about two times slower.
455
456TGraph::TGraph(const char *filename, const char *format, Option_t *option)
457 : TNamed("Graph", filename), TAttFill(0, 1000)
458{
459 Double_t x, y;
462
463 std::ifstream infile(fname.Data());
464 if (!infile.good()) {
465 MakeZombie();
466 Error("TGraph", "Cannot open file: %s, TGraph is Zombie", filename);
467 fNpoints = 0;
468 return;
469 } else {
470 fNpoints = 100; //initial number of points
471 }
472 if (!CtorAllocate()) return;
473 std::string line;
474 Int_t np = 0;
475
476 // No delimiters specified (standard constructor).
477 if (strcmp(option, "") == 0) {
478
479 while (std::getline(infile, line, '\n')) {
480 if (2 != sscanf(line.c_str(), format, &x, &y)) {
481 continue; //skip empty and ill-formed lines
482 }
483 SetPoint(np, x, y);
484 np++;
485 }
486 Set(np);
487
488 // A delimiter has been specified in "option"
489 } else {
490
491 // Checking format and creating its boolean counterpart
493 format_.ReplaceAll(" ", "") ;
494 format_.ReplaceAll("\t", "") ;
495 format_.ReplaceAll("lg", "") ;
496 format_.ReplaceAll("s", "") ;
497 format_.ReplaceAll("%*", "0") ;
498 format_.ReplaceAll("%", "1") ;
499 if (!format_.IsDigit()) {
500 Error("TGraph", "Incorrect input format! Allowed formats are {\"%%lg\",\"%%*lg\" or \"%%*s\"}");
501 return;
502 }
503 Int_t ntokens = format_.Length() ;
504 if (ntokens < 2) {
505 Error("TGraph", "Incorrect input format! Only %d tag(s) in format whereas 2 \"%%lg\" tags are expected!", ntokens);
506 return;
507 }
510 for (Int_t idx = 0; idx < ntokens; idx++) {
511 isTokenToBeSaved[idx] = TString::Format("%c", format_[idx]).Atoi() ; //atoi(&format_[idx]) does not work for some reason...
512 if (isTokenToBeSaved[idx] == 1) {
514 }
515 }
516 if (ntokens >= 2 && ntokensToBeSaved != 2) { //first condition not to repeat the previous error message
517 Error("TGraph", "Incorrect input format! There are %d \"%%lg\" tag(s) in format whereas 2 and only 2 are expected!", ntokensToBeSaved);
518 delete [] isTokenToBeSaved ;
519 return;
520 }
521
522 // Initializing loop variables
523 Bool_t isLineToBeSkipped = kFALSE ; //empty and ill-formed lines
524 char * token = nullptr ;
525 TString token_str = "" ;
526 Int_t token_idx = 0 ;
527 Double_t * value = new Double_t [2] ; //x,y buffers
528 Int_t value_idx = 0 ;
529
530 // Looping
531 char *rest;
532 while (std::getline(infile, line, '\n')) {
533 if (!line.empty()) {
534 if (line[line.size() - 1] == char(13)) { // removing DOS CR character
535 line.erase(line.end() - 1, line.end()) ;
536 }
537 //token = R__STRTOK_R(const_cast<char *>(line.c_str()), option, rest);
538 token = R__STRTOK_R(const_cast<char *>(line.c_str()), option, &rest);
539 while (token != nullptr && value_idx < 2) {
541 token_str = TString(token) ;
542 token_str.ReplaceAll("\t", "") ;
543 if (!token_str.IsFloat()) {
545 break ;
546 } else {
547 value[value_idx] = token_str.Atof() ;
548 value_idx++ ;
549 }
550 }
551 token = R__STRTOK_R(nullptr, option, &rest); // next token
552 token_idx++ ;
553 }
554 if (!isLineToBeSkipped && value_idx == 2) {
555 x = value[0] ;
556 y = value[1] ;
557 SetPoint(np, x, y) ;
558 np++ ;
559 }
560 }
562 token = nullptr ;
563 token_idx = 0 ;
564 value_idx = 0 ;
565 }
566 Set(np) ;
567
568 // Cleaning
569 delete [] isTokenToBeSaved ;
570 delete [] value ;
571 delete token ;
572 }
573 infile.close();
574 if (fNpoints == 0) {
575 Warning("TGraph", "No points were found in file %s with the specified input format %s", filename, format);
576 return;
577 }
578}
579
580////////////////////////////////////////////////////////////////////////////////
581/// Graph default destructor.
582
584{
585 delete [] fX;
586 delete [] fY;
587 if (fFunctions) {
589 //special logic to support the case where the same object is
590 //added multiple times in fFunctions.
591 //This case happens when the same object is added with different
592 //drawing modes
593 TObject *obj;
594 while ((obj = fFunctions->First())) {
595 while (fFunctions->Remove(obj)) { }
596 delete obj;
597 }
598 delete fFunctions;
599 fFunctions = nullptr; //to avoid accessing a deleted object in RecursiveRemove
600 }
601 delete fHistogram;
602}
603
604////////////////////////////////////////////////////////////////////////////////
605/// Allocate internal data structures for `newsize` points.
606
611
612////////////////////////////////////////////////////////////////////////////////
613/// Allocate arrays.
614
616{
617 if (arraySize < 0) {
618 arraySize = 0;
619 }
621 if (!arraySize) {
622 for (Int_t i = 0; i < Narrays; ++i)
623 newarrays[i] = nullptr;
624 } else {
625 for (Int_t i = 0; i < Narrays; ++i)
626 newarrays[i] = new Double_t[arraySize];
627 }
628 fMaxSize = arraySize;
629 return newarrays;
630}
631
632////////////////////////////////////////////////////////////////////////////////
633/// Performs the operation: `y = y + c1*f(x,y)`
634/// Errors are not recalculated.
635///
636/// \param f may be a 1-D function TF1 or 2-d function TF2
637/// \param c1 a scaling factor, 1 by default
638
640{
642
643 for (Int_t i = 0; i < fNpoints; i++) {
644 fY[i] += c1*f->Eval(fX[i], fY[i]);
645 }
646 if (gPad) gPad->Modified();
647}
648
649////////////////////////////////////////////////////////////////////////////////
650/// Apply function f to all the data points
651/// f may be a 1-D function TF1 or 2-d function TF2
652/// The Y values of the graph are replaced by the new values computed
653/// using the function
654
656{
658
659 for (Int_t i = 0; i < fNpoints; i++) {
660 fY[i] = f->Eval(fX[i], fY[i]);
661 }
662 if (gPad) gPad->Modified();
663}
664
665////////////////////////////////////////////////////////////////////////////////
666/// Browse
667
669{
670 TString opt = gEnv->GetValue("TGraph.BrowseOption", "");
671 if (opt.IsNull()) {
672 opt = b ? b->GetDrawOption() : "alp";
673 opt = (opt == "") ? "alp" : opt.Data();
674 }
675 Draw(opt.Data());
676 gPad->Update();
677}
678
679////////////////////////////////////////////////////////////////////////////////
680/// Return the chisquare of this graph with respect to f1.
681/// The chisquare is computed as the sum of the quantity below at each point:
682/// \f[
683/// \frac{(y-f1(x))^{2}}{ey^{2}+(\frac{1}{2}(exl+exh)f1'(x))^{2}}
684/// \f]
685/// where x and y are the graph point coordinates and f1'(x) is the derivative of function f1(x).
686/// This method to approximate the uncertainty in y because of the errors in x, is called
687/// "effective variance" method.
688/// In case of a pure TGraph, the denominator is 1.
689/// In case of a TGraphErrors or TGraphAsymmErrors the errors are taken
690/// into account.
691/// By default the range of the graph is used whatever function range.
692/// Use option "R" to use the function range
693
695{
696 if (!func) {
697 Error("Chisquare","Function pointer is Null - return -1");
698 return -1;
699 }
700
701 TString opt(option); opt.ToUpper();
702 bool useRange = opt.Contains("R");
703
704 return ROOT::Fit::Chisquare(*this, *func,useRange);
705}
706
707////////////////////////////////////////////////////////////////////////////////
708/// Return kTRUE if point number "left"'s argument (angle with respect to positive
709/// x-axis) is bigger than that of point number "right". Can be used by Sort.
710
712{
713 Double_t xl = 0, yl = 0, xr = 0, yr = 0;
714 gr->GetPoint(left, xl, yl);
715 gr->GetPoint(right, xr, yr);
716 return (TMath::ATan2(yl, xl) > TMath::ATan2(yr, xr));
717}
718
719////////////////////////////////////////////////////////////////////////////////
720/// Return kTRUE if fX[left] > fX[right]. Can be used by Sort.
721
723{
724 return gr->fX[left] > gr->fX[right];
725}
726
727////////////////////////////////////////////////////////////////////////////////
728/// Return kTRUE if fY[left] > fY[right]. Can be used by Sort.
729
731{
732 return gr->fY[left] > gr->fY[right];
733}
734
735////////////////////////////////////////////////////////////////////////////////
736/// Return kTRUE if point number "left"'s distance to origin is bigger than
737/// that of point number "right". Can be used by Sort.
738
740{
741 return gr->fX[left] * gr->fX[left] + gr->fY[left] * gr->fY[left]
742 > gr->fX[right] * gr->fX[right] + gr->fY[right] * gr->fY[right];
743}
744
745////////////////////////////////////////////////////////////////////////////////
746/// Compute the x/y range of the points in this graph
747
749{
750 if (fNpoints <= 0) {
751 xmin = xmax = ymin = ymax = 0;
752 return;
753 }
754 xmin = xmax = fX[0];
755 ymin = ymax = fY[0];
756
757 Double_t xminl = 0; // Positive minimum. Used in case of log scale along X axis.
758 Double_t yminl = 0; // Positive minimum. Used in case of log scale along Y axis.
759
760 for (Int_t i = 1; i < fNpoints; i++) {
761 if (fX[i] < xmin) xmin = fX[i];
762 if (fX[i] > xmax) xmax = fX[i];
763 if (fY[i] < ymin) ymin = fY[i];
764 if (fY[i] > ymax) ymax = fY[i];
765 if (ymin>0 && (yminl==0 || ymin<yminl)) yminl = ymin;
766 if (xmin>0 && (xminl==0 || xmin<xminl)) xminl = xmin;
767 }
768
769 if (gPad && gPad->GetLogy() && yminl>0) ymin = yminl;
770 if (gPad && gPad->GetLogx() && xminl>0) xmin = xminl;
771}
772
773////////////////////////////////////////////////////////////////////////////////
774/// Copy points from fX and fY to arrays[0] and arrays[1]
775/// or to fX and fY if arrays == 0 and ibegin != iend.
776/// If newarrays is non null, replace fX, fY with pointers from newarrays[0,1].
777/// Delete newarrays, old fX and fY
778
781{
783 if (newarrays) {
784 delete[] fX;
785 fX = newarrays[0];
786 delete[] fY;
787 fY = newarrays[1];
788 delete[] newarrays;
789 }
790}
791
792////////////////////////////////////////////////////////////////////////////////
793/// Copy points from fX and fY to arrays[0] and arrays[1]
794/// or to fX and fY if arrays == 0 and ibegin != iend.
795
798{
799 if (ibegin < 0 || iend <= ibegin || obegin < 0) { // Error;
800 return kFALSE;
801 }
802 if (!arrays && ibegin == obegin) { // No copying is needed
803 return kFALSE;
804 }
805 Int_t n = (iend - ibegin) * sizeof(Double_t);
806 if (arrays) {
807 memmove(&arrays[0][obegin], &fX[ibegin], n);
808 memmove(&arrays[1][obegin], &fY[ibegin], n);
809 } else {
810 memmove(&fX[obegin], &fX[ibegin], n);
811 memmove(&fY[obegin], &fY[ibegin], n);
812 }
813 return kTRUE;
814}
815
816////////////////////////////////////////////////////////////////////////////////
817/// In constructors set fNpoints than call this method.
818/// Return kFALSE if the graph will contain no points.
819///Note: This function should be called only from the constructor
820/// since it does not delete previously existing arrays
821
823{
824 fHistogram = nullptr;
825 fMaximum = -1111;
826 fMinimum = -1111;
828 fFunctions = new TList;
829 if (fNpoints <= 0) {
830 fNpoints = 0;
831 fMaxSize = 0;
832 fX = nullptr;
833 fY = nullptr;
834 return kFALSE;
835 } else {
837 fX = new Double_t[fMaxSize];
838 fY = new Double_t[fMaxSize];
839 }
840 return kTRUE;
841}
842
843////////////////////////////////////////////////////////////////////////////////
844/// Draw this graph with its current attributes.
845///
846/// The options to draw a graph are described in TGraphPainter class.
847
849{
850 TString opt = option;
851 opt.ToLower();
852
853 if (opt.Contains("same")) {
854 opt.ReplaceAll("same", "");
855 }
856
857 // in case of option *, set marker style to 3 (star) and replace
858 // * option by option P.
859 Ssiz_t pos;
860 if ((pos = opt.Index("*")) != kNPOS) {
862 opt.Replace(pos, 1, "p");
863 }
864
865 // If no option is specified, it is defined as "alp" in case there is
866 // no current pad or if the current pad has no axis defined and if there is
867 // no default option set using TGraph::SetOption. If fOption is set using
868 // TGraph::SetOption, it is used as default option.
869 if ((!option || !strlen(option))) {
870 Option_t *topt = (!fOption.IsNull()) ? fOption.Data() : "alp";
871 if (gPad) {
872 if (!gPad->GetListOfPrimitives()->FindObject("TFrame"))
873 opt = topt;
874 } else {
875 opt = topt;
876 }
877 }
878
879 if (gPad) {
880 if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
881 if (opt.Contains("a")) gPad->Clear();
882 }
883
884 AppendPad(opt);
885
886 gPad->IncrementPaletteColor(1, opt);
887
888}
889
890////////////////////////////////////////////////////////////////////////////////
891/// Compute distance from point px,py to a graph.
892///
893/// Compute the closest distance of approach from point px,py to this line.
894/// The distance is computed in pixels units.
895
897{
899 if (painter) return painter->DistancetoPrimitiveHelper(this, px, py);
900 else return 0;
901}
902
903////////////////////////////////////////////////////////////////////////////////
904/// Draw this graph with new attributes.
905
907{
908 TGraph *newgraph = new TGraph(n, x, y);
912 newgraph->SetBit(kCanDelete);
913 newgraph->AppendPad(option);
914}
915
916////////////////////////////////////////////////////////////////////////////////
917/// Draw this graph with new attributes.
918
920{
921 TGraph *newgraph = new TGraph(n, x, y);
925 newgraph->SetBit(kCanDelete);
926 newgraph->AppendPad(option);
927}
928
929////////////////////////////////////////////////////////////////////////////////
930/// Draw this graph with new attributes.
931
933{
934 const Double_t *xx = x;
935 const Double_t *yy = y;
936 if (!xx) xx = fX;
937 if (!yy) yy = fY;
938 TGraph *newgraph = new TGraph(n, xx, yy);
942 newgraph->SetBit(kCanDelete);
943 newgraph->AppendPad(option);
944}
945
946////////////////////////////////////////////////////////////////////////////////
947/// Display a panel with all graph drawing options.
948
950{
952 if (painter) painter->DrawPanelHelper(this);
953}
954
955////////////////////////////////////////////////////////////////////////////////
956/// Interpolate points in this graph at x using a TSpline.
957///
958/// - if spline==0 and option="" a linear interpolation between the two points
959/// close to x is computed. If x is outside the graph range, a linear
960/// extrapolation is computed.
961/// - if spline==0 and option="S" a TSpline3 object is created using this graph
962/// and the interpolated value from the spline is returned.
963/// the internally created spline is deleted on return.
964/// - if spline is specified, it is used to return the interpolated value.
965///
966/// If the points are sorted in X a binary search is used (significantly faster)
967/// One needs to set the bit TGraph::SetBit(TGraph::kIsSortedX) before calling
968/// TGraph::Eval to indicate that the graph is sorted in X.
969
971{
972
973 if (spline) {
974 //spline interpolation using the input spline
975 return spline->Eval(x);
976 }
977
978 if (fNpoints == 0) return 0;
979 if (fNpoints == 1) return fY[0];
980
981 if (option && *option) {
982 TString opt = option;
983 opt.ToLower();
984 // create a TSpline every time when using option "s" and no spline pointer is given
985 if (opt.Contains("s")) {
986
987 // points must be sorted before using a TSpline
988 std::vector<Double_t> xsort(fNpoints);
989 std::vector<Double_t> ysort(fNpoints);
990 std::vector<Int_t> indxsort(fNpoints);
991 TMath::Sort(fNpoints, fX, &indxsort[0], false);
992 for (Int_t i = 0; i < fNpoints; ++i) {
993 xsort[i] = fX[ indxsort[i] ];
994 ysort[i] = fY[ indxsort[i] ];
995 }
996
997 // spline interpolation creating a new spline
998 TSpline3 s("", &xsort[0], &ysort[0], fNpoints);
999 Double_t result = s.Eval(x);
1000 return result;
1001 }
1002 }
1003 //linear interpolation
1004 //In case x is < fX[0] or > fX[fNpoints-1] return the extrapolated point
1005
1006 //find points in graph around x assuming points are not sorted
1007 // (if point are sorted use a binary search)
1008 Int_t low = -1;
1009 Int_t up = -1;
1012 if (low == -1) {
1013 // use first two points for doing an extrapolation
1014 low = 0;
1015 }
1016 if (fX[low] == x) return fY[low];
1017 if (low == fNpoints-1) low--; // for extrapolating
1018 up = low+1;
1019 }
1020 else {
1021 // case TGraph is not sorted
1022
1023 // find neighbours simply looping all points
1024 // and find also the 2 adjacent points: (low2 < low < x < up < up2 )
1025 // needed in case x is outside the graph ascissa interval
1026 Int_t low2 = -1;
1027 Int_t up2 = -1;
1028
1029 for (Int_t i = 0; i < fNpoints; ++i) {
1030 if (fX[i] < x) {
1031 if (low == -1 || fX[i] > fX[low]) {
1032 low2 = low;
1033 low = i;
1034 } else if (low2 == -1) low2 = i;
1035 } else if (fX[i] > x) {
1036 if (up == -1 || fX[i] < fX[up]) {
1037 up2 = up;
1038 up = i;
1039 } else if (up2 == -1) up2 = i;
1040 } else // case x == fX[i]
1041 return fY[i]; // no interpolation needed
1042 }
1043
1044 // treat cases when x is outside graph min max abscissa
1045 if (up == -1) {
1046 up = low;
1047 low = low2;
1048 }
1049 if (low == -1) {
1050 low = up;
1051 up = up2;
1052 }
1053 }
1054 // do now the linear interpolation
1055 assert(low != -1 && up != -1);
1056
1057 if (fX[low] == fX[up]) return fY[low];
1058 Double_t yn = fY[up] + (x - fX[up]) * (fY[low] - fY[up]) / (fX[low] - fX[up]);
1059 return yn;
1060}
1061
1062////////////////////////////////////////////////////////////////////////////////
1063/// Execute action corresponding to one event.
1064///
1065/// This member function is called when a graph is clicked with the locator
1066///
1067/// If Left button clicked on one of the line end points, this point
1068/// follows the cursor until button is released.
1069///
1070/// if Middle button clicked, the line is moved parallel to itself
1071/// until the button is released.
1072
1074{
1076 if (painter) painter->ExecuteEventHelper(this, event, px, py);
1077}
1078
1079////////////////////////////////////////////////////////////////////////////////
1080/// If array sizes <= newsize, expand storage to 2*newsize.
1081
1083{
1085 CopyAndRelease(ps, 0, 0, 0);
1086}
1087
1088////////////////////////////////////////////////////////////////////////////////
1089/// If graph capacity is less than newsize points then make array sizes
1090/// equal to least multiple of step to contain newsize points.
1091
1093{
1094 if (newsize <= fMaxSize) {
1095 return;
1096 }
1097 Double_t **ps = Allocate(step * (newsize / step + (newsize % step ? 1 : 0)));
1098 CopyAndRelease(ps, 0, fNpoints, 0);
1099}
1100
1101////////////////////////////////////////////////////////////////////////////////
1102/// if size > fMaxSize allocate new arrays of 2*size points and copy iend first
1103/// points.
1104/// Return pointer to new arrays.
1105
1107{
1108 if (size <= fMaxSize)
1109 return nullptr;
1111 CopyPoints(newarrays, 0, iend, 0);
1112 return newarrays;
1113}
1114
1115////////////////////////////////////////////////////////////////////////////////
1116/// Set zero values for point arrays in the range [begin, end)
1117/// Should be redefined in descendant classes
1118
1120{
1121 memset(fX + begin, 0, (end - begin)*sizeof(Double_t));
1122 memset(fY + begin, 0, (end - begin)*sizeof(Double_t));
1123}
1124
1125////////////////////////////////////////////////////////////////////////////////
1126/// Search object named name in the list of functions
1127
1129{
1130 return fFunctions ? fFunctions->FindObject(name) : nullptr;
1131}
1132
1133////////////////////////////////////////////////////////////////////////////////
1134/// Search object obj in the list of functions
1135
1137{
1138 return fFunctions ? fFunctions->FindObject(obj) : nullptr;
1139}
1140
1141////////////////////////////////////////////////////////////////////////////////
1142/// Fit this graph with function f1.
1143///
1144/// \param[in] f1 pointer to the function object
1145/// \param[in] option string defining the fit options (see table below).
1146/// \param[in] goption specify a list of graphics options. See TGraph::Draw and TGraphPainter for a complete list of these possible options.
1147/// \param[in] rxmin lower fitting range
1148/// \param[in] rxmax upper fitting range
1149///
1150/// \anchor GFitOpt
1151/// ### Graph Fitting Options
1152/// The list of fit options is given in parameter option.
1153///
1154/// option | description
1155/// -------|------------
1156/// "S" | The full result of the fit is returned in the `TFitResultPtr`. This is needed to get the covariance matrix of the fit. See `TFitResult` and the base class `ROOT::Math::FitResult`.
1157/// "W" | Ignore all point errors when fitting a TGraphErrors or TGraphAsymmErrors
1158/// "F" | Uses the default minimizer (e.g. Minuit) when fitting a linear function (e.g. polN) instead of the linear fitter.
1159/// "U" | Uses a user specified objective function (e.g. user providedlikelihood function) defined using `TVirtualFitter::SetFCN`
1160/// "E" | Performs a better parameter errors estimation using the Minos technique for all fit parameters.
1161/// "M" | Uses the IMPROVE algorithm (available only in TMinuit). This algorithm attempts improve the found local minimum by searching for a better one.
1162/// "Q" | Quiet mode (minimum printing)
1163/// "V" | Verbose mode (default is between Q and V)
1164/// "+" | Adds this new fitted function to the list of fitted functions. By default, the previous function is deleted and only the last one is kept.
1165/// "N" | Does not store the graphics function, does not draw the histogram with the function after fitting.
1166/// "0" | Does not draw the histogram and the fitted function after fitting, but in contrast to option "N", it stores the fitted function in the histogram list of functions.
1167/// "R" | Fit using a fitting range specified in the function range with `TF1::SetRange`.
1168/// "B" | Use this option when you want to fix one or more parameters and the fitting function is a predefined one (e.g gaus, expo,..), otherwise in case of pre-defined functions, some default initial values and limits are set.
1169/// "C" | In case of linear fitting, do no calculate the chisquare (saves CPU time).
1170/// "G" | Uses the gradient implemented in `TF1::GradientPar` for the minimization. This allows to use Automatic Differentiation when it is supported by the provided TF1 function.
1171/// "EX0" | When fitting a TGraphErrors or TGraphAsymErrors do not consider errors in the X coordinates
1172/// "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
1173///
1174///
1175/// This function is used for fitting also the derived TGraph classes such as TGraphErrors or TGraphAsymmErrors.
1176/// See the note below on how the errors are used when fitting a TGraphErrors or TGraphAsymmErrors.
1177///
1178/// The fitting of the TGraph, i.e simple data points without any error associated, is performed using the
1179/// un-weighted least-square (chi-square) method.
1180///
1181///
1182///\anchor GFitErrors
1183/// ### TGraphErrors fit:
1184///
1185/// In case of a TGraphErrors or TGraphAsymmErrors object, when `x` errors are present, the error along x,
1186/// is projected along the y-direction by calculating the function at the points `x-ex_low` and
1187/// `x+ex_high`, where `ex_low` and `ex_high` are the corresponding lower and upper error in x.
1188/// The chi-square is then computed as the sum of the quantity below at each data point:
1189///
1190/// \f[
1191/// \frac{(y-f(x))^{2}}{ey^{2}+(\frac{1}{2}(exl+exh)f'(x))^{2}}
1192/// \f]
1193///
1194/// where `x` and `y` are the point coordinates, and `f'(x)` is the derivative of the
1195/// function `f(x)`.
1196///
1197/// In case of asymmetric errors, if the function lies below (above) the data point, `ey` is `ey_low` (`ey_high`).
1198///
1199/// The approach used to approximate the uncertainty in y because of the
1200/// errors in x is to make it equal the error in x times the slope of the line.
1201/// This approach is called "effective variance method" and
1202/// the implementation is provided in the function FitUtil::EvaluateChi2Effective
1203///
1204/// \anchor GFitLinear
1205/// ### Linear fitting:
1206/// When the fitting function is linear (contains the `++` sign) or the fitting
1207/// function is a polynomial, a linear fitter is initialised.
1208/// To create a linear function, use the following syntax: linear parts
1209/// separated by `++` sign.
1210/// Example: to fit the parameters of the function `p0*x + p1*sin(x)`, you can create a
1211/// TF1 object as
1212///
1213/// TF1 *f1 = new TF1("f1", "x++sin(x)", xmin, xmax);
1214///
1215/// For such a TF1 you don't have to set the initial conditions and the linear fitter is used.
1216/// Going via the linear fitter for functions, linear in parameters, gives a
1217/// considerable advantage in speed.
1218/// When using the linear fitting it is also possible to perform a robust fitting with the
1219/// Least Trimmed Square (LTS) regression algorithm, by using the fit option `ROB`.
1220/// See the tutorial `fitLinearRobust.C`.
1221///
1222/// ### Notes on TGraph/TGraphErrors Fitting:
1223///
1224/// 1. By using the "effective variance" method a simple linear regression
1225/// becomes a non-linear case, which takes several iterations
1226/// instead of 0 as in the linear case.
1227/// 2. The effective variance technique assumes that there is no correlation
1228/// between the x and y coordinate.
1229/// 3. The standard chi2 (least square) method without error in the coordinates (x) can
1230/// be forced by using option "EX0"
1231/// 4. The linear fitter doesn't take into account the errors in x. When fitting a
1232/// TGraphErrors with a linear functions the errors in x will not be considered.
1233/// If errors in x are important, use option "F" for linear function fitting.
1234/// 5. When fitting a TGraph (i.e. no errors associated with each point),
1235/// a correction is applied to the errors on the parameters with the following
1236/// formula:
1237/// `parameter_error *= sqrt(chisquare/(ndf-1))`
1238///
1239/// ### General Fitting documentation
1240///
1241/// See in TH1::Fit for the documentation of
1242/// - [Fit Result](\ref HFitRes)
1243/// - [Fit Status](\ref HFitStatus)
1244/// - [Fit Statistics Box](\ref HFitStatBox)
1245/// - [Fitting in a Range](\ref HFitRange)
1246/// - [Setting Initial Conditions](\ref HFitInitial)
1247
1257
1258////////////////////////////////////////////////////////////////////////////////
1259/// Fit this graph with function with name `fname`.
1260///
1261/// This is a different interface to TGraph fitting using TGraph::Fit(TF1 *f1,Option_t *, Option_t *, Axis_t, Axis_t)
1262/// See there for the details about fitting a TGraph.
1263///
1264/// The parameter `fname` is the name of an already predefined function created by TF1 or TF2
1265/// Predefined functions such as gaus, expo and poln are automatically
1266/// created by ROOT.
1267///
1268/// The parameter `fname` can also be a formula, accepted by the linear fitter (linear parts divided
1269/// by "++" sign), for example "x++sin(x)" for fitting "[0]*x+[1]*sin(x)"
1270
1272{
1273 const char *linear = fname ? strstr(fname, "++") : nullptr;
1274 if (linear) {
1275 TF1 f1(fname, fname, xmin, xmax);
1276 return Fit(&f1, option, "", xmin, xmax);
1277 }
1278 TF1 * f1 = (TF1*)gROOT->GetFunction(fname);
1279 if (!f1) {
1280 Printf("Unknown function: %s", fname);
1281 return -1;
1282 }
1283 return Fit(f1, option, "", xmin, xmax);
1284}
1285
1286////////////////////////////////////////////////////////////////////////////////
1287/// Display a GUI panel with all graph fit options.
1288///
1289/// See class TFitEditor for example
1290
1292{
1293 if (!gPad)
1294 gROOT->MakeDefCanvas();
1295
1296 if (!gPad) {
1297 Error("FitPanel", "Unable to create a default canvas");
1298 return;
1299 }
1300
1301 // use plugin manager to create instance of TFitEditor
1302 TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
1303 if (handler && handler->LoadPlugin() != -1) {
1304 if (handler->ExecPlugin(2, gPad, this) == 0)
1305 Error("FitPanel", "Unable to crate the FitPanel");
1306 } else
1307 Error("FitPanel", "Unable to find the FitPanel plug-in");
1308}
1309
1310////////////////////////////////////////////////////////////////////////////////
1311/// Return graph correlation factor
1312
1314{
1315 Double_t rms1 = GetRMS(1);
1316 if (rms1 == 0) return 0;
1317 Double_t rms2 = GetRMS(2);
1318 if (rms2 == 0) return 0;
1319 return GetCovariance() / rms1 / rms2;
1320}
1321
1322////////////////////////////////////////////////////////////////////////////////
1323/// Return covariance of vectors x,y
1324
1326{
1327 if (fNpoints <= 0) return 0;
1328 Double_t sum = fNpoints, sumx = 0, sumy = 0, sumxy = 0;
1329
1330 for (Int_t i = 0; i < fNpoints; i++) {
1331 sumx += fX[i];
1332 sumy += fY[i];
1333 sumxy += fX[i] * fY[i];
1334 }
1335 return sumxy / sum - sumx / sum * sumy / sum;
1336}
1337
1338////////////////////////////////////////////////////////////////////////////////
1339/// Return mean value of X (axis=1) or Y (axis=2)
1340
1342{
1343 if (axis < 1 || axis > 2) return 0;
1344 if (fNpoints <= 0) return 0;
1345 Double_t sumx = 0;
1346 for (Int_t i = 0; i < fNpoints; i++) {
1347 if (axis == 1) sumx += fX[i];
1348 else sumx += fY[i];
1349 }
1350 return sumx / fNpoints;
1351}
1352
1353////////////////////////////////////////////////////////////////////////////////
1354/// Return RMS of X (axis=1) or Y (axis=2)
1355
1357{
1358 if (axis < 1 || axis > 2) return 0;
1359 if (fNpoints <= 0) return 0;
1360 Double_t sumx = 0, sumx2 = 0;
1361 for (Int_t i = 0; i < fNpoints; i++) {
1362 if (axis == 1) {
1363 sumx += fX[i];
1364 sumx2 += fX[i] * fX[i];
1365 } else {
1366 sumx += fY[i];
1367 sumx2 += fY[i] * fY[i];
1368 }
1369 }
1370 Double_t x = sumx / fNpoints;
1372 return TMath::Sqrt(rms2);
1373}
1374
1375////////////////////////////////////////////////////////////////////////////////
1376/// It always returns a negative value. Real implementation in TGraphErrors
1377
1379{
1380 return -1;
1381}
1382
1383////////////////////////////////////////////////////////////////////////////////
1384/// It always returns a negative value. Real implementation in TGraphErrors
1385
1387{
1388 return -1;
1389}
1390
1391////////////////////////////////////////////////////////////////////////////////
1392/// It always returns a negative value. Real implementation in TGraphErrors
1393/// and TGraphAsymmErrors
1394
1396{
1397 return -1;
1398}
1399
1400////////////////////////////////////////////////////////////////////////////////
1401/// It always returns a negative value. Real implementation in TGraphErrors
1402/// and TGraphAsymmErrors
1403
1405{
1406 return -1;
1407}
1408
1409////////////////////////////////////////////////////////////////////////////////
1410/// It always returns a negative value. Real implementation in TGraphErrors
1411/// and TGraphAsymmErrors
1412
1414{
1415 return -1;
1416}
1417
1418////////////////////////////////////////////////////////////////////////////////
1419/// It always returns a negative value. Real implementation in TGraphErrors
1420/// and TGraphAsymmErrors
1421
1423{
1424 return -1;
1425}
1426
1427////////////////////////////////////////////////////////////////////////////////
1428/// Return pointer to function with name.
1429///
1430/// Functions such as TGraph::Fit store the fitted function in the list of
1431/// functions of this graph.
1432
1433TF1 *TGraph::GetFunction(const char *name) const
1434{
1435 return dynamic_cast<TF1*>(FindObject(name));
1436}
1437
1438////////////////////////////////////////////////////////////////////////////////
1439/// Returns a pointer to the histogram used to draw the axis
1440/// Takes into account the two following cases.
1441/// 1. option 'A' was specified in TGraph::Draw. Return fHistogram
1442/// 2. user had called TPad::DrawFrame. return pointer to hframe histogram
1443
1445{
1448
1449 ComputeRange(rwxmin, rwymin, rwxmax, rwymax); //this is redefined in TGraphErrors
1450
1451 // (if fHistogram exist) && (if the log scale is on) &&
1452 // (if the computed range minimum is > 0) && (if the fHistogram minimum is zero)
1453 // then it means fHistogram limits have been computed in linear scale
1454 // therefore they might be too strict and cut some points. In that case the
1455 // fHistogram limits should be recomputed ie: the existing fHistogram
1456 // should not be returned.
1457 TH1F *histogr = nullptr;
1458 if (fHistogram) {
1459 if (!TestBit(kResetHisto)) {
1460 if (gPad && gPad->GetLogx()) {
1461 if (rwxmin <= 0 || fHistogram->GetXaxis()->GetXmin() != 0) return fHistogram;
1462 } else if (gPad && gPad->GetLogy()) {
1464 } else {
1465 return fHistogram;
1466 }
1467 } else {
1468 const_cast <TGraph*>(this)->ResetBit(kResetHisto);
1469 }
1471 }
1472
1473 if (rwxmin == rwxmax) rwxmax += 1.;
1474 if (rwymin == rwymax) rwymax += 1.;
1475 dx = 0.1 * (rwxmax - rwxmin);
1476 dy = 0.1 * (rwymax - rwymin);
1477 uxmin = rwxmin - dx;
1478 uxmax = rwxmax + dx;
1479 minimum = rwymin - dy;
1480 maximum = rwymax + dy;
1481
1482 if (fMinimum != -1111) minimum = fMinimum;
1483 if (fMaximum != -1111) maximum = fMaximum;
1484
1485 // the graph is created with at least as many channels as there are points
1486 // to permit zooming on the full range
1487 if (uxmin < 0 && rwxmin >= 0) {
1488 if (gPad && gPad->GetLogx()) uxmin = 0.9 * rwxmin;
1489 else uxmin = 0;
1490 }
1491 if (uxmax > 0 && rwxmax <= 0) {
1492 if (gPad && gPad->GetLogx()) uxmax = 1.1 * rwxmax;
1493 else uxmax = 0;
1494 }
1495
1496 if (minimum < 0 && rwymin >= 0) minimum = 0.9 * rwymin;
1497
1498 if (minimum <= 0 && gPad && gPad->GetLogy()) minimum = 0.001 * maximum;
1499 if (uxmin <= 0 && gPad && gPad->GetLogx()) {
1500 if (uxmax > 1000) uxmin = 1;
1501 else uxmin = 0.001 * uxmax;
1502 }
1503
1504 rwxmin = uxmin;
1505 rwxmax = uxmax;
1506 Int_t npt = 100;
1507 if (fNpoints > npt) npt = fNpoints;
1508 const char *gname = GetName();
1509 if (!gname[0]) gname = "Graph";
1510 // do not add the histogram to gDirectory
1511 // use local TDirectory::TContext that will set temporarly gDirectory to a nullptr and
1512 // will avoid that histogram is added in the global directory
1513 {
1514 TDirectory::TContext ctx(nullptr);
1515 ((TGraph*)this)->fHistogram = new TH1F(gname, GetTitle(), npt, rwxmin, rwxmax);
1516 }
1517 if (!fHistogram) return nullptr;
1522 // Restore the axis attributes if needed
1523 if (histogr) {
1524 fHistogram->GetXaxis()->SetTitle(histogr->GetXaxis()->GetTitle());
1525 fHistogram->GetXaxis()->CenterTitle(histogr->GetXaxis()->GetCenterTitle());
1526 fHistogram->GetXaxis()->RotateTitle(histogr->GetXaxis()->GetRotateTitle());
1527 fHistogram->GetXaxis()->SetNoExponent(histogr->GetXaxis()->GetNoExponent());
1528 fHistogram->GetXaxis()->SetTimeDisplay(histogr->GetXaxis()->GetTimeDisplay());
1529 fHistogram->GetXaxis()->SetTimeFormat(histogr->GetXaxis()->GetTimeFormat());
1530 histogr->GetXaxis()->TAttAxis::Copy(*(fHistogram->GetXaxis()));
1531
1532 fHistogram->GetYaxis()->SetTitle(histogr->GetYaxis()->GetTitle());
1533 fHistogram->GetYaxis()->CenterTitle(histogr->GetYaxis()->GetCenterTitle());
1534 fHistogram->GetYaxis()->RotateTitle(histogr->GetYaxis()->GetRotateTitle());
1535 fHistogram->GetYaxis()->SetNoExponent(histogr->GetYaxis()->GetNoExponent());
1536 fHistogram->GetYaxis()->SetTimeDisplay(histogr->GetYaxis()->GetTimeDisplay());
1537 fHistogram->GetYaxis()->SetTimeFormat(histogr->GetYaxis()->GetTimeFormat());
1538 histogr->GetYaxis()->TAttAxis::Copy(*(fHistogram->GetYaxis()));
1539
1540 delete histogr;
1541 }
1542 return fHistogram;
1543}
1544
1545////////////////////////////////////////////////////////////////////////////////
1546/// Get x and y values for point number i.
1547/// The function returns -1 in case of an invalid request or the point number otherwise
1548
1550{
1551 if (i < 0 || i >= fNpoints || !fX || !fY) return -1;
1552 x = fX[i];
1553 y = fY[i];
1554 return i;
1555}
1556
1557////////////////////////////////////////////////////////////////////////////////
1558/// Get x value for point i.
1559
1561{
1562 if (i < 0 || i >= fNpoints || !fX)
1563 return -1.;
1564
1565 return fX[i];
1566}
1567
1568////////////////////////////////////////////////////////////////////////////////
1569/// Get y value for point i.
1570
1572{
1573 if (i < 0 || i >= fNpoints || !fY)
1574 return -1.;
1575
1576 return fY[i];
1577}
1578
1579////////////////////////////////////////////////////////////////////////////////
1580/// Get x axis of the graph.
1581
1583{
1584 auto h = GetHistogram();
1585 return h ? h->GetXaxis() : nullptr;
1586}
1587
1588////////////////////////////////////////////////////////////////////////////////
1589/// Get y axis of the graph.
1590
1592{
1593 auto h = GetHistogram();
1594 return h ? h->GetYaxis() : nullptr;
1595}
1596
1597////////////////////////////////////////////////////////////////////////////////
1598/// Implementation to get information on point of graph at cursor position
1599/// Adapted from class TH1
1600
1602{
1603 if (!gPad) {
1604 Error("GetObjectInfo", "Cannot be used without gPad");
1605 return nullptr;
1606 }
1607
1608 // localize point
1609 Int_t ipoint = -2;
1610 // start with a small window (in case the mouse is very close to one point)
1611 for (Int_t i = 0; i < fNpoints; i++) {
1612 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
1613 Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
1614
1615 if (dpx * dpx + dpy * dpy < 25) {
1616 ipoint = i;
1617 break;
1618 }
1619 }
1620
1621 Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
1622 Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py));
1623
1624 if (ipoint == -2)
1625 return Form("x=%g, y=%g", x, y);
1626
1629
1630 return Form("x=%g, y=%g, point=%d, xval=%g, yval=%g", x, y, ipoint, xval, yval);
1631}
1632
1633////////////////////////////////////////////////////////////////////////////////
1634/// Compute Initial values of parameters for a gaussian.
1635
1637{
1638 Double_t allcha, sumx, sumx2, x, val, rms, mean;
1639 Int_t bin;
1640 const Double_t sqrtpi = 2.506628;
1641
1642 // Compute mean value and RMS of the graph in the given range
1643 if (xmax <= xmin) {
1644 xmin = fX[0];
1645 xmax = fX[fNpoints-1];
1646 }
1647 Int_t np = 0;
1648 allcha = sumx = sumx2 = 0;
1649 for (bin = 0; bin < fNpoints; bin++) {
1650 x = fX[bin];
1651 if (x < xmin || x > xmax) continue;
1652 np++;
1653 val = fY[bin];
1654 sumx += val * x;
1655 sumx2 += val * x * x;
1656 allcha += val;
1657 }
1658 if (np == 0 || allcha == 0) return;
1659 mean = sumx / allcha;
1660 rms = TMath::Sqrt(sumx2 / allcha - mean * mean);
1662 if (rms == 0) rms = 1;
1664 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1665 f1->SetParameter(0, binwidx * allcha / (sqrtpi * rms));
1666 f1->SetParameter(1, mean);
1667 f1->SetParameter(2, rms);
1668 f1->SetParLimits(2, 0, 10 * rms);
1669}
1670
1671////////////////////////////////////////////////////////////////////////////////
1672/// Compute Initial values of parameters for an exponential.
1673
1675{
1677 Int_t ifail;
1678 if (xmax <= xmin) {
1679 xmin = fX[0];
1680 xmax = fX[fNpoints-1];
1681 }
1683
1685
1687 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1689 f1->SetParameter(1, slope);
1690}
1691
1692////////////////////////////////////////////////////////////////////////////////
1693/// Compute Initial values of parameters for a polynom.
1694
1696{
1697 Double_t fitpar[25];
1698
1700 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1701 Int_t npar = f1->GetNpar();
1702 if (xmax <= xmin) {
1703 xmin = fX[0];
1704 xmax = fX[fNpoints-1];
1705 }
1706
1708
1709 for (Int_t i = 0; i < npar; i++) f1->SetParameter(i, fitpar[i]);
1710}
1711
1712////////////////////////////////////////////////////////////////////////////////
1713/// Insert a new point at the mouse position
1714
1716{
1717 if (!gPad) {
1718 Error("InsertPoint", "Cannot be used without gPad, requires last mouse position");
1719 return -1;
1720 }
1721
1722 Int_t px = gPad->GetEventX();
1723 Int_t py = gPad->GetEventY();
1724
1725 //localize point where to insert
1726 Int_t ipoint = -2;
1727 Int_t i, d = 0;
1728 // start with a small window (in case the mouse is very close to one point)
1729 for (i = 0; i < fNpoints - 1; i++) {
1730 d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
1731 if (d < 5) {
1732 ipoint = i + 1;
1733 break;
1734 }
1735 }
1736 if (ipoint == -2) {
1737 //may be we are far from one point, try again with a larger window
1738 for (i = 0; i < fNpoints - 1; i++) {
1739 d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
1740 if (d < 10) {
1741 ipoint = i + 1;
1742 break;
1743 }
1744 }
1745 }
1746 if (ipoint == -2) {
1747 //distinguish between first and last point
1748 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[0]));
1749 Int_t dpy = py - gPad->YtoAbsPixel(gPad->XtoPad(fY[0]));
1750 if (dpx * dpx + dpy * dpy < 25) ipoint = 0;
1751 else ipoint = fNpoints;
1752 }
1753
1754
1755 InsertPointBefore(ipoint, gPad->AbsPixeltoX(px), gPad->AbsPixeltoY(py));
1756
1757 gPad->Modified();
1758 return ipoint;
1759}
1760
1761
1762////////////////////////////////////////////////////////////////////////////////
1763/// Insert a new point with coordinates (x,y) before the point number `ipoint`.
1764
1766{
1767 if (ipoint < 0) {
1768 Error("TGraph", "Inserted point index should be >= 0");
1769 return;
1770 }
1771
1772 if (ipoint > fNpoints) {
1773 Error("TGraph", "Inserted point index should be <= %d", fNpoints);
1774 return;
1775 }
1776
1777 if (ipoint == fNpoints) {
1778 SetPoint(ipoint, x, y);
1779 return;
1780 }
1781
1783 CopyAndRelease(ps, ipoint, fNpoints++, ipoint + 1);
1784
1785 // To avoid redefinitions in descendant classes
1786 FillZero(ipoint, ipoint + 1);
1787
1788 fX[ipoint] = x;
1789 fY[ipoint] = y;
1790}
1791
1792
1793////////////////////////////////////////////////////////////////////////////////
1794/// Integrate the TGraph data within a given (index) range.
1795/// Note that this function computes the area of the polygon enclosed by the points of the TGraph.
1796/// The polygon segments, which are defined by the points of the TGraph, do not need to form a closed polygon,
1797/// since the last polygon segment, which closes the polygon, is taken as the line connecting the last TGraph point
1798/// with the first one. It is clear that the order of the point is essential in defining the polygon.
1799/// Also note that the segments should not intersect.
1800///
1801/// NB:
1802/// - if last=-1 (default) last is set to the last point.
1803/// - if (first <0) the first point (0) is taken.
1804///
1805/// ### Method:
1806///
1807/// There are many ways to calculate the surface of a polygon. It all depends on what kind of data
1808/// you have to deal with. The most evident solution would be to divide the polygon in triangles and
1809/// calculate the surface of them. But this can quickly become complicated as you will have to test
1810/// every segments of every triangles and check if they are intersecting with a current polygon's
1811/// segment or if it goes outside the polygon. Many calculations that would lead to many problems...
1812///
1813/// ### The solution (implemented by R.Brun)
1814/// Fortunately for us, there is a simple way to solve this problem, as long as the polygon's
1815/// segments don't intersect.
1816/// It takes the x coordinate of the current vertex and multiply it by the y coordinate of the next
1817/// vertex. Then it subtracts from it the result of the y coordinate of the current vertex multiplied
1818/// by the x coordinate of the next vertex. Then divide the result by 2 to get the surface/area.
1819///
1820/// ### Sources
1821/// - http://forums.wolfram.com/mathgroup/archive/1998/Mar/msg00462.html
1822/// - http://stackoverflow.com/questions/451426/how-do-i-calculate-the-surface-area-of-a-2d-polygon
1823
1825{
1826 if (first < 0) first = 0;
1827 if (last < 0) last = fNpoints - 1;
1828 if (last >= fNpoints) last = fNpoints - 1;
1829 if (first >= last) return 0;
1830 Int_t np = last - first + 1;
1831 Double_t sum = 0.0;
1832 //for(Int_t i=first;i<=last;i++) {
1833 // Int_t j = first + (i-first+1)%np;
1834 // sum += TMath::Abs(fX[i]*fY[j]);
1835 // sum -= TMath::Abs(fY[i]*fX[j]);
1836 //}
1837 for (Int_t i = first; i <= last; i++) {
1838 Int_t j = first + (i - first + 1) % np;
1839 sum += (fY[i] + fY[j]) * (fX[j] - fX[i]);
1840 }
1841 return 0.5 * TMath::Abs(sum);
1842}
1843
1844////////////////////////////////////////////////////////////////////////////////
1845/// Return 1 if the point (x,y) is inside the polygon defined by
1846/// the graph vertices 0 otherwise.
1847///
1848/// Algorithm:
1849///
1850/// The loop is executed with the end-point coordinates of a line segment
1851/// (X1,Y1)-(X2,Y2) and the Y-coordinate of a horizontal line.
1852/// The counter inter is incremented if the line (X1,Y1)-(X2,Y2) intersects
1853/// the horizontal line. In this case XINT is set to the X-coordinate of the
1854/// intersection point. If inter is an odd number, then the point x,y is within
1855/// the polygon.
1856
1858{
1859 return (Int_t)TMath::IsInside(x, y, fNpoints, fX, fY);
1860}
1861
1862////////////////////////////////////////////////////////////////////////////////
1863/// Least squares polynomial fitting without weights.
1864///
1865/// \param [in] m number of parameters
1866/// \param [in] a array of parameters
1867/// \param [in] xmin 1st point number to fit (default =0)
1868/// \param [in] xmax last point number to fit (default=fNpoints-1)
1869///
1870/// based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
1871
1873{
1874 const Double_t zero = 0.;
1875 const Double_t one = 1.;
1876 const Int_t idim = 20;
1877
1878 Double_t b[400] /* was [20][20] */;
1879 Int_t i, k, l, ifail;
1881 Double_t da[20], xk, yk;
1882 Int_t n = fNpoints;
1883 if (xmax <= xmin) {
1884 xmin = fX[0];
1885 xmax = fX[fNpoints-1];
1886 }
1887
1888 if (m <= 2) {
1889 LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
1890 return;
1891 }
1892 if (m > idim || m > n) return;
1893 da[0] = zero;
1894 for (l = 2; l <= m; ++l) {
1895 b[l-1] = zero;
1896 b[m + l*20 - 21] = zero;
1897 da[l-1] = zero;
1898 }
1899 Int_t np = 0;
1900 for (k = 0; k < fNpoints; ++k) {
1901 xk = fX[k];
1902 if (xk < xmin || xk > xmax) continue;
1903 np++;
1904 yk = fY[k];
1905 power = one;
1906 da[0] += yk;
1907 for (l = 2; l <= m; ++l) {
1908 power *= xk;
1909 b[l-1] += power;
1910 da[l-1] += power * yk;
1911 }
1912 for (l = 2; l <= m; ++l) {
1913 power *= xk;
1914 b[m + l*20 - 21] += power;
1915 }
1916 }
1917 b[0] = Double_t(np);
1918 for (i = 3; i <= m; ++i) {
1919 for (k = i; k <= m; ++k) {
1920 b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
1921 }
1922 }
1924
1925 if (ifail < 0) {
1926 a[0] = fY[0];
1927 for (i = 1; i < m; ++i) a[i] = 0;
1928 return;
1929 }
1930 for (i = 0; i < m; ++i) a[i] = da[i];
1931}
1932
1933////////////////////////////////////////////////////////////////////////////////
1934/// Least square linear fit without weights.
1935///
1936/// Fit a straight line (a0 + a1*x) to the data in this graph.
1937///
1938/// \param [in] ndata if ndata<0, fits the logarithm of the graph (used in InitExpo() to set
1939/// the initial parameter values for a fit with exponential function.
1940/// \param [in] a0 constant
1941/// \param [in] a1 slope
1942/// \param [in] ifail return parameter indicating the status of the fit (ifail=0, fit is OK)
1943/// \param [in] xmin, xmax fitting range
1944///
1945/// extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
1946
1948{
1950 Int_t i;
1952 Double_t fn, xk, yk;
1953 Double_t det;
1954 if (xmax <= xmin) {
1955 xmin = fX[0];
1956 xmax = fX[fNpoints-1];
1957 }
1958
1959 ifail = -2;
1960 xbar = ybar = x2bar = xybar = 0;
1961 Int_t np = 0;
1962 for (i = 0; i < fNpoints; ++i) {
1963 xk = fX[i];
1964 if (xk < xmin || xk > xmax) continue;
1965 np++;
1966 yk = fY[i];
1967 if (ndata < 0) {
1968 if (yk <= 0) yk = 1e-9;
1969 yk = TMath::Log(yk);
1970 }
1971 xbar += xk;
1972 ybar += yk;
1973 x2bar += xk * xk;
1974 xybar += xk * yk;
1975 }
1976 fn = Double_t(np);
1977 det = fn * x2bar - xbar * xbar;
1978 ifail = -1;
1979 if (det <= 0) {
1980 if (fn > 0) a0 = ybar / fn;
1981 else a0 = 0;
1982 a1 = 0;
1983 return;
1984 }
1985 ifail = 0;
1986 a0 = (x2bar * ybar - xbar * xybar) / det;
1987 a1 = (fn * xybar - xbar * ybar) / det;
1988}
1989
1990////////////////////////////////////////////////////////////////////////////////
1991/// Draw this graph with its current attributes.
1992
1998
1999////////////////////////////////////////////////////////////////////////////////
2000/// Draw the (x,y) as a graph.
2001
2007
2008////////////////////////////////////////////////////////////////////////////////
2009/// Draw the (x,y) as a histogram.
2010
2016
2017////////////////////////////////////////////////////////////////////////////////
2018/// Draw the stats
2019
2021{
2023 if (painter) painter->PaintStats(this, fit);
2024}
2025
2026////////////////////////////////////////////////////////////////////////////////
2027/// Print graph values.
2028
2030{
2031 for (Int_t i = 0; i < fNpoints; i++) {
2032 printf("x[%d]=%g, y[%d]=%g\n", i, fX[i], i, fY[i]);
2033 }
2034}
2035
2036////////////////////////////////////////////////////////////////////////////////
2037/// Recursively remove object from the list of functions
2038
2040{
2041 if (fFunctions) {
2044 }
2045 if (fHistogram == obj)
2046 fHistogram = nullptr;
2047}
2048
2049////////////////////////////////////////////////////////////////////////////////
2050/// Delete point close to the mouse position
2051/// Returns index of removed point (or -1 if nothing was changed)
2052
2054{
2055 if (!gPad) {
2056 Error("RemovePoint", "Cannot be used without gPad, requires last mouse position");
2057 return -1;
2058 }
2059
2060 Int_t px = gPad->GetEventX();
2061 Int_t py = gPad->GetEventY();
2062
2063 //localize point to be deleted
2064 Int_t ipoint = -2;
2065 // start with a small window (in case the mouse is very close to one point)
2066 for (Int_t i = 0; i < fNpoints; i++) {
2067 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
2068 Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
2069 if (dpx * dpx + dpy * dpy < 100) {
2070 ipoint = i;
2071 break;
2072 }
2073 }
2074 return RemovePoint(ipoint);
2075}
2076
2077////////////////////////////////////////////////////////////////////////////////
2078/// Delete point number ipoint
2079/// Returns index of removed point (or -1 if nothing was changed)
2080
2082{
2083 if ((ipoint < 0) || (ipoint >= fNpoints))
2084 return -1;
2085
2087 CopyAndRelease(ps, ipoint + 1, fNpoints--, ipoint);
2088 if (gPad) gPad->Modified();
2089 return ipoint;
2090}
2091
2092////////////////////////////////////////////////////////////////////////////////
2093/// Save the graph as .csv, .tsv or .txt. In case of any other extension, fall
2094/// back to TObject::SaveAs
2095///
2096/// The result can be immediately imported into Excel, gnuplot, Python or whatever,
2097/// without the needing to install pyroot, etc.
2098///
2099/// \param filename the name of the file where to store the graph
2100/// \param option some tuning options
2101///
2102/// The file extension defines the delimiter used:
2103/// - `.csv` : comma
2104/// - `.tsv` : tab
2105/// - `.txt` : space
2106///
2107/// By default file contains lines with (X, Y) coordinates. If errors are present,
2108/// (X, EX, Y, EY) are written. With asymmetric errors (X, EXL, EXH, Y, EYL, EYH) are stored.
2109/// If option contains "asroot" string, order of values will match such order in TGraph constructors.
2110/// So one will get (X, Y, EX, EY) or (X, Y, EXL, EXH, EYL, EYH).
2111///
2112/// Also one can directly select that kind of errors are stored:
2113/// - "errors" - (X, Y, EX, EY) will be stored
2114/// - "asymmerrors" - (X, Y, EXL, EXH, EYL, EYH) will be stored
2115/// - "noerrors" - just (X, Y) will be stored disregard of graph kind
2116///
2117/// If option contains "title" a title line is generated with the axis titles.
2118
2119
2120void TGraph::SaveAs(const char *filename, Option_t *option) const
2121{
2122 char del = '\0';
2123 TString ext = "";
2125 TString opt = option;
2126 opt.ToLower();
2127
2128 if (filename) {
2129 if (fname.EndsWith(".csv")) {del = ','; ext = "csv";}
2130 else if (fname.EndsWith(".tsv")) {del = '\t'; ext = "tsv";}
2131 else if (fname.EndsWith(".txt")) {del = ' '; ext = "txt";}
2132 }
2133 if (del) {
2134 std::ofstream out;
2135 out.open(filename, std::ios::out);
2136 if (!out.good ()) {
2137 Error("SaveAs", "cannot open file: %s", filename);
2138 return;
2139 }
2140 Bool_t store_title = opt.Contains("title");
2142 Bool_t as_root = opt.Contains("asroot") || opt.Contains("native");
2143 if (opt.Contains("noerrors"))
2144 no_errors = kTRUE;
2145 else if (opt.Contains("asymmerrors"))
2147 else if (opt.Contains("errors"))
2149 else if (InheritsFrom("TGraphErrors"))
2151 else if (InheritsFrom("TGraphAsymmErrors") || InheritsFrom("TGraphBentErrors"))
2153 else
2154 no_errors = kTRUE;
2155
2157 if (fHistogram) {
2160 }
2161 if (xtitle.IsNull())
2162 xtitle = "x";
2163 if (ytitle.IsNull())
2164 ytitle = "y";
2165
2166 if (plain_errors) {
2167 if(store_title) {
2168 if (as_root)
2169 out << "# " << xtitle << "\t" << ytitle << "\tex\tey\n";
2170 else
2171 out << "# " << xtitle << "\tex\t" << ytitle << "\tey\n";
2172 }
2173 for(int i = 0; i < fNpoints ; i++) {
2174 Double_t x = GetPointX(i);
2175 Double_t y = GetPointY(i);
2176 Double_t ex = GetErrorX(i);
2177 Double_t ey = GetErrorY(i);
2178 if (as_root)
2179 out << x << del << y << del << ex << del << ey << "\n";
2180 else
2181 out << x << del << ex << del << y << del << ey << "\n";
2182 }
2183 } else if (asymm_erros) {
2184 if(store_title) {
2185 if (as_root)
2186 out << "# " << xtitle << "\t" << ytitle << "\texl\texh\teyl\teyh\n";
2187 else
2188 out << "# " << xtitle << "\texl\texh\t" << ytitle << "\teyl\teyh\n";
2189 }
2190 for(int i = 0; i < GetN(); i++) {
2191 Double_t x = GetPointX(i);
2192 Double_t y = GetPointY(i);
2197 if (as_root)
2198 out << x << del << y << del << exl << del << exh << del << eyl << del << eyh << "\n";
2199 else
2200 out << x << del << exl << del << exh << del << y << del << eyl << del << eyh << "\n";
2201 }
2202 } else if (no_errors) {
2203 if(store_title)
2204 out << "# " << xtitle << "\t" << ytitle << "\n";
2205 for (Int_t i = 0 ; i < GetN(); i++)
2206 out << GetPointX(i) << del << GetPointY(i) << "\n";
2207 }
2208 out.close();
2209 Info("SaveAs", "%s file: %s has been generated", ext.Data(), filename);
2210 } else {
2212 }
2213}
2214
2215
2216////////////////////////////////////////////////////////////////////////////////
2217/// Save primitive as a C++ statement(s) on output stream out
2218
2219void TGraph::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
2220{
2221 TString args;
2222 if (fNpoints >= 1) {
2223 TString xname = SavePrimitiveVector(out, "graph_x", fNpoints, fX, kTRUE);
2224 TString yname = SavePrimitiveVector(out, "graph_y", fNpoints, fY);
2225 args.Form("%d, %s.data(), %s.data()", fNpoints, xname.Data(), yname.Data());
2226 }
2227
2228 SavePrimitiveConstructor(out, Class(), "graph", args, fNpoints < 1);
2229
2230 SaveHistogramAndFunctions(out, "graph", option);
2231}
2232
2233////////////////////////////////////////////////////////////////////////////////
2234/// Save histogram and list of functions of TGraph as C++ statement
2235/// Used in all TGraph-derived classes
2236
2237void TGraph::SaveHistogramAndFunctions(std::ostream &out, const char *varname, Option_t *option)
2238{
2239 thread_local Int_t frameNumber = 0;
2240
2241 TString ref = "Graph";
2242 if ((ref != GetName()) || (ref != GetTitle()))
2244
2245 SaveFillAttributes(out, varname, 0, 1000);
2246 SaveLineAttributes(out, varname, 1, 1, 1);
2247 SaveMarkerAttributes(out, varname, 1, 1, 1);
2248
2249 if (fHistogram) {
2251 fHistogram->SetName(TString::Format("Graph_histogram%d", ++frameNumber).Data());
2252 fHistogram->SavePrimitive(out, "nodraw");
2253 out << " " <<varname << "->SetHistogram(" << fHistogram->GetName() << ");\n";
2254 out << " \n";
2255 fHistogram->SetName(hname.Data());
2256 }
2257
2259
2260 if (!option)
2261 option = "";
2262 const char *l = strstr(option, "multigraph");
2263 if (l) {
2264 out << " multigraph->Add(" << varname << ",\"" << l + 10 << "\");\n";
2265 return;
2266 }
2267 l = strstr(option, "th2poly");
2268 if (l) {
2269 out << " " << l + 7 << "->AddBin(" << varname << ");\n";
2270 return;
2271 }
2272
2274}
2275
2276////////////////////////////////////////////////////////////////////////////////
2277/// Multiply the values of a TGraph by a constant c1.
2278///
2279/// If option contains "x" the x values are scaled
2280/// If option contains "y" the y values are scaled
2281/// If option contains "xy" both x and y values are scaled
2282
2284{
2285 TString opt = option; opt.ToLower();
2286 if (opt.Contains("x")) {
2287 for (Int_t i=0; i<GetN(); i++)
2288 GetX()[i] *= c1;
2289 }
2290 if (opt.Contains("y")) {
2291 for (Int_t i=0; i<GetN(); i++)
2292 GetY()[i] *= c1;
2293 }
2294}
2295
2296////////////////////////////////////////////////////////////////////////////////
2297/// Set number of points in the graph
2298/// Existing coordinates are preserved
2299/// New coordinates above fNpoints are preset to 0.
2300
2302{
2303 if (n < 0) n = 0;
2304 if (n == fNpoints) return;
2305 Double_t **ps = Allocate(n);
2306 CopyAndRelease(ps, 0, TMath::Min(fNpoints, n), 0);
2307 if (n > fNpoints) {
2309 }
2310 fNpoints = n;
2311}
2312
2313////////////////////////////////////////////////////////////////////////////////
2314/// Return kTRUE if kNotEditable bit is not set, kFALSE otherwise.
2315
2317{
2318 return TestBit(kNotEditable) ? kFALSE : kTRUE;
2319}
2320
2321////////////////////////////////////////////////////////////////////////////////
2322/// if editable=kFALSE, the graph cannot be modified with the mouse
2323/// by default a TGraph is editable
2324
2330
2331////////////////////////////////////////////////////////////////////////////////
2332/// Set highlight (enable/disable) mode for the graph
2333/// by default highlight mode is disable
2334
2336{
2337 if (IsHighlight() == set) return;
2338
2340 if (!painter) return;
2341 SetBit(kIsHighlight, set);
2342 painter->SetHighlight(this);
2343}
2344
2345/// Set the histogram underlying the TGraph. This transfers the ownership of h to the TGraph. The preexisting fHistogram will be deleted.
2347{
2348 delete fHistogram;
2349 fHistogram = h;
2350}
2351
2352////////////////////////////////////////////////////////////////////////////////
2353/// Set the maximum of the graph.
2354
2356{
2357 fMaximum = maximum;
2358 GetHistogram()->SetMaximum(maximum);
2359}
2360
2361////////////////////////////////////////////////////////////////////////////////
2362/// Set the minimum of the graph.
2363
2365{
2366 fMinimum = minimum;
2367 GetHistogram()->SetMinimum(minimum);
2368}
2369
2370////////////////////////////////////////////////////////////////////////////////
2371/// Set x and y values for point number i.
2372
2374{
2375 if (i < 0) return;
2377
2378 if (i >= fMaxSize) {
2379 Double_t **ps = ExpandAndCopy(i + 1, fNpoints);
2380 CopyAndRelease(ps, 0, 0, 0);
2381 }
2382 if (i >= fNpoints) {
2383 // points above i can be not initialized
2384 // set zero up to i-th point to avoid redefinition
2385 // of this method in descendant classes
2386 FillZero(fNpoints, i + 1);
2387 fNpoints = i + 1;
2388 }
2389 fX[i] = x;
2390 fY[i] = y;
2391 if (gPad) gPad->Modified();
2392}
2393
2394////////////////////////////////////////////////////////////////////////////////
2395/// Set x value for point i.
2396
2398{
2399 SetPoint(i, x, GetPointY(i));
2400}
2401
2402////////////////////////////////////////////////////////////////////////////////
2403/// Set y value for point i.
2404
2406{
2407 SetPoint(i, GetPointX(i), y);
2408}
2409
2410////////////////////////////////////////////////////////////////////////////////
2411/// Set graph name.
2412void TGraph::SetName(const char *name)
2413{
2414 fName = name;
2416}
2417
2418////////////////////////////////////////////////////////////////////////////////
2419/// Change (i.e. set) the title
2420///
2421/// if title is in the form `stringt;stringx;stringy;stringz`
2422/// the graph title is set to `stringt`, the x axis title to `stringx`,
2423/// the y axis title to `stringy`, and the z axis title to `stringz`.
2424///
2425/// To insert the character `;` in one of the titles, one should use `#;`
2426/// or `#semicolon`.
2427
2428void TGraph::SetTitle(const char* title)
2429{
2430 fTitle = title;
2431 fTitle.ReplaceAll("#;",2,"#semicolon",10);
2432 Int_t p = fTitle.Index(";");
2433
2434 if (p>0) {
2435 if (!fHistogram) GetHistogram();
2436 fHistogram->SetTitle(title);
2437 Int_t n = fTitle.Length()-p;
2438 if (p>0) fTitle.Remove(p,n);
2439 fTitle.ReplaceAll("#semicolon",10,"#;",2);
2440 } else {
2441 if (fHistogram) fHistogram->SetTitle(title);
2442 }
2443}
2444
2445////////////////////////////////////////////////////////////////////////////////
2446/// Set graph name and title
2447
2448void TGraph::SetNameTitle(const char *name, const char *title)
2449{
2450 SetName(name);
2451 SetTitle(title);
2452}
2453
2454////////////////////////////////////////////////////////////////////////////////
2455/// Set statistics option on/off.
2456///
2457/// By default, the statistics box is drawn.
2458/// The paint options can be selected via gStyle->SetOptStat.
2459/// This function sets/resets the kNoStats bit in the graph object.
2460/// It has priority over the Style option.
2461
2463{
2465 if (!stats) {
2467 //remove the "stats" object from the list of functions
2468 if (fFunctions) {
2469 TObject *obj = fFunctions->FindObject("stats");
2470 if (obj) {
2471 fFunctions->Remove(obj);
2472 delete obj;
2473 }
2474 }
2475 }
2476}
2477
2478////////////////////////////////////////////////////////////////////////////////
2479/// if size*2 <= fMaxSize allocate new arrays of size points,
2480/// copy points [0,oend).
2481/// Return newarray (passed or new instance if it was zero
2482/// and allocations are needed)
2483
2485{
2486 if (size * 2 > fMaxSize || !fMaxSize)
2487 return nullptr;
2488
2490 CopyPoints(newarrays, 0, oend, 0);
2491 return newarrays;
2492}
2493
2494////////////////////////////////////////////////////////////////////////////////
2495/// Sorts the points of this TGraph using in-place quicksort (see e.g. older glibc).
2496/// To compare two points the function parameter greaterfunc is used (see TGraph::CompareX for an
2497/// example of such a method, which is also the default comparison function for Sort). After
2498/// the sort, greaterfunc(this, i, j) will return kTRUE for all i>j if ascending == kTRUE, and
2499/// kFALSE otherwise.
2500///
2501/// The last two parameters are used for the recursive quick sort, stating the range to be sorted
2502///
2503/// Examples:
2504/// ~~~ {.cpp}
2505/// // sort points along x axis
2506/// graph->Sort();
2507/// // sort points along their distance to origin
2508/// graph->Sort(&TGraph::CompareRadius);
2509///
2510/// Bool_t CompareErrors(const TGraph* gr, Int_t i, Int_t j) {
2511/// const TGraphErrors* ge=(const TGraphErrors*)gr;
2512/// return (ge->GetEY()[i]>ge->GetEY()[j]); }
2513/// // sort using the above comparison function, largest errors first
2514/// graph->Sort(&CompareErrors, kFALSE);
2515/// ~~~
2516
2517void TGraph::Sort(Bool_t (*greaterfunc)(const TGraph *, Int_t, Int_t) /*=TGraph::CompareX()*/,
2518 Bool_t ascending /*=kTRUE*/, Int_t low /*=0*/, Int_t high /*=-1111*/)
2519{
2520 // set the bit in case of an ascending =sort in X
2521 if (greaterfunc == TGraph::CompareX && ascending && low == 0 && high == -1111)
2523
2524 if (high == -1111)
2525 high = fNpoints - 1;
2526
2527 // Create a vector to store the indices of the graph data points.
2528 // We use std::vector<Int_t> instead of std::vector<ULong64_t> to match the input type
2529 // required by the comparison operator's signature provided as `greaterfunc`
2530 std::vector<int> sorting_indices(fNpoints);
2531 std::iota(sorting_indices.begin(), sorting_indices.end(), 0);
2532
2533 // Sort the indices using the provided comparison function
2534 // We use std::stable_sort here because the libc++ implementation of std::sort
2535 // is not standard-compliant until LLVM 14 which caused errors on the mac nodes
2536 // of our CI, related issue: https://github.com/llvm/llvm-project/issues/21211
2537 std::stable_sort(sorting_indices.begin() + low, sorting_indices.begin() + high + 1,
2538 [&](int left, int right) { return left != right && greaterfunc(this, left, right) != ascending; });
2539
2540 Int_t numSortedPoints = high - low + 1;
2542}
2543
2544////////////////////////////////////////////////////////////////////////////////
2545/// Stream an object of class TGraph.
2546
2548{
2549 if (b.IsReading()) {
2550 UInt_t R__s, R__c;
2551 Version_t R__v = b.ReadVersion(&R__s, &R__c);
2552 if (R__v > 2) {
2553 b.ReadClassBuffer(TGraph::Class(), this, R__v, R__s, R__c);
2554 if (fHistogram) fHistogram->SetDirectory(nullptr);
2555 TIter next(fFunctions);
2556 TObject *obj;
2557 while ((obj = next())) {
2558 if (obj->InheritsFrom(TF1::Class())) {
2559 TF1 *f1 = (TF1*)obj;
2560 f1->SetParent(this);
2561 }
2562 }
2564 return;
2565 }
2566 //====process old versions before automatic schema evolution
2571 b >> fNpoints;
2573 fX = new Double_t[fNpoints];
2574 fY = new Double_t[fNpoints];
2575 if (R__v < 2) {
2576 Float_t *x = new Float_t[fNpoints];
2577 Float_t *y = new Float_t[fNpoints];
2578 b.ReadFastArray(x, fNpoints);
2579 b.ReadFastArray(y, fNpoints);
2580 for (Int_t i = 0; i < fNpoints; i++) {
2581 fX[i] = x[i];
2582 fY[i] = y[i];
2583 }
2584 delete [] y;
2585 delete [] x;
2586 } else {
2587 b.ReadFastArray(fX, fNpoints);
2588 b.ReadFastArray(fY, fNpoints);
2589 }
2590 b >> fFunctions;
2591 b >> fHistogram;
2592 if (fHistogram) fHistogram->SetDirectory(nullptr);
2593 if (R__v < 2) {
2594 Float_t mi, ma;
2595 b >> mi;
2596 b >> ma;
2597 fMinimum = mi;
2598 fMaximum = ma;
2599 } else {
2600 b >> fMinimum;
2601 b >> fMaximum;
2602 }
2603 b.CheckByteCount(R__s, R__c, TGraph::IsA());
2604 //====end of old versions
2605
2606 } else {
2607 b.WriteClassBuffer(TGraph::Class(), this);
2608 }
2609}
2610
2611////////////////////////////////////////////////////////////////////////////////
2612/// Swap points.
2613
2619
2620////////////////////////////////////////////////////////////////////////////////
2621/// Update the fX and fY arrays with the sorted values.
2622
2624{
2625 std::vector<Double_t> fXSorted(numSortedPoints);
2626 std::vector<Double_t> fYSorted(numSortedPoints);
2627
2628 // Fill the sorted X and Y values based on the sorted indices
2629 std::generate(fXSorted.begin(), fXSorted.end(),
2630 [begin = low, &sorting_indices, this]() mutable { return fX[sorting_indices[begin++]]; });
2631 std::generate(fYSorted.begin(), fYSorted.end(),
2632 [begin = low, &sorting_indices, this]() mutable { return fY[sorting_indices[begin++]]; });
2633
2634 // Copy the sorted X and Y values back to the original arrays
2635 std::copy(fXSorted.begin(), fXSorted.end(), fX + low);
2636 std::copy(fYSorted.begin(), fYSorted.end(), fY + low);
2637}
2638
2639////////////////////////////////////////////////////////////////////////////////
2640/// Swap values.
2641
2643{
2644 Double_t tmp = arr[pos1];
2645 arr[pos1] = arr[pos2];
2646 arr[pos2] = tmp;
2647}
2648
2649////////////////////////////////////////////////////////////////////////////////
2650/// Set current style settings in this graph
2651/// This function is called when either TCanvas::UseCurrentStyle
2652/// or TROOT::ForceStyle have been invoked.
2653
2684
2685////////////////////////////////////////////////////////////////////////////////
2686/// Adds all graphs from the collection to this graph.
2687/// Returns the total number of points in the result or -1 in case of an error.
2688
2690{
2691 TIter next(li);
2692 while (TObject* o = next()) {
2693 TGraph *g = dynamic_cast<TGraph*>(o);
2694 if (!g) {
2695 Error("Merge",
2696 "Cannot merge - an object which doesn't inherit from TGraph found in the list");
2697 return -1;
2698 }
2699 DoMerge(g);
2700 }
2701 return GetN();
2702}
2703
2704////////////////////////////////////////////////////////////////////////////////
2705/// protected function to perform the merge operation of a graph
2706
2708{
2709 Double_t x = 0, y = 0;
2710 for (Int_t i = 0 ; i < g->GetN(); i++) {
2711 g->GetPoint(i, x, y);
2712 SetPoint(GetN(), x, y);
2713 }
2714 return kTRUE;
2715}
2716
2717////////////////////////////////////////////////////////////////////////////////
2718/// Move all graph points on specified values dx,dy
2719/// If log argument specified, calculation done in logarithmic scale like:
2720/// new_value = exp( log(old_value) + delta );
2721
2723{
2724 Double_t x = 0, y = 0;
2725 for (Int_t i = 0 ; i < GetN(); i++) {
2726 GetPoint(i, x, y);
2727 if (!logx) {
2728 x += dx;
2729 } else if (x > 0) {
2730 x = TMath::Exp(TMath::Log(x) + dx);
2731 }
2732 if (!logy) {
2733 y += dy;
2734 } else if (y > 0) {
2735 y = TMath::Exp(TMath::Log(y) + dy);
2736 }
2737 SetPoint(i, x, y);
2738 }
2739}
2740
2741
2742////////////////////////////////////////////////////////////////////////////////
2743/// Find zero of a continuous function.
2744/// This function finds a real zero of the continuous real
2745/// function Y(X) in a given interval (A,B). See accompanying
2746/// notes for details of the argument list and calling sequence
2747
2750{
2751 static Double_t a, b, ya, ytest, y1, x1, h;
2752 static Int_t j1, it, j3, j2;
2753 Double_t yb, x2;
2754 yb = 0;
2755
2756 // Calculate Y(X) at X=AZ.
2757 if (k <= 0) {
2758 a = AZ;
2759 b = BZ;
2760 X = a;
2761 j1 = 1;
2762 it = 1;
2763 k = j1;
2764 return;
2765 }
2766
2767 // Test whether Y(X) is sufficiently small.
2768
2769 if (TMath::Abs(Y) <= E2) {
2770 k = 2;
2771 return;
2772 }
2773
2774 // Calculate Y(X) at X=BZ.
2775
2776 if (j1 == 1) {
2777 ya = Y;
2778 X = b;
2779 j1 = 2;
2780 return;
2781 }
2782 // Test whether the signs of Y(AZ) and Y(BZ) are different.
2783 // if not, begin the binary subdivision.
2784
2785 if (j1 != 2) goto L100;
2786 if (ya * Y < 0) goto L120;
2787 x1 = a;
2788 y1 = ya;
2789 j1 = 3;
2790 h = b - a;
2791 j2 = 1;
2792 x2 = a + 0.5 * h;
2793 j3 = 1;
2794 it++; //*-*- Check whether (maxiterations) function values have been calculated.
2795 if (it >= maxiterations) k = j1;
2796 else X = x2;
2797 return;
2798
2799 // Test whether a bracket has been found .
2800 // If not,continue the search
2801
2802L100:
2803 if (j1 > 3) goto L170;
2804 if (ya*Y >= 0) {
2805 if (j3 >= j2) {
2806 h = 0.5 * h;
2807 j2 = 2 * j2;
2808 a = x1;
2809 ya = y1;
2810 x2 = a + 0.5 * h;
2811 j3 = 1;
2812 } else {
2813 a = X;
2814 ya = Y;
2815 x2 = X + h;
2816 j3++;
2817 }
2818 it++;
2819 if (it >= maxiterations) k = j1;
2820 else X = x2;
2821 return;
2822 }
2823
2824 // The first bracket has been found.calculate the next X by the
2825 // secant method based on the bracket.
2826
2827L120:
2828 b = X;
2829 yb = Y;
2830 j1 = 4;
2831L130:
2832 if (TMath::Abs(ya) > TMath::Abs(yb)) {
2833 x1 = a;
2834 y1 = ya;
2835 X = b;
2836 Y = yb;
2837 } else {
2838 x1 = b;
2839 y1 = yb;
2840 X = a;
2841 Y = ya;
2842 }
2843
2844 // Use the secant method based on the function values y1 and Y.
2845 // check that x2 is inside the interval (a,b).
2846
2847L150:
2848 x2 = X - Y * (X - x1) / (Y - y1);
2849 x1 = X;
2850 y1 = Y;
2852 if ((x2 - a)*(x2 - b) < 0) {
2853 it++;
2854 if (it >= maxiterations) k = j1;
2855 else X = x2;
2856 return;
2857 }
2858
2859 // Calculate the next value of X by bisection . Check whether
2860 // the maximum accuracy has been achieved.
2861
2862L160:
2863 x2 = 0.5 * (a + b);
2864 ytest = 0;
2865 if ((x2 - a)*(x2 - b) >= 0) {
2866 k = 2;
2867 return;
2868 }
2869 it++;
2870 if (it >= maxiterations) k = j1;
2871 else X = x2;
2872 return;
2873
2874
2875 // Revise the bracket (a,b).
2876
2877L170:
2878 if (j1 != 4) return;
2879 if (ya * Y < 0) {
2880 b = X;
2881 yb = Y;
2882 } else {
2883 a = X;
2884 ya = Y;
2885 }
2886
2887 // Use ytest to decide the method for the next value of X.
2888
2889 if (ytest <= 0) goto L130;
2890 if (TMath::Abs(Y) - ytest <= 0) goto L150;
2891 goto L160;
2892}
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define g(i)
Definition RSha256.hxx:105
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
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 Ssiz_t kNPOS
The equivalent of std::string::npos for the ROOT class TString.
Definition RtypesCore.h:131
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
#define X(type, name)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t SetLineWidth
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 SetFillStyle
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 del
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 char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t SetLineColor
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 x2
Option_t Option_t TPoint TPoint const char x1
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
Option_t Option_t SetFillColor
Option_t Option_t SetMarkerStyle
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b)
Extracted from CERN Program library routine DSEQN.
Definition TH1.cxx:4885
float xmin
float ymin
float xmax
float ymax
#define gROOT
Definition TROOT.h:411
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
void Printf(const char *fmt,...)
Formats a string in a circular formatting buffer and prints the string.
Definition TString.cxx:2509
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
R__EXTERN TSystem * gSystem
Definition TSystem.h:572
#define gPad
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition DataRange.h:35
const_iterator begin() const
const_iterator end() const
Fill Area Attributes class.
Definition TAttFill.h:20
virtual void Streamer(TBuffer &)
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:31
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition TAttFill.cxx:206
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:32
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:238
Line Attributes class.
Definition TAttLine.h:20
virtual void Streamer(TBuffer &)
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:37
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:36
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:176
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition TAttLine.cxx:210
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:274
Marker Attributes class.
Definition TAttMarker.h:20
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:33
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:32
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition TAttMarker.h:34
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
virtual void Streamer(TBuffer &)
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:46
Class to manage histogram axis.
Definition TAxis.h:32
const char * GetTitle() const override
Returns title of object.
Definition TAxis.h:137
virtual void SetTimeDisplay(Int_t value)
Definition TAxis.h:173
void RotateTitle(Bool_t rotate=kTRUE)
Rotate title by 180 degrees.
Definition TAxis.h:205
void CenterTitle(Bool_t center=kTRUE)
Center axis title.
Definition TAxis.h:196
void SetNoExponent(Bool_t noExponent=kTRUE)
Set the NoExponent flag By default, an exponent of the form 10^N is used when the label value are eit...
Definition TAxis.h:235
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition TAxis.h:166
virtual void SetTimeFormat(const char *format="")
Change the format used for time plotting.
Definition TAxis.cxx:1149
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
Buffer base class used for serializing objects.
Definition TBuffer.h:43
Collection abstract base class.
Definition TCollection.h:65
virtual Bool_t IsEmpty() const
TObject * Clone(const char *newname="") const override
Make a clone of an collection using the Streamer facility.
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:505
1-Dim function class
Definition TF1.h:182
static TClass * Class()
virtual Int_t GetNpar() const
Definition TF1.h:461
virtual void SetParent(TObject *p=nullptr)
Definition TF1.h:662
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set lower and upper limits for parameter ipar.
Definition TF1.cxx:3538
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:623
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual Double_t GetPointX(Int_t i) const
Get x value for point i.
Definition TGraph.cxx:1560
static TClass * Class()
virtual Double_t Integral(Int_t first=0, Int_t last=-1) const
Integrate the TGraph data within a given (index) range.
Definition TGraph.cxx:1824
Int_t fNpoints
Number of points <= fMaxSize.
Definition TGraph.h:46
virtual Int_t IsInside(Double_t x, Double_t y) const
Return 1 if the point (x,y) is inside the polygon defined by the graph vertices 0 otherwise.
Definition TGraph.cxx:1857
virtual void LeastSquareFit(Int_t m, Double_t *a, Double_t xmin=0, Double_t xmax=0)
Least squares polynomial fitting without weights.
Definition TGraph.cxx:1872
virtual Double_t Chisquare(TF1 *f1, Option_t *option="") const
Return the chisquare of this graph with respect to f1.
Definition TGraph.cxx:694
void UseCurrentStyle() override
Set current style settings in this graph This function is called when either TCanvas::UseCurrentStyle...
Definition TGraph.cxx:2654
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2373
Double_t * GetY() const
Definition TGraph.h:139
virtual Int_t Merge(TCollection *list)
Adds all graphs from the collection to this graph.
Definition TGraph.cxx:2689
Int_t fMaxSize
!Current dimension of arrays fX and fY
Definition TGraph.h:45
TString fOption
Options used for drawing the graph.
Definition TGraph.h:53
~TGraph() override
Graph default destructor.
Definition TGraph.cxx:583
Double_t ** ShrinkAndCopy(Int_t size, Int_t iend)
if size*2 <= fMaxSize allocate new arrays of size points, copy points [0,oend).
Definition TGraph.cxx:2484
virtual Double_t GetRMS(Int_t axis=1) const
Return RMS of X (axis=1) or Y (axis=2)
Definition TGraph.cxx:1356
TH1F * fHistogram
Pointer to histogram used for drawing axis.
Definition TGraph.h:50
void Paint(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:1993
@ kNotEditable
Bit set if graph is non editable.
Definition TGraph.h:77
@ kIsHighlight
Bit set if graph is highlight.
Definition TGraph.h:79
@ kIsSortedX
Graph is sorted in X points.
Definition TGraph.h:78
@ kClipFrame
Clip to the frame boundary.
Definition TGraph.h:75
@ kResetHisto
fHistogram must be reset in GetHistogram
Definition TGraph.h:76
@ kNoStats
Don't draw stats box.
Definition TGraph.h:74
virtual Double_t GetErrorXlow(Int_t bin) const
It always returns a negative value.
Definition TGraph.cxx:1404
virtual void MovePoints(Double_t dx, Double_t dy, Bool_t logx=kFALSE, Bool_t logy=kFALSE)
Move all graph points on specified values dx,dy If log argument specified, calculation done in logari...
Definition TGraph.cxx:2722
virtual Double_t GetErrorYlow(Int_t bin) const
It always returns a negative value.
Definition TGraph.cxx:1422
virtual void UpdateArrays(const std::vector< Int_t > &sorting_indices, Int_t numSortedPoints, Int_t low)
Update the fX and fY arrays with the sorted values.
Definition TGraph.cxx:2623
virtual void CopyAndRelease(Double_t **newarrays, Int_t ibegin, Int_t iend, Int_t obegin)
Copy points from fX and fY to arrays[0] and arrays[1] or to fX and fY if arrays == 0 and ibegin !...
Definition TGraph.cxx:779
Double_t GetMinimum() const
Definition TGraph.h:151
void Print(Option_t *chopt="") const override
Print graph values.
Definition TGraph.cxx:2029
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum of the graph.
Definition TGraph.cxx:2355
static Bool_t CompareY(const TGraph *gr, Int_t left, Int_t right)
Return kTRUE if fY[left] > fY[right]. Can be used by Sort.
Definition TGraph.cxx:730
static Bool_t CompareRadius(const TGraph *gr, Int_t left, Int_t right)
Return kTRUE if point number "left"'s distance to origin is bigger than that of point number "right".
Definition TGraph.cxx:739
virtual Double_t GetErrorYhigh(Int_t bin) const
It always returns a negative value.
Definition TGraph.cxx:1413
TClass * IsA() const override
Definition TGraph.h:202
static Bool_t CompareX(const TGraph *gr, Int_t left, Int_t right)
Return kTRUE if fX[left] > fX[right]. Can be used by Sort.
Definition TGraph.cxx:722
Int_t GetN() const
Definition TGraph.h:131
TF1 * GetFunction(const char *name) const
Return pointer to function with name.
Definition TGraph.cxx:1433
virtual void LeastSquareLinearFit(Int_t n, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin=0, Double_t xmax=0)
Least square linear fit without weights.
Definition TGraph.cxx:1947
Double_t * fY
[fNpoints] array of Y points
Definition TGraph.h:48
Bool_t CtorAllocate()
In constructors set fNpoints than call this method.
Definition TGraph.cxx:822
virtual void DrawGraph(Int_t n, const Int_t *x, const Int_t *y, Option_t *option="")
Draw this graph with new attributes.
Definition TGraph.cxx:906
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0)
Fit this graph with function with name fname.
Definition TGraph.cxx:1271
virtual void Sort(Bool_t(*greater)(const TGraph *, Int_t, Int_t)=&TGraph::CompareX, Bool_t ascending=kTRUE, Int_t low=0, Int_t high=-1111)
Sorts the points of this TGraph using in-place quicksort (see e.g.
Definition TGraph.cxx:2517
static Bool_t CompareArg(const TGraph *gr, Int_t left, Int_t right)
Return kTRUE if point number "left"'s argument (angle with respect to positive x-axis) is bigger than...
Definition TGraph.cxx:711
virtual void ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
Compute the x/y range of the points in this graph.
Definition TGraph.cxx:748
char * GetObjectInfo(Int_t px, Int_t py) const override
Implementation to get information on point of graph at cursor position Adapted from class TH1.
Definition TGraph.cxx:1601
Double_t ** AllocateArrays(Int_t Narrays, Int_t arraySize)
Allocate arrays.
Definition TGraph.cxx:615
virtual void Scale(Double_t c1=1., Option_t *option="y")
Multiply the values of a TGraph by a constant c1.
Definition TGraph.cxx:2283
TList * fFunctions
Pointer to list of functions (fits and user)
Definition TGraph.h:49
virtual Double_t GetCovariance() const
Return covariance of vectors x,y.
Definition TGraph.cxx:1325
static void SwapValues(Double_t *arr, Int_t pos1, Int_t pos2)
Swap values.
Definition TGraph.cxx:2642
void Streamer(TBuffer &) override
Stream an object of class TGraph.
Definition TGraph.cxx:2547
void Zero(Int_t &k, Double_t AZ, Double_t BZ, Double_t E2, Double_t &X, Double_t &Y, Int_t maxiterations)
Find zero of a continuous function.
Definition TGraph.cxx:2748
virtual Double_t ** Allocate(Int_t newsize)
Allocate internal data structures for newsize points.
Definition TGraph.cxx:607
virtual void FitPanel()
Display a GUI panel with all graph fit options.
Definition TGraph.cxx:1291
void Browse(TBrowser *b) override
Browse.
Definition TGraph.cxx:668
virtual Bool_t DoMerge(const TGraph *g)
protected function to perform the merge operation of a graph
Definition TGraph.cxx:2707
virtual void InsertPointBefore(Int_t ipoint, Double_t x, Double_t y)
Insert a new point with coordinates (x,y) before the point number ipoint.
Definition TGraph.cxx:1765
TList * GetListOfFunctions() const
Definition TGraph.h:125
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TGraph.cxx:1073
Double_t * GetX() const
Definition TGraph.h:138
void SaveAs(const char *filename="graph", Option_t *option="") const override
Save the graph as .csv, .tsv or .txt.
Definition TGraph.cxx:2120
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
Interpolate points in this graph at x using a TSpline.
Definition TGraph.cxx:970
virtual void InitExpo(Double_t xmin=0, Double_t xmax=0)
Compute Initial values of parameters for an exponential.
Definition TGraph.cxx:1674
virtual Int_t RemovePoint()
Delete point close to the mouse position Returns index of removed point (or -1 if nothing was changed...
Definition TGraph.cxx:2053
virtual void InitGaus(Double_t xmin=0, Double_t xmax=0)
Compute Initial values of parameters for a gaussian.
Definition TGraph.cxx:1636
virtual Bool_t IsHighlight() const
Definition TGraph.h:166
virtual void Add(TF1 *f, Double_t c1=1)
Performs the operation: y = y + c1*f(x,y) Errors are not recalculated.
Definition TGraph.cxx:639
virtual void Apply(TF1 *f)
Apply function f to all the data points f may be a 1-D function TF1 or 2-d function TF2 The Y values ...
Definition TGraph.cxx:655
void SetName(const char *name="") override
Set graph name.
Definition TGraph.cxx:2412
virtual void SetHighlight(Bool_t set=kTRUE)
Set highlight (enable/disable) mode for the graph by default highlight mode is disable.
Definition TGraph.cxx:2335
virtual void SwapPoints(Int_t pos1, Int_t pos2)
Swap points.
Definition TGraph.cxx:2614
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:848
TAxis * GetXaxis() const
Get x axis of the graph.
Definition TGraph.cxx:1582
Bool_t GetEditable() const
Return kTRUE if kNotEditable bit is not set, kFALSE otherwise.
Definition TGraph.cxx:2316
virtual Double_t GetCorrelationFactor() const
Return graph correlation factor.
Definition TGraph.cxx:1313
virtual void FillZero(Int_t begin, Int_t end, Bool_t from_ctor=kTRUE)
Set zero values for point arrays in the range [begin, end) Should be redefined in descendant classes.
Definition TGraph.cxx:1119
Double_t ** ExpandAndCopy(Int_t size, Int_t iend)
if size > fMaxSize allocate new arrays of 2*size points and copy iend first points.
Definition TGraph.cxx:1106
virtual void Expand(Int_t newsize)
If array sizes <= newsize, expand storage to 2*newsize.
Definition TGraph.cxx:1082
void SaveHistogramAndFunctions(std::ostream &out, const char *varname, Option_t *option)
Save histogram and list of functions of TGraph as C++ statement Used in all TGraph-derived classes.
Definition TGraph.cxx:2237
virtual Double_t GetMean(Int_t axis=1) const
Return mean value of X (axis=1) or Y (axis=2)
Definition TGraph.cxx:1341
Double_t * fX
[fNpoints] array of X points
Definition TGraph.h:47
virtual void PaintStats(TF1 *fit)
Draw the stats.
Definition TGraph.cxx:2020
TAxis * GetYaxis() const
Get y axis of the graph.
Definition TGraph.cxx:1591
TObject * FindObject(const char *name) const override
Search object named name in the list of functions.
Definition TGraph.cxx:1128
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TGraph.cxx:2462
virtual TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1444
virtual Double_t GetErrorY(Int_t bin) const
It always returns a negative value. Real implementation in TGraphErrors.
Definition TGraph.cxx:1386
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TGraph.cxx:2219
virtual Double_t GetPointY(Int_t i) const
Get y value for point i.
Definition TGraph.cxx:1571
Double_t fMinimum
Minimum value for plotting along y.
Definition TGraph.h:51
void PaintGraph(Int_t npoints, const Double_t *x, const Double_t *y, Option_t *chopt)
Draw the (x,y) as a graph.
Definition TGraph.cxx:2002
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2428
virtual Int_t InsertPoint()
Insert a new point at the mouse position.
Definition TGraph.cxx:1715
void RecursiveRemove(TObject *obj) override
Recursively remove object from the list of functions.
Definition TGraph.cxx:2039
virtual void SetPointY(Int_t i, Double_t y)
Set y value for point i.
Definition TGraph.cxx:2405
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a graph.
Definition TGraph.cxx:896
virtual void SetHistogram(TH1F *h)
Set the histogram underlying the TGraph. This transfers the ownership of h to the TGraph....
Definition TGraph.cxx:2346
virtual void DrawPanel()
Display a panel with all graph drawing options.
Definition TGraph.cxx:949
void PaintGrapHist(Int_t npoints, const Double_t *x, const Double_t *y, Option_t *chopt)
Draw the (x,y) as a histogram.
Definition TGraph.cxx:2011
void SetNameTitle(const char *name="", const char *title="") override
Set graph name and title.
Definition TGraph.cxx:2448
virtual void SetPointX(Int_t i, Double_t x)
Set x value for point i.
Definition TGraph.cxx:2397
virtual void Set(Int_t n)
Set number of points in the graph Existing coordinates are preserved New coordinates above fNpoints a...
Definition TGraph.cxx:2301
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition TGraph.cxx:1549
virtual void SetEditable(Bool_t editable=kTRUE)
if editable=kFALSE, the graph cannot be modified with the mouse by default a TGraph is editable
Definition TGraph.cxx:2325
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum of the graph.
Definition TGraph.cxx:2364
TGraph()
Graph default constructor.
Definition TGraph.cxx:125
virtual Bool_t CopyPoints(Double_t **newarrays, Int_t ibegin, Int_t iend, Int_t obegin)
Copy points from fX and fY to arrays[0] and arrays[1] or to fX and fY if arrays == 0 and ibegin !...
Definition TGraph.cxx:796
Double_t fMaximum
Maximum value for plotting along y.
Definition TGraph.h:52
virtual Double_t GetErrorXhigh(Int_t bin) const
It always returns a negative value.
Definition TGraph.cxx:1395
TGraph & operator=(const TGraph &)
Equal operator for this graph.
Definition TGraph.cxx:249
virtual void InitPolynom(Double_t xmin=0, Double_t xmax=0)
Compute Initial values of parameters for a polynom.
Definition TGraph.cxx:1695
virtual Double_t GetErrorX(Int_t bin) const
It always returns a negative value. Real implementation in TGraphErrors.
Definition TGraph.cxx:1378
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:878
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
virtual void SetDirectory(TDirectory *dir)
By default, when a histogram is created, it is added to the list of histogram objects in the current ...
Definition TH1.cxx:8988
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6769
void UseCurrentStyle() override
Copy current attributes from/to current style.
Definition TH1.cxx:7515
@ kNoStats
Don't draw stats box.
Definition TH1.h:403
TAxis * GetXaxis()
Definition TH1.h:571
virtual void SetMaximum(Double_t maximum=-1111)
Definition TH1.h:652
TAxis * GetYaxis()
Definition TH1.h:572
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:653
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TH1.cxx:7322
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:9011
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
Definition TH1.cxx:2734
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:7481
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:575
void RecursiveRemove(TObject *obj) override
Remove object from this collection and recursively remove the object from all other objects (and coll...
Definition TList.cxx:761
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:819
TObject * First() const override
Return the first object in the list. Returns 0 when list is empty.
Definition TList.cxx:656
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
void SavePrimitiveNameTitle(std::ostream &out, const char *variable_name)
Save object name and title into the output stream "out".
Definition TNamed.cxx:135
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
void Streamer(TBuffer &) override
Stream an object of class TObject.
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
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:50
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:202
virtual void UseCurrentStyle()
Set current style settings in this object This function is called when either TCanvas::UseCurrentStyl...
Definition TObject.cxx:902
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1074
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:203
virtual void SaveAs(const char *filename="", Option_t *option="") const
Save this object in the file specified by filename.
Definition TObject.cxx:705
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:881
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:543
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1088
void MakeZombie()
Definition TObject.h:53
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:839
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:771
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:790
void ResetBit(UInt_t f)
Definition TObject.h:201
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:68
@ kInvalidObject
if object ctor succeeded but object should not be used
Definition TObject.h:78
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1062
Longptr_t ExecPlugin(int nargs)
Int_t LoadPlugin()
Load the plugin library for this handler.
Class to create third splines to interpolate knots Arbitrary conditions can be introduced for first a...
Definition TSpline.h:182
Double_t Eval(Double_t x) const override
Eval this spline at x.
Definition TSpline.cxx:780
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
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
void Clear()
Clear string without changing its capacity.
Definition TString.cxx:1241
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
Definition TString.h:703
const char * Data() const
Definition TString.h:384
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:713
void ToUpper()
Change string to upper case.
Definition TString.cxx:1202
Bool_t IsNull() const
Definition TString.h:422
TString & Remove(Ssiz_t pos)
Definition TString.h:694
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
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2362
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:641
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:660
void SetHistFillColor(Color_t color=1)
Definition TStyle.h:383
Color_t GetHistLineColor() const
Definition TStyle.h:235
Bool_t IsReading() const
Definition TStyle.h:300
void SetHistLineStyle(Style_t styl=0)
Definition TStyle.h:386
Style_t GetHistFillStyle() const
Definition TStyle.h:236
Color_t GetHistFillColor() const
Definition TStyle.h:234
void SetHistLineColor(Color_t color=1)
Definition TStyle.h:384
Style_t GetHistLineStyle() const
Definition TStyle.h:237
void SetHistFillStyle(Style_t styl=0)
Definition TStyle.h:385
Width_t GetHistLineWidth() const
Definition TStyle.h:238
void SetHistLineWidth(Width_t width=1)
Definition TStyle.h:387
virtual Bool_t ExpandPathName(TString &path)
Expand a pathname getting rid of special shell characters like ~.
Definition TSystem.cxx:1285
TVectorT.
Definition TVectorT.h:29
Abstract Base Class for Fitting.
static TVirtualFitter * GetFitter()
static: return the current Fitter
Abstract interface to a histogram painter.
static TVirtualGraphPainter * GetPainter()
Static function returning a pointer to the current graph painter.
TLine * line
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
Double_t ey[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TGraphErrors * gr
Definition legend1.C:25
Double_t ex[n]
Definition legend1.C:17
TF1 * f1
Definition legend1.C:11
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
double Chisquare(const TH1 &h1, TF1 &f1, bool useRange, EChisquareType type)
compute the chi2 value for an histogram given a function (see TH1::Chisquare for the documentation)
void FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption)
Decode list of options into fitOption.
Definition HFitImpl.cxx:685
Bool_t IsInside(T xp, T yp, Int_t np, T *x, T *y)
Function which returns kTRUE if point xp,yp lies inside the polygon defined by the np points in array...
Definition TMath.h:1320
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:720
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:657
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:767
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:197
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
Definition TMathBase.h:413
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:329
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
TLine l
Definition textangle.C:4
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2339