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 = h->GetXaxis()->GetNbins();
362 }
363
364 if (!CtorAllocate())
365 return;
366
367 auto xaxis = h->GetXaxis();
368 for (Int_t i = 0; i < fNpoints; i++) {
369 fX[i] = xaxis->GetBinCenter(i + 1);
370 fY[i] = h->GetBinContent(i + 1);
371 }
372 h->TAttLine::Copy(*this);
373 h->TAttFill::Copy(*this);
374 h->TAttMarker::Copy(*this);
375
376 std::string gname = "Graph_from_" + std::string(h->GetName());
377 SetName(gname.c_str());
378 SetTitle(h->GetTitle());
379}
380
381////////////////////////////////////////////////////////////////////////////////
382/// Graph constructor importing its parameters from the TF1 object passed as argument
383/// - if option =="" (default), a TGraph is created with points computed
384/// at the fNpx points of f.
385/// - if option =="d", a TGraph is created with points computed with the derivatives
386/// at the fNpx points of f.
387/// - if option =="i", a TGraph is created with points computed with the integral
388/// at the fNpx points of f.
389/// - if option =="I", a TGraph is created with points computed with the integral
390/// at the fNpx+1 points of f and the integral is normalized to 1.
391
393 : TNamed("Graph", "Graph"), TAttFill(0, 1000)
394{
395 char coption = ' ';
396 if (!f) {
397 Error("TGraph", "Pointer to function is null");
398 fNpoints = 0;
399 } else {
400 fNpoints = f->GetNpx();
401 if (option) coption = *option;
402 if (coption == 'i' || coption == 'I') fNpoints++;
403 }
404 if (!CtorAllocate()) return;
405
406 Double_t xmin = f->GetXmin();
407 Double_t xmax = f->GetXmax();
408 Double_t dx = (xmax - xmin) / fNpoints;
409 Double_t integ = 0;
410 Int_t i;
411 for (i = 0; i < fNpoints; i++) {
412 if (coption == 'i' || coption == 'I') {
413 fX[i] = xmin + i * dx;
414 if (i == 0) fY[i] = 0;
415 else fY[i] = integ + ((TF1*)f)->Integral(fX[i] - dx, fX[i]);
416 integ = fY[i];
417 } else if (coption == 'd' || coption == 'D') {
418 fX[i] = xmin + (i + 0.5) * dx;
419 fY[i] = ((TF1*)f)->Derivative(fX[i]);
420 } else {
421 fX[i] = xmin + (i + 0.5) * dx;
422 fY[i] = ((TF1*)f)->Eval(fX[i]);
423 }
424 }
425 if (integ != 0 && coption == 'I') {
426 for (i = 1; i < fNpoints; i++) fY[i] /= integ;
427 }
428
429 f->TAttLine::Copy(*this);
430 f->TAttFill::Copy(*this);
431 f->TAttMarker::Copy(*this);
432
433 SetName(f->GetName());
434 SetTitle(f->GetTitle());
435}
436
437////////////////////////////////////////////////////////////////////////////////
438/// Graph constructor reading input from filename.
439///
440/// `filename` is assumed to contain at least two columns of numbers.
441/// The string format is by default `"%lg %lg"`.
442/// This is a standard c formatting for `scanf()`.
443/// For example, set format to `"%lg,%lg"` for a comma-separated file.
444/// If format string is empty, suitable value will be provided based on file extension
445///
446/// If columns of numbers should be skipped, a `"%*lg"` or `"%*s"` for each column
447/// can be added, e.g. `"%lg %*lg %lg"` would read x-values from the first and
448/// y-values from the third column.
449///
450/// For files separated by a specific delimiter different from ' ' and '\\t' (e.g.
451/// ';' in csv files) you can avoid using `%*s` to bypass this delimiter by explicitly
452/// specify the `option` argument,
453/// e.g. option=`" \\t,;"` for columns of figures separated by any of these characters
454/// (' ', '\\t', ',', ';')
455/// used once (e.g. `"1;1"`) or in a combined way (`" 1;,;; 1"`).
456/// Note in that case, the instantiation is about two times slower.
457
458TGraph::TGraph(const char *filename, const char *format, Option_t *option)
459 : TNamed("Graph", filename), TAttFill(0, 1000)
460{
461 Double_t x, y;
464
465 std::ifstream infile(fname.Data());
466 if (!infile.good()) {
467 MakeZombie();
468 Error("TGraph", "Cannot open file: %s, TGraph is Zombie", filename);
469 fNpoints = 0;
470 return;
471 } else {
472 fNpoints = 100; //initial number of points
473 }
474 if (!CtorAllocate()) return;
475 std::string line;
476 Int_t np = 0;
477
479
480 if (!option || !*option) { // No delimiters specified (standard constructor).
481 // is empty format string specified - try to guess format from the file extension
482 if (format_.IsNull()) {
483 if (fname.EndsWith(".txt", TString::kIgnoreCase))
484 format_ = "%lg %lg";
485 else if (fname.EndsWith(".tsv", TString::kIgnoreCase))
486 format_ = "%lg\t%lg";
487 else
488 format_ = "%lg,%lg";
489 }
490
491 while (std::getline(infile, line, '\n')) {
492 if (2 != sscanf(line.c_str(), format_.Data(), &x, &y)) {
493 continue; //skip empty and ill-formed lines
494 }
495 SetPoint(np, x, y);
496 np++;
497 }
498 Set(np);
499
500 // A delimiter has been specified in "option"
501 } else {
502
503 // Checking format and creating its boolean counterpart
504 format_.ReplaceAll(" ", "") ;
505 format_.ReplaceAll("\t", "") ;
506 format_.ReplaceAll("lg", "") ;
507 format_.ReplaceAll("s", "") ;
508 format_.ReplaceAll("%*", "0") ;
509 format_.ReplaceAll("%", "1") ;
510 if (!format_.IsDigit()) {
511 Error("TGraph", "Incorrect input format! Allowed formats are {\"%%lg\",\"%%*lg\" or \"%%*s\"}");
512 return;
513 }
514 Int_t ntokens = format_.Length() ;
515 if (ntokens < 2) {
516 Error("TGraph", "Incorrect input format! Only %d tag(s) in format whereas 2 \"%%lg\" tags are expected!", ntokens);
517 return;
518 }
521 for (Int_t idx = 0; idx < ntokens; idx++) {
522 isTokenToBeSaved[idx] = TString::Format("%c", format_[idx]).Atoi() ; //atoi(&format_[idx]) does not work for some reason...
523 if (isTokenToBeSaved[idx] == 1) {
525 }
526 }
527 if (ntokens >= 2 && ntokensToBeSaved != 2) { //first condition not to repeat the previous error message
528 Error("TGraph", "Incorrect input format! There are %d \"%%lg\" tag(s) in format whereas 2 and only 2 are expected!", ntokensToBeSaved);
529 delete [] isTokenToBeSaved ;
530 return;
531 }
532
533 // Initializing loop variables
534 Bool_t isLineToBeSkipped = kFALSE ; //empty and ill-formed lines
535 char * token = nullptr ;
536 TString token_str = "" ;
537 Int_t token_idx = 0 ;
538 Double_t * value = new Double_t [2] ; //x,y buffers
539 Int_t value_idx = 0 ;
540
541 // Looping
542 char *rest;
543 while (std::getline(infile, line, '\n')) {
544 if (!line.empty()) {
545 if (line[line.size() - 1] == char(13)) { // removing DOS CR character
546 line.erase(line.end() - 1, line.end()) ;
547 }
548 //token = R__STRTOK_R(const_cast<char *>(line.c_str()), option, rest);
549 token = R__STRTOK_R(const_cast<char *>(line.c_str()), option, &rest);
550 while (token != nullptr && value_idx < 2) {
552 token_str = TString(token) ;
553 token_str.ReplaceAll("\t", "") ;
554 if (!token_str.IsFloat()) {
556 break ;
557 } else {
558 value[value_idx] = token_str.Atof() ;
559 value_idx++ ;
560 }
561 }
562 token = R__STRTOK_R(nullptr, option, &rest); // next token
563 token_idx++ ;
564 }
565 if (!isLineToBeSkipped && value_idx == 2) {
566 x = value[0] ;
567 y = value[1] ;
568 SetPoint(np, x, y) ;
569 np++ ;
570 }
571 }
573 token = nullptr ;
574 token_idx = 0 ;
575 value_idx = 0 ;
576 }
577 Set(np) ;
578
579 // Cleaning
580 delete [] isTokenToBeSaved ;
581 delete [] value ;
582 delete token ;
583 }
584 infile.close();
585 if (fNpoints == 0) {
586 Warning("TGraph", "No points were found in file %s with the specified input format %s", filename, format);
587 return;
588 }
589}
590
591////////////////////////////////////////////////////////////////////////////////
592/// Graph default destructor.
593
595{
596 delete [] fX;
597 delete [] fY;
598 if (fFunctions) {
600 //special logic to support the case where the same object is
601 //added multiple times in fFunctions.
602 //This case happens when the same object is added with different
603 //drawing modes
604 TObject *obj;
605 while ((obj = fFunctions->First())) {
606 while (fFunctions->Remove(obj)) { }
607 delete obj;
608 }
609 delete fFunctions;
610 fFunctions = nullptr; //to avoid accessing a deleted object in RecursiveRemove
611 }
612 delete fHistogram;
613}
614
615////////////////////////////////////////////////////////////////////////////////
616/// Allocate internal data structures for `newsize` points.
617
622
623////////////////////////////////////////////////////////////////////////////////
624/// Allocate arrays.
625
627{
628 if (arraySize < 0) {
629 arraySize = 0;
630 }
632 if (!arraySize) {
633 for (Int_t i = 0; i < Narrays; ++i)
634 newarrays[i] = nullptr;
635 } else {
636 for (Int_t i = 0; i < Narrays; ++i)
637 newarrays[i] = new Double_t[arraySize];
638 }
639 fMaxSize = arraySize;
640 return newarrays;
641}
642
643////////////////////////////////////////////////////////////////////////////////
644/// Performs the operation: `y = y + c1*f(x,y)`
645/// Errors are not recalculated.
646///
647/// \param f may be a 1-D function TF1 or 2-d function TF2
648/// \param c1 a scaling factor, 1 by default
649
651{
653
654 for (Int_t i = 0; i < fNpoints; i++) {
655 fY[i] += c1*f->Eval(fX[i], fY[i]);
656 }
657 if (gPad) gPad->Modified();
658}
659
660////////////////////////////////////////////////////////////////////////////////
661/// Apply function f to all the data points
662/// f may be a 1-D function TF1 or 2-d function TF2
663/// The Y values of the graph are replaced by the new values computed
664/// using the function
665
667{
669
670 for (Int_t i = 0; i < fNpoints; i++) {
671 fY[i] = f->Eval(fX[i], fY[i]);
672 }
673 if (gPad) gPad->Modified();
674}
675
676////////////////////////////////////////////////////////////////////////////////
677/// Browse
678
680{
681 TString opt = gEnv->GetValue("TGraph.BrowseOption", "");
682 if (opt.IsNull()) {
683 opt = b ? b->GetDrawOption() : "alp";
684 opt = (opt == "") ? "alp" : opt.Data();
685 }
686 Draw(opt.Data());
687 gPad->Update();
688}
689
690////////////////////////////////////////////////////////////////////////////////
691/// Return the chisquare of this graph with respect to f1.
692/// The chisquare is computed as the sum of the quantity below at each point:
693/// \f[
694/// \frac{(y-f1(x))^{2}}{ey^{2}+(\frac{1}{2}(exl+exh)f1'(x))^{2}}
695/// \f]
696/// where x and y are the graph point coordinates and f1'(x) is the derivative of function f1(x).
697/// This method to approximate the uncertainty in y because of the errors in x, is called
698/// "effective variance" method.
699/// In case of a pure TGraph, the denominator is 1.
700/// In case of a TGraphErrors or TGraphAsymmErrors the errors are taken
701/// into account.
702/// By default the range of the graph is used whatever function range.
703/// Use option "R" to use the function range
704
706{
707 if (!func) {
708 Error("Chisquare","Function pointer is Null - return -1");
709 return -1;
710 }
711
712 TString opt(option); opt.ToUpper();
713 bool useRange = opt.Contains("R");
714
715 return ROOT::Fit::Chisquare(*this, *func,useRange);
716}
717
718////////////////////////////////////////////////////////////////////////////////
719/// Return kTRUE if point number "left"'s argument (angle with respect to positive
720/// x-axis) is bigger than that of point number "right". Can be used by Sort.
721
723{
724 Double_t xl = 0, yl = 0, xr = 0, yr = 0;
725 gr->GetPoint(left, xl, yl);
726 gr->GetPoint(right, xr, yr);
727 return (TMath::ATan2(yl, xl) > TMath::ATan2(yr, xr));
728}
729
730////////////////////////////////////////////////////////////////////////////////
731/// Return kTRUE if fX[left] > fX[right]. Can be used by Sort.
732
734{
735 return gr->fX[left] > gr->fX[right];
736}
737
738////////////////////////////////////////////////////////////////////////////////
739/// Return kTRUE if fY[left] > fY[right]. Can be used by Sort.
740
742{
743 return gr->fY[left] > gr->fY[right];
744}
745
746////////////////////////////////////////////////////////////////////////////////
747/// Return kTRUE if point number "left"'s distance to origin is bigger than
748/// that of point number "right". Can be used by Sort.
749
751{
752 return gr->fX[left] * gr->fX[left] + gr->fY[left] * gr->fY[left]
753 > gr->fX[right] * gr->fX[right] + gr->fY[right] * gr->fY[right];
754}
755
756////////////////////////////////////////////////////////////////////////////////
757/// Compute the x/y range of the points in this graph
758
760{
761 if (fNpoints <= 0) {
762 xmin = xmax = ymin = ymax = 0;
763 return;
764 }
765 xmin = xmax = fX[0];
766 ymin = ymax = fY[0];
767
768 Double_t xminl = 0; // Positive minimum. Used in case of log scale along X axis.
769 Double_t yminl = 0; // Positive minimum. Used in case of log scale along Y axis.
770
771 for (Int_t i = 1; i < fNpoints; i++) {
772 if (fX[i] < xmin) xmin = fX[i];
773 if (fX[i] > xmax) xmax = fX[i];
774 if (fY[i] < ymin) ymin = fY[i];
775 if (fY[i] > ymax) ymax = fY[i];
776 if (ymin>0 && (yminl==0 || ymin<yminl)) yminl = ymin;
777 if (xmin>0 && (xminl==0 || xmin<xminl)) xminl = xmin;
778 }
779
780 if (gPad && gPad->GetLogy() && yminl>0) ymin = yminl;
781 if (gPad && gPad->GetLogx() && xminl>0) xmin = xminl;
782}
783
784////////////////////////////////////////////////////////////////////////////////
785/// Copy points from fX and fY to arrays[0] and arrays[1]
786/// or to fX and fY if arrays == 0 and ibegin != iend.
787/// If newarrays is non null, replace fX, fY with pointers from newarrays[0,1].
788/// Delete newarrays, old fX and fY
789
792{
794 if (newarrays) {
795 delete[] fX;
796 fX = newarrays[0];
797 delete[] fY;
798 fY = newarrays[1];
799 delete[] newarrays;
800 }
801}
802
803////////////////////////////////////////////////////////////////////////////////
804/// Copy points from fX and fY to arrays[0] and arrays[1]
805/// or to fX and fY if arrays == 0 and ibegin != iend.
806
809{
810 if (ibegin < 0 || iend <= ibegin || obegin < 0) { // Error;
811 return kFALSE;
812 }
813 if (!arrays && ibegin == obegin) { // No copying is needed
814 return kFALSE;
815 }
816 Int_t n = (iend - ibegin) * sizeof(Double_t);
817 if (arrays) {
818 memmove(&arrays[0][obegin], &fX[ibegin], n);
819 memmove(&arrays[1][obegin], &fY[ibegin], n);
820 } else {
821 memmove(&fX[obegin], &fX[ibegin], n);
822 memmove(&fY[obegin], &fY[ibegin], n);
823 }
824 return kTRUE;
825}
826
827////////////////////////////////////////////////////////////////////////////////
828/// In constructors set fNpoints than call this method.
829/// Return kFALSE if the graph will contain no points.
830///Note: This function should be called only from the constructor
831/// since it does not delete previously existing arrays
832
834{
835 fHistogram = nullptr;
836 fMaximum = -1111;
837 fMinimum = -1111;
839 fFunctions = new TList;
840 if (fNpoints <= 0) {
841 fNpoints = 0;
842 fMaxSize = 0;
843 fX = nullptr;
844 fY = nullptr;
845 return kFALSE;
846 } else {
848 fX = new Double_t[fMaxSize];
849 fY = new Double_t[fMaxSize];
850 }
851 return kTRUE;
852}
853
854////////////////////////////////////////////////////////////////////////////////
855/// Draw this graph with its current attributes.
856///
857/// The options to draw a graph are described in TGraphPainter class.
858
860{
861 TString opt = option;
862 opt.ToLower();
863
864 if (opt.Contains("same")) {
865 opt.ReplaceAll("same", "");
866 }
867
868 // in case of option *, set marker style to 3 (star) and replace
869 // * option by option P.
870 Ssiz_t pos;
871 if ((pos = opt.Index("*")) != kNPOS) {
873 opt.Replace(pos, 1, "p");
874 }
875
876 // If no option is specified, it is defined as "alp" in case there is
877 // no current pad or if the current pad has no axis defined and if there is
878 // no default option set using TGraph::SetOption. If fOption is set using
879 // TGraph::SetOption, it is used as default option.
880 if ((!option || !strlen(option))) {
881 Option_t *topt = (!fOption.IsNull()) ? fOption.Data() : "alp";
882 if (gPad) {
883 if (!gPad->GetListOfPrimitives()->FindObject("TFrame"))
884 opt = topt;
885 } else {
886 opt = topt;
887 }
888 }
889
890 if (gPad) {
891 if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
892 if (opt.Contains("a")) gPad->Clear();
893 }
894
895 AppendPad(opt);
896
897 gPad->IncrementPaletteColor(1, opt);
898
899}
900
901////////////////////////////////////////////////////////////////////////////////
902/// Compute distance from point px,py to a graph.
903///
904/// Compute the closest distance of approach from point px,py to this line.
905/// The distance is computed in pixels units.
906
908{
910 if (painter) return painter->DistancetoPrimitiveHelper(this, px, py);
911 else return 0;
912}
913
914////////////////////////////////////////////////////////////////////////////////
915/// Draw this graph with new attributes.
916
918{
919 TGraph *newgraph = new TGraph(n, x, y);
923 newgraph->SetBit(kCanDelete);
924 newgraph->AppendPad(option);
925}
926
927////////////////////////////////////////////////////////////////////////////////
928/// Draw this graph with new attributes.
929
931{
932 TGraph *newgraph = new TGraph(n, x, y);
936 newgraph->SetBit(kCanDelete);
937 newgraph->AppendPad(option);
938}
939
940////////////////////////////////////////////////////////////////////////////////
941/// Draw this graph with new attributes.
942
944{
945 const Double_t *xx = x;
946 const Double_t *yy = y;
947 if (!xx) xx = fX;
948 if (!yy) yy = fY;
949 TGraph *newgraph = new TGraph(n, xx, yy);
953 newgraph->SetBit(kCanDelete);
954 newgraph->AppendPad(option);
955}
956
957////////////////////////////////////////////////////////////////////////////////
958/// Display a panel with all graph drawing options.
959
961{
963 if (painter) painter->DrawPanelHelper(this);
964}
965
966////////////////////////////////////////////////////////////////////////////////
967/// Interpolate points in this graph at x using a TSpline.
968///
969/// - if spline==0 and option="" a linear interpolation between the two points
970/// close to x is computed. If x is outside the graph range, a linear
971/// extrapolation is computed.
972/// - if spline==0 and option="S" a TSpline3 object is created using this graph
973/// and the interpolated value from the spline is returned.
974/// the internally created spline is deleted on return.
975/// - if spline is specified, it is used to return the interpolated value.
976///
977/// If the points are sorted in X a binary search is used (significantly faster)
978/// One needs to set the bit TGraph::SetBit(TGraph::kIsSortedX) before calling
979/// TGraph::Eval to indicate that the graph is sorted in X.
980
982{
983
984 if (spline) {
985 //spline interpolation using the input spline
986 return spline->Eval(x);
987 }
988
989 if (fNpoints == 0) return 0;
990 if (fNpoints == 1) return fY[0];
991
992 if (option && *option) {
993 TString opt = option;
994 opt.ToLower();
995 // create a TSpline every time when using option "s" and no spline pointer is given
996 if (opt.Contains("s")) {
997
998 // points must be sorted before using a TSpline
999 std::vector<Double_t> xsort(fNpoints);
1000 std::vector<Double_t> ysort(fNpoints);
1001 std::vector<Int_t> indxsort(fNpoints);
1002 TMath::Sort(fNpoints, fX, &indxsort[0], false);
1003 for (Int_t i = 0; i < fNpoints; ++i) {
1004 xsort[i] = fX[ indxsort[i] ];
1005 ysort[i] = fY[ indxsort[i] ];
1006 }
1007
1008 // spline interpolation creating a new spline
1009 TSpline3 s("", &xsort[0], &ysort[0], fNpoints);
1010 Double_t result = s.Eval(x);
1011 return result;
1012 }
1013 }
1014 //linear interpolation
1015 //In case x is < fX[0] or > fX[fNpoints-1] return the extrapolated point
1016
1017 //find points in graph around x assuming points are not sorted
1018 // (if point are sorted use a binary search)
1019 Int_t low = -1;
1020 Int_t up = -1;
1023 if (low == -1) {
1024 // use first two points for doing an extrapolation
1025 low = 0;
1026 }
1027 if (fX[low] == x) return fY[low];
1028 if (low == fNpoints-1) low--; // for extrapolating
1029 up = low+1;
1030 }
1031 else {
1032 // case TGraph is not sorted
1033
1034 // find neighbours simply looping all points
1035 // and find also the 2 adjacent points: (low2 < low < x < up < up2 )
1036 // needed in case x is outside the graph ascissa interval
1037 Int_t low2 = -1;
1038 Int_t up2 = -1;
1039
1040 for (Int_t i = 0; i < fNpoints; ++i) {
1041 if (fX[i] < x) {
1042 if (low == -1 || fX[i] > fX[low]) {
1043 low2 = low;
1044 low = i;
1045 } else if (low2 == -1) low2 = i;
1046 } else if (fX[i] > x) {
1047 if (up == -1 || fX[i] < fX[up]) {
1048 up2 = up;
1049 up = i;
1050 } else if (up2 == -1) up2 = i;
1051 } else // case x == fX[i]
1052 return fY[i]; // no interpolation needed
1053 }
1054
1055 // treat cases when x is outside graph min max abscissa
1056 if (up == -1) {
1057 up = low;
1058 low = low2;
1059 }
1060 if (low == -1) {
1061 low = up;
1062 up = up2;
1063 }
1064 }
1065 // do now the linear interpolation
1066 assert(low != -1 && up != -1);
1067
1068 if (fX[low] == fX[up]) return fY[low];
1069 Double_t yn = fY[up] + (x - fX[up]) * (fY[low] - fY[up]) / (fX[low] - fX[up]);
1070 return yn;
1071}
1072
1073////////////////////////////////////////////////////////////////////////////////
1074/// Execute action corresponding to one event.
1075///
1076/// This member function is called when a graph is clicked with the locator
1077///
1078/// If Left button clicked on one of the line end points, this point
1079/// follows the cursor until button is released.
1080///
1081/// if Middle button clicked, the line is moved parallel to itself
1082/// until the button is released.
1083
1085{
1087 if (painter) painter->ExecuteEventHelper(this, event, px, py);
1088}
1089
1090////////////////////////////////////////////////////////////////////////////////
1091/// If array sizes <= newsize, expand storage to 2*newsize.
1092
1094{
1096 CopyAndRelease(ps, 0, 0, 0);
1097}
1098
1099////////////////////////////////////////////////////////////////////////////////
1100/// If graph capacity is less than newsize points then make array sizes
1101/// equal to least multiple of step to contain newsize points.
1102
1104{
1105 if (newsize <= fMaxSize) {
1106 return;
1107 }
1108 Double_t **ps = Allocate(step * (newsize / step + (newsize % step ? 1 : 0)));
1109 CopyAndRelease(ps, 0, fNpoints, 0);
1110}
1111
1112////////////////////////////////////////////////////////////////////////////////
1113/// if size > fMaxSize allocate new arrays of 2*size points and copy iend first
1114/// points.
1115/// Return pointer to new arrays.
1116
1118{
1119 if (size <= fMaxSize)
1120 return nullptr;
1122 CopyPoints(newarrays, 0, iend, 0);
1123 return newarrays;
1124}
1125
1126////////////////////////////////////////////////////////////////////////////////
1127/// Set zero values for point arrays in the range [begin, end)
1128/// Should be redefined in descendant classes
1129
1131{
1132 memset(fX + begin, 0, (end - begin)*sizeof(Double_t));
1133 memset(fY + begin, 0, (end - begin)*sizeof(Double_t));
1134}
1135
1136////////////////////////////////////////////////////////////////////////////////
1137/// Search object named name in the list of functions
1138
1140{
1141 return fFunctions ? fFunctions->FindObject(name) : nullptr;
1142}
1143
1144////////////////////////////////////////////////////////////////////////////////
1145/// Search object obj in the list of functions
1146
1148{
1149 return fFunctions ? fFunctions->FindObject(obj) : nullptr;
1150}
1151
1152////////////////////////////////////////////////////////////////////////////////
1153/// Fit this graph with function f1.
1154///
1155/// \param[in] f1 pointer to the function object
1156/// \param[in] option string defining the fit options (see table below).
1157/// \param[in] goption specify a list of graphics options. See TGraph::Draw and TGraphPainter for a complete list of these possible options.
1158/// \param[in] rxmin lower fitting range
1159/// \param[in] rxmax upper fitting range
1160///
1161/// \anchor GFitOpt
1162/// ### Graph Fitting Options
1163/// The list of fit options is given in parameter option.
1164///
1165/// option | description
1166/// -------|------------
1167/// "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`.
1168/// "W" | Ignore all point errors when fitting a TGraphErrors or TGraphAsymmErrors
1169/// "F" | Uses the default minimizer (e.g. Minuit) when fitting a linear function (e.g. polN) instead of the linear fitter.
1170/// "U" | Uses a user specified objective function (e.g. user providedlikelihood function) defined using `TVirtualFitter::SetFCN`
1171/// "E" | Performs a better parameter errors estimation using the Minos technique for all fit parameters.
1172/// "M" | Uses the IMPROVE algorithm (available only in TMinuit). This algorithm attempts improve the found local minimum by searching for a better one.
1173/// "Q" | Quiet mode (minimum printing)
1174/// "V" | Verbose mode (default is between Q and V)
1175/// "+" | 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.
1176/// "N" | Does not store the graphics function, does not draw the histogram with the function after fitting.
1177/// "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.
1178/// "R" | Fit using a fitting range specified in the function range with `TF1::SetRange`.
1179/// "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.
1180/// "C" | In case of linear fitting, do no calculate the chisquare (saves CPU time).
1181/// "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.
1182/// "EX0" | When fitting a TGraphErrors or TGraphAsymErrors do not consider errors in the X coordinates
1183/// "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
1184///
1185///
1186/// This function is used for fitting also the derived TGraph classes such as TGraphErrors or TGraphAsymmErrors.
1187/// See the note below on how the errors are used when fitting a TGraphErrors or TGraphAsymmErrors.
1188///
1189/// The fitting of the TGraph, i.e simple data points without any error associated, is performed using the
1190/// un-weighted least-square (chi-square) method.
1191///
1192///
1193///\anchor GFitErrors
1194/// ### TGraphErrors fit:
1195///
1196/// In case of a TGraphErrors or TGraphAsymmErrors object, when `x` errors are present, the error along x,
1197/// is projected along the y-direction by calculating the function at the points `x-ex_low` and
1198/// `x+ex_high`, where `ex_low` and `ex_high` are the corresponding lower and upper error in x.
1199/// The chi-square is then computed as the sum of the quantity below at each data point:
1200///
1201/// \f[
1202/// \frac{(y-f(x))^{2}}{ey^{2}+(\frac{1}{2}(exl+exh)f'(x))^{2}}
1203/// \f]
1204///
1205/// where `x` and `y` are the point coordinates, and `f'(x)` is the derivative of the
1206/// function `f(x)`.
1207///
1208/// In case of asymmetric errors, if the function lies below (above) the data point, `ey` is `ey_low` (`ey_high`).
1209///
1210/// The approach used to approximate the uncertainty in y because of the
1211/// errors in x is to make it equal the error in x times the slope of the line.
1212/// This approach is called "effective variance method" and
1213/// the implementation is provided in the function FitUtil::EvaluateChi2Effective
1214///
1215/// \anchor GFitLinear
1216/// ### Linear fitting:
1217/// When the fitting function is linear (contains the `++` sign) or the fitting
1218/// function is a polynomial, a linear fitter is initialised.
1219/// To create a linear function, use the following syntax: linear parts
1220/// separated by `++` sign.
1221/// Example: to fit the parameters of the function `p0*x + p1*sin(x)`, you can create a
1222/// TF1 object as
1223///
1224/// TF1 *f1 = new TF1("f1", "x++sin(x)", xmin, xmax);
1225///
1226/// For such a TF1 you don't have to set the initial conditions and the linear fitter is used.
1227/// Going via the linear fitter for functions, linear in parameters, gives a
1228/// considerable advantage in speed.
1229/// When using the linear fitting it is also possible to perform a robust fitting with the
1230/// Least Trimmed Square (LTS) regression algorithm, by using the fit option `ROB`.
1231/// See the tutorial `fitLinearRobust.C`.
1232///
1233/// ### Notes on TGraph/TGraphErrors Fitting:
1234///
1235/// 1. By using the "effective variance" method a simple linear regression
1236/// becomes a non-linear case, which takes several iterations
1237/// instead of 0 as in the linear case.
1238/// 2. The effective variance technique assumes that there is no correlation
1239/// between the x and y coordinate.
1240/// 3. The standard chi2 (least square) method without error in the coordinates (x) can
1241/// be forced by using option "EX0"
1242/// 4. The linear fitter doesn't take into account the errors in x. When fitting a
1243/// TGraphErrors with a linear functions the errors in x will not be considered.
1244/// If errors in x are important, use option "F" for linear function fitting.
1245/// 5. When fitting a TGraph (i.e. no errors associated with each point),
1246/// a correction is applied to the errors on the parameters with the following
1247/// formula:
1248/// `parameter_error *= sqrt(chisquare/(ndf-1))`
1249///
1250/// ### General Fitting documentation
1251///
1252/// See in TH1::Fit for the documentation of
1253/// - [Fit Result](\ref HFitRes)
1254/// - [Fit Status](\ref HFitStatus)
1255/// - [Fit Statistics Box](\ref HFitStatBox)
1256/// - [Fitting in a Range](\ref HFitRange)
1257/// - [Setting Initial Conditions](\ref HFitInitial)
1258
1268
1269////////////////////////////////////////////////////////////////////////////////
1270/// Fit this graph with function with name `fname`.
1271///
1272/// This is a different interface to TGraph fitting using TGraph::Fit(TF1 *f1,Option_t *, Option_t *, Axis_t, Axis_t)
1273/// See there for the details about fitting a TGraph.
1274///
1275/// The parameter `fname` is the name of an already predefined function created by TF1 or TF2
1276/// Predefined functions such as gaus, expo and poln are automatically
1277/// created by ROOT.
1278///
1279/// The parameter `fname` can also be a formula, accepted by the linear fitter (linear parts divided
1280/// by "++" sign), for example "x++sin(x)" for fitting "[0]*x+[1]*sin(x)"
1281
1283{
1284 const char *linear = fname ? strstr(fname, "++") : nullptr;
1285 if (linear) {
1286 TF1 f1(fname, fname, xmin, xmax);
1287 return Fit(&f1, option, "", xmin, xmax);
1288 }
1289 TF1 * f1 = (TF1*)gROOT->GetFunction(fname);
1290 if (!f1) {
1291 Printf("Unknown function: %s", fname);
1292 return -1;
1293 }
1294 return Fit(f1, option, "", xmin, xmax);
1295}
1296
1297////////////////////////////////////////////////////////////////////////////////
1298/// Display a GUI panel with all graph fit options.
1299///
1300/// See class TFitEditor for example
1301
1303{
1304 if (!gPad)
1305 gROOT->MakeDefCanvas();
1306
1307 if (!gPad) {
1308 Error("FitPanel", "Unable to create a default canvas");
1309 return;
1310 }
1311
1312 // use plugin manager to create instance of TFitEditor
1313 TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
1314 if (handler && handler->LoadPlugin() != -1) {
1315 if (handler->ExecPlugin(2, gPad, this) == 0)
1316 Error("FitPanel", "Unable to crate the FitPanel");
1317 } else
1318 Error("FitPanel", "Unable to find the FitPanel plug-in");
1319}
1320
1321////////////////////////////////////////////////////////////////////////////////
1322/// Return graph correlation factor
1323
1325{
1326 Double_t rms1 = GetRMS(1);
1327 if (rms1 == 0) return 0;
1328 Double_t rms2 = GetRMS(2);
1329 if (rms2 == 0) return 0;
1330 return GetCovariance() / rms1 / rms2;
1331}
1332
1333////////////////////////////////////////////////////////////////////////////////
1334/// Return covariance of vectors x,y
1335
1337{
1338 if (fNpoints <= 0) return 0;
1339 Double_t sum = fNpoints, sumx = 0, sumy = 0, sumxy = 0;
1340
1341 for (Int_t i = 0; i < fNpoints; i++) {
1342 sumx += fX[i];
1343 sumy += fY[i];
1344 sumxy += fX[i] * fY[i];
1345 }
1346 return sumxy / sum - sumx / sum * sumy / sum;
1347}
1348
1349////////////////////////////////////////////////////////////////////////////////
1350/// Return mean value of X (axis=1) or Y (axis=2)
1351
1353{
1354 if (axis < 1 || axis > 2) return 0;
1355 if (fNpoints <= 0) return 0;
1356 Double_t sumx = 0;
1357 for (Int_t i = 0; i < fNpoints; i++) {
1358 if (axis == 1) sumx += fX[i];
1359 else sumx += fY[i];
1360 }
1361 return sumx / fNpoints;
1362}
1363
1364////////////////////////////////////////////////////////////////////////////////
1365/// Return RMS of X (axis=1) or Y (axis=2)
1366
1368{
1369 if (axis < 1 || axis > 2) return 0;
1370 if (fNpoints <= 0) return 0;
1371 Double_t sumx = 0, sumx2 = 0;
1372 for (Int_t i = 0; i < fNpoints; i++) {
1373 if (axis == 1) {
1374 sumx += fX[i];
1375 sumx2 += fX[i] * fX[i];
1376 } else {
1377 sumx += fY[i];
1378 sumx2 += fY[i] * fY[i];
1379 }
1380 }
1381 Double_t x = sumx / fNpoints;
1383 return TMath::Sqrt(rms2);
1384}
1385
1386////////////////////////////////////////////////////////////////////////////////
1387/// It always returns a negative value. Real implementation in TGraphErrors
1388
1390{
1391 return -1;
1392}
1393
1394////////////////////////////////////////////////////////////////////////////////
1395/// It always returns a negative value. Real implementation in TGraphErrors
1396
1398{
1399 return -1;
1400}
1401
1402////////////////////////////////////////////////////////////////////////////////
1403/// It always returns a negative value. Real implementation in TGraphErrors
1404/// and TGraphAsymmErrors
1405
1407{
1408 return -1;
1409}
1410
1411////////////////////////////////////////////////////////////////////////////////
1412/// It always returns a negative value. Real implementation in TGraphErrors
1413/// and TGraphAsymmErrors
1414
1416{
1417 return -1;
1418}
1419
1420////////////////////////////////////////////////////////////////////////////////
1421/// It always returns a negative value. Real implementation in TGraphErrors
1422/// and TGraphAsymmErrors
1423
1425{
1426 return -1;
1427}
1428
1429////////////////////////////////////////////////////////////////////////////////
1430/// It always returns a negative value. Real implementation in TGraphErrors
1431/// and TGraphAsymmErrors
1432
1434{
1435 return -1;
1436}
1437
1438////////////////////////////////////////////////////////////////////////////////
1439/// Return pointer to function with name.
1440///
1441/// Functions such as TGraph::Fit store the fitted function in the list of
1442/// functions of this graph.
1443
1444TF1 *TGraph::GetFunction(const char *name) const
1445{
1446 return dynamic_cast<TF1*>(FindObject(name));
1447}
1448
1449////////////////////////////////////////////////////////////////////////////////
1450/// Returns a pointer to the histogram used to draw the axis
1451/// Takes into account the two following cases.
1452/// 1. option 'A' was specified in TGraph::Draw. Return fHistogram
1453/// 2. user had called TPad::DrawFrame. return pointer to hframe histogram
1454
1456{
1459
1460 ComputeRange(rwxmin, rwymin, rwxmax, rwymax); //this is redefined in TGraphErrors
1461
1462 // (if fHistogram exist) && (if the log scale is on) &&
1463 // (if the computed range minimum is > 0) && (if the fHistogram minimum is zero)
1464 // then it means fHistogram limits have been computed in linear scale
1465 // therefore they might be too strict and cut some points. In that case the
1466 // fHistogram limits should be recomputed ie: the existing fHistogram
1467 // should not be returned.
1468 TH1F *histogr = nullptr;
1469 if (fHistogram) {
1470 if (!TestBit(kResetHisto)) {
1471 if (gPad && gPad->GetLogx()) {
1472 if (rwxmin <= 0 || fHistogram->GetXaxis()->GetXmin() != 0) return fHistogram;
1473 } else if (gPad && gPad->GetLogy()) {
1475 } else {
1476 return fHistogram;
1477 }
1478 } else {
1479 const_cast <TGraph*>(this)->ResetBit(kResetHisto);
1480 }
1482 }
1483
1484 if (rwxmin == rwxmax) rwxmax += 1.;
1485 if (rwymin == rwymax) rwymax += 1.;
1486 dx = 0.1 * (rwxmax - rwxmin);
1487 dy = 0.1 * (rwymax - rwymin);
1488 uxmin = rwxmin - dx;
1489 uxmax = rwxmax + dx;
1490 minimum = rwymin - dy;
1491 maximum = rwymax + dy;
1492
1493 if (fMinimum != -1111) minimum = fMinimum;
1494 if (fMaximum != -1111) maximum = fMaximum;
1495
1496 // the graph is created with at least as many channels as there are points
1497 // to permit zooming on the full range
1498 if (uxmin < 0 && rwxmin >= 0) {
1499 if (gPad && gPad->GetLogx()) uxmin = 0.9 * rwxmin;
1500 else uxmin = 0;
1501 }
1502 if (uxmax > 0 && rwxmax <= 0) {
1503 if (gPad && gPad->GetLogx()) uxmax = 1.1 * rwxmax;
1504 else uxmax = 0;
1505 }
1506
1507 if (minimum < 0 && rwymin >= 0) minimum = 0.9 * rwymin;
1508
1509 if (minimum <= 0 && gPad && gPad->GetLogy()) minimum = 0.001 * maximum;
1510 if (uxmin <= 0 && gPad && gPad->GetLogx()) {
1511 if (uxmax > 1000) uxmin = 1;
1512 else uxmin = 0.001 * uxmax;
1513 }
1514
1515 rwxmin = uxmin;
1516 rwxmax = uxmax;
1517 Int_t npt = 100;
1518 if (fNpoints > npt) npt = fNpoints;
1519 const char *gname = GetName();
1520 if (!gname[0]) gname = "Graph";
1521 // do not add the histogram to gDirectory
1522 // use local TDirectory::TContext that will set temporarly gDirectory to a nullptr and
1523 // will avoid that histogram is added in the global directory
1524 {
1525 TDirectory::TContext ctx(nullptr);
1526 ((TGraph*)this)->fHistogram = new TH1F(gname, GetTitle(), npt, rwxmin, rwxmax);
1527 }
1528 if (!fHistogram) return nullptr;
1533 // Restore the axis attributes if needed
1534 if (histogr) {
1535 fHistogram->GetXaxis()->SetTitle(histogr->GetXaxis()->GetTitle());
1536 fHistogram->GetXaxis()->CenterTitle(histogr->GetXaxis()->GetCenterTitle());
1537 fHistogram->GetXaxis()->RotateTitle(histogr->GetXaxis()->GetRotateTitle());
1538 fHistogram->GetXaxis()->SetNoExponent(histogr->GetXaxis()->GetNoExponent());
1539 fHistogram->GetXaxis()->SetTimeDisplay(histogr->GetXaxis()->GetTimeDisplay());
1540 fHistogram->GetXaxis()->SetTimeFormat(histogr->GetXaxis()->GetTimeFormat());
1541 histogr->GetXaxis()->TAttAxis::Copy(*(fHistogram->GetXaxis()));
1542
1543 fHistogram->GetYaxis()->SetTitle(histogr->GetYaxis()->GetTitle());
1544 fHistogram->GetYaxis()->CenterTitle(histogr->GetYaxis()->GetCenterTitle());
1545 fHistogram->GetYaxis()->RotateTitle(histogr->GetYaxis()->GetRotateTitle());
1546 fHistogram->GetYaxis()->SetNoExponent(histogr->GetYaxis()->GetNoExponent());
1547 fHistogram->GetYaxis()->SetTimeDisplay(histogr->GetYaxis()->GetTimeDisplay());
1548 fHistogram->GetYaxis()->SetTimeFormat(histogr->GetYaxis()->GetTimeFormat());
1549 histogr->GetYaxis()->TAttAxis::Copy(*(fHistogram->GetYaxis()));
1550
1551 delete histogr;
1552 }
1553 return fHistogram;
1554}
1555
1556////////////////////////////////////////////////////////////////////////////////
1557/// Get x and y values for point number i.
1558/// The function returns -1 in case of an invalid request or the point number otherwise
1559
1561{
1562 if (i < 0 || i >= fNpoints || !fX || !fY) return -1;
1563 x = fX[i];
1564 y = fY[i];
1565 return i;
1566}
1567
1568////////////////////////////////////////////////////////////////////////////////
1569/// Get x value for point i.
1570
1572{
1573 if (i < 0 || i >= fNpoints || !fX)
1574 return -1.;
1575
1576 return fX[i];
1577}
1578
1579////////////////////////////////////////////////////////////////////////////////
1580/// Get y value for point i.
1581
1583{
1584 if (i < 0 || i >= fNpoints || !fY)
1585 return -1.;
1586
1587 return fY[i];
1588}
1589
1590////////////////////////////////////////////////////////////////////////////////
1591/// Get x axis of the graph.
1592
1594{
1595 auto h = GetHistogram();
1596 return h ? h->GetXaxis() : nullptr;
1597}
1598
1599////////////////////////////////////////////////////////////////////////////////
1600/// Get y axis of the graph.
1601
1603{
1604 auto h = GetHistogram();
1605 return h ? h->GetYaxis() : nullptr;
1606}
1607
1608////////////////////////////////////////////////////////////////////////////////
1609/// Implementation to get information on point of graph at cursor position
1610/// Adapted from class TH1
1611
1613{
1614 if (!gPad) {
1615 Error("GetObjectInfo", "Cannot be used without gPad");
1616 return nullptr;
1617 }
1618
1619 // localize point
1620 Int_t ipoint = -2;
1621 // start with a small window (in case the mouse is very close to one point)
1622 for (Int_t i = 0; i < fNpoints; i++) {
1623 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
1624 Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
1625
1626 if (dpx * dpx + dpy * dpy < 25) {
1627 ipoint = i;
1628 break;
1629 }
1630 }
1631
1632 Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
1633 Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py));
1634
1635 if (ipoint == -2)
1636 return Form("x=%g, y=%g", x, y);
1637
1640
1641 return Form("x=%g, y=%g, point=%d, xval=%g, yval=%g", x, y, ipoint, xval, yval);
1642}
1643
1644////////////////////////////////////////////////////////////////////////////////
1645/// Compute Initial values of parameters for a gaussian.
1646
1648{
1649 Double_t allcha, sumx, sumx2, x, val, rms, mean;
1650 Int_t bin;
1651 const Double_t sqrtpi = 2.506628;
1652
1653 // Compute mean value and RMS of the graph in the given range
1654 if (xmax <= xmin) {
1655 xmin = fX[0];
1656 xmax = fX[fNpoints-1];
1657 }
1658 Int_t np = 0;
1659 allcha = sumx = sumx2 = 0;
1660 for (bin = 0; bin < fNpoints; bin++) {
1661 x = fX[bin];
1662 if (x < xmin || x > xmax) continue;
1663 np++;
1664 val = fY[bin];
1665 sumx += val * x;
1666 sumx2 += val * x * x;
1667 allcha += val;
1668 }
1669 if (np == 0 || allcha == 0) return;
1670 mean = sumx / allcha;
1671 rms = TMath::Sqrt(sumx2 / allcha - mean * mean);
1673 if (rms == 0) rms = 1;
1675 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1676 f1->SetParameter(0, binwidx * allcha / (sqrtpi * rms));
1677 f1->SetParameter(1, mean);
1678 f1->SetParameter(2, rms);
1679 f1->SetParLimits(2, 0, 10 * rms);
1680}
1681
1682////////////////////////////////////////////////////////////////////////////////
1683/// Compute Initial values of parameters for an exponential.
1684
1686{
1688 Int_t ifail;
1689 if (xmax <= xmin) {
1690 xmin = fX[0];
1691 xmax = fX[fNpoints-1];
1692 }
1694
1696
1698 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1700 f1->SetParameter(1, slope);
1701}
1702
1703////////////////////////////////////////////////////////////////////////////////
1704/// Compute Initial values of parameters for a polynom.
1705
1707{
1708 Double_t fitpar[25];
1709
1711 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
1712 Int_t npar = f1->GetNpar();
1713 if (xmax <= xmin) {
1714 xmin = fX[0];
1715 xmax = fX[fNpoints-1];
1716 }
1717
1719
1720 for (Int_t i = 0; i < npar; i++) f1->SetParameter(i, fitpar[i]);
1721}
1722
1723////////////////////////////////////////////////////////////////////////////////
1724/// Insert a new point at the mouse position
1725
1727{
1728 if (!gPad) {
1729 Error("InsertPoint", "Cannot be used without gPad, requires last mouse position");
1730 return -1;
1731 }
1732
1733 Int_t px = gPad->GetEventX();
1734 Int_t py = gPad->GetEventY();
1735
1736 //localize point where to insert
1737 Int_t ipoint = -2;
1738 Int_t i, d = 0;
1739 // start with a small window (in case the mouse is very close to one point)
1740 for (i = 0; i < fNpoints - 1; i++) {
1741 d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
1742 if (d < 5) {
1743 ipoint = i + 1;
1744 break;
1745 }
1746 }
1747 if (ipoint == -2) {
1748 //may be we are far from one point, try again with a larger window
1749 for (i = 0; i < fNpoints - 1; i++) {
1750 d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
1751 if (d < 10) {
1752 ipoint = i + 1;
1753 break;
1754 }
1755 }
1756 }
1757 if (ipoint == -2) {
1758 //distinguish between first and last point
1759 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[0]));
1760 Int_t dpy = py - gPad->YtoAbsPixel(gPad->XtoPad(fY[0]));
1761 if (dpx * dpx + dpy * dpy < 25) ipoint = 0;
1762 else ipoint = fNpoints;
1763 }
1764
1765
1766 InsertPointBefore(ipoint, gPad->AbsPixeltoX(px), gPad->AbsPixeltoY(py));
1767
1768 gPad->Modified();
1769 return ipoint;
1770}
1771
1772
1773////////////////////////////////////////////////////////////////////////////////
1774/// Insert a new point with coordinates (x,y) before the point number `ipoint`.
1775
1777{
1778 if (ipoint < 0) {
1779 Error("TGraph", "Inserted point index should be >= 0");
1780 return;
1781 }
1782
1783 if (ipoint > fNpoints) {
1784 Error("TGraph", "Inserted point index should be <= %d", fNpoints);
1785 return;
1786 }
1787
1788 if (ipoint == fNpoints) {
1789 SetPoint(ipoint, x, y);
1790 return;
1791 }
1792
1794 CopyAndRelease(ps, ipoint, fNpoints++, ipoint + 1);
1795
1796 // To avoid redefinitions in descendant classes
1797 FillZero(ipoint, ipoint + 1);
1798
1799 fX[ipoint] = x;
1800 fY[ipoint] = y;
1801}
1802
1803
1804////////////////////////////////////////////////////////////////////////////////
1805/// Integrate the TGraph data within a given (index) range.
1806/// Note that this function computes the area of the polygon enclosed by the points of the TGraph.
1807/// The polygon segments, which are defined by the points of the TGraph, do not need to form a closed polygon,
1808/// since the last polygon segment, which closes the polygon, is taken as the line connecting the last TGraph point
1809/// with the first one. It is clear that the order of the point is essential in defining the polygon.
1810/// Also note that the segments should not intersect.
1811///
1812/// NB:
1813/// - if last=-1 (default) last is set to the last point.
1814/// - if (first <0) the first point (0) is taken.
1815///
1816/// ### Method:
1817///
1818/// There are many ways to calculate the surface of a polygon. It all depends on what kind of data
1819/// you have to deal with. The most evident solution would be to divide the polygon in triangles and
1820/// calculate the surface of them. But this can quickly become complicated as you will have to test
1821/// every segments of every triangles and check if they are intersecting with a current polygon's
1822/// segment or if it goes outside the polygon. Many calculations that would lead to many problems...
1823///
1824/// ### The solution (implemented by R.Brun)
1825/// Fortunately for us, there is a simple way to solve this problem, as long as the polygon's
1826/// segments don't intersect.
1827/// It takes the x coordinate of the current vertex and multiply it by the y coordinate of the next
1828/// vertex. Then it subtracts from it the result of the y coordinate of the current vertex multiplied
1829/// by the x coordinate of the next vertex. Then divide the result by 2 to get the surface/area.
1830///
1831/// ### Sources
1832/// - http://forums.wolfram.com/mathgroup/archive/1998/Mar/msg00462.html
1833/// - http://stackoverflow.com/questions/451426/how-do-i-calculate-the-surface-area-of-a-2d-polygon
1834
1836{
1837 if (first < 0) first = 0;
1838 if (last < 0) last = fNpoints - 1;
1839 if (last >= fNpoints) last = fNpoints - 1;
1840 if (first >= last) return 0;
1841 Int_t np = last - first + 1;
1842 Double_t sum = 0.0;
1843 //for(Int_t i=first;i<=last;i++) {
1844 // Int_t j = first + (i-first+1)%np;
1845 // sum += TMath::Abs(fX[i]*fY[j]);
1846 // sum -= TMath::Abs(fY[i]*fX[j]);
1847 //}
1848 for (Int_t i = first; i <= last; i++) {
1849 Int_t j = first + (i - first + 1) % np;
1850 sum += (fY[i] + fY[j]) * (fX[j] - fX[i]);
1851 }
1852 return 0.5 * TMath::Abs(sum);
1853}
1854
1855////////////////////////////////////////////////////////////////////////////////
1856/// Return 1 if the point (x,y) is inside the polygon defined by
1857/// the graph vertices 0 otherwise.
1858///
1859/// Algorithm:
1860///
1861/// The loop is executed with the end-point coordinates of a line segment
1862/// (X1,Y1)-(X2,Y2) and the Y-coordinate of a horizontal line.
1863/// The counter inter is incremented if the line (X1,Y1)-(X2,Y2) intersects
1864/// the horizontal line. In this case XINT is set to the X-coordinate of the
1865/// intersection point. If inter is an odd number, then the point x,y is within
1866/// the polygon.
1867
1869{
1870 return (Int_t)TMath::IsInside(x, y, fNpoints, fX, fY);
1871}
1872
1873////////////////////////////////////////////////////////////////////////////////
1874/// Least squares polynomial fitting without weights.
1875///
1876/// \param [in] m number of parameters
1877/// \param [in] a array of parameters
1878/// \param [in] xmin 1st point number to fit (default =0)
1879/// \param [in] xmax last point number to fit (default=fNpoints-1)
1880///
1881/// based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
1882
1884{
1885 const Double_t zero = 0.;
1886 const Double_t one = 1.;
1887 const Int_t idim = 20;
1888
1889 Double_t b[400] /* was [20][20] */;
1890 Int_t i, k, l, ifail;
1892 Double_t da[20], xk, yk;
1893 Int_t n = fNpoints;
1894 if (xmax <= xmin) {
1895 xmin = fX[0];
1896 xmax = fX[fNpoints-1];
1897 }
1898
1899 if (m <= 2) {
1900 LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
1901 return;
1902 }
1903 if (m > idim || m > n) return;
1904 da[0] = zero;
1905 for (l = 2; l <= m; ++l) {
1906 b[l-1] = zero;
1907 b[m + l*20 - 21] = zero;
1908 da[l-1] = zero;
1909 }
1910 Int_t np = 0;
1911 for (k = 0; k < fNpoints; ++k) {
1912 xk = fX[k];
1913 if (xk < xmin || xk > xmax) continue;
1914 np++;
1915 yk = fY[k];
1916 power = one;
1917 da[0] += yk;
1918 for (l = 2; l <= m; ++l) {
1919 power *= xk;
1920 b[l-1] += power;
1921 da[l-1] += power * yk;
1922 }
1923 for (l = 2; l <= m; ++l) {
1924 power *= xk;
1925 b[m + l*20 - 21] += power;
1926 }
1927 }
1928 b[0] = Double_t(np);
1929 for (i = 3; i <= m; ++i) {
1930 for (k = i; k <= m; ++k) {
1931 b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
1932 }
1933 }
1935
1936 if (ifail < 0) {
1937 a[0] = fY[0];
1938 for (i = 1; i < m; ++i) a[i] = 0;
1939 return;
1940 }
1941 for (i = 0; i < m; ++i) a[i] = da[i];
1942}
1943
1944////////////////////////////////////////////////////////////////////////////////
1945/// Least square linear fit without weights.
1946///
1947/// Fit a straight line (a0 + a1*x) to the data in this graph.
1948///
1949/// \param [in] ndata if ndata<0, fits the logarithm of the graph (used in InitExpo() to set
1950/// the initial parameter values for a fit with exponential function.
1951/// \param [in] a0 constant
1952/// \param [in] a1 slope
1953/// \param [in] ifail return parameter indicating the status of the fit (ifail=0, fit is OK)
1954/// \param [in] xmin, xmax fitting range
1955///
1956/// extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
1957
1959{
1961 Int_t i;
1963 Double_t fn, xk, yk;
1964 Double_t det;
1965 if (xmax <= xmin) {
1966 xmin = fX[0];
1967 xmax = fX[fNpoints-1];
1968 }
1969
1970 ifail = -2;
1971 xbar = ybar = x2bar = xybar = 0;
1972 Int_t np = 0;
1973 for (i = 0; i < fNpoints; ++i) {
1974 xk = fX[i];
1975 if (xk < xmin || xk > xmax) continue;
1976 np++;
1977 yk = fY[i];
1978 if (ndata < 0) {
1979 if (yk <= 0) yk = 1e-9;
1980 yk = TMath::Log(yk);
1981 }
1982 xbar += xk;
1983 ybar += yk;
1984 x2bar += xk * xk;
1985 xybar += xk * yk;
1986 }
1987 fn = Double_t(np);
1988 det = fn * x2bar - xbar * xbar;
1989 ifail = -1;
1990 if (det <= 0) {
1991 if (fn > 0) a0 = ybar / fn;
1992 else a0 = 0;
1993 a1 = 0;
1994 return;
1995 }
1996 ifail = 0;
1997 a0 = (x2bar * ybar - xbar * xybar) / det;
1998 a1 = (fn * xybar - xbar * ybar) / det;
1999}
2000
2001////////////////////////////////////////////////////////////////////////////////
2002/// Draw this graph with its current attributes.
2003
2009
2010////////////////////////////////////////////////////////////////////////////////
2011/// Draw the (x,y) as a graph.
2012
2018
2019////////////////////////////////////////////////////////////////////////////////
2020/// Draw the (x,y) as a histogram.
2021
2027
2028////////////////////////////////////////////////////////////////////////////////
2029/// Draw the stats
2030
2032{
2034 if (painter) painter->PaintStats(this, fit);
2035}
2036
2037////////////////////////////////////////////////////////////////////////////////
2038/// Print graph values.
2039
2041{
2042 for (Int_t i = 0; i < fNpoints; i++) {
2043 printf("x[%d]=%g, y[%d]=%g\n", i, fX[i], i, fY[i]);
2044 }
2045}
2046
2047////////////////////////////////////////////////////////////////////////////////
2048/// Recursively remove object from the list of functions
2049
2051{
2052 if (fFunctions) {
2055 }
2056 if (fHistogram == obj)
2057 fHistogram = nullptr;
2058}
2059
2060////////////////////////////////////////////////////////////////////////////////
2061/// Delete point close to the mouse position
2062/// Returns index of removed point (or -1 if nothing was changed)
2063
2065{
2066 if (!gPad) {
2067 Error("RemovePoint", "Cannot be used without gPad, requires last mouse position");
2068 return -1;
2069 }
2070
2071 Int_t px = gPad->GetEventX();
2072 Int_t py = gPad->GetEventY();
2073
2074 //localize point to be deleted
2075 Int_t ipoint = -2;
2076 // start with a small window (in case the mouse is very close to one point)
2077 for (Int_t i = 0; i < fNpoints; i++) {
2078 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
2079 Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
2080 if (dpx * dpx + dpy * dpy < 100) {
2081 ipoint = i;
2082 break;
2083 }
2084 }
2085 return RemovePoint(ipoint);
2086}
2087
2088////////////////////////////////////////////////////////////////////////////////
2089/// Delete point number ipoint
2090/// Returns index of removed point (or -1 if nothing was changed)
2091
2093{
2094 if ((ipoint < 0) || (ipoint >= fNpoints))
2095 return -1;
2096
2098 CopyAndRelease(ps, ipoint + 1, fNpoints--, ipoint);
2099 if (gPad) gPad->Modified();
2100 return ipoint;
2101}
2102
2103////////////////////////////////////////////////////////////////////////////////
2104/// Save the graph as .csv, .tsv or .txt. In case of any other extension, fall
2105/// back to TObject::SaveAs
2106///
2107/// The result can be immediately imported into Excel, gnuplot, Python or whatever,
2108/// without the needing to install pyroot, etc.
2109///
2110/// \param filename the name of the file where to store the graph
2111/// \param option some tuning options
2112///
2113/// The file extension defines the delimiter used:
2114/// - `.csv` : comma
2115/// - `.tsv` : tab
2116/// - `.txt` : space
2117///
2118/// By default file contains lines with (X, Y) coordinates. If errors are present,
2119/// (X, EX, Y, EY) are written. With asymmetric errors (X, EXL, EXH, Y, EYL, EYH) are stored.
2120/// If option contains "asroot" string, order of values will match such order in TGraph constructors.
2121/// So one will get (X, Y, EX, EY) or (X, Y, EXL, EXH, EYL, EYH).
2122///
2123/// Also one can directly select that kind of errors are stored:
2124/// - "errors" - (X, Y, EX, EY) will be stored
2125/// - "asymmerrors" - (X, Y, EXL, EXH, EYL, EYH) will be stored
2126/// - "noerrors" - just (X, Y) will be stored disregard of graph kind
2127///
2128/// If option contains "title" a title line is generated with the axis titles.
2129
2130
2131void TGraph::SaveAs(const char *filename, Option_t *option) const
2132{
2133 char del = '\0';
2134 TString ext = "";
2136 TString opt = option;
2137 opt.ToLower();
2138
2139 if (filename) {
2140 if (fname.EndsWith(".csv")) {del = ','; ext = "csv";}
2141 else if (fname.EndsWith(".tsv")) {del = '\t'; ext = "tsv";}
2142 else if (fname.EndsWith(".txt")) {del = ' '; ext = "txt";}
2143 }
2144 if (del) {
2145 std::ofstream out;
2146 out.open(filename, std::ios::out);
2147 if (!out.good ()) {
2148 Error("SaveAs", "cannot open file: %s", filename);
2149 return;
2150 }
2151 Bool_t store_title = opt.Contains("title");
2153 Bool_t as_root = opt.Contains("asroot") || opt.Contains("native");
2154 if (opt.Contains("noerrors"))
2155 no_errors = kTRUE;
2156 else if (opt.Contains("asymmerrors"))
2158 else if (opt.Contains("errors"))
2160 else if (InheritsFrom("TGraphErrors"))
2162 else if (InheritsFrom("TGraphAsymmErrors") || InheritsFrom("TGraphBentErrors"))
2164 else
2165 no_errors = kTRUE;
2166
2168 if (fHistogram) {
2171 }
2172 if (xtitle.IsNull())
2173 xtitle = "x";
2174 if (ytitle.IsNull())
2175 ytitle = "y";
2176
2177 if (plain_errors) {
2178 if(store_title) {
2179 if (as_root)
2180 out << "# " << xtitle << "\t" << ytitle << "\tex\tey\n";
2181 else
2182 out << "# " << xtitle << "\tex\t" << ytitle << "\tey\n";
2183 }
2184 for(int i = 0; i < fNpoints ; i++) {
2185 Double_t x = GetPointX(i);
2186 Double_t y = GetPointY(i);
2187 Double_t ex = GetErrorX(i);
2188 Double_t ey = GetErrorY(i);
2189 if (as_root)
2190 out << x << del << y << del << ex << del << ey << "\n";
2191 else
2192 out << x << del << ex << del << y << del << ey << "\n";
2193 }
2194 } else if (asymm_erros) {
2195 if(store_title) {
2196 if (as_root)
2197 out << "# " << xtitle << "\t" << ytitle << "\texl\texh\teyl\teyh\n";
2198 else
2199 out << "# " << xtitle << "\texl\texh\t" << ytitle << "\teyl\teyh\n";
2200 }
2201 for(int i = 0; i < GetN(); i++) {
2202 Double_t x = GetPointX(i);
2203 Double_t y = GetPointY(i);
2208 if (as_root)
2209 out << x << del << y << del << exl << del << exh << del << eyl << del << eyh << "\n";
2210 else
2211 out << x << del << exl << del << exh << del << y << del << eyl << del << eyh << "\n";
2212 }
2213 } else if (no_errors) {
2214 if(store_title)
2215 out << "# " << xtitle << "\t" << ytitle << "\n";
2216 for (Int_t i = 0 ; i < GetN(); i++)
2217 out << GetPointX(i) << del << GetPointY(i) << "\n";
2218 }
2219 out.close();
2220 Info("SaveAs", "%s file: %s has been generated", ext.Data(), filename);
2221 } else {
2223 }
2224}
2225
2226
2227////////////////////////////////////////////////////////////////////////////////
2228/// Save primitive as a C++ statement(s) on output stream out
2229
2230void TGraph::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
2231{
2232 TString args;
2233 if (fNpoints >= 1) {
2234 TString xname = SavePrimitiveVector(out, "graph_x", fNpoints, fX, kTRUE);
2235 TString yname = SavePrimitiveVector(out, "graph_y", fNpoints, fY);
2236 args.Form("%d, %s.data(), %s.data()", fNpoints, xname.Data(), yname.Data());
2237 }
2238
2239 SavePrimitiveConstructor(out, Class(), "graph", args, fNpoints < 1);
2240
2241 SaveHistogramAndFunctions(out, "graph", option);
2242}
2243
2244////////////////////////////////////////////////////////////////////////////////
2245/// Save histogram and list of functions of TGraph as C++ statement
2246/// Used in all TGraph-derived classes
2247
2248void TGraph::SaveHistogramAndFunctions(std::ostream &out, const char *varname, Option_t *option)
2249{
2250 thread_local Int_t frameNumber = 0;
2251
2252 TString ref = "Graph";
2253 if ((ref != GetName()) || (ref != GetTitle()))
2255
2256 SaveFillAttributes(out, varname, 0, 1000);
2257 SaveLineAttributes(out, varname, 1, 1, 1);
2258 SaveMarkerAttributes(out, varname, 1, 1, 1);
2259
2260 if (fHistogram) {
2262 fHistogram->SetName(TString::Format("Graph_histogram%d", ++frameNumber).Data());
2263 fHistogram->SavePrimitive(out, "nodraw");
2264 out << " " <<varname << "->SetHistogram(" << fHistogram->GetName() << ");\n";
2265 out << " \n";
2266 fHistogram->SetName(hname.Data());
2267 }
2268
2270
2271 if (!option)
2272 option = "";
2273 const char *l = strstr(option, "multigraph");
2274 if (l) {
2275 out << " multigraph->Add(" << varname << ",\"" << l + 10 << "\");\n";
2276 return;
2277 }
2278 l = strstr(option, "th2poly");
2279 if (l) {
2280 out << " " << l + 7 << "->AddBin(" << varname << ");\n";
2281 return;
2282 }
2283
2285}
2286
2287////////////////////////////////////////////////////////////////////////////////
2288/// Multiply the values of a TGraph by a constant c1.
2289///
2290/// If option contains "x" the x values are scaled
2291/// If option contains "y" the y values are scaled
2292/// If option contains "xy" both x and y values are scaled
2293
2295{
2296 TString opt = option; opt.ToLower();
2297 if (opt.Contains("x")) {
2298 for (Int_t i=0; i<GetN(); i++)
2299 GetX()[i] *= c1;
2300 }
2301 if (opt.Contains("y")) {
2302 for (Int_t i=0; i<GetN(); i++)
2303 GetY()[i] *= c1;
2304 }
2305}
2306
2307////////////////////////////////////////////////////////////////////////////////
2308/// Set number of points in the graph
2309/// Existing coordinates are preserved
2310/// New coordinates above fNpoints are preset to 0.
2311
2313{
2314 if (n < 0) n = 0;
2315 if (n == fNpoints) return;
2316 Double_t **ps = Allocate(n);
2317 CopyAndRelease(ps, 0, TMath::Min(fNpoints, n), 0);
2318 if (n > fNpoints) {
2320 }
2321 fNpoints = n;
2322}
2323
2324////////////////////////////////////////////////////////////////////////////////
2325/// Return kTRUE if kNotEditable bit is not set, kFALSE otherwise.
2326
2328{
2329 return TestBit(kNotEditable) ? kFALSE : kTRUE;
2330}
2331
2332////////////////////////////////////////////////////////////////////////////////
2333/// if editable=kFALSE, the graph cannot be modified with the mouse
2334/// by default a TGraph is editable
2335
2341
2342////////////////////////////////////////////////////////////////////////////////
2343/// Set highlight (enable/disable) mode for the graph
2344/// by default highlight mode is disable
2345
2347{
2348 if (IsHighlight() == set) return;
2349
2351 if (!painter) return;
2352 SetBit(kIsHighlight, set);
2353 painter->SetHighlight(this);
2354}
2355
2356/// Set the histogram underlying the TGraph. This transfers the ownership of h to the TGraph. The preexisting fHistogram will be deleted.
2358{
2359 delete fHistogram;
2360 fHistogram = h;
2361}
2362
2363////////////////////////////////////////////////////////////////////////////////
2364/// Set the maximum of the graph.
2365
2367{
2368 fMaximum = maximum;
2369 GetHistogram()->SetMaximum(maximum);
2370}
2371
2372////////////////////////////////////////////////////////////////////////////////
2373/// Set the minimum of the graph.
2374
2376{
2377 fMinimum = minimum;
2378 GetHistogram()->SetMinimum(minimum);
2379}
2380
2381////////////////////////////////////////////////////////////////////////////////
2382/// Set x and y values for point number i.
2383
2385{
2386 if (i < 0) return;
2388
2389 if (i >= fMaxSize) {
2390 Double_t **ps = ExpandAndCopy(i + 1, fNpoints);
2391 CopyAndRelease(ps, 0, 0, 0);
2392 }
2393 if (i >= fNpoints) {
2394 // points above i can be not initialized
2395 // set zero up to i-th point to avoid redefinition
2396 // of this method in descendant classes
2397 FillZero(fNpoints, i + 1);
2398 fNpoints = i + 1;
2399 }
2400 fX[i] = x;
2401 fY[i] = y;
2402 if (gPad) gPad->Modified();
2403}
2404
2405////////////////////////////////////////////////////////////////////////////////
2406/// Set x value for point i.
2407
2409{
2410 SetPoint(i, x, GetPointY(i));
2411}
2412
2413////////////////////////////////////////////////////////////////////////////////
2414/// Set y value for point i.
2415
2417{
2418 SetPoint(i, GetPointX(i), y);
2419}
2420
2421////////////////////////////////////////////////////////////////////////////////
2422/// Set graph name.
2423void TGraph::SetName(const char *name)
2424{
2425 fName = name;
2427}
2428
2429////////////////////////////////////////////////////////////////////////////////
2430/// Change (i.e. set) the title
2431///
2432/// if title is in the form `stringt;stringx;stringy;stringz`
2433/// the graph title is set to `stringt`, the x axis title to `stringx`,
2434/// the y axis title to `stringy`, and the z axis title to `stringz`.
2435///
2436/// To insert the character `;` in one of the titles, one should use `#;`
2437/// or `#semicolon`.
2438
2439void TGraph::SetTitle(const char* title)
2440{
2441 fTitle = title;
2442 fTitle.ReplaceAll("#;",2,"#semicolon",10);
2443 Int_t p = fTitle.Index(";");
2444
2445 if (p>0) {
2446 if (!fHistogram) GetHistogram();
2447 fHistogram->SetTitle(title);
2448 Int_t n = fTitle.Length()-p;
2449 if (p>0) fTitle.Remove(p,n);
2450 fTitle.ReplaceAll("#semicolon",10,"#;",2);
2451 } else {
2452 if (fHistogram) fHistogram->SetTitle(title);
2453 }
2454}
2455
2456////////////////////////////////////////////////////////////////////////////////
2457/// Set graph name and title
2458
2459void TGraph::SetNameTitle(const char *name, const char *title)
2460{
2461 SetName(name);
2462 SetTitle(title);
2463}
2464
2465////////////////////////////////////////////////////////////////////////////////
2466/// Set statistics option on/off.
2467///
2468/// By default, the statistics box is drawn.
2469/// The paint options can be selected via gStyle->SetOptStat.
2470/// This function sets/resets the kNoStats bit in the graph object.
2471/// It has priority over the Style option.
2472
2474{
2476 if (!stats) {
2478 //remove the "stats" object from the list of functions
2479 if (fFunctions) {
2480 TObject *obj = fFunctions->FindObject("stats");
2481 if (obj) {
2482 fFunctions->Remove(obj);
2483 delete obj;
2484 }
2485 }
2486 }
2487}
2488
2489////////////////////////////////////////////////////////////////////////////////
2490/// if size*2 <= fMaxSize allocate new arrays of size points,
2491/// copy points [0,oend).
2492/// Return newarray (passed or new instance if it was zero
2493/// and allocations are needed)
2494
2496{
2497 if (size * 2 > fMaxSize || !fMaxSize)
2498 return nullptr;
2499
2501 CopyPoints(newarrays, 0, oend, 0);
2502 return newarrays;
2503}
2504
2505////////////////////////////////////////////////////////////////////////////////
2506/// Sorts the points of this TGraph using in-place quicksort (see e.g. older glibc).
2507/// To compare two points the function parameter greaterfunc is used (see TGraph::CompareX for an
2508/// example of such a method, which is also the default comparison function for Sort). After
2509/// the sort, greaterfunc(this, i, j) will return kTRUE for all i>j if ascending == kTRUE, and
2510/// kFALSE otherwise.
2511///
2512/// The last two parameters are used for the recursive quick sort, stating the range to be sorted
2513///
2514/// Examples:
2515/// ~~~ {.cpp}
2516/// // sort points along x axis
2517/// graph->Sort();
2518/// // sort points along their distance to origin
2519/// graph->Sort(&TGraph::CompareRadius);
2520///
2521/// Bool_t CompareErrors(const TGraph* gr, Int_t i, Int_t j) {
2522/// const TGraphErrors* ge=(const TGraphErrors*)gr;
2523/// return (ge->GetEY()[i]>ge->GetEY()[j]); }
2524/// // sort using the above comparison function, largest errors first
2525/// graph->Sort(&CompareErrors, kFALSE);
2526/// ~~~
2527
2528void TGraph::Sort(Bool_t (*greaterfunc)(const TGraph *, Int_t, Int_t) /*=TGraph::CompareX()*/,
2529 Bool_t ascending /*=kTRUE*/, Int_t low /*=0*/, Int_t high /*=-1111*/)
2530{
2531 // set the bit in case of an ascending =sort in X
2532 if (greaterfunc == TGraph::CompareX && ascending && low == 0 && high == -1111)
2534
2535 if (high == -1111)
2536 high = fNpoints - 1;
2537
2538 // Create a vector to store the indices of the graph data points.
2539 // We use std::vector<Int_t> instead of std::vector<ULong64_t> to match the input type
2540 // required by the comparison operator's signature provided as `greaterfunc`
2541 std::vector<int> sorting_indices(fNpoints);
2542 std::iota(sorting_indices.begin(), sorting_indices.end(), 0);
2543
2544 // Sort the indices using the provided comparison function
2545 // We use std::stable_sort here because the libc++ implementation of std::sort
2546 // is not standard-compliant until LLVM 14 which caused errors on the mac nodes
2547 // of our CI, related issue: https://github.com/llvm/llvm-project/issues/21211
2548 std::stable_sort(sorting_indices.begin() + low, sorting_indices.begin() + high + 1,
2549 [&](int left, int right) { return left != right && greaterfunc(this, left, right) != ascending; });
2550
2551 Int_t numSortedPoints = high - low + 1;
2553}
2554
2555////////////////////////////////////////////////////////////////////////////////
2556/// Stream an object of class TGraph.
2557
2559{
2560 if (b.IsReading()) {
2561 UInt_t R__s, R__c;
2562 Version_t R__v = b.ReadVersion(&R__s, &R__c);
2563 if (R__v > 2) {
2564 b.ReadClassBuffer(TGraph::Class(), this, R__v, R__s, R__c);
2565 if (fHistogram) fHistogram->SetDirectory(nullptr);
2566 TIter next(fFunctions);
2567 TObject *obj;
2568 while ((obj = next())) {
2569 if (obj->InheritsFrom(TF1::Class())) {
2570 TF1 *f1 = (TF1*)obj;
2571 f1->SetParent(this);
2572 }
2573 }
2575 return;
2576 }
2577 //====process old versions before automatic schema evolution
2582 b >> fNpoints;
2584 fX = new Double_t[fNpoints];
2585 fY = new Double_t[fNpoints];
2586 if (R__v < 2) {
2587 Float_t *x = new Float_t[fNpoints];
2588 Float_t *y = new Float_t[fNpoints];
2589 b.ReadFastArray(x, fNpoints);
2590 b.ReadFastArray(y, fNpoints);
2591 for (Int_t i = 0; i < fNpoints; i++) {
2592 fX[i] = x[i];
2593 fY[i] = y[i];
2594 }
2595 delete [] y;
2596 delete [] x;
2597 } else {
2598 b.ReadFastArray(fX, fNpoints);
2599 b.ReadFastArray(fY, fNpoints);
2600 }
2601 b >> fFunctions;
2602 b >> fHistogram;
2603 if (fHistogram) fHistogram->SetDirectory(nullptr);
2604 if (R__v < 2) {
2605 Float_t mi, ma;
2606 b >> mi;
2607 b >> ma;
2608 fMinimum = mi;
2609 fMaximum = ma;
2610 } else {
2611 b >> fMinimum;
2612 b >> fMaximum;
2613 }
2614 b.CheckByteCount(R__s, R__c, TGraph::IsA());
2615 //====end of old versions
2616
2617 } else {
2618 b.WriteClassBuffer(TGraph::Class(), this);
2619 }
2620}
2621
2622////////////////////////////////////////////////////////////////////////////////
2623/// Swap points.
2624
2630
2631////////////////////////////////////////////////////////////////////////////////
2632/// Update the fX and fY arrays with the sorted values.
2633
2635{
2636 std::vector<Double_t> fXSorted(numSortedPoints);
2637 std::vector<Double_t> fYSorted(numSortedPoints);
2638
2639 // Fill the sorted X and Y values based on the sorted indices
2640 std::generate(fXSorted.begin(), fXSorted.end(),
2641 [begin = low, &sorting_indices, this]() mutable { return fX[sorting_indices[begin++]]; });
2642 std::generate(fYSorted.begin(), fYSorted.end(),
2643 [begin = low, &sorting_indices, this]() mutable { return fY[sorting_indices[begin++]]; });
2644
2645 // Copy the sorted X and Y values back to the original arrays
2646 std::copy(fXSorted.begin(), fXSorted.end(), fX + low);
2647 std::copy(fYSorted.begin(), fYSorted.end(), fY + low);
2648}
2649
2650////////////////////////////////////////////////////////////////////////////////
2651/// Swap values.
2652
2654{
2655 Double_t tmp = arr[pos1];
2656 arr[pos1] = arr[pos2];
2657 arr[pos2] = tmp;
2658}
2659
2660////////////////////////////////////////////////////////////////////////////////
2661/// Set current style settings in this graph
2662/// This function is called when either TCanvas::UseCurrentStyle
2663/// or TROOT::ForceStyle have been invoked.
2664
2695
2696////////////////////////////////////////////////////////////////////////////////
2697/// Adds all graphs from the collection to this graph.
2698/// Returns the total number of points in the result or -1 in case of an error.
2699
2701{
2702 TIter next(li);
2703 while (TObject* o = next()) {
2704 TGraph *g = dynamic_cast<TGraph*>(o);
2705 if (!g) {
2706 Error("Merge",
2707 "Cannot merge - an object which doesn't inherit from TGraph found in the list");
2708 return -1;
2709 }
2710 DoMerge(g);
2711 }
2712 return GetN();
2713}
2714
2715////////////////////////////////////////////////////////////////////////////////
2716/// protected function to perform the merge operation of a graph
2717
2719{
2720 Double_t x = 0, y = 0;
2721 for (Int_t i = 0 ; i < g->GetN(); i++) {
2722 g->GetPoint(i, x, y);
2723 SetPoint(GetN(), x, y);
2724 }
2725 return kTRUE;
2726}
2727
2728////////////////////////////////////////////////////////////////////////////////
2729/// Move all graph points on specified values dx,dy
2730/// If log argument specified, calculation done in logarithmic scale like:
2731/// new_value = exp( log(old_value) + delta );
2732
2734{
2735 Double_t x = 0, y = 0;
2736 for (Int_t i = 0 ; i < GetN(); i++) {
2737 GetPoint(i, x, y);
2738 if (!logx) {
2739 x += dx;
2740 } else if (x > 0) {
2741 x = TMath::Exp(TMath::Log(x) + dx);
2742 }
2743 if (!logy) {
2744 y += dy;
2745 } else if (y > 0) {
2746 y = TMath::Exp(TMath::Log(y) + dy);
2747 }
2748 SetPoint(i, x, y);
2749 }
2750}
2751
2752
2753////////////////////////////////////////////////////////////////////////////////
2754/// Find zero of a continuous function.
2755/// This function finds a real zero of the continuous real
2756/// function Y(X) in a given interval (A,B). See accompanying
2757/// notes for details of the argument list and calling sequence
2758
2761{
2762 static Double_t a, b, ya, ytest, y1, x1, h;
2763 static Int_t j1, it, j3, j2;
2764 Double_t yb, x2;
2765 yb = 0;
2766
2767 // Calculate Y(X) at X=AZ.
2768 if (k <= 0) {
2769 a = AZ;
2770 b = BZ;
2771 X = a;
2772 j1 = 1;
2773 it = 1;
2774 k = j1;
2775 return;
2776 }
2777
2778 // Test whether Y(X) is sufficiently small.
2779
2780 if (TMath::Abs(Y) <= E2) {
2781 k = 2;
2782 return;
2783 }
2784
2785 // Calculate Y(X) at X=BZ.
2786
2787 if (j1 == 1) {
2788 ya = Y;
2789 X = b;
2790 j1 = 2;
2791 return;
2792 }
2793 // Test whether the signs of Y(AZ) and Y(BZ) are different.
2794 // if not, begin the binary subdivision.
2795
2796 if (j1 != 2) goto L100;
2797 if (ya * Y < 0) goto L120;
2798 x1 = a;
2799 y1 = ya;
2800 j1 = 3;
2801 h = b - a;
2802 j2 = 1;
2803 x2 = a + 0.5 * h;
2804 j3 = 1;
2805 it++; //*-*- Check whether (maxiterations) function values have been calculated.
2806 if (it >= maxiterations) k = j1;
2807 else X = x2;
2808 return;
2809
2810 // Test whether a bracket has been found .
2811 // If not,continue the search
2812
2813L100:
2814 if (j1 > 3) goto L170;
2815 if (ya*Y >= 0) {
2816 if (j3 >= j2) {
2817 h = 0.5 * h;
2818 j2 = 2 * j2;
2819 a = x1;
2820 ya = y1;
2821 x2 = a + 0.5 * h;
2822 j3 = 1;
2823 } else {
2824 a = X;
2825 ya = Y;
2826 x2 = X + h;
2827 j3++;
2828 }
2829 it++;
2830 if (it >= maxiterations) k = j1;
2831 else X = x2;
2832 return;
2833 }
2834
2835 // The first bracket has been found.calculate the next X by the
2836 // secant method based on the bracket.
2837
2838L120:
2839 b = X;
2840 yb = Y;
2841 j1 = 4;
2842L130:
2843 if (TMath::Abs(ya) > TMath::Abs(yb)) {
2844 x1 = a;
2845 y1 = ya;
2846 X = b;
2847 Y = yb;
2848 } else {
2849 x1 = b;
2850 y1 = yb;
2851 X = a;
2852 Y = ya;
2853 }
2854
2855 // Use the secant method based on the function values y1 and Y.
2856 // check that x2 is inside the interval (a,b).
2857
2858L150:
2859 x2 = X - Y * (X - x1) / (Y - y1);
2860 x1 = X;
2861 y1 = Y;
2863 if ((x2 - a)*(x2 - b) < 0) {
2864 it++;
2865 if (it >= maxiterations) k = j1;
2866 else X = x2;
2867 return;
2868 }
2869
2870 // Calculate the next value of X by bisection . Check whether
2871 // the maximum accuracy has been achieved.
2872
2873L160:
2874 x2 = 0.5 * (a + b);
2875 ytest = 0;
2876 if ((x2 - a)*(x2 - b) >= 0) {
2877 k = 2;
2878 return;
2879 }
2880 it++;
2881 if (it >= maxiterations) k = j1;
2882 else X = x2;
2883 return;
2884
2885
2886 // Revise the bracket (a,b).
2887
2888L170:
2889 if (j1 != 4) return;
2890 if (ya * Y < 0) {
2891 b = X;
2892 yb = Y;
2893 } else {
2894 a = X;
2895 ya = Y;
2896 }
2897
2898 // Use ytest to decide the method for the next value of X.
2899
2900 if (ytest <= 0) goto L130;
2901 if (TMath::Abs(Y) - ytest <= 0) goto L150;
2902 goto L160;
2903}
#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:1571
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:1835
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:1868
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:1883
virtual Double_t Chisquare(TF1 *f1, Option_t *option="") const
Return the chisquare of this graph with respect to f1.
Definition TGraph.cxx:705
void UseCurrentStyle() override
Set current style settings in this graph This function is called when either TCanvas::UseCurrentStyle...
Definition TGraph.cxx:2665
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2384
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:2700
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:594
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:2495
virtual Double_t GetRMS(Int_t axis=1) const
Return RMS of X (axis=1) or Y (axis=2)
Definition TGraph.cxx:1367
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:2004
@ 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:1415
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:2733
virtual Double_t GetErrorYlow(Int_t bin) const
It always returns a negative value.
Definition TGraph.cxx:1433
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:2634
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:790
Double_t GetMinimum() const
Definition TGraph.h:151
void Print(Option_t *chopt="") const override
Print graph values.
Definition TGraph.cxx:2040
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum of the graph.
Definition TGraph.cxx:2366
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:741
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:750
virtual Double_t GetErrorYhigh(Int_t bin) const
It always returns a negative value.
Definition TGraph.cxx:1424
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:733
Int_t GetN() const
Definition TGraph.h:131
TF1 * GetFunction(const char *name) const
Return pointer to function with name.
Definition TGraph.cxx:1444
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:1958
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:833
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:917
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:1282
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:2528
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:722
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:759
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:1612
Double_t ** AllocateArrays(Int_t Narrays, Int_t arraySize)
Allocate arrays.
Definition TGraph.cxx:626
virtual void Scale(Double_t c1=1., Option_t *option="y")
Multiply the values of a TGraph by a constant c1.
Definition TGraph.cxx:2294
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:1336
static void SwapValues(Double_t *arr, Int_t pos1, Int_t pos2)
Swap values.
Definition TGraph.cxx:2653
void Streamer(TBuffer &) override
Stream an object of class TGraph.
Definition TGraph.cxx:2558
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:2759
virtual Double_t ** Allocate(Int_t newsize)
Allocate internal data structures for newsize points.
Definition TGraph.cxx:618
virtual void FitPanel()
Display a GUI panel with all graph fit options.
Definition TGraph.cxx:1302
void Browse(TBrowser *b) override
Browse.
Definition TGraph.cxx:679
virtual Bool_t DoMerge(const TGraph *g)
protected function to perform the merge operation of a graph
Definition TGraph.cxx:2718
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:1776
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:1084
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:2131
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:981
virtual void InitExpo(Double_t xmin=0, Double_t xmax=0)
Compute Initial values of parameters for an exponential.
Definition TGraph.cxx:1685
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:2064
virtual void InitGaus(Double_t xmin=0, Double_t xmax=0)
Compute Initial values of parameters for a gaussian.
Definition TGraph.cxx:1647
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:650
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:666
void SetName(const char *name="") override
Set graph name.
Definition TGraph.cxx:2423
virtual void SetHighlight(Bool_t set=kTRUE)
Set highlight (enable/disable) mode for the graph by default highlight mode is disable.
Definition TGraph.cxx:2346
virtual void SwapPoints(Int_t pos1, Int_t pos2)
Swap points.
Definition TGraph.cxx:2625
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:859
TAxis * GetXaxis() const
Get x axis of the graph.
Definition TGraph.cxx:1593
Bool_t GetEditable() const
Return kTRUE if kNotEditable bit is not set, kFALSE otherwise.
Definition TGraph.cxx:2327
virtual Double_t GetCorrelationFactor() const
Return graph correlation factor.
Definition TGraph.cxx:1324
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:1130
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:1117
virtual void Expand(Int_t newsize)
If array sizes <= newsize, expand storage to 2*newsize.
Definition TGraph.cxx:1093
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:2248
virtual Double_t GetMean(Int_t axis=1) const
Return mean value of X (axis=1) or Y (axis=2)
Definition TGraph.cxx:1352
Double_t * fX
[fNpoints] array of X points
Definition TGraph.h:47
virtual void PaintStats(TF1 *fit)
Draw the stats.
Definition TGraph.cxx:2031
TAxis * GetYaxis() const
Get y axis of the graph.
Definition TGraph.cxx:1602
TObject * FindObject(const char *name) const override
Search object named name in the list of functions.
Definition TGraph.cxx:1139
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TGraph.cxx:2473
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:1455
virtual Double_t GetErrorY(Int_t bin) const
It always returns a negative value. Real implementation in TGraphErrors.
Definition TGraph.cxx:1397
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TGraph.cxx:2230
virtual Double_t GetPointY(Int_t i) const
Get y value for point i.
Definition TGraph.cxx:1582
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:2013
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2439
virtual Int_t InsertPoint()
Insert a new point at the mouse position.
Definition TGraph.cxx:1726
void RecursiveRemove(TObject *obj) override
Recursively remove object from the list of functions.
Definition TGraph.cxx:2050
virtual void SetPointY(Int_t i, Double_t y)
Set y value for point i.
Definition TGraph.cxx:2416
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a graph.
Definition TGraph.cxx:907
virtual void SetHistogram(TH1F *h)
Set the histogram underlying the TGraph. This transfers the ownership of h to the TGraph....
Definition TGraph.cxx:2357
virtual void DrawPanel()
Display a panel with all graph drawing options.
Definition TGraph.cxx:960
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:2022
void SetNameTitle(const char *name="", const char *title="") override
Set graph name and title.
Definition TGraph.cxx:2459
virtual void SetPointX(Int_t i, Double_t x)
Set x value for point i.
Definition TGraph.cxx:2408
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:2312
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:1560
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:2336
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum of the graph.
Definition TGraph.cxx:2375
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:807
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:1406
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:1706
virtual Double_t GetErrorX(Int_t bin) const
It always returns a negative value. Real implementation in TGraphErrors.
Definition TGraph.cxx:1389
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
@ kIgnoreCase
Definition TString.h:285
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