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