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