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