Logo ROOT  
Reference Guide
TMultiGraph.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Rene Brun 12/10/2000
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 "TROOT.h"
13#include "TEnv.h"
14#include "TBrowser.h"
15#include "TMultiGraph.h"
16#include "TGraph.h"
17#include "TH1.h"
18#include "TH2.h"
19#include "TPolyLine3D.h"
20#include "TVirtualPad.h"
21#include "Riostream.h"
22#include "TVirtualFitter.h"
23#include "TPluginManager.h"
24#include "TClass.h"
25#include "TMath.h"
26#include "TSystem.h"
27#include <stdlib.h>
28
29#include "HFitInterface.h"
30#include "Fit/DataRange.h"
32
33#include <ctype.h>
34
35extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
36
38
39
40////////////////////////////////////////////////////////////////////////////////
41
42/** \class TMultiGraph
43 \ingroup Hist
44 \brief A TMultiGraph is a collection of TGraph (or derived) objects.
45
46- [Introduction](#MG00)
47- [MultiGraphs' drawing](#MG01)
48 - [Setting drawing options](#MG01a)
49 - [Titles setting](#MG01b)
50 - [The option "3D"](#MG01c)
51 - [Legend drawing](#MG01d)
52 - [Automatic coloring](#MG01e)
53 - [Reverse axis](#MG01f)
54- [MultiGraphs' fitting](#MG02)
55 - [Fit box position](#MG02a)
56- [Axis' limits setting](#MG03)
57
58
59### <a name="MG00"></a> Introduction
60
61A TMultiGraph allows to manipulate a set of graphs as a single entity. In particular,
62when drawn, the X and Y axis ranges are automatically computed such as all the graphs
63will be visible.
64
65`TMultiGraph::Add` should be used to add a new graph to the list.
66
67The TMultiGraph owns the objects in the list.
68
69The number of graphs in a multigraph can be retrieve with:
70~~~ {.cpp}
71mg->GetListOfGraphs()->GetEntries();
72~~~
73
74### <a name="MG00"></a> MultiGraphs' Drawing
75
76The drawing options are the same as for TGraph.
77Like for TGraph, the painting is performed thanks to the TGraphPainter
78class. All details about the various painting options are given in this class.
79
80Example:
81~~~ {.cpp}
82 TGraph *gr1 = new TGraph(...
83 TGraphErrors *gr2 = new TGraphErrors(...
84 TMultiGraph *mg = new TMultiGraph();
85 mg->Add(gr1,"lp");
86 mg->Add(gr2,"cp");
87 mg->Draw("a");
88~~~
89
90#### <a name="MG01a"></a> Setting drawing options
91
92The drawing option for each TGraph may be specified as an optional
93second argument of the `Add` function.
94
95If a draw option is specified, it will be used to draw the graph,
96otherwise the graph will be drawn with the option specified in
97`TMultiGraph::Draw`
98
99#### <a name="MG01b"></a> Titles setting
100
101The global title and the axis titles can be modified the following way:
102
103~~~ {.cpp}
104 [...]
105 auto mg = new TMultiGraph;
106 mg->SetTitle("title;xaxis title; yaxis title");
107 mg->Add(g1);
108 mg->Add(g2);
109 mg->Draw("apl");
110~~~
111
112#### <a name="MG01c"></a> The option "3D"
113
114A special option `3D` allows to draw the graphs in a 3D space. See the
115following example:
116
117Begin_Macro(source)
118{
119 auto c0 = new TCanvas("c1","multigraph L3",200,10,700,500);
120
121 auto mg = new TMultiGraph();
122
123 auto gr1 = new TGraph(); gr1->SetLineColor(kBlue);
124 auto gr2 = new TGraph(); gr2->SetLineColor(kRed);
125 auto gr3 = new TGraph(); gr3->SetLineColor(kGreen);
126 auto gr4 = new TGraph(); gr4->SetLineColor(kOrange);
127
128 Double_t dx = 6.28/1000;
129 Double_t x = -3.14;
130
131 for (int i=0; i<=1000; i++) {
132 x = x+dx;
133 gr1->SetPoint(i,x,2.*TMath::Sin(x));
134 gr2->SetPoint(i,x,TMath::Cos(x));
135 gr3->SetPoint(i,x,TMath::Cos(x*x));
136 gr4->SetPoint(i,x,TMath::Cos(x*x*x));
137 }
138
139 mg->Add(gr4); gr4->SetTitle("Cos(x*x*x)"); gr4->SetLineWidth(3);
140 mg->Add(gr3); gr3->SetTitle("Cos(x*x)") ; gr3->SetLineWidth(3);
141 mg->Add(gr2); gr2->SetTitle("Cos(x)") ; gr2->SetLineWidth(3);
142 mg->Add(gr1); gr1->SetTitle("2*Sin(x)") ; gr1->SetLineWidth(3);
143
144 mg->SetTitle("Multi-graph Title; X-axis Title; Y-axis Title");
145
146 mg->Draw("a fb l3d");
147
148 mg->GetHistogram()->GetXaxis()->SetRangeUser(0.,2.5);
149 gPad->Modified();
150 gPad->Update();
151}
152End_Macro
153
154#### <a name="MG01d"></a> Legend drawing
155
156The method TPad::BuildLegend is able to extract the graphs inside a
157multigraph. The following example demonstrate this.
158
159Begin_Macro(source)
160{
161 auto c3 = new TCanvas("c3","c3",600, 400);
162
163 auto mg = new TMultiGraph("mg","mg");
164
165 const Int_t size = 10;
166
167 double px[size];
168 double py1[size];
169 double py2[size];
170 double py3[size];
171
172 for ( int i = 0; i < size ; ++i ) {
173 px[i] = i;
174 py1[i] = size - i;
175 py2[i] = size - 0.5 * i;
176 py3[i] = size - 0.6 * i;
177 }
178
179 auto gr1 = new TGraph( size, px, py1 );
180 gr1->SetName("gr1");
181 gr1->SetTitle("graph 1");
182 gr1->SetMarkerStyle(21);
183 gr1->SetDrawOption("AP");
184 gr1->SetLineColor(2);
185 gr1->SetLineWidth(4);
186 gr1->SetFillStyle(0);
187
188 auto gr2 = new TGraph( size, px, py2 );
189 gr2->SetName("gr2");
190 gr2->SetTitle("graph 2");
191 gr2->SetMarkerStyle(22);
192 gr2->SetMarkerColor(2);
193 gr2->SetDrawOption("P");
194 gr2->SetLineColor(3);
195 gr2->SetLineWidth(4);
196 gr2->SetFillStyle(0);
197
198 auto gr3 = new TGraph( size, px, py3 );
199 gr3->SetName("gr3");
200 gr3->SetTitle("graph 3");
201 gr3->SetMarkerStyle(23);
202 gr3->SetLineColor(4);
203 gr3->SetLineWidth(4);
204 gr3->SetFillStyle(0);
205
206 mg->Add( gr1 );
207 mg->Add( gr2 );
208
209 gr3->Draw("ALP");
210 mg->Draw("LP");
211 c3->BuildLegend();
212}
213End_Macro
214
215#### <a name="MG01e"></a> Automatic coloring
216
217Automatic coloring according to the current palette is available as shown in the
218following example:
219
220Begin_Macro(source)
221../../../tutorials/graphs/multigraphpalettecolor.C
222End_Macro
223
224#### <a name="MG01f"></a> Reverse axis
225
226\since **ROOT version 6.19/02**
227
228When a TMultiGraph is drawn, the X-axis is drawn with increasing values from left to
229right and the Y-axis from bottom to top. The two options RX and RY allow to change
230this order. The option RX allows to draw the X-axis with increasing values from
231right to left and the RY option allows to draw the Y-axis with increasing values
232from top to bottom. The following example illustrate how to use these options.
233
234Begin_Macro(source)
235{
236 auto *c = new TCanvas();
237 c->Divide(2,1);
238
239 auto *g1 = new TGraphErrors();
240 g1->SetPoint(0,-4,-3);
241 g1->SetPoint(1,1,1);
242 g1->SetPoint(2,2,1);
243 g1->SetPoint(3,3,4);
244 g1->SetPoint(4,5,5);
245 g1->SetPointError(0,1.,2.);
246 g1->SetPointError(1,2,1);
247 g1->SetPointError(2,2,3);
248 g1->SetPointError(3,3,2);
249 g1->SetPointError(4,4,5);
250 g1->SetMarkerStyle(21);
251
252 auto *g2 = new TGraph();
253 g2->SetPoint(0,4,8);
254 g2->SetPoint(1,5,9);
255 g2->SetPoint(2,6,10);
256 g2->SetPoint(3,10,11);
257 g2->SetPoint(4,15,12);
258 g2->SetLineColor(kRed);
259 g2->SetLineWidth(5);
260
261 auto mg = new TMultiGraph();
262 mg->Add(g1,"P");
263 mg->Add(g2,"L");
264
265 c->cd(1); gPad->SetGrid(1,1);
266 mg->Draw("A");
267
268 c->cd(2); gPad->SetGrid(1,1);
269 mg->Draw("A RX RY");
270}
271End_Macro
272
273### <a name="MG02"></a> MultiGraphs' fitting
274
275The following example shows how to fit a TMultiGraph.
276
277Begin_Macro(source)
278{
279 auto c1 = new TCanvas("c1","c1",600,400);
280
281 Double_t px1[2] = {2.,4.};
282 Double_t dx1[2] = {0.1,0.1};
283 Double_t py1[2] = {2.1,4.0};
284 Double_t dy1[2] = {0.3,0.2};
285
286 Double_t px2[2] = {3.,5.};
287 Double_t dx2[2] = {0.1,0.1};
288 Double_t py2[2] = {3.2,4.8};
289 Double_t dy2[2] = {0.3,0.2};
290
291 gStyle->SetOptFit(0001);
292
293 auto g1 = new TGraphErrors(2,px1,py1,dx1,dy1);
294 g1->SetMarkerStyle(21);
295 g1->SetMarkerColor(2);
296
297 auto g2 = new TGraphErrors(2,px2,py2,dx2,dy2);
298 g2->SetMarkerStyle(22);
299 g2->SetMarkerColor(3);
300
301 auto g = new TMultiGraph();
302 g->Add(g1);
303 g->Add(g2);
304
305 g->Draw("AP");
306
307 g->Fit("pol1","FQ");
308}
309End_Macro
310
311#### <a name="MG02a"></a> Fit box position
312
313When the graphs in a TMultiGraph are fitted, the fit parameters boxes
314overlap. The following example shows how to make them all visible.
315
316
317Begin_Macro(source)
318../../../tutorials/graphs/multigraph.C
319End_Macro
320
321### <a name="MG03"></a> Axis' limits setting
322
323The axis limits can be changed the like for TGraph. The same methods apply on
324the multigraph.
325Note the two differents ways to change limits on X and Y axis.
326
327Begin_Macro(source)
328{
329 auto c2 = new TCanvas("c2","c2",600,400);
330
331 TGraph *g[3];
332 Double_t x[10] = {0,1,2,3,4,5,6,7,8,9};
333 Double_t y[10] = {1,2,3,4,5,5,4,3,2,1};
334 auto mg = new TMultiGraph();
335 for (int i=0; i<3; i++) {
336 g[i] = new TGraph(10, x, y);
337 g[i]->SetMarkerStyle(20);
338 g[i]->SetMarkerColor(i+2);
339 for (int j=0; j<10; j++) y[j] = y[j]-1;
340 mg->Add(g[i]);
341 }
342 mg->Draw("APL");
343 mg->GetXaxis()->SetTitle("E_{#gamma} (GeV)");
344 mg->GetYaxis()->SetTitle("Coefficients");
345
346 // Change the axis limits
347 gPad->Modified();
348 mg->GetXaxis()->SetLimits(1.5,7.5);
349 mg->SetMinimum(0.);
350 mg->SetMaximum(10.);
351}
352End_Macro
353*/
354
355
356////////////////////////////////////////////////////////////////////////////////
357/// TMultiGraph default constructor.
358
360{
361 fGraphs = 0;
362 fFunctions = 0;
363 fHistogram = 0;
364 fMaximum = -1111;
365 fMinimum = -1111;
366}
367
368
369////////////////////////////////////////////////////////////////////////////////
370/// Constructor with name and title.
371
372TMultiGraph::TMultiGraph(const char *name, const char *title)
373 : TNamed(name,title)
374{
375 fGraphs = 0;
376 fFunctions = 0;
377 fHistogram = 0;
378 fMaximum = -1111;
379 fMinimum = -1111;
380}
381
382
383////////////////////////////////////////////////////////////////////////////////
384/// Copy constructor.
385
387 TNamed (mg),
388 fGraphs(mg.fGraphs),
389 fFunctions(mg.fFunctions),
390 fHistogram(mg.fHistogram),
391 fMaximum(mg.fMaximum),
392 fMinimum(mg.fMinimum)
393{
394}
395
396
397////////////////////////////////////////////////////////////////////////////////
398/// Assignment operator.
399
401{
402 if (this!=&mg) {
404 fGraphs=mg.fGraphs;
405 fFunctions=mg.fFunctions;
406 fHistogram=mg.fHistogram;
407 fMaximum=mg.fMaximum;
408 fMinimum=mg.fMinimum;
409 }
410 return *this;
411}
412
413
414////////////////////////////////////////////////////////////////////////////////
415/// TMultiGraph destructor.
416
418{
419 if (!fGraphs) return;
420 TGraph *g;
421 TIter next(fGraphs);
422 while ((g = (TGraph*) next())) {
423 g->ResetBit(kMustCleanup);
424 }
425 fGraphs->Delete();
426 delete fGraphs;
427 fGraphs = 0;
428 delete fHistogram;
429 fHistogram = 0;
430 if (fFunctions) {
432 //special logic to support the case where the same object is
433 //added multiple times in fFunctions.
434 //This case happens when the same object is added with different
435 //drawing modes
436 TObject *obj;
437 while ((obj = fFunctions->First())) {
438 while (fFunctions->Remove(obj)) { }
439 delete obj;
440 }
441 delete fFunctions;
442 }
443}
444
445
446////////////////////////////////////////////////////////////////////////////////
447/// Add a new graph to the list of graphs.
448/// Note that the graph is now owned by the TMultigraph.
449/// Deleting the TMultiGraph object will automatically delete the graphs.
450/// You should not delete the graphs when the TMultigraph is still active.
451
453{
454 if (!fGraphs) fGraphs = new TList();
455 graph->SetBit(kMustCleanup);
456 fGraphs->Add(graph,chopt);
457}
458
459
460////////////////////////////////////////////////////////////////////////////////
461/// Add all the graphs in "multigraph" to the list of graphs.
462///
463/// - If "chopt" is defined all the graphs in "multigraph" will be added with
464/// the "chopt" option.
465/// - If "chopt" is undefined each graph will be added with the option it had
466/// in "multigraph".
467
468void TMultiGraph::Add(TMultiGraph *multigraph, Option_t *chopt)
469{
470 TList *graphlist = multigraph->GetListOfGraphs();
471 if (!graphlist) return;
472
473 if (!fGraphs) fGraphs = new TList();
474
475 TObjOptLink *lnk = (TObjOptLink*)graphlist->FirstLink();
476 TObject *obj = 0;
477
478 while (lnk) {
479 obj = lnk->GetObject();
480 if (!strlen(chopt)) fGraphs->Add(obj,lnk->GetOption());
481 else fGraphs->Add(obj,chopt);
482 lnk = (TObjOptLink*)lnk->Next();
483 }
484}
485
486
487////////////////////////////////////////////////////////////////////////////////
488/// Browse multigraph.
489
491{
492 TString opt = gEnv->GetValue("TGraph.BrowseOption", "");
493 if (opt.IsNull()) {
494 opt = b ? b->GetDrawOption() : "alp";
495 opt = (opt == "") ? "alp" : opt.Data();
496 }
497 Draw(opt.Data());
498 gPad->Update();
499}
500
501
502////////////////////////////////////////////////////////////////////////////////
503/// Compute distance from point px,py to each graph.
504
506{
507 // Are we on the axis?
508 const Int_t kMaxDiff = 10;
509 Int_t distance = 9999;
510 if (fHistogram) {
511 distance = fHistogram->DistancetoPrimitive(px,py);
512 if (distance <= 0) return distance;
513 }
514
515 // Loop on the list of graphs
516 if (!fGraphs) return distance;
517 TGraph *g;
518 TIter next(fGraphs);
519 while ((g = (TGraph*) next())) {
520 Int_t dist = g->DistancetoPrimitive(px,py);
521 if (dist <= 0) return 0;
522 if (dist < kMaxDiff) {gPad->SetSelected(g); return dist;}
523 }
524 return distance;
525}
526
527
528////////////////////////////////////////////////////////////////////////////////
529/// Draw this multigraph with its current attributes.
530///
531/// Options to draw a graph are described in TGraphPainter.
532///
533/// The drawing option for each TGraph may be specified as an optional
534/// second argument of the Add function. You can use GetGraphDrawOption
535/// to return this option.
536///
537/// If a draw option is specified, it will be used to draw the graph,
538/// otherwise the graph will be drawn with the option specified in
539/// TMultiGraph::Draw. Use GetDrawOption to return the option specified
540/// when drawing the TMultiGraph.
541
543{
544 TString opt = option;
545 opt.ToLower();
546
547 if (gPad) {
548 if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
549 if (opt.Contains("a")) gPad->Clear();
550 }
551 AppendPad(option);
552}
553
554
555////////////////////////////////////////////////////////////////////////////////
556/// Fit this graph with function with name fname.
557///
558/// interface to TF1::Fit(TF1 *f1...
559
561{
562 char *linear;
563 linear= (char*)strstr(fname, "++");
564 if (linear) {
565 TF1 f1(fname, fname, xmin, xmax);
566 return Fit(&f1,option,"",xmin,xmax);
567 }
568 TF1 * f1 = (TF1*)gROOT->GetFunction(fname);
569 if (!f1) { Printf("Unknown function: %s",fname); return -1; }
570
571 return Fit(f1,option,"",xmin,xmax);
572}
573
574
575////////////////////////////////////////////////////////////////////////////////
576/// Fit this multigraph with function f1.
577///
578/// In this function all graphs of the multigraph are fitted simultaneously
579///
580/// f1 is an already predefined function created by TF1.
581/// Predefined functions such as gaus, expo and poln are automatically
582/// created by ROOT.
583///
584/// The list of fit options is given in parameter `option`which may takes the
585/// following values:
586///
587/// - "W" Set all errors to 1
588/// - "U" Use a User specified fitting algorithm (via SetFCN)
589/// - "Q" Quiet mode (minimum printing)
590/// - "V" Verbose mode (default is between Q and V)
591/// - "B" Use this option when you want to fix one or more parameters
592/// and the fitting function is like "gaus","expo","poln","landau".
593/// - "R" Use the Range specified in the function range
594/// - "N" Do not store the graphics function, do not draw
595/// - "0" Do not plot the result of the fit. By default the fitted function
596/// is drawn unless the option"N" above is specified.
597/// - "+" Add this new fitted function to the list of fitted functions
598/// (by default, any previous function is deleted)
599/// - "C" In case of linear fitting, not calculate the chisquare (saves time)
600/// - "F" If fitting a polN, switch to minuit fitter
601/// - "ROB" In case of linear fitting, compute the LTS regression
602/// coefficients (robust(resistant) regression), using
603/// the default fraction of good points
604/// - "ROB=0.x" - compute the LTS regression coefficients, using
605/// 0.x as a fraction of good points
606///
607/// When the fit is drawn (by default), the parameter goption may be used
608/// to specify a list of graphics options. See TGraph::Paint for a complete
609/// list of these options.
610///
611/// In order to use the Range option, one must first create a function
612/// with the expression to be fitted. For example, if your graph
613/// has a defined range between -4 and 4 and you want to fit a gaussian
614/// only in the interval 1 to 3, you can do:
615/// ~~~ {.cpp}
616/// TF1 *f1 = new TF1("f1","gaus",1,3);
617/// graph->Fit("f1","R");
618/// ~~~
619///
620/// ### Who is calling this function ?
621///
622/// Note that this function is called when calling TGraphErrors::Fit
623/// or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
624/// see the discussion below on the errors calculation.
625///
626/// ### Setting initial conditions
627///
628/// Parameters must be initialized before invoking the Fit function.
629/// The setting of the parameter initial values is automatic for the
630/// predefined functions : poln, expo, gaus, landau. One can however disable
631/// this automatic computation by specifying the option "B".
632/// You can specify boundary limits for some or all parameters via
633/// ~~~ {.cpp}
634/// f1->SetParLimits(p_number, parmin, parmax);
635/// ~~~
636/// if `parmin>=parmax`, the parameter is fixed
637/// Note that you are not forced to fix the limits for all parameters.
638/// For example, if you fit a function with 6 parameters, you can do:
639/// ~~~ {.cpp}
640/// func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
641/// func->SetParLimits(4,-10,-4);
642/// func->SetParLimits(5, 1,1);
643/// ~~~
644/// With this setup, parameters 0->3 can vary freely
645/// Parameter 4 has boundaries [-10,-4] with initial value -8
646/// Parameter 5 is fixed to 100.
647///
648/// ### Fit range
649///
650/// The fit range can be specified in two ways:
651///
652/// - specify rxmax > rxmin (default is rxmin=rxmax=0)
653/// - specify the option "R". In this case, the function will be taken
654/// instead of the full graph range.
655///
656/// ### Changing the fitting function
657///
658/// By default a chi2 fitting function is used for fitting the TGraphs's.
659/// The function is implemented in `FitUtil::EvaluateChi2`.
660/// In case of TGraphErrors an effective chi2 is used
661/// (see TGraphErrors fit in TGraph::Fit) and is implemented in
662/// `FitUtil::EvaluateChi2Effective`
663/// To specify a User defined fitting function, specify option "U" and
664/// call the following function:
665/// ~~~ {.cpp}
666/// TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
667/// ~~~
668/// where MyFittingFunction is of type:
669/// ~~~ {.cpp}
670/// extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
671/// ~~~
672///
673/// ### Access to the fit result
674///
675/// The function returns a TFitResultPtr which can hold a pointer to a TFitResult object.
676/// By default the TFitResultPtr contains only the status of the fit and it converts
677/// automatically to an integer. If the option "S" is instead used, TFitResultPtr contains
678/// the TFitResult and behaves as a smart pointer to it. For example one can do:
679/// ~~~ {.cpp}
680/// TFitResultPtr r = graph->Fit("myFunc","S");
681/// TMatrixDSym cov = r->GetCovarianceMatrix(); // to access the covariance matrix
682/// Double_t par0 = r->Parameter(0); // retrieve the value for the parameter 0
683/// Double_t err0 = r->ParError(0); // retrieve the error for the parameter 0
684/// r->Print("V"); // print full information of fit including covariance matrix
685/// r->Write(); // store the result in a file
686/// ~~~
687///
688/// The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also
689/// from the fitted function.
690///
691/// ### Associated functions
692///
693/// One or more object (typically a TF1*) can be added to the list
694/// of functions (fFunctions) associated to each graph.
695/// When TGraph::Fit is invoked, the fitted function is added to this list.
696/// Given a graph gr, one can retrieve an associated function
697/// with:
698/// ~~~ {.cpp}
699/// TF1 *myfunc = gr->GetFunction("myfunc");
700/// ~~~
701///
702/// If the graph is made persistent, the list of
703/// associated functions is also persistent. Given a pointer (see above)
704/// to an associated function myfunc, one can retrieve the function/fit
705/// parameters with calls such as:
706/// ~~~ {.cpp}
707/// Double_t chi2 = myfunc->GetChisquare();
708/// Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
709/// Double_t err0 = myfunc->GetParError(0); //error on first parameter
710/// ~~~
711///
712/// ### Fit Statistics
713///
714/// You can change the statistics box to display the fit parameters with
715/// the TStyle::SetOptFit(mode) method. This mode has four digits.
716/// mode = pcev (default = 0111)
717///
718/// - v = 1; print name/values of parameters
719/// - e = 1; print errors (if e=1, v must be 1)
720/// - c = 1; print Chisquare/Number of degrees of freedom
721/// - p = 1; print Probability
722///
723/// For example: `gStyle->SetOptFit(1011);`
724/// prints the fit probability, parameter names/values, and errors.
725/// You can change the position of the statistics box with these lines
726/// (where g is a pointer to the TGraph):
727///
728/// ~~~ {.cpp}
729/// Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
730/// Root > st->SetX1NDC(newx1); //new x start position
731/// Root > st->SetX2NDC(newx2); //new x end position
732/// ~~~
733
735{
736 // internal multigraph fitting methods
737 Foption_t fitOption;
739
740 // create range and minimizer options with default values
741 ROOT::Fit::DataRange range(rxmin,rxmax);
743 return ROOT::Fit::FitObject(this, f1 , fitOption , minOption, goption, range);
744
745}
746
747////////////////////////////////////////////////////////////////////////////////
748/// Display a panel with all histogram fit options.
749/// See class TFitPanel for example
750
752{
753 if (!gPad)
754 gROOT->MakeDefCanvas();
755
756 if (!gPad) {
757 Error("FitPanel", "Unable to create a default canvas");
758 return;
759 }
760
761 // use plugin manager to create instance of TFitEditor
762 TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
763 if (handler && handler->LoadPlugin() != -1) {
764 if (handler->ExecPlugin(2, gPad, this) == 0)
765 Error("FitPanel", "Unable to crate the FitPanel");
766 }
767 else
768 Error("FitPanel", "Unable to find the FitPanel plug-in");
769}
770
771////////////////////////////////////////////////////////////////////////////////
772/// Return the draw option for the TGraph `gr` in this TMultiGraph.
773/// The return option is the one specified when calling TMultiGraph::Add(gr,option).
774
776{
777 if (!fGraphs || !gr) return "";
778 TListIter next(fGraphs);
779 TObject *obj;
780 while ((obj = next())) {
781 if (obj == (TObject*)gr) return next.GetOption();
782 }
783 return "";
784}
785
786
787////////////////////////////////////////////////////////////////////////////////
788/// Compute Initial values of parameters for a gaussian.
789
791{
792 Double_t allcha, sumx, sumx2, x, val, rms, mean;
793 Int_t bin;
794 const Double_t sqrtpi = 2.506628;
795
796 // Compute mean value and RMS of the graph in the given range
797 Int_t np = 0;
798 allcha = sumx = sumx2 = 0;
799 TGraph *g;
800 TIter next(fGraphs);
801 Double_t *px, *py;
802 Int_t npp; //number of points in each graph
803 while ((g = (TGraph*) next())) {
804 px=g->GetX();
805 py=g->GetY();
806 npp=g->GetN();
807 for (bin=0; bin<npp; bin++) {
808 x=px[bin];
809 if (x<xmin || x>xmax) continue;
810 np++;
811 val=py[bin];
812 sumx+=val*x;
813 sumx2+=val*x*x;
814 allcha+=val;
815 }
816 }
817 if (np == 0 || allcha == 0) return;
818 mean = sumx/allcha;
819 rms = TMath::Sqrt(sumx2/allcha - mean*mean);
820
821 Double_t binwidx = TMath::Abs((xmax-xmin)/np);
822 if (rms == 0) rms = 1;
824 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
825 f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
826 f1->SetParameter(1,mean);
827 f1->SetParameter(2,rms);
828 f1->SetParLimits(2,0,10*rms);
829}
830
831
832////////////////////////////////////////////////////////////////////////////////
833/// Compute Initial values of parameters for an exponential.
834
836{
837 Double_t constant, slope;
838 Int_t ifail;
839
840 LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
841
843 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
844 f1->SetParameter(0,constant);
845 f1->SetParameter(1,slope);
846}
847
848
849////////////////////////////////////////////////////////////////////////////////
850/// Compute Initial values of parameters for a polynom.
851
853{
854 Double_t fitpar[25];
855
857 TF1 *f1 = (TF1*)grFitter->GetUserFunc();
858 Int_t npar = f1->GetNpar();
859
860 LeastSquareFit(npar, fitpar, xmin, xmax);
861
862 for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
863}
864
865
866////////////////////////////////////////////////////////////////////////////////
867/// Least squares lpolynomial fitting without weights.
868///
869/// - m number of parameters
870/// - a array of parameters
871/// - first 1st point number to fit (default =0)
872/// - last last point number to fit (default=fNpoints-1)
873///
874/// based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
875
877{
878 const Double_t zero = 0.;
879 const Double_t one = 1.;
880 const Int_t idim = 20;
881
882 Double_t b[400] /* was [20][20] */;
883 Int_t i, k, l, ifail, bin;
884 Double_t power;
885 Double_t da[20], xk, yk;
886
887
888 //count the total number of points to fit
889 TGraph *g;
890 TIter next(fGraphs);
891 Double_t *px, *py;
892 Int_t n=0;
893 Int_t npp;
894 while ((g = (TGraph*) next())) {
895 px=g->GetX();
896 npp=g->GetN();
897 for (bin=0; bin<npp; bin++) {
898 xk=px[bin];
899 if (xk < xmin || xk > xmax) continue;
900 n++;
901 }
902 }
903 if (m <= 2) {
904 LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
905 return;
906 }
907 if (m > idim || m > n) return;
908 da[0] = zero;
909 for (l = 2; l <= m; ++l) {
910 b[l-1] = zero;
911 b[m + l*20 - 21] = zero;
912 da[l-1] = zero;
913 }
914 Int_t np = 0;
915
916 next.Reset();
917 while ((g = (TGraph*) next())) {
918 px=g->GetX();
919 py=g->GetY();
920 npp=g->GetN();
921
922 for (k = 0; k <= npp; ++k) {
923 xk = px[k];
924 if (xk < xmin || xk > xmax) continue;
925 np++;
926 yk = py[k];
927 power = one;
928 da[0] += yk;
929 for (l = 2; l <= m; ++l) {
930 power *= xk;
931 b[l-1] += power;
932 da[l-1] += power*yk;
933 }
934 for (l = 2; l <= m; ++l) {
935 power *= xk;
936 b[m + l*20 - 21] += power;
937 }
938 }
939 }
940 b[0] = Double_t(np);
941 for (i = 3; i <= m; ++i) {
942 for (k = i; k <= m; ++k) {
943 b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
944 }
945 }
946 H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
947
948 if (ifail < 0) {
949 //a[0] = fY[0];
950 py=((TGraph *)fGraphs->First())->GetY();
951 a[0]=py[0];
952 for (i=1; i<m; ++i) a[i] = 0;
953 return;
954 }
955 for (i=0; i<m; ++i) a[i] = da[i];
956}
957
958
959////////////////////////////////////////////////////////////////////////////////
960/// Least square linear fit without weights.
961///
962/// Fit a straight line (a0 + a1*x) to the data in this graph.
963///
964/// - ndata: number of points to fit
965/// - first: first point number to fit
966/// - last: last point to fit O(ndata should be last-first
967/// - ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
968///
969/// extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
970
973{
974 Double_t xbar, ybar, x2bar;
975 Int_t i;
976 Double_t xybar;
977 Double_t fn, xk, yk;
978 Double_t det;
979
980 ifail = -2;
981 xbar = ybar = x2bar = xybar = 0;
982 Int_t np = 0;
983 TGraph *g;
984 TIter next(fGraphs);
985 Double_t *px, *py;
986 Int_t npp;
987 while ((g = (TGraph*) next())) {
988 px=g->GetX();
989 py=g->GetY();
990 npp=g->GetN();
991 for (i = 0; i < npp; ++i) {
992 xk = px[i];
993 if (xk < xmin || xk > xmax) continue;
994 np++;
995 yk = py[i];
996 if (ndata < 0) {
997 if (yk <= 0) yk = 1e-9;
998 yk = TMath::Log(yk);
999 }
1000 xbar += xk;
1001 ybar += yk;
1002 x2bar += xk*xk;
1003 xybar += xk*yk;
1004 }
1005 }
1006 fn = Double_t(np);
1007 det = fn*x2bar - xbar*xbar;
1008 ifail = -1;
1009 if (det <= 0) {
1010 if (fn > 0) a0 = ybar/fn;
1011 else a0 = 0;
1012 a1 = 0;
1013 return;
1014 }
1015 ifail = 0;
1016 a0 = (x2bar*ybar - xbar*xybar) / det;
1017 a1 = (fn*xybar - xbar*ybar) / det;
1018}
1019
1020
1021////////////////////////////////////////////////////////////////////////////////
1022/// Return 1 if the point (x,y) is inside one of the graphs 0 otherwise.
1023
1025{
1026 Int_t in = 0;
1027 if (!fGraphs) return in;
1028 TGraph *g;
1029 TIter next(fGraphs);
1030 while ((g = (TGraph*) next())) {
1031 in = g->IsInside(x, y);
1032 if (in) return in;
1033 }
1034 return in;
1035}
1036
1037
1038////////////////////////////////////////////////////////////////////////////////
1039/// Returns a pointer to the histogram used to draw the axis.
1040/// Takes into account following cases.
1041///
1042/// 1. if `fHistogram` exists it is returned
1043/// 2. if `fHistogram` doesn't exists and `gPad` exists `gPad` is updated. That
1044/// may trigger the creation of `fHistogram`. If `fHistogram` still does not
1045/// exit but `hframe` does (if user called `TPad::DrawFrame`) the pointer to
1046/// `hframe` histogram is returned
1047/// 3. after the two previous steps, if `fHistogram` still doesn't exist, then
1048/// it is created.
1049
1051{
1052 if (fHistogram) return fHistogram;
1053
1054 if (gPad) {
1055 gPad->Modified();
1056 gPad->Update();
1057 if (fHistogram) return fHistogram;
1058 TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
1059 if (h1) return h1;
1060 }
1061
1062 Bool_t initialrangeset = kFALSE;
1063 Double_t rwxmin = 0.,rwxmax = 0.,rwymin = 0.,rwymax = 0.;
1064 TGraph *g;
1065 Int_t npt = 100 ;
1066 TIter next(fGraphs);
1067 while ((g = (TGraph*) next())) {
1068 if (g->GetN() <= 0) continue;
1069 if (initialrangeset) {
1070 Double_t rx1,ry1,rx2,ry2;
1071 g->ComputeRange(rx1, ry1, rx2, ry2);
1072 if (rx1 < rwxmin) rwxmin = rx1;
1073 if (ry1 < rwymin) rwymin = ry1;
1074 if (rx2 > rwxmax) rwxmax = rx2;
1075 if (ry2 > rwymax) rwymax = ry2;
1076 } else {
1077 g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1078 initialrangeset = kTRUE;
1079 }
1080 if (g->GetN() > npt) npt = g->GetN();
1081 }
1082 fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1083 if (!fHistogram) return 0;
1084 fHistogram->SetMinimum(rwymin);
1086 fHistogram->SetMaximum(rwymax);
1087 fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1089 return fHistogram;
1090}
1091
1092
1093////////////////////////////////////////////////////////////////////////////////
1094/// Return pointer to function with name.
1095///
1096/// Functions such as TGraph::Fit store the fitted function in the list of
1097/// functions of this graph.
1098
1100{
1101 if (!fFunctions) return 0;
1102 return (TF1*)fFunctions->FindObject(name);
1103}
1104
1105////////////////////////////////////////////////////////////////////////////////
1106/// Return pointer to list of functions.
1107/// If pointer is null create the list
1108
1110{
1111 if (!fFunctions) fFunctions = new TList();
1112 return fFunctions;
1113}
1114
1115
1116////////////////////////////////////////////////////////////////////////////////
1117/// Get x axis of the graph.
1118/// This method returns a valid axis only after the TMultigraph has been drawn.
1119
1121{
1122 TH1 *h = GetHistogram();
1123 if (!h) return 0;
1124 return h->GetXaxis();
1125}
1126
1127
1128////////////////////////////////////////////////////////////////////////////////
1129/// Get y axis of the graph.
1130/// This method returns a valid axis only after the TMultigraph has been drawn.
1131
1133{
1134 TH1 *h = GetHistogram();
1135 if (!h) return 0;
1136 return h->GetYaxis();
1137}
1138
1139
1140////////////////////////////////////////////////////////////////////////////////
1141/// Paint all the graphs of this multigraph.
1142
1144{
1145 const TPickerStackGuard pushGuard(this);
1146
1147 if (!fGraphs) return;
1148 if (fGraphs->GetSize() == 0) return;
1149
1150 char option[128];
1151 strlcpy(option,choptin,128);
1152 Int_t nch = strlen(choptin);
1153 for (Int_t i=0;i<nch;i++) option[i] = toupper(option[i]);
1154
1155 // Automatic color
1156 char *l1 = strstr(option,"PFC"); // Automatic Fill Color
1157 char *l2 = strstr(option,"PLC"); // Automatic Line Color
1158 char *l3 = strstr(option,"PMC"); // Automatic Marker Color
1159 if (l1 || l2 || l3) {
1160 TString opt1 = option; opt1.ToLower();
1161 if (l1) memcpy(l1," ",3);
1162 if (l2) memcpy(l2," ",3);
1163 if (l3) memcpy(l3," ",3);
1165 TGraph* gAti;
1166 Int_t ngraphs = fGraphs->GetSize();
1167 Int_t ic;
1168 gPad->IncrementPaletteColor(ngraphs, opt1);
1169 for (Int_t i=0;i<ngraphs;i++) {
1170 ic = gPad->NextPaletteColor();
1171 gAti = (TGraph*)(fGraphs->At(i));
1172 if (l1) gAti->SetFillColor(ic);
1173 if (l2) gAti->SetLineColor(ic);
1174 if (l3) gAti->SetMarkerColor(ic);
1175 lnk = (TObjOptLink*)lnk->Next();
1176 }
1177 }
1178
1179 char *l;
1180
1181 TString chopt = option;
1182
1183 l = (char*)strstr(chopt.Data(),"3D");
1184 if (l) {
1185 l = (char*)strstr(chopt.Data(),"L");
1186 if (l) PaintPolyLine3D(chopt.Data());
1187 return;
1188 }
1189
1190 l = (char*)strstr(chopt.Data(),"PADS");
1191 if (l) {
1192 chopt.ReplaceAll("PADS","");
1193 PaintPads(chopt.Data());
1194 return;
1195 }
1196
1197 char *lrx = (char *)strstr(chopt.Data(), "RX"); // Reverse graphs along X axis
1198 char *lry = (char *)strstr(chopt.Data(), "RY"); // Reverse graphs along Y axis
1199 if (lrx || lry) {
1200 PaintReverse(chopt.Data());
1201 return;
1202 }
1203
1204 TGraph *g;
1205
1206 l = (char*)strstr(chopt.Data(),"A");
1207 if (l) {
1208 *l = ' ';
1209 TIter next(fGraphs);
1210 Int_t npt = 100;
1211 Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
1212 rwxmin = gPad->GetUxmin();
1213 rwxmax = gPad->GetUxmax();
1214 rwymin = gPad->GetUymin();
1215 rwymax = gPad->GetUymax();
1216 char *xtitle = 0;
1217 char *ytitle = 0;
1218 Int_t firstx = 0;
1219 Int_t lastx = 0;
1220 Bool_t timedisplay = kFALSE;
1221 char *timeformat = 0;
1222
1223 if (fHistogram) {
1224 //cleanup in case of a previous unzoom and in case one of the TGraph has changed
1226 TGraph* gAti;
1227 Int_t ngraphs = fGraphs->GetSize();
1228 Bool_t reset_hist = kFALSE;
1229 for (Int_t i=0;i<ngraphs;i++) {
1230 gAti = (TGraph*)(fGraphs->At(i));
1231 if(gAti->TestBit(TGraph::kResetHisto)) {reset_hist = kTRUE; break;}
1232 lnk = (TObjOptLink*)lnk->Next();
1233 }
1234 if (fHistogram->GetMinimum() >= fHistogram->GetMaximum() || reset_hist) {
1235 nch = strlen(fHistogram->GetXaxis()->GetTitle());
1236 firstx = fHistogram->GetXaxis()->GetFirst();
1237 lastx = fHistogram->GetXaxis()->GetLast();
1238 timedisplay = fHistogram->GetXaxis()->GetTimeDisplay();
1239 if (nch) {
1240 xtitle = new char[nch+1];
1241 strlcpy(xtitle,fHistogram->GetXaxis()->GetTitle(),nch+1);
1242 }
1243 nch = strlen(fHistogram->GetYaxis()->GetTitle());
1244 if (nch) {
1245 ytitle = new char[nch+1];
1246 strlcpy(ytitle,fHistogram->GetYaxis()->GetTitle(),nch+1);
1247 }
1248 nch = strlen(fHistogram->GetXaxis()->GetTimeFormat());
1249 if (nch) {
1250 timeformat = new char[nch+1];
1251 strlcpy(timeformat,fHistogram->GetXaxis()->GetTimeFormat(),nch+1);
1252 }
1253 delete fHistogram;
1254 fHistogram = 0;
1255 }
1256 }
1257 if (fHistogram) {
1258 minimum = fHistogram->GetYaxis()->GetXmin();
1259 maximum = fHistogram->GetYaxis()->GetXmax();
1260 uxmin = gPad->PadtoX(rwxmin);
1261 uxmax = gPad->PadtoX(rwxmax);
1262 } else {
1263 Bool_t initialrangeset = kFALSE;
1264 while ((g = (TGraph*) next())) {
1265 if (g->GetN() <= 0) continue;
1266 if (initialrangeset) {
1267 Double_t rx1,ry1,rx2,ry2;
1268 g->ComputeRange(rx1, ry1, rx2, ry2);
1269 if (rx1 < rwxmin) rwxmin = rx1;
1270 if (ry1 < rwymin) rwymin = ry1;
1271 if (rx2 > rwxmax) rwxmax = rx2;
1272 if (ry2 > rwymax) rwymax = ry2;
1273 } else {
1274 g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1275 initialrangeset = kTRUE;
1276 }
1277 if (g->GetN() > npt) npt = g->GetN();
1278 }
1279 if (rwxmin == rwxmax) rwxmax += 1.;
1280 if (rwymin == rwymax) rwymax += 1.;
1281 dx = 0.05*(rwxmax-rwxmin);
1282 dy = 0.05*(rwymax-rwymin);
1283 uxmin = rwxmin - dx;
1284 uxmax = rwxmax + dx;
1285 if (gPad->GetLogy()) {
1286 if (rwymin <= 0) rwymin = 0.001*rwymax;
1287 minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
1288 maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
1289 } else {
1290 minimum = rwymin - dy;
1291 maximum = rwymax + dy;
1292 }
1293 if (minimum < 0 && rwymin >= 0) minimum = 0;
1294 if (maximum > 0 && rwymax <= 0) maximum = 0;
1295 }
1296
1297 if (fMinimum != -1111) rwymin = minimum = fMinimum;
1298 if (fMaximum != -1111) rwymax = maximum = fMaximum;
1299 if (uxmin < 0 && rwxmin >= 0) {
1300 if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
1301 //else uxmin = 0;
1302 }
1303 if (uxmax > 0 && rwxmax <= 0) {
1304 if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
1305 //else uxmax = 0;
1306 }
1307 if (minimum < 0 && rwymin >= 0) {
1308 if (gPad->GetLogy()) minimum = 0.9*rwymin;
1309 //else minimum = 0;
1310 }
1311 if (maximum > 0 && rwymax <= 0) {
1312 if (gPad->GetLogy()) maximum = 1.1*rwymax;
1313 //else maximum = 0;
1314 }
1315 if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
1316 if (uxmin <= 0 && gPad->GetLogx()) {
1317 if (uxmax > 1000) uxmin = 1;
1318 else uxmin = 0.001*uxmax;
1319 }
1320 rwymin = minimum;
1321 rwymax = maximum;
1322 if (fHistogram) {
1323 fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1324 }
1325
1326 // Create a temporary histogram to draw the axis
1327 if (!fHistogram) {
1328 // the graph is created with at least as many channels as there are points
1329 // to permit zooming on the full range
1330 rwxmin = uxmin;
1331 rwxmax = uxmax;
1332 fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1333 if (!fHistogram) return;
1334 fHistogram->SetMinimum(rwymin);
1336 fHistogram->SetMaximum(rwymax);
1337 fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1339 if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
1340 if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
1341 if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
1342 if (timedisplay) {fHistogram->GetXaxis()->SetTimeDisplay(timedisplay);}
1343 if (timeformat) {fHistogram->GetXaxis()->SetTimeFormat(timeformat); delete [] timeformat;}
1344 }
1345 fHistogram->Paint("0");
1346 }
1347
1348 TGraph *gfit = nullptr;
1349 if (fGraphs) {
1351 TObject *obj = 0;
1352
1353 chopt.ReplaceAll("A","");
1354
1355 while (lnk) {
1356
1357 obj = lnk->GetObject();
1358
1359 gPad->PushSelectableObject(obj);
1360
1361 if (!gPad->PadInHighlightMode() || (gPad->PadInHighlightMode() && obj == gPad->GetSelected())) {
1362 TString opt = lnk->GetOption();
1363 if (!opt.IsWhitespace())
1364 obj->Paint(opt.ReplaceAll("A","").Data());
1365 else {
1366 if (!chopt.IsWhitespace()) obj->Paint(chopt.Data());
1367 else obj->Paint("L");
1368 }
1369 }
1370
1371 lnk = (TObjOptLink*)lnk->Next();
1372 }
1373
1374 gfit = (TGraph*)obj; // pick one TGraph in the list to paint the fit parameters.
1375 }
1376
1377 TObject *f;
1378 TF1 *fit = nullptr;
1379 if (fFunctions) {
1380 TIter next(fFunctions);
1381 while ((f = (TObject*) next())) {
1382 if (f->InheritsFrom(TF1::Class())) {
1383 if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
1384 fit = (TF1*)f;
1385 } else {
1386 f->Paint();
1387 }
1388 }
1389 }
1390
1391 if (gfit && fit) gfit->PaintStats(fit);
1392}
1393
1394
1395////////////////////////////////////////////////////////////////////////////////
1396/// Divides the active pad and draws all Graphs in the Multigraph separately.
1397
1399{
1400 TIter next(fGraphs);
1401 Int_t neededPads = fGraphs->GetSize();
1402 Int_t existingPads = 0;
1403 TString opt = (TString)option;
1404
1405 TVirtualPad *curPad = gPad;
1406 TObject *obj;
1407 TIter nextPad(curPad->GetListOfPrimitives());
1408
1409 while ((obj = nextPad())) {
1410 if (obj->InheritsFrom(TVirtualPad::Class())) existingPads++;
1411 }
1412 if (existingPads < neededPads) {
1413 curPad->Clear();
1414 Int_t nx = (Int_t)TMath::Sqrt((Double_t)neededPads);
1415 if (nx*nx < neededPads) nx++;
1416 Int_t ny = nx;
1417 if (((nx*ny)-nx) >= neededPads) ny--;
1418 curPad->Divide(nx,ny);
1419 }
1420 Int_t i = 0;
1421 TGraph *g;
1422
1424 obj = 0;
1425
1426 while (lnk) {
1427 g = (TGraph*)lnk->GetObject();
1428 i++;
1429 curPad->cd(i);
1430 TString apopt = lnk->GetOption();
1431 if (strlen(apopt)) {
1432 g->Draw((apopt.Append("A")).Data());
1433 } else {
1434 if (strlen(opt)) g->Draw(opt.Append("A"));
1435 else g->Draw("LA");
1436 }
1437 lnk = (TObjOptLink*)lnk->Next();
1438 }
1439
1440 curPad->cd();
1441}
1442
1443
1444////////////////////////////////////////////////////////////////////////////////
1445/// Paint all the graphs of this multigraph as 3D lines.
1446
1448{
1449 Int_t i, npt=0;
1450 char *l;
1451 Double_t rwxmin=0., rwxmax=0., rwymin=0., rwymax=0.;
1452 TIter next(fGraphs);
1453 TGraph *g;
1454
1455 g = (TGraph*) next();
1456 if (g) {
1457 g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1458 npt = g->GetN();
1459 }
1460
1461 if (!fHistogram) {
1462 fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1463 }
1464
1465 while ((g = (TGraph*) next())) {
1466 Double_t rx1,ry1,rx2,ry2;
1467 g->ComputeRange(rx1, ry1, rx2, ry2);
1468 if (rx1 < rwxmin) rwxmin = rx1;
1469 if (ry1 < rwymin) rwymin = ry1;
1470 if (rx2 > rwxmax) rwxmax = rx2;
1471 if (ry2 > rwymax) rwymax = ry2;
1472 if (g->GetN() > npt) npt = g->GetN();
1473 }
1474
1475 Int_t ndiv = fGraphs->GetSize();
1476
1477 TH2F* frame = new TH2F("frame","", ndiv, 0., (Double_t)(ndiv), npt, rwxmin, rwxmax);
1478 if (fHistogram) {
1479 frame->SetTitle(fHistogram->GetTitle());
1480 frame->GetYaxis()->SetTitle(fHistogram->GetXaxis()->GetTitle());
1482 frame->GetZaxis()->SetTitle(fHistogram->GetYaxis()->GetTitle());
1483 }
1484
1485 TAxis *Xaxis = frame->GetXaxis();
1486 Xaxis->SetNdivisions(-ndiv);
1487 next.Reset();
1488 for (i=ndiv; i>=1; i--) {
1489 g = (TGraph*) next();
1490 Xaxis->SetBinLabel(i, g->GetTitle());
1491 }
1492
1493 frame->SetStats(kFALSE);
1494 if (fMinimum != -1111) frame->SetMinimum(fMinimum);
1495 else frame->SetMinimum(rwymin);
1496 if (fMaximum != -1111) frame->SetMaximum(fMaximum);
1497 else frame->SetMaximum(rwymax);
1498
1499 l = (char*)strstr(option,"A");
1500 if (l) frame->Paint("lego9,fb,bb");
1501 l = (char*)strstr(option,"BB");
1502 if (!l) frame->Paint("lego9,fb,a,same");
1503
1504 Double_t *x, *y;
1505 Double_t xyz1[3], xyz2[3];
1506
1507 Double_t xl = frame->GetYaxis()->GetBinLowEdge(frame->GetYaxis()->GetFirst());
1508 Double_t xu = frame->GetYaxis()->GetBinUpEdge(frame->GetYaxis()->GetLast());
1509 Double_t yl = frame->GetMinimum();
1510 Double_t yu = frame->GetMaximum();
1511 Double_t xc[2],yc[2];
1512 next.Reset();
1513 Int_t j = ndiv;
1514
1515 while ((g = (TGraph*) next())) {
1516 npt = g->GetN();
1517 x = g->GetX();
1518 y = g->GetY();
1519 gPad->SetLineColor(g->GetLineColor());
1520 gPad->SetLineWidth(g->GetLineWidth());
1521 gPad->SetLineStyle(g->GetLineStyle());
1522 gPad->TAttLine::Modify();
1523 for (i=0; i<npt-1; i++) {
1524 xc[0] = x[i];
1525 xc[1] = x[i+1];
1526 yc[0] = y[i];
1527 yc[1] = y[i+1];
1528 if (gPad->Clip(&xc[0], &yc[0], xl, yl, xu, yu)<2) {
1529 xyz1[0] = j-0.5;
1530 xyz1[1] = xc[0];
1531 xyz1[2] = yc[0];
1532 xyz2[0] = j-0.5;
1533 xyz2[1] = xc[1];
1534 xyz2[2] = yc[1];
1535 gPad->PaintLine3D(xyz1, xyz2);
1536 }
1537 }
1538 j--;
1539 }
1540
1541 l = (char*)strstr(option,"FB");
1542 if (!l) frame->Paint("lego9,bb,a,same");
1543 delete frame;
1544}
1545
1546
1547////////////////////////////////////////////////////////////////////////////////
1548/// Paint all the graphs of this multigraph reverting values along X and/or Y axis.
1549/// New graphs are created.
1550
1552{
1553 auto *h = GetHistogram();
1554 TH1F *hg = nullptr;
1555 TGraph *fg = nullptr;
1556 if (!h)
1557 return;
1558 TString mgopt = option;
1559 mgopt.ToLower();
1560
1561 TIter next(fGraphs);
1562 TGraph *g;
1563 Bool_t first = kTRUE;
1564 TString gopt;
1565 while ((g = (TGraph *)next())) {
1566 gopt = GetGraphDrawOption(g);
1567 gopt.Append(mgopt);
1568 if (first) {
1569 fg = g;
1570 hg = fg->GetHistogram();
1571 fg->SetHistogram(h);
1572 fg->Paint(gopt.Data());
1573 first = kFALSE;
1574 } else {
1575 g->Paint(gopt.ReplaceAll("a", "").Data());
1576 }
1577 }
1578 if (fg)
1579 fg->SetHistogram(hg);
1580}
1581
1582
1583////////////////////////////////////////////////////////////////////////////////
1584/// Print the list of graphs.
1585
1586void TMultiGraph::Print(Option_t *option) const
1587{
1588 TGraph *g;
1589 if (fGraphs) {
1590 TIter next(fGraphs);
1591 while ((g = (TGraph*) next())) {
1592 g->Print(option);
1593 }
1594 }
1595}
1596
1597
1598////////////////////////////////////////////////////////////////////////////////
1599/// Recursively remove this object from a list. Typically implemented
1600/// by classes that can contain multiple references to a same object.
1601
1603{
1604 if (!fGraphs) return;
1605 TObject *objr = fGraphs->Remove(obj);
1606 if (!objr) return;
1607 delete fHistogram; fHistogram = 0;
1608 if (gPad) gPad->Modified();
1609}
1610
1611
1612////////////////////////////////////////////////////////////////////////////////
1613/// Save primitive as a C++ statement(s) on output stream out.
1614
1615void TMultiGraph::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1616{
1617 char quote = '"';
1618 out<<" "<<std::endl;
1619 if (gROOT->ClassSaved(TMultiGraph::Class())) {
1620 out<<" ";
1621 } else {
1622 out<<" TMultiGraph *";
1623 }
1624 out<<"multigraph = new TMultiGraph();"<<std::endl;
1625 out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
1626 out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<std::endl;
1627
1628 if (fGraphs) {
1630 TObject *g;
1631
1632 while (lnk) {
1633 g = lnk->GetObject();
1634 g->SavePrimitive(out, Form("multigraph%s",lnk->GetOption()));
1635 lnk = (TObjOptLink*)lnk->Next();
1636 }
1637 }
1638 const char *l = strstr(option,"th2poly");
1639 if (l) {
1640 out<<" "<<l+7<<"->AddBin(multigraph);"<<std::endl;
1641 } else {
1642 out<<" multigraph->Draw(" <<quote<<option<<quote<<");"<<std::endl;
1643 }
1644 TAxis *xaxis = GetXaxis();
1645 TAxis *yaxis = GetYaxis();
1646
1647 if (xaxis) {
1648 out<<" multigraph->GetXaxis()->SetLimits("<<xaxis->GetXmin()<<", "<<xaxis->GetXmax()<<");"<<std::endl;
1649 xaxis->SaveAttributes(out, "multigraph","->GetXaxis()");
1650 }
1651 if (yaxis) yaxis->SaveAttributes(out, "multigraph","->GetYaxis()");
1652 if (fMinimum != -1111) out<<" multigraph->SetMinimum("<<fMinimum<<");"<<std::endl;
1653 if (fMaximum != -1111) out<<" multigraph->SetMaximum("<<fMaximum<<");"<<std::endl;
1654}
1655
1656
1657////////////////////////////////////////////////////////////////////////////////
1658/// Set multigraph maximum.
1659
1661{
1662 fMaximum = maximum;
1663 if (fHistogram) fHistogram->SetMaximum(maximum);
1664}
1665
1666
1667////////////////////////////////////////////////////////////////////////////////
1668/// Set multigraph minimum.
1669
1671{
1672 fMinimum = minimum;
1673 if (fHistogram) fHistogram->SetMinimum(minimum);
1674}
1675
1676
1677////////////////////////////////////////////////////////////////////////////////
1678/// Get iterator over internal graphs list.
1679
1681{
1682 return TIter(fGraphs);
1683}
void Class()
Definition: Class.C:29
#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
int Int_t
Definition: RtypesCore.h:41
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
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
R__EXTERN TEnv * gEnv
Definition: TEnv.h:171
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
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
#define gROOT
Definition: TROOT.h:415
char * Form(const char *fmt,...)
void Printf(const char *fmt,...)
#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 void SetNdivisions(Int_t n=510, Bool_t optim=kTRUE)
Set the number of divisions for this axis.
Definition: TAttAxis.cxx:229
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
Class to manage histogram axis.
Definition: TAxis.h:30
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition: TAxis.cxx:809
virtual Bool_t GetTimeDisplay() const
Definition: TAxis.h:126
Double_t GetXmax() const
Definition: TAxis.h:134
virtual void SaveAttributes(std::ostream &out, const char *name, const char *subname)
Save axis attributes as C++ statement(s) on output stream out.
Definition: TAxis.cxx:647
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:504
virtual void SetTimeDisplay(Int_t value)
Definition: TAxis.h:161
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition: TAxis.h:154
Double_t GetXmin() const
Definition: TAxis.h:133
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
virtual void SetRange(Int_t first=0, Int_t last=0)
Set the viewing range for the axis from bin first to last.
Definition: TAxis.cxx:903
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:514
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
Definition: TCollection.h:182
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 Int_t GetNpar() const
Definition: TF1.h:475
@ kNotDraw
Definition: TF1.h:321
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
@ kResetHisto
fHistogram must be reset in GetHistogram
Definition: TGraph.h:71
virtual void Paint(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:2018
virtual void PaintStats(TF1 *fit)
Draw the stats.
Definition: TGraph.cxx:2045
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 void SetHistogram(TH1F *h)
Definition: TGraph.h:175
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 SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6333
TAxis * GetZaxis()
Definition: TH1.h:318
@ 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
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition: TH1.cxx:7994
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 Paint(Option_t *option="")
Control routine to paint any kind of histograms.
Definition: TH1.cxx:5825
virtual Double_t GetMinimum(Double_t minval=-FLT_MAX) const
Return minimum value larger than minval of bins in the range, unless the value has been overridden by...
Definition: TH1.cxx:8079
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a line.
Definition: TH1.cxx:2736
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8434
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
void Reset()
Definition: TCollection.h:252
Iterator of linked list.
Definition: TList.h:200
Option_t * GetOption() const
Returns the object option stored in the list.
Definition: TList.cxx:1141
A doubly linked list.
Definition: TList.h:44
virtual void Add(TObject *obj)
Definition: TList.h:87
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 TObjLink * FirstLink() const
Definition: TList.h:108
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:354
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:656
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:35
Double_t fMinimum
Definition: TMultiGraph.h:42
TH1F * fHistogram
Definition: TMultiGraph.h:40
TMultiGraph()
TMultiGraph default constructor.
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.
TList * GetListOfGraphs() const
Definition: TMultiGraph.h:69
TList * fGraphs
Definition: TMultiGraph.h:38
virtual void Add(TGraph *graph, Option_t *chopt="")
Add a new graph to the list of graphs.
Double_t fMaximum
Definition: TMultiGraph.h:41
virtual void SetMinimum(Double_t minimum=-1111)
Set multigraph minimum.
TH1F * GetHistogram()
Returns a pointer to the histogram used to draw the axis.
TF1 * GetFunction(const char *name) const
Return pointer to function with name.
virtual void LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
Least squares lpolynomial fitting without weights.
virtual void InitPolynom(Double_t xmin, Double_t xmax)
Compute Initial values of parameters for a polynom.
void PaintPolyLine3D(Option_t *chopt="")
Paint all the graphs of this multigraph as 3D lines.
virtual void InitExpo(Double_t xmin, Double_t xmax)
Compute Initial values of parameters for an exponential.
virtual void FitPanel()
Display a panel with all histogram fit options.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to each graph.
virtual void Paint(Option_t *chopt="")
Paint all the graphs of this multigraph.
virtual void Draw(Option_t *chopt="")
Draw this multigraph with its current attributes.
TIter begin() const
Get iterator over internal graphs list.
TMultiGraph & operator=(const TMultiGraph &)
Assignment operator.
virtual void LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
Least square linear fit without weights.
void PaintPads(Option_t *chopt="")
Divides the active pad and draws all Graphs in the Multigraph separately.
virtual Option_t * GetGraphDrawOption(const TGraph *gr) const
Return the draw option for the TGraph gr in this TMultiGraph.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
TAxis * GetYaxis()
Get y axis of the graph.
virtual void InitGaus(Double_t xmin, Double_t xmax)
Compute Initial values of parameters for a gaussian.
virtual ~TMultiGraph()
TMultiGraph destructor.
virtual Int_t IsInside(Double_t x, Double_t y) const
Return 1 if the point (x,y) is inside one of the graphs 0 otherwise.
TList * GetListOfFunctions()
Return pointer to list of functions.
virtual void Print(Option_t *chopt="") const
Print the list of graphs.
virtual void SetMaximum(Double_t maximum=-1111)
Set multigraph maximum.
virtual void Browse(TBrowser *b)
Browse multigraph.
void PaintReverse(Option_t *chopt="")
Paint all the graphs of this multigraph reverting values along X and/or Y axis.
TList * fFunctions
Definition: TMultiGraph.h:39
TAxis * GetXaxis()
Get x axis of the graph.
virtual void RecursiveRemove(TObject *obj)
Recursively remove this object from a list.
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
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
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
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
virtual void Paint(Option_t *option="")
This method must be overridden if a class wants to paint itself.
Definition: TObject.cxx:519
@ kInvalidObject
if object ctor succeeded but object should not be used
Definition: TObject.h:68
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition: TObject.h:60
Long_t ExecPlugin(int nargs, const T &... params)
Int_t LoadPlugin()
Load the plugin library for this handler.
Basic string class.
Definition: TString.h:131
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
const char * Data() const
Definition: TString.h:364
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
Bool_t IsNull() const
Definition: TString.h:402
TString & Append(const char *cs)
Definition: TString.h:559
Bool_t IsWhitespace() const
Definition: TString.h:403
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
Abstract Base Class for Fitting.
static TVirtualFitter * GetFitter()
static: return the current Fitter
virtual TObject * GetUserFunc() const
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:50
virtual TList * GetListOfPrimitives() const =0
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)=0
virtual void Clear(Option_t *option="")=0
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
TH1F * h1
Definition: legend1.C:5
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 dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
static constexpr double mg
Double_t Log(Double_t x)
Definition: TMath.h:750
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
Double_t Log10(Double_t x)
Definition: TMath.h:754
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
Definition: first.py:1
Definition: graph.py:1
auto * m
Definition: textangle.C:8
auto * l
Definition: textangle.C:4
auto * a
Definition: textangle.C:12