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