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