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