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