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