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