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