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