Logo ROOT   6.12/07
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"
31 #include "Math/MinimizerOptions.h"
32 
33 #include <ctype.h>
34 
35 extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
36 
38 
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 
42 /** \class TMultiGraph
43  \ingroup Hist
44 A TMultiGraph is a collection of TGraph (or derived) objects. It allows to
45 manipulate a set of graphs as a single entity. In particular, when drawn,
46 the X and Y axis ranges are automatically computed such as all the graphs
47 will be visible.
48 
49 `TMultiGraph::Add` should be used to add a new graph to the list.
50 
51 The TMultiGraph owns the objects in the list.
52 
53 The drawing options are the same as for TGraph.
54 Like for TGraph, the painting is performed thanks to the TGraphPainter
55 class. All details about the various painting options are given in this class.
56 
57 Example:
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 ~~~
66 A special option `3D` allows to draw the graphs in a 3D space. See the
67 following example:
68 
69 Begin_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 }
104 End_Macro
105 
106 
107 The number of graphs in a multigraph can be retrieve with:
108 ~~~ {.cpp}
109 mg->GetListOfGraphs()->GetSize();
110 ~~~
111 
112 The drawing option for each TGraph may be specified as an optional
113 second argument of the `Add` function.
114 
115 If a draw option is specified, it will be used to draw the graph,
116 otherwise the graph will be drawn with the option specified in
117 `TMultiGraph::Draw`.
118 
119 The following example shows how to fit a TMultiGraph.
120 
121 Begin_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 }
153 End_Macro
154 
155 
156 The 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 
167 When the graphs in a TMultiGraph are fitted, the fit parameters boxes
168 overlap. The following example shows how to make them all visible.
169 
170 
171 Begin_Macro(source)
172 ../../../tutorials/graphs/multigraph.C
173 End_Macro
174 
175 
176 The axis limits can be changed the like for TGraph. The same methods apply on
177 the multigraph.
178 Note the two differents ways to change limits on X and Y axis.
179 
180 Begin_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 }
205 End_Macro
206 
207 
208 The method TPad::BuildLegend is able to extract the graphs inside a
209 multigraph. The following example demonstrate this.
210 
211 Begin_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 }
265 End_Macro
266 
267 Automatic coloring according to the current palette is available as shown in the
268 following example:
269 
270 Begin_Macro(source)
271 ../../../tutorials/graphs/multigraphpalettecolor.C
272 End_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 
292 TMultiGraph::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) {
323  TNamed::operator=(mg);
324  fGraphs=mg.fGraphs;
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())) {
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 
372 void TMultiGraph::Add(TGraph *graph, Option_t *chopt)
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 
388 void 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 
462 void TMultiGraph::Draw(Option_t *option)
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 
480 TFitResultPtr TMultiGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
481 {
482  char *linear;
483  linear= (char*)strstr(fname, "++");
484  TF1 *f1=0;
485  if (linear)
486  f1=new TF1(fname, fname, xmin, xmax);
487  else {
488  f1 = (TF1*)gROOT->GetFunction(fname);
489  if (!f1) { Printf("Unknown function: %s",fname); return -1; }
490  }
491 
492  return Fit(f1,option,"",xmin,xmax);
493 }
494 
495 
496 ////////////////////////////////////////////////////////////////////////////////
497 /// Fit this multigraph with function f1.
498 ///
499 /// In this function all graphs of the multigraph are fitted simultaneously
500 ///
501 /// f1 is an already predefined function created by TF1.
502 /// Predefined functions such as gaus, expo and poln are automatically
503 /// created by ROOT.
504 ///
505 /// The list of fit options is given in parameter `option`which may takes the
506 /// following values:
507 ///
508 /// - "W" Set all errors to 1
509 /// - "U" Use a User specified fitting algorithm (via SetFCN)
510 /// - "Q" Quiet mode (minimum printing)
511 /// - "V" Verbose mode (default is between Q and V)
512 /// - "B" Use this option when you want to fix one or more parameters
513 /// and the fitting function is like "gaus","expo","poln","landau".
514 /// - "R" Use the Range specified in the function range
515 /// - "N" Do not store the graphics function, do not draw
516 /// - "0" Do not plot the result of the fit. By default the fitted function
517 /// is drawn unless the option"N" above is specified.
518 /// - "+" Add this new fitted function to the list of fitted functions
519 /// (by default, any previous function is deleted)
520 /// - "C" In case of linear fitting, not calculate the chisquare (saves time)
521 /// - "F" If fitting a polN, switch to minuit fitter
522 /// - "ROB" In case of linear fitting, compute the LTS regression
523 /// coefficients (robust(resistant) regression), using
524 /// the default fraction of good points
525 /// - "ROB=0.x" - compute the LTS regression coefficients, using
526 /// 0.x as a fraction of good points
527 ///
528 /// When the fit is drawn (by default), the parameter goption may be used
529 /// to specify a list of graphics options. See TGraph::Paint for a complete
530 /// list of these options.
531 ///
532 /// In order to use the Range option, one must first create a function
533 /// with the expression to be fitted. For example, if your graph
534 /// has a defined range between -4 and 4 and you want to fit a gaussian
535 /// only in the interval 1 to 3, you can do:
536 /// ~~~ {.cpp}
537 /// TF1 *f1 = new TF1("f1","gaus",1,3);
538 /// graph->Fit("f1","R");
539 /// ~~~
540 ///
541 /// ### Who is calling this function ?
542 ///
543 /// Note that this function is called when calling TGraphErrors::Fit
544 /// or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
545 /// see the discussion below on the errors calculation.
546 ///
547 /// ### Setting initial conditions
548 ///
549 /// Parameters must be initialized before invoking the Fit function.
550 /// The setting of the parameter initial values is automatic for the
551 /// predefined functions : poln, expo, gaus, landau. One can however disable
552 /// this automatic computation by specifying the option "B".
553 /// You can specify boundary limits for some or all parameters via
554 /// ~~~ {.cpp}
555 /// f1->SetParLimits(p_number, parmin, parmax);
556 /// ~~~
557 /// if `parmin>=parmax`, the parameter is fixed
558 /// Note that you are not forced to fix the limits for all parameters.
559 /// For example, if you fit a function with 6 parameters, you can do:
560 /// ~~~ {.cpp}
561 /// func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
562 /// func->SetParLimits(4,-10,-4);
563 /// func->SetParLimits(5, 1,1);
564 /// ~~~
565 /// With this setup, parameters 0->3 can vary freely
566 /// Parameter 4 has boundaries [-10,-4] with initial value -8
567 /// Parameter 5 is fixed to 100.
568 ///
569 /// ### Fit range
570 ///
571 /// The fit range can be specified in two ways:
572 ///
573 /// - specify rxmax > rxmin (default is rxmin=rxmax=0)
574 /// - specify the option "R". In this case, the function will be taken
575 /// instead of the full graph range.
576 ///
577 /// ### Changing the fitting function
578 ///
579 /// By default a chi2 fitting function is used for fitting the TGraphs's.
580 /// The function is implemented in `FitUtil::EvaluateChi2`.
581 /// In case of TGraphErrors an effective chi2 is used
582 /// (see TGraphErrors fit in TGraph::Fit) and is implemented in
583 /// `FitUtil::EvaluateChi2Effective`
584 /// To specify a User defined fitting function, specify option "U" and
585 /// call the following function:
586 /// ~~~ {.cpp}
587 /// TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
588 /// ~~~
589 /// where MyFittingFunction is of type:
590 /// ~~~ {.cpp}
591 /// extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
592 /// ~~~
593 ///
594 /// ### Access to the fit result
595 ///
596 /// The function returns a TFitResultPtr which can hold a pointer to a TFitResult object.
597 /// By default the TFitResultPtr contains only the status of the fit and it converts
598 /// automatically to an integer. If the option "S" is instead used, TFitResultPtr contains
599 /// the TFitResult and behaves as a smart pointer to it. For example one can do:
600 /// ~~~ {.cpp}
601 /// TFitResultPtr r = graph->Fit("myFunc","S");
602 /// TMatrixDSym cov = r->GetCovarianceMatrix(); // to access the covariance matrix
603 /// Double_t par0 = r->Parameter(0); // retrieve the value for the parameter 0
604 /// Double_t err0 = r->ParError(0); // retrieve the error for the parameter 0
605 /// r->Print("V"); // print full information of fit including covariance matrix
606 /// r->Write(); // store the result in a file
607 /// ~~~
608 ///
609 /// The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also
610 /// from the fitted function.
611 ///
612 /// ### Associated functions
613 ///
614 /// One or more object (typically a TF1*) can be added to the list
615 /// of functions (fFunctions) associated to each graph.
616 /// When TGraph::Fit is invoked, the fitted function is added to this list.
617 /// Given a graph gr, one can retrieve an associated function
618 /// with:
619 /// ~~~ {.cpp}
620 /// TF1 *myfunc = gr->GetFunction("myfunc");
621 /// ~~~
622 ///
623 /// If the graph is made persistent, the list of
624 /// associated functions is also persistent. Given a pointer (see above)
625 /// to an associated function myfunc, one can retrieve the function/fit
626 /// parameters with calls such as:
627 /// ~~~ {.cpp}
628 /// Double_t chi2 = myfunc->GetChisquare();
629 /// Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
630 /// Double_t err0 = myfunc->GetParError(0); //error on first parameter
631 /// ~~~
632 ///
633 /// ### Fit Statistics
634 ///
635 /// You can change the statistics box to display the fit parameters with
636 /// the TStyle::SetOptFit(mode) method. This mode has four digits.
637 /// mode = pcev (default = 0111)
638 ///
639 /// - v = 1; print name/values of parameters
640 /// - e = 1; print errors (if e=1, v must be 1)
641 /// - c = 1; print Chisquare/Number of degrees of freedom
642 /// - p = 1; print Probability
643 ///
644 /// For example: `gStyle->SetOptFit(1011);`
645 /// prints the fit probability, parameter names/values, and errors.
646 /// You can change the position of the statistics box with these lines
647 /// (where g is a pointer to the TGraph):
648 ///
649 /// ~~~ {.cpp}
650 /// Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
651 /// Root > st->SetX1NDC(newx1); //new x start position
652 /// Root > st->SetX2NDC(newx2); //new x end position
653 /// ~~~
654 
655 TFitResultPtr TMultiGraph::Fit(TF1 *f1, Option_t *option, Option_t *goption, Axis_t rxmin, Axis_t rxmax)
656 {
657  // internal multigraph fitting methods
658  Foption_t fitOption;
660 
661  // create range and minimizer options with default values
662  ROOT::Fit::DataRange range(rxmin,rxmax);
664  return ROOT::Fit::FitObject(this, f1 , fitOption , minOption, goption, range);
665 
666 }
667 
668 ////////////////////////////////////////////////////////////////////////////////
669 /// Display a panel with all histogram fit options.
670 /// See class TFitPanel for example
671 
673 {
674  if (!gPad)
675  gROOT->MakeDefCanvas();
676 
677  if (!gPad) {
678  Error("FitPanel", "Unable to create a default canvas");
679  return;
680  }
681 
682  // use plugin manager to create instance of TFitEditor
683  TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
684  if (handler && handler->LoadPlugin() != -1) {
685  if (handler->ExecPlugin(2, gPad, this) == 0)
686  Error("FitPanel", "Unable to crate the FitPanel");
687  }
688  else
689  Error("FitPanel", "Unable to find the FitPanel plug-in");
690 }
691 
692 ////////////////////////////////////////////////////////////////////////////////
693 /// Return the draw option for the TGraph `gr` in this TMultiGraph.
694 /// The return option is the one specified when calling TMultiGraph::Add(gr,option).
695 
697 {
698  if (!fGraphs || !gr) return "";
699  TListIter next(fGraphs);
700  TObject *obj;
701  while ((obj = next())) {
702  if (obj == (TObject*)gr) return next.GetOption();
703  }
704  return "";
705 }
706 
707 
708 ////////////////////////////////////////////////////////////////////////////////
709 /// Compute Initial values of parameters for a gaussian.
710 
711 void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
712 {
713  Double_t allcha, sumx, sumx2, x, val, rms, mean;
714  Int_t bin;
715  const Double_t sqrtpi = 2.506628;
716 
717  // Compute mean value and RMS of the graph in the given range
718  Int_t np = 0;
719  allcha = sumx = sumx2 = 0;
720  TGraph *g;
721  TIter next(fGraphs);
722  Double_t *px, *py;
723  Int_t npp; //number of points in each graph
724  while ((g = (TGraph*) next())) {
725  px=g->GetX();
726  py=g->GetY();
727  npp=g->GetN();
728  for (bin=0; bin<npp; bin++) {
729  x=px[bin];
730  if (x<xmin || x>xmax) continue;
731  np++;
732  val=py[bin];
733  sumx+=val*x;
734  sumx2+=val*x*x;
735  allcha+=val;
736  }
737  }
738  if (np == 0 || allcha == 0) return;
739  mean = sumx/allcha;
740  rms = TMath::Sqrt(sumx2/allcha - mean*mean);
741 
742  Double_t binwidx = TMath::Abs((xmax-xmin)/np);
743  if (rms == 0) rms = 1;
745  TF1 *f1 = (TF1*)grFitter->GetUserFunc();
746  f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
747  f1->SetParameter(1,mean);
748  f1->SetParameter(2,rms);
749  f1->SetParLimits(2,0,10*rms);
750 }
751 
752 
753 ////////////////////////////////////////////////////////////////////////////////
754 /// Compute Initial values of parameters for an exponential.
755 
756 void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
757 {
758  Double_t constant, slope;
759  Int_t ifail;
760 
761  LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
762 
764  TF1 *f1 = (TF1*)grFitter->GetUserFunc();
765  f1->SetParameter(0,constant);
766  f1->SetParameter(1,slope);
767 }
768 
769 
770 ////////////////////////////////////////////////////////////////////////////////
771 /// Compute Initial values of parameters for a polynom.
772 
774 {
775  Double_t fitpar[25];
776 
778  TF1 *f1 = (TF1*)grFitter->GetUserFunc();
779  Int_t npar = f1->GetNpar();
780 
781  LeastSquareFit(npar, fitpar, xmin, xmax);
782 
783  for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
784 }
785 
786 
787 ////////////////////////////////////////////////////////////////////////////////
788 /// Least squares lpolynomial fitting without weights.
789 ///
790 /// - m number of parameters
791 /// - a array of parameters
792 /// - first 1st point number to fit (default =0)
793 /// - last last point number to fit (default=fNpoints-1)
794 ///
795 /// based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
796 
798 {
799  const Double_t zero = 0.;
800  const Double_t one = 1.;
801  const Int_t idim = 20;
802 
803  Double_t b[400] /* was [20][20] */;
804  Int_t i, k, l, ifail, bin;
805  Double_t power;
806  Double_t da[20], xk, yk;
807 
808 
809  //count the total number of points to fit
810  TGraph *g;
811  TIter next(fGraphs);
812  Double_t *px, *py;
813  Int_t n=0;
814  Int_t npp;
815  while ((g = (TGraph*) next())) {
816  px=g->GetX();
817  npp=g->GetN();
818  for (bin=0; bin<npp; bin++) {
819  xk=px[bin];
820  if (xk < xmin || xk > xmax) continue;
821  n++;
822  }
823  }
824  if (m <= 2) {
825  LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
826  return;
827  }
828  if (m > idim || m > n) return;
829  da[0] = zero;
830  for (l = 2; l <= m; ++l) {
831  b[l-1] = zero;
832  b[m + l*20 - 21] = zero;
833  da[l-1] = zero;
834  }
835  Int_t np = 0;
836 
837  next.Reset();
838  while ((g = (TGraph*) next())) {
839  px=g->GetX();
840  py=g->GetY();
841  npp=g->GetN();
842 
843  for (k = 0; k <= npp; ++k) {
844  xk = px[k];
845  if (xk < xmin || xk > xmax) continue;
846  np++;
847  yk = py[k];
848  power = one;
849  da[0] += yk;
850  for (l = 2; l <= m; ++l) {
851  power *= xk;
852  b[l-1] += power;
853  da[l-1] += power*yk;
854  }
855  for (l = 2; l <= m; ++l) {
856  power *= xk;
857  b[m + l*20 - 21] += power;
858  }
859  }
860  }
861  b[0] = Double_t(np);
862  for (i = 3; i <= m; ++i) {
863  for (k = i; k <= m; ++k) {
864  b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
865  }
866  }
867  H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
868 
869  if (ifail < 0) {
870  //a[0] = fY[0];
871  py=((TGraph *)fGraphs->First())->GetY();
872  a[0]=py[0];
873  for (i=1; i<m; ++i) a[i] = 0;
874  return;
875  }
876  for (i=0; i<m; ++i) a[i] = da[i];
877 }
878 
879 
880 ////////////////////////////////////////////////////////////////////////////////
881 /// Least square linear fit without weights.
882 ///
883 /// Fit a straight line (a0 + a1*x) to the data in this graph.
884 ///
885 /// - ndata: number of points to fit
886 /// - first: first point number to fit
887 /// - last: last point to fit O(ndata should be last-first
888 /// - ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
889 ///
890 /// extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
891 
893  Int_t &ifail, Double_t xmin, Double_t xmax)
894 {
895  Double_t xbar, ybar, x2bar;
896  Int_t i;
897  Double_t xybar;
898  Double_t fn, xk, yk;
899  Double_t det;
900 
901  ifail = -2;
902  xbar = ybar = x2bar = xybar = 0;
903  Int_t np = 0;
904  TGraph *g;
905  TIter next(fGraphs);
906  Double_t *px, *py;
907  Int_t npp;
908  while ((g = (TGraph*) next())) {
909  px=g->GetX();
910  py=g->GetY();
911  npp=g->GetN();
912  for (i = 0; i < npp; ++i) {
913  xk = px[i];
914  if (xk < xmin || xk > xmax) continue;
915  np++;
916  yk = py[i];
917  if (ndata < 0) {
918  if (yk <= 0) yk = 1e-9;
919  yk = TMath::Log(yk);
920  }
921  xbar += xk;
922  ybar += yk;
923  x2bar += xk*xk;
924  xybar += xk*yk;
925  }
926  }
927  fn = Double_t(np);
928  det = fn*x2bar - xbar*xbar;
929  ifail = -1;
930  if (det <= 0) {
931  if (fn > 0) a0 = ybar/fn;
932  else a0 = 0;
933  a1 = 0;
934  return;
935  }
936  ifail = 0;
937  a0 = (x2bar*ybar - xbar*xybar) / det;
938  a1 = (fn*xybar - xbar*ybar) / det;
939 }
940 
941 
942 ////////////////////////////////////////////////////////////////////////////////
943 /// Return 1 if the point (x,y) is inside one of the graphs 0 otherwise.
944 
946 {
947  Int_t in = 0;
948  if (!fGraphs) return in;
949  TGraph *g;
950  TIter next(fGraphs);
951  while ((g = (TGraph*) next())) {
952  in = g->IsInside(x, y);
953  if (in) return in;
954  }
955  return in;
956 }
957 
958 
959 ////////////////////////////////////////////////////////////////////////////////
960 /// Returns a pointer to the histogram used to draw the axis.
961 /// Takes into account the two following cases.
962 ///
963 /// 1. option 'A' was specified in TMultiGraph::Draw. Return fHistogram
964 /// 2. user had called TPad::DrawFrame. return pointer to hframe histogram
965 
967 {
968  if (fHistogram) return fHistogram;
969 
970  if (gPad) {
971  gPad->Modified();
972  gPad->Update();
973  if (fHistogram) return fHistogram;
974  TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
975  return h1;
976  } else {
977  Bool_t initialrangeset = kFALSE;
978  Double_t rwxmin = 0.,rwxmax = 0.,rwymin = 0.,rwymax = 0.;
979  TGraph *g;
980  Int_t npt = 100 ;
981  TIter next(fGraphs);
982  while ((g = (TGraph*) next())) {
983  if (g->GetN() <= 0) continue;
984  if (initialrangeset) {
985  Double_t rx1,ry1,rx2,ry2;
986  g->ComputeRange(rx1, ry1, rx2, ry2);
987  if (rx1 < rwxmin) rwxmin = rx1;
988  if (ry1 < rwymin) rwymin = ry1;
989  if (rx2 > rwxmax) rwxmax = rx2;
990  if (ry2 > rwymax) rwymax = ry2;
991  } else {
992  g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
993  initialrangeset = kTRUE;
994  }
995  if (g->GetN() > npt) npt = g->GetN();
996  }
997  fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
998  if (!fHistogram) return 0;
999  fHistogram->SetMinimum(rwymin);
1001  fHistogram->SetMaximum(rwymax);
1002  fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1004  return fHistogram;
1005  }
1006 }
1007 
1008 
1009 ////////////////////////////////////////////////////////////////////////////////
1010 /// Return pointer to function with name.
1011 ///
1012 /// Functions such as TGraph::Fit store the fitted function in the list of
1013 /// functions of this graph.
1014 
1015 TF1 *TMultiGraph::GetFunction(const char *name) const
1016 {
1017  if (!fFunctions) return 0;
1018  return (TF1*)fFunctions->FindObject(name);
1019 }
1020 
1021 ////////////////////////////////////////////////////////////////////////////////
1022 /// Return pointer to list of functions.
1023 /// If pointer is null create the list
1024 
1026 {
1027  if (!fFunctions) fFunctions = new TList();
1028  return fFunctions;
1029 }
1030 
1031 
1032 ////////////////////////////////////////////////////////////////////////////////
1033 /// Get x axis of the graph.
1034 /// This method returns a valid axis only after the TMultigraph has been drawn.
1035 
1037 {
1039  if (!h) return 0;
1040  return h->GetXaxis();
1041 }
1042 
1043 
1044 ////////////////////////////////////////////////////////////////////////////////
1045 /// Get y axis of the graph.
1046 /// This method returns a valid axis only after the TMultigraph has been drawn.
1047 
1049 {
1051  if (!h) return 0;
1052  return h->GetYaxis();
1053 }
1054 
1055 
1056 ////////////////////////////////////////////////////////////////////////////////
1057 /// Paint all the graphs of this multigraph.
1058 
1059 void TMultiGraph::Paint(Option_t *choptin)
1060 {
1061  const TPickerStackGuard pushGuard(this);
1062 
1063  if (!fGraphs) return;
1064  if (fGraphs->GetSize() == 0) return;
1065 
1066  char option[128];
1067  strlcpy(option,choptin,128);
1068  Int_t nch = strlen(choptin);
1069  for (Int_t i=0;i<nch;i++) option[i] = toupper(option[i]);
1070 
1071  // Automatic color
1072  char *l1 = strstr(option,"PFC"); // Automatic Fill Color
1073  char *l2 = strstr(option,"PLC"); // Automatic Line Color
1074  char *l3 = strstr(option,"PMC"); // Automatic Marker Color
1075  if (l1 || l2 || l3) {
1076  TString opt1 = option; opt1.ToLower();
1077  if (l1) strncpy(l1," ",3);
1078  if (l2) strncpy(l2," ",3);
1079  if (l3) strncpy(l3," ",3);
1081  TGraph* gAti;
1082  Int_t ngraphs = fGraphs->GetSize();
1083  Int_t ic;
1084  gPad->IncrementPaletteColor(ngraphs, opt1);
1085  for (Int_t i=0;i<ngraphs;i++) {
1086  ic = gPad->NextPaletteColor();
1087  gAti = (TGraph*)(fGraphs->At(i));
1088  if (l1) gAti->SetFillColor(ic);
1089  if (l2) gAti->SetLineColor(ic);
1090  if (l3) gAti->SetMarkerColor(ic);
1091  lnk = (TObjOptLink*)lnk->Next();
1092  }
1093  }
1094 
1095  char *l;
1096 
1097  TString chopt = option;
1098 
1099  l = (char*)strstr(chopt.Data(),"3D");
1100  if (l) {
1101  l = (char*)strstr(chopt.Data(),"L");
1102  if (l) PaintPolyLine3D(chopt.Data());
1103  return;
1104  }
1105 
1106  l = (char*)strstr(chopt.Data(),"PADS");
1107  if (l) {
1108  chopt.ReplaceAll("PADS","");
1109  PaintPads(chopt.Data());
1110  return;
1111  }
1112 
1113  TGraph *g;
1114 
1115  l = (char*)strstr(chopt.Data(),"A");
1116  if (l) {
1117  *l = ' ';
1118  TIter next(fGraphs);
1119  Int_t npt = 100;
1120  Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
1121  rwxmin = gPad->GetUxmin();
1122  rwxmax = gPad->GetUxmax();
1123  rwymin = gPad->GetUymin();
1124  rwymax = gPad->GetUymax();
1125  char *xtitle = 0;
1126  char *ytitle = 0;
1127  Int_t firstx = 0;
1128  Int_t lastx = 0;
1129  Bool_t timedisplay = kFALSE;
1130  char *timeformat = 0;
1131 
1132  if (fHistogram) {
1133  //cleanup in case of a previous unzoom and in case one of the TGraph has changed
1135  TGraph* gAti;
1136  Int_t ngraphs = fGraphs->GetSize();
1137  Bool_t reset_hist = kFALSE;
1138  for (Int_t i=0;i<ngraphs;i++) {
1139  gAti = (TGraph*)(fGraphs->At(i));
1140  if(gAti->TestBit(TGraph::kResetHisto)) {reset_hist = kTRUE; break;}
1141  lnk = (TObjOptLink*)lnk->Next();
1142  }
1143  if (fHistogram->GetMinimum() >= fHistogram->GetMaximum() || reset_hist) {
1144  nch = strlen(fHistogram->GetXaxis()->GetTitle());
1145  firstx = fHistogram->GetXaxis()->GetFirst();
1146  lastx = fHistogram->GetXaxis()->GetLast();
1147  timedisplay = fHistogram->GetXaxis()->GetTimeDisplay();
1148  if (nch) {
1149  xtitle = new char[nch+1];
1150  strlcpy(xtitle,fHistogram->GetXaxis()->GetTitle(),nch+1);
1151  }
1152  nch = strlen(fHistogram->GetYaxis()->GetTitle());
1153  if (nch) {
1154  ytitle = new char[nch+1];
1155  strlcpy(ytitle,fHistogram->GetYaxis()->GetTitle(),nch+1);
1156  }
1157  nch = strlen(fHistogram->GetXaxis()->GetTimeFormat());
1158  if (nch) {
1159  timeformat = new char[nch+1];
1160  strlcpy(timeformat,fHistogram->GetXaxis()->GetTimeFormat(),nch+1);
1161  }
1162  delete fHistogram;
1163  fHistogram = 0;
1164  }
1165  }
1166  if (fHistogram) {
1167  minimum = fHistogram->GetYaxis()->GetXmin();
1168  maximum = fHistogram->GetYaxis()->GetXmax();
1169  uxmin = gPad->PadtoX(rwxmin);
1170  uxmax = gPad->PadtoX(rwxmax);
1171  } else {
1172  Bool_t initialrangeset = kFALSE;
1173  while ((g = (TGraph*) next())) {
1174  if (g->GetN() <= 0) continue;
1175  if (initialrangeset) {
1176  Double_t rx1,ry1,rx2,ry2;
1177  g->ComputeRange(rx1, ry1, rx2, ry2);
1178  if (rx1 < rwxmin) rwxmin = rx1;
1179  if (ry1 < rwymin) rwymin = ry1;
1180  if (rx2 > rwxmax) rwxmax = rx2;
1181  if (ry2 > rwymax) rwymax = ry2;
1182  } else {
1183  g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1184  initialrangeset = kTRUE;
1185  }
1186  if (g->GetN() > npt) npt = g->GetN();
1187  }
1188  if (rwxmin == rwxmax) rwxmax += 1.;
1189  if (rwymin == rwymax) rwymax += 1.;
1190  dx = 0.05*(rwxmax-rwxmin);
1191  dy = 0.05*(rwymax-rwymin);
1192  uxmin = rwxmin - dx;
1193  uxmax = rwxmax + dx;
1194  if (gPad->GetLogy()) {
1195  if (rwymin <= 0) rwymin = 0.001*rwymax;
1196  minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
1197  maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
1198  } else {
1199  minimum = rwymin - dy;
1200  maximum = rwymax + dy;
1201  }
1202  if (minimum < 0 && rwymin >= 0) minimum = 0;
1203  if (maximum > 0 && rwymax <= 0) maximum = 0;
1204  }
1205 
1206  if (fMinimum != -1111) rwymin = minimum = fMinimum;
1207  if (fMaximum != -1111) rwymax = maximum = fMaximum;
1208  if (uxmin < 0 && rwxmin >= 0) {
1209  if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
1210  //else uxmin = 0;
1211  }
1212  if (uxmax > 0 && rwxmax <= 0) {
1213  if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
1214  //else uxmax = 0;
1215  }
1216  if (minimum < 0 && rwymin >= 0) {
1217  if (gPad->GetLogy()) minimum = 0.9*rwymin;
1218  //else minimum = 0;
1219  }
1220  if (maximum > 0 && rwymax <= 0) {
1221  if (gPad->GetLogy()) maximum = 1.1*rwymax;
1222  //else maximum = 0;
1223  }
1224  if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
1225  if (uxmin <= 0 && gPad->GetLogx()) {
1226  if (uxmax > 1000) uxmin = 1;
1227  else uxmin = 0.001*uxmax;
1228  }
1229  rwymin = minimum;
1230  rwymax = maximum;
1231  if (fHistogram) {
1232  fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1233  }
1234 
1235  // Create a temporary histogram to draw the axis
1236  if (!fHistogram) {
1237  // the graph is created with at least as many channels as there are points
1238  // to permit zooming on the full range
1239  rwxmin = uxmin;
1240  rwxmax = uxmax;
1241  fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1242  if (!fHistogram) return;
1243  fHistogram->SetMinimum(rwymin);
1245  fHistogram->SetMaximum(rwymax);
1246  fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1248  if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
1249  if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
1250  if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
1251  if (timedisplay) {fHistogram->GetXaxis()->SetTimeDisplay(timedisplay);}
1252  if (timeformat) {fHistogram->GetXaxis()->SetTimeFormat(timeformat); delete [] timeformat;}
1253  }
1254  fHistogram->Paint("0");
1255  }
1256 
1257  TGraph *gfit = 0;
1258  if (fGraphs) {
1260  TObject *obj = 0;
1261 
1262  chopt.ReplaceAll("A","");
1263 
1264  while (lnk) {
1265 
1266  obj = lnk->GetObject();
1267 
1268  gPad->PushSelectableObject(obj);
1269 
1270  if (!gPad->PadInHighlightMode() || (gPad->PadInHighlightMode() && obj == gPad->GetSelected())) {
1271  TString opt = lnk->GetOption();
1272  if (!opt.IsWhitespace())
1273  obj->Paint(opt.ReplaceAll("A","").Data());
1274  else {
1275  if (!chopt.IsWhitespace()) obj->Paint(chopt.Data());
1276  else obj->Paint("L");
1277  }
1278  }
1279 
1280  lnk = (TObjOptLink*)lnk->Next();
1281  }
1282 
1283  gfit = (TGraph*)obj; // pick one TGraph in the list to paint the fit parameters.
1284  }
1285 
1286  TObject *f;
1287  TF1 *fit = 0;
1288  if (fFunctions) {
1289  TIter next(fFunctions);
1290  while ((f = (TObject*) next())) {
1291  if (f->InheritsFrom(TF1::Class())) {
1292  if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
1293  fit = (TF1*)f;
1294  } else {
1295  f->Paint();
1296  }
1297  }
1298  }
1299 
1300  if (fit) gfit->PaintStats(fit);
1301 }
1302 
1303 
1304 ////////////////////////////////////////////////////////////////////////////////
1305 /// Divides the active pad and draws all Graphs in the Multigraph separately.
1306 
1307 void TMultiGraph::PaintPads(Option_t *option)
1308 {
1310  Int_t neededPads = fGraphs->GetSize();
1311  Int_t existingPads = 0;
1312  TString opt = (TString)option;
1313 
1314  TVirtualPad *curPad = gPad;
1315  TObject *obj;
1316  TIter nextPad(curPad->GetListOfPrimitives());
1317 
1318  while ((obj = nextPad())) {
1319  if (obj->InheritsFrom(TVirtualPad::Class())) existingPads++;
1320  }
1321  if (existingPads < neededPads) {
1322  curPad->Clear();
1323  Int_t nx = (Int_t)TMath::Sqrt((Double_t)neededPads);
1324  if (nx*nx < neededPads) nx++;
1325  Int_t ny = nx;
1326  if (((nx*ny)-nx) >= neededPads) ny--;
1327  curPad->Divide(nx,ny);
1328  }
1329  Int_t i = 0;
1330  TGraph *g;
1331 
1333  obj = 0;
1334 
1335  while (lnk) {
1336  g = (TGraph*)lnk->GetObject();
1337  i++;
1338  curPad->cd(i);
1339  TString apopt = lnk->GetOption();
1340  if (strlen(apopt)) {
1341  g->Draw((apopt.Append("A")).Data());
1342  } else {
1343  if (strlen(opt)) g->Draw(opt.Append("A"));
1344  else g->Draw("LA");
1345  }
1346  lnk = (TObjOptLink*)lnk->Next();
1347  }
1348 
1349  curPad->cd();
1350 }
1351 
1352 
1353 ////////////////////////////////////////////////////////////////////////////////
1354 /// Paint all the graphs of this multigraph as 3D lines.
1355 
1357 {
1358  Int_t i, npt=0;
1359  char *l;
1360  Double_t rwxmin=0., rwxmax=0., rwymin=0., rwymax=0.;
1361  TIter next(fGraphs);
1362  TGraph *g;
1363 
1364  g = (TGraph*) next();
1365  if (g) {
1366  g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1367  npt = g->GetN();
1368  }
1369 
1370  if (!fHistogram) {
1371  fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1372  }
1373 
1374  while ((g = (TGraph*) next())) {
1375  Double_t rx1,ry1,rx2,ry2;
1376  g->ComputeRange(rx1, ry1, rx2, ry2);
1377  if (rx1 < rwxmin) rwxmin = rx1;
1378  if (ry1 < rwymin) rwymin = ry1;
1379  if (rx2 > rwxmax) rwxmax = rx2;
1380  if (ry2 > rwymax) rwymax = ry2;
1381  if (g->GetN() > npt) npt = g->GetN();
1382  }
1383 
1384  Int_t ndiv = fGraphs->GetSize();
1385 
1386  TH2F* frame = new TH2F("frame","", ndiv, 0., (Double_t)(ndiv), npt, rwxmin, rwxmax);
1387  if (fHistogram) {
1388  frame->SetTitle(fHistogram->GetTitle());
1389  frame->GetYaxis()->SetTitle(fHistogram->GetXaxis()->GetTitle());
1391  frame->GetZaxis()->SetTitle(fHistogram->GetYaxis()->GetTitle());
1392  }
1393 
1394  TAxis *Xaxis = frame->GetXaxis();
1395  Xaxis->SetNdivisions(-ndiv);
1396  next.Reset();
1397  for (i=ndiv; i>=1; i--) {
1398  g = (TGraph*) next();
1399  Xaxis->SetBinLabel(i, g->GetTitle());
1400  }
1401 
1402  frame->SetStats(kFALSE);
1403  if (fMinimum != -1111) frame->SetMinimum(fMinimum);
1404  else frame->SetMinimum(rwymin);
1405  if (fMaximum != -1111) frame->SetMaximum(fMaximum);
1406  else frame->SetMaximum(rwymax);
1407 
1408  l = (char*)strstr(option,"A");
1409  if (l) frame->Paint("lego9,fb,bb");
1410  l = (char*)strstr(option,"BB");
1411  if (!l) frame->Paint("lego9,fb,a,same");
1412 
1413  Double_t *x, *y;
1414  Double_t xyz1[3], xyz2[3];
1415 
1416  Double_t xl = frame->GetYaxis()->GetBinLowEdge(frame->GetYaxis()->GetFirst());
1417  Double_t xu = frame->GetYaxis()->GetBinUpEdge(frame->GetYaxis()->GetLast());
1418  Double_t yl = frame->GetMinimum();
1419  Double_t yu = frame->GetMaximum();
1420  Double_t xc[2],yc[2];
1421  next.Reset();
1422  Int_t j = ndiv;
1423 
1424  while ((g = (TGraph*) next())) {
1425  npt = g->GetN();
1426  x = g->GetX();
1427  y = g->GetY();
1428  gPad->SetLineColor(g->GetLineColor());
1429  gPad->SetLineWidth(g->GetLineWidth());
1430  gPad->SetLineStyle(g->GetLineStyle());
1431  gPad->TAttLine::Modify();
1432  for (i=0; i<npt-1; i++) {
1433  xc[0] = x[i];
1434  xc[1] = x[i+1];
1435  yc[0] = y[i];
1436  yc[1] = y[i+1];
1437  if (gPad->Clip(&xc[0], &yc[0], xl, yl, xu, yu)<2) {
1438  xyz1[0] = j-0.5;
1439  xyz1[1] = xc[0];
1440  xyz1[2] = yc[0];
1441  xyz2[0] = j-0.5;
1442  xyz2[1] = xc[1];
1443  xyz2[2] = yc[1];
1444  gPad->PaintLine3D(xyz1, xyz2);
1445  }
1446  }
1447  j--;
1448  }
1449 
1450  l = (char*)strstr(option,"FB");
1451  if (!l) frame->Paint("lego9,bb,a,same");
1452  delete frame;
1453 }
1454 
1455 
1456 ////////////////////////////////////////////////////////////////////////////////
1457 /// Print the list of graphs.
1458 
1459 void TMultiGraph::Print(Option_t *option) const
1460 {
1462  if (fGraphs) {
1463  TIter next(fGraphs);
1464  while ((g = (TGraph*) next())) {
1465  g->Print(option);
1466  }
1467  }
1468 }
1469 
1470 
1471 ////////////////////////////////////////////////////////////////////////////////
1472 /// Recursively remove this object from a list. Typically implemented
1473 /// by classes that can contain multiple references to a same object.
1474 
1476 {
1477  if (!fGraphs) return;
1478  TObject *objr = fGraphs->Remove(obj);
1479  if (!objr) return;
1480  delete fHistogram; fHistogram = 0;
1481  if (gPad) gPad->Modified();
1482 }
1483 
1484 
1485 ////////////////////////////////////////////////////////////////////////////////
1486 /// Save primitive as a C++ statement(s) on output stream out.
1487 
1488 void TMultiGraph::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1489 {
1490  char quote = '"';
1491  out<<" "<<std::endl;
1492  if (gROOT->ClassSaved(TMultiGraph::Class())) {
1493  out<<" ";
1494  } else {
1495  out<<" TMultiGraph *";
1496  }
1497  out<<"multigraph = new TMultiGraph();"<<std::endl;
1498  out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
1499  out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<std::endl;
1500 
1501  if (fGraphs) {
1503  TObject *g;
1504 
1505  while (lnk) {
1506  g = lnk->GetObject();
1507  g->SavePrimitive(out, Form("multigraph%s",lnk->GetOption()));
1508  lnk = (TObjOptLink*)lnk->Next();
1509  }
1510  }
1511  const char *l = strstr(option,"th2poly");
1512  if (l) {
1513  out<<" "<<l+7<<"->AddBin(multigraph);"<<std::endl;
1514  } else {
1515  out<<" multigraph->Draw(" <<quote<<option<<quote<<");"<<std::endl;
1516  }
1517  TAxis *xaxis = GetXaxis();
1518  TAxis *yaxis = GetYaxis();
1519 
1520  if (xaxis) xaxis->SaveAttributes(out, "multigraph","->GetXaxis()");
1521  if (yaxis) yaxis->SaveAttributes(out, "multigraph","->GetYaxis()");
1522 }
1523 
1524 
1525 ////////////////////////////////////////////////////////////////////////////////
1526 /// Set multigraph maximum.
1527 
1528 void TMultiGraph::SetMaximum(Double_t maximum)
1529 {
1530  fMaximum = maximum;
1531  if (fHistogram) fHistogram->SetMaximum(maximum);
1532 }
1533 
1534 
1535 ////////////////////////////////////////////////////////////////////////////////
1536 /// Set multigraph minimum.
1537 
1538 void TMultiGraph::SetMinimum(Double_t minimum)
1539 {
1540  fMinimum = minimum;
1541  if (fHistogram) fHistogram->SetMinimum(minimum);
1542 }
1543 
1544 
1545 ////////////////////////////////////////////////////////////////////////////////
1546 /// Get iterator over internal graphs list.
1547 TIter TMultiGraph::begin() const
1548 {
1549  return TIter(fGraphs);
1550 }
void PaintPads(Option_t *chopt="")
Divides the active pad and draws all Graphs in the Multigraph separately.
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
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:4576
virtual void SetMinimum(Double_t minimum=-1111)
Set multigraph minimum.
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.
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
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:7807
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 void Paint(Option_t *option="")
Control routine to paint any kind of histograms.
Definition: TH1.cxx:5663
float xmin
Definition: THbookFile.cxx:93
TF1 * GetFunction(const char *name) const
Return pointer to function with name.
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
virtual void FitPanel()
Display a panel with all histogram fit options.
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:390
auto * m
Definition: textangle.C:8
virtual void SetLimits(Double_t xmin, Double_t xmax)
Definition: TAxis.h:154
Double_t Log(Double_t x)
Definition: TMath.h:648
virtual void SetTimeFormat(const char *format="")
Change the format used for time plotting.
Definition: TAxis.cxx:1002
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:8194
const char Option_t
Definition: RtypesCore.h:62
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:504
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.
double Axis_t
Definition: RtypesCore.h:72
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:638
THist< 1, float, THistStatContent, THistStatUncertainty > TH1F
Definition: THist.hxx:285
TH1 * h
Definition: legend2.C:5
virtual void InitPolynom(Double_t xmin, Double_t xmax)
Compute Initial values of parameters for a polynom.
virtual void SetNdivisions(Int_t n=510, Bool_t optim=kTRUE)
Set the number of divisions for this axis.
Definition: TAttAxis.cxx:229
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:35
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:391
#define gROOT
Definition: TROOT.h:402
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
Int_t LoadPlugin()
Load the plugin library for this handler.
Basic string class.
Definition: TString.h:125
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:567
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual void Paint(Option_t *chopt="")
Paint all the graphs of this multigraph.
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.
virtual void RecursiveRemove(TObject *obj)
Recursively remove this object from a list.
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:745
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
don&#39;t draw stats box
Definition: TH1.h:160
Iterator of linked list.
Definition: TList.h:197
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
static constexpr double mg
virtual TObject * FindObject(const char *name) const
Delete a TObjLink object.
Definition: TList.cxx:574
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:514
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
Double_t GetXmin() const
Definition: TAxis.h:133
Double_t x[n]
Definition: legend1.C:17
TAxis * GetXaxis()
Get x axis of the graph.
void Class()
Definition: Class.C:29
virtual Bool_t GetTimeDisplay() const
Definition: TAxis.h:126
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual void Print(Option_t *chopt="") const
Print graph values.
Definition: TGraph.cxx:1995
Double_t Log10(Double_t x)
Definition: TMath.h:651
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
TString & Append(const char *cs)
Definition: TString.h:495
virtual ~TMultiGraph()
TMultiGraph destructor.
TList * GetListOfFunctions()
Return pointer to list of functions.
TH1F * GetHistogram()
Returns a pointer to the histogram used to draw the axis.
virtual void Draw(Option_t *chopt="")
Draw this multigraph with its current attributes.
TH1F * h1
Definition: legend1.C:5
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition: TVirtualPad.h:49
TMultiGraph & operator=(const TMultiGraph &)
Assignment operator.
virtual const char * GetTimeFormat() const
Definition: TAxis.h:127
virtual void LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
Least squares lpolynomial fitting without weights.
virtual void SetTimeDisplay(Int_t value)
Definition: TAxis.h:161
A doubly linked list.
Definition: TList.h:44
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a line.
Definition: TH1.cxx:2707
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
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 void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set limits for parameter ipar.
Definition: TF1.cxx:3385
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition: TNamed.cxx:51
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:655
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
const char * GetTitle() const
Returns title of object.
Definition: TAxis.h:129
Class to manage histogram axis.
Definition: TAxis.h:30
if object ctor succeeded but object should not be used
Definition: TObject.h:68
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
auto * a
Definition: textangle.C:12
virtual void Print(Option_t *chopt="") const
Print the list of graphs.
virtual void Browse(TBrowser *b)
Browse multigraph.
Long_t ExecPlugin(int nargs, const T &... params)
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:249
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:818
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Definition: TFitResultPtr.h:31
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
fHistogram must be reset in GetHistogram
Definition: TGraph.h:71
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
char * Form(const char *fmt,...)
TList * fFunctions
Definition: TMultiGraph.h:39
virtual void Clear(Option_t *option="")=0
static TVirtualFitter * GetFitter()
static: return the current Fitter
TAxis * GetYaxis()
Get y axis of the graph.
Int_t GetN() const
Definition: TGraph.h:122
TAxis * GetYaxis()
Definition: TH1.h:316
virtual TList * GetListOfPrimitives() const =0
float xmax
Definition: THbookFile.cxx:93
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 InitExpo(Double_t xmin, Double_t xmax)
Compute Initial values of parameters for an exponential.
Definition: graph.py:1
if object destructor must call RecursiveRemove()
Definition: TObject.h:60
TGraphErrors * gr
Definition: legend1.C:25
virtual TObjLink * FirstLink() const
Definition: TList.h:108
#define Printf
Definition: TGeoToOCC.h:18
TMultiGraph()
TMultiGraph default constructor.
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:7892
const Bool_t kFALSE
Definition: RtypesCore.h:88
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to each graph.
Double_t * GetX() const
Definition: TGraph.h:129
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual void Paint(Option_t *option="")
This method must be overridden if a class wants to paint itself.
Definition: TObject.cxx:519
#define ClassImp(name)
Definition: Rtypes.h:359
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition: DataRange.h:34
double Double_t
Definition: RtypesCore.h:55
virtual void PaintStats(TF1 *fit)
Draw the stats.
Definition: TGraph.cxx:1986
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a graph.
Definition: TGraph.cxx:789
R__EXTERN TEnv * gEnv
Definition: TEnv.h:171
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Double_t y[n]
Definition: legend1.C:17
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:570
virtual Option_t * GetOption() const
Definition: TObject.h:120
The TH1 histogram class.
Definition: TH1.h:56
Double_t fMaximum
Definition: TMultiGraph.h:41
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition: TAxis.cxx:809
Bool_t IsNull() const
Definition: TString.h:383
void FitOptionsMake(EFitObjectType type, const char *option, Foption_t &fitOption)
Decode list of options into fitOption.
Definition: HFitImpl.cxx:681
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
Abstract Base Class for Fitting.
TAxis * GetZaxis()
Definition: TH1.h:317
void PaintPolyLine3D(Option_t *chopt="")
Paint all the graphs of this multigraph as 3D lines.
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Int_t GetNpar() const
Definition: TF1.h:465
TList * fGraphs
Definition: TMultiGraph.h:38
virtual void SetMaximum(Double_t maximum=-1111)
Set multigraph maximum.
TList * GetListOfGraphs() const
Definition: TMultiGraph.h:69
Double_t * GetY() const
Definition: TGraph.h:130
virtual void Add(TObject *obj)
Definition: TList.h:87
auto * l
Definition: textangle.C:4
TH1F * fHistogram
Definition: TMultiGraph.h:40
1-Dim function class
Definition: TF1.h:211
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:41
TF1 * f1
Definition: legend1.C:11
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
#define gPad
Definition: TVirtualPad.h:285
void ResetBit(UInt_t f)
Definition: TObject.h:171
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:618
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6154
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
Option_t * GetDrawOption() const
Get option used by the graphics system to draw this object.
Definition: TBrowser.h:104
virtual Int_t IsInside(Double_t x, Double_t y) const
Return 1 if the point (x,y) is inside the polygon defined by the graph vertices 0 otherwise...
Definition: TGraph.cxx:1823
virtual void Add(TGraph *graph, Option_t *chopt="")
Add a new graph to the list of graphs.
Double_t Sqrt(Double_t x)
Definition: TMath.h:590
virtual void ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
Compute the x/y range of the points in this graph.
Definition: TGraph.cxx:645
virtual Option_t * GetGraphDrawOption(const TGraph *gr) const
Return the draw option for the TGraph gr in this TMultiGraph.
virtual Int_t GetSize() const
Definition: TCollection.h:180
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual TObject * GetUserFunc() const
THist< 2, float, THistStatContent, THistStatUncertainty > TH2F
Definition: THist.hxx:291
Bool_t IsWhitespace() const
Definition: TString.h:384
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition: TEnv.cxx:491
const Bool_t kTRUE
Definition: RtypesCore.h:87
Double_t GetXmax() const
Definition: TAxis.h:134
virtual void InitGaus(Double_t xmin, Double_t xmax)
Compute Initial values of parameters for a gaussian.
const Int_t n
Definition: legend1.C:16
TIter begin() const
Get iterator over internal graphs list.
Double_t fMinimum
Definition: TMultiGraph.h:42
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8247
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:315
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
const char * Data() const
Definition: TString.h:345
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition: TObject.cxx:664
static constexpr double g