[root] / trunk / hist / hist / src / TMultiGraph.cxx Repository:
ViewVC logotype

Annotation of /trunk/hist/hist/src/TMultiGraph.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 48282 - (view) (download) (as text)

1 : brun 24702 // @(#)root/hist:$Id$
2 : brun 754 // 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 "TMultiGraph.h"
14 :     #include "TGraph.h"
15 :     #include "TH1.h"
16 : couet 43622 #include "TH2.h"
17 :     #include "TPolyLine3D.h"
18 : brun 754 #include "TVirtualPad.h"
19 : rdm 3748 #include "Riostream.h"
20 : brun 11226 #include "TVirtualFitter.h"
21 : moneta 25487 #include "TPluginManager.h"
22 : pcanal 14336 #include "TClass.h"
23 : brun 17338 #include "TMath.h"
24 : moneta 25487 #include "TSystem.h"
25 : rdm 22419 #include <stdlib.h>
26 : brun 754
27 : moneta 31207 #include "HFitInterface.h"
28 :     #include "Fit/DataRange.h"
29 :     #include "Math/MinimizerOptions.h"
30 :    
31 : brun 754 #include <ctype.h>
32 :    
33 : brun 11226 extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
34 : brun 754
35 :     ClassImp(TMultiGraph)
36 :    
37 : couet 13304
38 : brun 754 //______________________________________________________________________________
39 : couet 20297 /* Begin_Html
40 :     <center><h2>TMultiGraph class</h2></center>
41 : couet 37721
42 :     A TMultiGraph is a collection of TGraph (or derived) objects. It allows to
43 :     manipulate a set of graphs as a single entity. In particular, when drawn,
44 :     the X and Y axis ranges are automatically computed such as all the graphs
45 :     will be visible.
46 :     <p>
47 :     <tt>TMultiGraph::Add</tt> should be used to add a new graph to the list.
48 :     <p>
49 : couet 20297 The TMultiGraph owns the objects in the list.
50 :     <p>
51 : couet 37721 The drawing options are the same as for TGraph.
52 :     Like for TGraph, the painting is performed thanks to the
53 :     <a href="http://root.cern.ch/root/html/TGraphPainter.html">TGraphPainter</a>
54 :     class. All details about the various painting options are given in
55 :     <a href="http://root.cern.ch/root/html/TGraphPainter.html">this class</a>.
56 : couet 20297 Example:
57 : couet 20433 <pre>
58 : couet 20297 TGraph *gr1 = new TGraph(...
59 :     TGraphErrors *gr2 = new TGraphErrors(...
60 :     TMultiGraph *mg = new TMultiGraph();
61 :     mg->Add(gr1,"lp");
62 :     mg->Add(gr2,"cp");
63 :     mg->Draw("a");
64 :     </pre>
65 : couet 43622 A special option <tt>3D</tt> allows to draw the graphs in a 3D space. See the
66 :     following example:
67 :     End_Html
68 :     Begin_Macro(source)
69 :     {
70 :     c0 = new TCanvas("c1","multigraph L3",200,10,700,500);
71 :     c0->SetFrameFillColor(30);
72 :    
73 :     TMultiGraph *mg = new TMultiGraph();
74 :    
75 :     TGraph *gr1 = new TGraph(); gr1->SetLineColor(kBlue);
76 :     TGraph *gr2 = new TGraph(); gr2->SetLineColor(kRed);
77 :     TGraph *gr3 = new TGraph(); gr3->SetLineColor(kGreen);
78 :     TGraph *gr4 = new TGraph(); gr4->SetLineColor(kOrange);
79 :    
80 :     Double_t dx = 6.28/100;
81 :     Double_t x = -3.14;
82 :    
83 :     for (int i=0; i<=100; 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->Draw("a fb l3d");
97 :     return c0;
98 :     }
99 :     End_Macro
100 :     Begin_Html
101 : couet 37731 <p>
102 :     The number of graphs in a multigraph can be retrieve with:
103 :     <pre>
104 :     mg->GetListOfGraphs()->GetSize();
105 :     </pre>
106 :     <p>
107 : couet 20297 The drawing option for each TGraph may be specified as an optional
108 : couet 37721 second argument of the <tt>Add</tt> function.
109 :     <p>
110 : couet 20297 If a draw option is specified, it will be used to draw the graph,
111 :     otherwise the graph will be drawn with the option specified in
112 :     <tt>TMultiGraph::Draw</tt>.
113 : couet 37721 <p>
114 :     The following example shows how to fit a TMultiGraph.
115 :     End_Html
116 :     Begin_Macro(source)
117 :     {
118 :     TCanvas *c1 = new TCanvas("c1","c1",600,400);
119 :    
120 :     Double_t x1[2] = {2.,4.};
121 :     Double_t dx1[2] = {0.1,0.1};
122 :     Double_t y1[2] = {2.1,4.0};
123 :     Double_t dy1[2] = {0.3,0.2};
124 :    
125 :     Double_t x2[2] = {3.,5.};
126 :     Double_t dx2[2] = {0.1,0.1};
127 :     Double_t y2[2] = {3.2,4.8};
128 :     Double_t dy2[2] = {0.3,0.2};
129 :    
130 :     gStyle->SetOptFit(0001);
131 :    
132 :     TGraphErrors *g1 = new TGraphErrors(2,x1,y1,dx1,dy1);
133 :     g1->SetMarkerStyle(21);
134 :     g1->SetMarkerColor(2);
135 :    
136 :     TGraphErrors *g2 = new TGraphErrors(2,x2,y2,dx2,dy2);
137 :     g2->SetMarkerStyle(22);
138 :     g2->SetMarkerColor(3);
139 :    
140 :     TMultiGraph *g = new TMultiGraph();
141 :     g->Add(g1);
142 :     g->Add(g2);
143 :    
144 :     g->Draw("AP");
145 :    
146 :     g->Fit("pol1","FQ");
147 :     return c1;
148 :     }
149 :     End_Macro
150 :     Begin_Html
151 :     <p>
152 :     The axis titles can be modified the following way:
153 :     <p>
154 :     <pre>
155 :     [...]
156 :     TMultiGraph *mg = new TMultiGraph;
157 :     mg->SetTitle("title;xaxis title; yaxis title");
158 :     mg->Add(g1);
159 :     mg->Add(g2);
160 :     mg->Draw("apl");
161 :     </pre>
162 : couet 37731 <p>
163 :     When the graphs in a TMultiGraph are fitted, the fit parameters boxes
164 :     overlap. The following example shows how to make them all visible.
165 : couet 37721
166 : couet 37731 End_Html
167 :     Begin_Macro(source)
168 :     ../../../tutorials/graphs/multigraph.C
169 :     End_Macro
170 :     Begin_Html
171 :    
172 :     <p>
173 :     The axis limits can be changed the like for TGraph. The same methods apply on
174 :     the multigraph.
175 :     Note the two differents ways to change limits on X and Y axis.
176 :    
177 :     End_Html
178 :     Begin_Macro(source)
179 :     {
180 :     TCanvas *c2 = new TCanvas("c2","c2",600,400);
181 :    
182 :     TGraph *g[3];
183 :     Double_t x[10] = {0,1,2,3,4,5,6,7,8,9};
184 :     Double_t y[10] = {1,2,3,4,5,5,4,3,2,1};
185 :     TMultiGraph *mg = new TMultiGraph();
186 :     for (int i=0; i<3; i++) {
187 :     g[i] = new TGraph(10, x, y);
188 :     g[i]->SetMarkerStyle(20);
189 :     g[i]->SetMarkerColor(i+2);
190 :     for (int j=0; j<10; j++) y[j] = y[j]-1;
191 :     mg->Add(g[i]);
192 :     }
193 :     mg->Draw("APL");
194 :     mg->GetXaxis()->SetTitle("E_{#gamma} (GeV)");
195 :     mg->GetYaxis()->SetTitle("Coefficients");
196 :    
197 :     // Change the axis limits
198 :     gPad->Modified();
199 :     mg->GetXaxis()->SetLimits(1.5,7.5);
200 :     mg->SetMinimum(0.);
201 :     mg->SetMaximum(10.);
202 :    
203 :     return c2;
204 :     }
205 :     End_Macro
206 :     Begin_Html
207 :     <p>
208 :     The method <a href="http://root.cern.ch/root/html/TPad.html#TPad:BuildLegend">
209 :     <tt>TPad::BuildLegend</tt></a> is able to extract the graphs inside a
210 :     multigraph. The following example demonstrate this.
211 :    
212 :     End_Html
213 :     Begin_Macro(source)
214 :     {
215 :     TCanvas *c3 = new TCanvas("c3","c3",600, 400);
216 :    
217 :     TMultiGraph * mg = new TMultiGraph("mg","mg");
218 :    
219 :     const Int_t size = 10;
220 :    
221 :     double x[size];
222 :     double y1[size];
223 :     double y2[size];
224 :     double y3[size];
225 :    
226 :     for ( int i = 0; i < size ; ++i ) {
227 :     x[i] = i;
228 :     y1[i] = size - i;
229 :     y2[i] = size - 0.5 * i;
230 :     y3[i] = size - 0.6 * i;
231 :     }
232 :    
233 :     TGraph * gr1 = new TGraph( size, x, y1 );
234 :     gr1->SetName("gr1");
235 :     gr1->SetTitle("graph 1");
236 :     gr1->SetMarkerStyle(21);
237 :     gr1->SetDrawOption("AP");
238 :     gr1->SetLineColor(2);
239 :     gr1->SetLineWidth(4);
240 :     gr1->SetFillStyle(0);
241 :    
242 :     TGraph * gr2 = new TGraph( size, x, y2 );
243 :     gr2->SetName("gr2");
244 :     gr2->SetTitle("graph 2");
245 :     gr2->SetMarkerStyle(22);
246 :     gr2->SetMarkerColor(2);
247 :     gr2->SetDrawOption("P");
248 :     gr2->SetLineColor(3);
249 :     gr2->SetLineWidth(4);
250 :     gr2->SetFillStyle(0);
251 :    
252 :     TGraph * gr3 = new TGraph( size, x, y3 );
253 :     gr3->SetName("gr3");
254 :     gr3->SetTitle("graph 3");
255 :     gr3->SetMarkerStyle(23);
256 :     gr3->SetLineColor(4);
257 :     gr3->SetLineWidth(4);
258 :     gr3->SetFillStyle(0);
259 :    
260 :     mg->Add( gr1 );
261 :     mg->Add( gr2 );
262 :    
263 :     gr3->Draw("ALP");
264 :     mg->Draw("LP");
265 :     c3->BuildLegend();
266 :    
267 :     return c3;
268 :     }
269 :     End_Macro
270 :     Begin_Html
271 :    
272 :    
273 : couet 20297 End_Html */
274 : brun 754
275 : couet 13304
276 : brun 754 //______________________________________________________________________________
277 :     TMultiGraph::TMultiGraph(): TNamed()
278 :     {
279 : couet 13304 // TMultiGraph default constructor
280 : brun 754
281 :     fGraphs = 0;
282 : brun 11226 fFunctions = 0;
283 : brun 754 fHistogram = 0;
284 :     fMaximum = -1111;
285 :     fMinimum = -1111;
286 :     }
287 :    
288 : couet 13304
289 : brun 754 //______________________________________________________________________________
290 :     TMultiGraph::TMultiGraph(const char *name, const char *title)
291 :     : TNamed(name,title)
292 :     {
293 : couet 43622 // Constructor with name and title
294 : couet 13304
295 : brun 754 fGraphs = 0;
296 : brun 11226 fFunctions = 0;
297 : brun 754 fHistogram = 0;
298 :     fMaximum = -1111;
299 :     fMinimum = -1111;
300 :     }
301 :    
302 : couet 20433
303 : brun 15134 //______________________________________________________________________________
304 :     TMultiGraph::TMultiGraph(const TMultiGraph& mg) :
305 :     TNamed (mg),
306 :     fGraphs(mg.fGraphs),
307 :     fFunctions(mg.fFunctions),
308 :     fHistogram(mg.fHistogram),
309 :     fMaximum(mg.fMaximum),
310 :     fMinimum(mg.fMinimum)
311 : couet 20433 {
312 : couet 43622 // Copy constructor
313 : brun 15171 }
314 : couet 13304
315 : couet 20433
316 : brun 754 //______________________________________________________________________________
317 : brun 15134 TMultiGraph& TMultiGraph::operator=(const TMultiGraph& mg)
318 :     {
319 : couet 43622 // Assignement operator
320 :    
321 : brun 15171 if(this!=&mg) {
322 :     TNamed::operator=(mg);
323 :     fGraphs=mg.fGraphs;
324 :     fFunctions=mg.fFunctions;
325 :     fHistogram=mg.fHistogram;
326 :     fMaximum=mg.fMaximum;
327 :     fMinimum=mg.fMinimum;
328 : couet 20433 }
329 : brun 15171 return *this;
330 : brun 15134 }
331 :    
332 : couet 20433
333 : brun 15134 //______________________________________________________________________________
334 : brun 754 TMultiGraph::~TMultiGraph()
335 :     {
336 : couet 13304 // TMultiGraph destructor
337 : brun 754
338 :     if (!fGraphs) return;
339 : brun 5552 TGraph *g;
340 :     TIter next(fGraphs);
341 :     while ((g = (TGraph*) next())) {
342 : couet 13304 g->ResetBit(kMustCleanup);
343 : brun 5552 }
344 : brun 754 fGraphs->Delete();
345 :     delete fGraphs;
346 :     fGraphs = 0;
347 :     delete fHistogram;
348 :     fHistogram = 0;
349 : brun 11226 if (fFunctions) {
350 :     fFunctions->SetBit(kInvalidObject);
351 :     //special logic to support the case where the same object is
352 :     //added multiple times in fFunctions.
353 :     //This case happens when the same object is added with different
354 :     //drawing modes
355 :     TObject *obj;
356 :     while ((obj = fFunctions->First())) {
357 : rdm 22546 while(fFunctions->Remove(obj)) { }
358 : brun 11226 delete obj;
359 :     }
360 :     delete fFunctions;
361 :     }
362 : brun 754 }
363 :    
364 : couet 13304
365 : brun 754 //______________________________________________________________________________
366 : brun 4037 void TMultiGraph::Add(TGraph *graph, Option_t *chopt)
367 : brun 754 {
368 : couet 43622 // Add a new graph to the list of graphs.
369 :     // Note that the graph is now owned by the TMultigraph.
370 : brun 6639 // Deleting the TMultiGraph object will automatically delete the graphs.
371 :     // You should not delete the graphs when the TMultigraph is still active.
372 : rdm 3742
373 : brun 754 if (!fGraphs) fGraphs = new TList();
374 : brun 5552 graph->SetBit(kMustCleanup);
375 : brun 4037 fGraphs->Add(graph,chopt);
376 : brun 754 }
377 :    
378 : couet 13304
379 : brun 754 //______________________________________________________________________________
380 : couet 16436 void TMultiGraph::Add(TMultiGraph *multigraph, Option_t *chopt)
381 :     {
382 : couet 43622 // Add all the graphs in "multigraph" to the list of graphs.
383 : couet 16436
384 :     TList *graphlist = multigraph->GetListOfGraphs();
385 :     if (!graphlist) return;
386 :    
387 :     if (!fGraphs) fGraphs = new TList();
388 : couet 20433
389 : couet 16436 TGraph *gr;
390 :     gr = (TGraph*)graphlist->First();
391 :     fGraphs->Add(gr,chopt);
392 : couet 22379 for(Int_t i = 1; i < graphlist->GetSize(); i++){
393 : couet 16436 gr = (TGraph*)graphlist->After(gr);
394 :     fGraphs->Add(gr,chopt);
395 :     }
396 :     }
397 :    
398 :    
399 :     //______________________________________________________________________________
400 : brun 754 void TMultiGraph::Browse(TBrowser *)
401 :     {
402 : couet 13304 // Browse multigraph.
403 :    
404 :     Draw("alp");
405 :     gPad->Update();
406 : brun 754 }
407 :    
408 : couet 13304
409 : brun 754 //______________________________________________________________________________
410 :     Int_t TMultiGraph::DistancetoPrimitive(Int_t px, Int_t py)
411 :     {
412 : couet 13304 // Compute distance from point px,py to each graph
413 : brun 754
414 : couet 13304 // Are we on the axis?
415 : brun 754 const Int_t kMaxDiff = 10;
416 :     Int_t distance = 9999;
417 :     if (fHistogram) {
418 :     distance = fHistogram->DistancetoPrimitive(px,py);
419 :     if (distance <= 0) return distance;
420 :     }
421 :    
422 : couet 13304 // Loop on the list of graphs
423 : brun 754 if (!fGraphs) return distance;
424 :     TGraph *g;
425 :     TIter next(fGraphs);
426 :     while ((g = (TGraph*) next())) {
427 :     Int_t dist = g->DistancetoPrimitive(px,py);
428 : brun 4915 if (dist <= 0) return 0;
429 : brun 754 if (dist < kMaxDiff) {gPad->SetSelected(g); return dist;}
430 :     }
431 :     return distance;
432 :     }
433 :    
434 : couet 13304
435 : brun 754 //______________________________________________________________________________
436 :     void TMultiGraph::Draw(Option_t *option)
437 :     {
438 : couet 13304 // Draw this multigraph with its current attributes.
439 :     //
440 : couet 37735 // Options to draw a graph are described in TGraph::PaintGraph
441 : couet 13304 //
442 :     // The drawing option for each TGraph may be specified as an optional
443 :     // second argument of the Add function. You can use GetGraphDrawOption
444 :     // to return this option.
445 :     // If a draw option is specified, it will be used to draw the graph,
446 :     // otherwise the graph will be drawn with the option specified in
447 :     // TMultiGraph::Draw. Use GetDrawOption to return the option specified
448 : couet 37735 // when drawing the TMultiGraph.
449 : brun 754
450 : couet 13304 AppendPad(option);
451 : brun 754 }
452 :    
453 : couet 13304
454 : brun 754 //______________________________________________________________________________
455 : moneta 31491 TFitResultPtr TMultiGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
456 : brun 11226 {
457 : couet 13304 // Fit this graph with function with name fname.
458 :     //
459 :     // interface to TF1::Fit(TF1 *f1...
460 : brun 11226
461 :     char *linear;
462 : brun 11231 linear= (char*)strstr(fname, "++");
463 : brun 11226 TF1 *f1=0;
464 :     if (linear)
465 : brun 12643 f1=new TF1(fname, fname, xmin, xmax);
466 : brun 11226 else {
467 :     f1 = (TF1*)gROOT->GetFunction(fname);
468 :     if (!f1) { Printf("Unknown function: %s",fname); return -1; }
469 :     }
470 :    
471 :     return Fit(f1,option,"",xmin,xmax);
472 :     }
473 :    
474 : couet 13304
475 : brun 11226 //______________________________________________________________________________
476 : moneta 31491 TFitResultPtr TMultiGraph::Fit(TF1 *f1, Option_t *option, Option_t *goption, Axis_t rxmin, Axis_t rxmax)
477 : brun 11226 {
478 : couet 13304 // Fit this multigraph with function f1.
479 :     //
480 :     // In this function all graphs of the multigraph are fitted simultaneously
481 :     //
482 :     // f1 is an already predefined function created by TF1.
483 :     // Predefined functions such as gaus, expo and poln are automatically
484 :     // created by ROOT.
485 :     //
486 :     // The list of fit options is given in parameter option.
487 :     // option = "W" Set all errors to 1
488 :     // = "U" Use a User specified fitting algorithm (via SetFCN)
489 :     // = "Q" Quiet mode (minimum printing)
490 :     // = "V" Verbose mode (default is between Q and V)
491 :     // = "B" Use this option when you want to fix one or more parameters
492 :     // and the fitting function is like "gaus","expo","poln","landau".
493 :     // = "R" Use the Range specified in the function range
494 :     // = "N" Do not store the graphics function, do not draw
495 :     // = "0" Do not plot the result of the fit. By default the fitted function
496 :     // is drawn unless the option"N" above is specified.
497 :     // = "+" Add this new fitted function to the list of fitted functions
498 :     // (by default, any previous function is deleted)
499 :     // = "C" In case of linear fitting, not calculate the chisquare
500 :     // (saves time)
501 :     // = "F" If fitting a polN, switch to minuit fitter
502 :     // = "ROB" In case of linear fitting, compute the LTS regression
503 : couet 20433 // coefficients (robust(resistant) regression), using
504 : couet 13304 // the default fraction of good points
505 :     // "ROB=0.x" - compute the LTS regression coefficients, using
506 :     // 0.x as a fraction of good points
507 :     //
508 :     // When the fit is drawn (by default), the parameter goption may be used
509 :     // to specify a list of graphics options. See TGraph::Paint for a complete
510 :     // list of these options.
511 :     //
512 :     // In order to use the Range option, one must first create a function
513 :     // with the expression to be fitted. For example, if your graph
514 :     // has a defined range between -4 and 4 and you want to fit a gaussian
515 :     // only in the interval 1 to 3, you can do:
516 :     // TF1 *f1 = new TF1("f1","gaus",1,3);
517 :     // graph->Fit("f1","R");
518 :     //
519 :     //
520 :     // who is calling this function
521 :     // ============================
522 :     // Note that this function is called when calling TGraphErrors::Fit
523 :     // or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
524 :     // see the discussion below on the errors calulation.
525 :     //
526 :     // Setting initial conditions
527 :     // ==========================
528 :     // Parameters must be initialized before invoking the Fit function.
529 :     // The setting of the parameter initial values is automatic for the
530 :     // predefined functions : poln, expo, gaus, landau. One can however disable
531 :     // this automatic computation by specifying the option "B".
532 :     // You can specify boundary limits for some or all parameters via
533 :     // f1->SetParLimits(p_number, parmin, parmax);
534 :     // if parmin>=parmax, the parameter is fixed
535 :     // Note that you are not forced to fix the limits for all parameters.
536 :     // For example, if you fit a function with 6 parameters, you can do:
537 :     // func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
538 :     // func->SetParLimits(4,-10,-4);
539 :     // func->SetParLimits(5, 1,1);
540 :     // With this setup, parameters 0->3 can vary freely
541 :     // Parameter 4 has boundaries [-10,-4] with initial value -8
542 :     // Parameter 5 is fixed to 100.
543 :     //
544 :     // Fit range
545 :     // =========
546 :     // The fit range can be specified in two ways:
547 :     // - specify rxmax > rxmin (default is rxmin=rxmax=0)
548 :     // - specify the option "R". In this case, the function will be taken
549 :     // instead of the full graph range.
550 :     //
551 : moneta 34699 // Changing the fitting function
552 :     // =============================
553 :     // By default a chi2 fitting function is used for fitting the TGraphs's.
554 :     // The function is implemented in FitUtil::EvaluateChi2.
555 : couet 37721 // In case of TGraphErrors an effective chi2 is used
556 :     // (see TGraphErrors fit in TGraph::Fit) and is implemented in
557 : moneta 34699 // FitUtil::EvaluateChi2Effective
558 :     // To specify a User defined fitting function, specify option "U" and
559 :     // call the following functions:
560 :     // TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
561 :     // where MyFittingFunction is of type:
562 :     // extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
563 : couet 13304 //
564 : moneta 34699 // Access to the fit result
565 :     // ========================
566 :     // The function returns a TFitResultPtr which can hold a pointer to a TFitResult object.
567 :     // By default the TFitResultPtr contains only the status of the fit and it converts
568 :     // automatically to an integer. If the option "S" is instead used, TFitResultPtr contains
569 :     // the TFitResult and behaves as a smart pointer to it. For example one can do:
570 :     // TFitResultPtr r = graph->Fit("myFunc","S");
571 :     // TMatrixDSym cov = r->GetCovarianceMatrix(); // to access the covariance matrix
572 :     // Double_t par0 = r->Parameter(0); // retrieve the value for the parameter 0
573 :     // Double_t err0 = r->ParError(0); // retrieve the error for the parameter 0
574 :     // r->Print("V"); // print full information of fit including covariance matrix
575 :     // r->Write(); // store the result in a file
576 : couet 13304 //
577 : moneta 34699 // The fit parameters, error and chi2 (but not covariance matrix) can be retrieved also
578 :     // from the fitted function.
579 : couet 13304 //
580 :     //
581 :     // Associated functions
582 :     // ====================
583 :     // One or more object (typically a TF1*) can be added to the list
584 :     // of functions (fFunctions) associated to each graph.
585 :     // When TGraph::Fit is invoked, the fitted function is added to this list.
586 :     // Given a graph gr, one can retrieve an associated function
587 :     // with: TF1 *myfunc = gr->GetFunction("myfunc");
588 :     //
589 :     // If the graph is made persistent, the list of
590 :     // associated functions is also persistent. Given a pointer (see above)
591 :     // to an associated function myfunc, one can retrieve the function/fit
592 :     // parameters with calls such as:
593 :     // Double_t chi2 = myfunc->GetChisquare();
594 :     // Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
595 :     // Double_t err0 = myfunc->GetParError(0); //error on first parameter
596 :     //
597 :     // Fit Statistics
598 :     // ==============
599 :     // You can change the statistics box to display the fit parameters with
600 :     // the TStyle::SetOptFit(mode) method. This mode has four digits.
601 :     // mode = pcev (default = 0111)
602 :     // v = 1; print name/values of parameters
603 :     // e = 1; print errors (if e=1, v must be 1)
604 :     // c = 1; print Chisquare/Number of degress of freedom
605 :     // p = 1; print Probability
606 :     //
607 :     // For example: gStyle->SetOptFit(1011);
608 :     // prints the fit probability, parameter names/values, and errors.
609 :     // You can change the position of the statistics box with these lines
610 :     // (where g is a pointer to the TGraph):
611 :     //
612 :     // Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
613 :     // Root > st->SetX1NDC(newx1); //new x start position
614 :     // Root > st->SetX2NDC(newx2); //new x end position
615 : brun 11226
616 : moneta 31207 // internal multigraph fitting methods
617 :     Foption_t fitOption;
618 : moneta 48282 ROOT::Fit::FitOptionsMake(ROOT::Fit::kGraph,option,fitOption);
619 : rdm 11750
620 : couet 37721 // create range and minimizer options with default values
621 :     ROOT::Fit::DataRange range(rxmin,rxmax);
622 :     ROOT::Math::MinimizerOptions minOption;
623 :     return ROOT::Fit::FitObject(this, f1 , fitOption , minOption, goption, range);
624 : moneta 31207
625 : moneta 25487 }
626 : rdm 11750
627 : moneta 25487 //______________________________________________________________________________
628 :     void TMultiGraph::FitPanel()
629 :     {
630 : couet 43622 // Display a panel with all histogram fit options
631 :     // See class TFitPanel for example
632 : rdm 11750
633 : moneta 25487 if (!gPad)
634 :     gROOT->MakeDefCanvas();
635 : brun 11226
636 : moneta 25487 if (!gPad) {
637 :     Error("FitPanel", "Unable to create a default canvas");
638 :     return;
639 : brun 12136 }
640 : brun 11226
641 : moneta 25487 // use plugin manager to create instance of TFitEditor
642 :     TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
643 :     if (handler && handler->LoadPlugin() != -1) {
644 :     if (handler->ExecPlugin(2, gPad, this) == 0)
645 :     Error("FitPanel", "Unable to crate the FitPanel");
646 : brun 11226 }
647 : couet 26029 else
648 : moneta 25487 Error("FitPanel", "Unable to find the FitPanel plug-in");
649 : brun 11226 }
650 :    
651 :     //______________________________________________________________________________
652 : brun 11584 Option_t *TMultiGraph::GetGraphDrawOption(const TGraph *gr) const
653 :     {
654 : couet 13304 // Return the draw option for the TGraph gr in this TMultiGraph
655 :     // The return option is the one specified when calling TMultiGraph::Add(gr,option).
656 : rdm 11750
657 : brun 11584 if (!fGraphs || !gr) return "";
658 :     TListIter next(fGraphs);
659 :     TObject *obj;
660 :     while ((obj = next())) {
661 : brun 11600 if (obj == (TObject*)gr) return next.GetOption();
662 : brun 11584 }
663 :     return "";
664 :     }
665 :    
666 : couet 13304
667 : brun 11584 //______________________________________________________________________________
668 : brun 11226 void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
669 :     {
670 : couet 13304 // Compute Initial values of parameters for a gaussian.
671 : brun 11226
672 :     Double_t allcha, sumx, sumx2, x, val, rms, mean;
673 :     Int_t bin;
674 :     const Double_t sqrtpi = 2.506628;
675 :    
676 : couet 13304 // Compute mean value and RMS of the graph in the given range
677 : brun 11226 Int_t np = 0;
678 :     allcha = sumx = sumx2 = 0;
679 :     TGraph *g;
680 :     TIter next(fGraphs);
681 :     Double_t *px, *py;
682 :     Int_t npp; //number of points in each graph
683 :     while ((g = (TGraph*) next())) {
684 :     px=g->GetX();
685 :     py=g->GetY();
686 :     npp=g->GetN();
687 :     for (bin=0; bin<npp; bin++){
688 : brun 12643 x=px[bin];
689 :     if (x<xmin || x>xmax) continue;
690 :     np++;
691 :     val=py[bin];
692 :     sumx+=val*x;
693 :     sumx2+=val*x*x;
694 :     allcha+=val;
695 : brun 11226 }
696 :     }
697 :     if (np == 0 || allcha == 0) return;
698 :     mean = sumx/allcha;
699 :     rms = TMath::Sqrt(sumx2/allcha - mean*mean);
700 :    
701 :     Double_t binwidx = TMath::Abs((xmax-xmin)/np);
702 :     if (rms == 0) rms = 1;
703 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
704 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
705 :     f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
706 :     f1->SetParameter(1,mean);
707 :     f1->SetParameter(2,rms);
708 :     f1->SetParLimits(2,0,10*rms);
709 :     }
710 :    
711 : couet 13304
712 : brun 11226 //______________________________________________________________________________
713 :     void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
714 :     {
715 : couet 13304 // Compute Initial values of parameters for an exponential.
716 : brun 11226
717 :     Double_t constant, slope;
718 :     Int_t ifail;
719 :    
720 :     LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
721 :    
722 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
723 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
724 :     f1->SetParameter(0,constant);
725 :     f1->SetParameter(1,slope);
726 :     }
727 :    
728 : couet 13304
729 : brun 11226 //______________________________________________________________________________
730 :     void TMultiGraph::InitPolynom(Double_t xmin, Double_t xmax)
731 :     {
732 : couet 13304 // Compute Initial values of parameters for a polynom.
733 : brun 11226
734 :     Double_t fitpar[25];
735 :    
736 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
737 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
738 :     Int_t npar = f1->GetNpar();
739 :    
740 :     LeastSquareFit(npar, fitpar, xmin, xmax);
741 :    
742 :     for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
743 :     }
744 :    
745 : couet 13304
746 : brun 11226 //______________________________________________________________________________
747 :     void TMultiGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
748 :     {
749 : couet 13304 // Least squares lpolynomial fitting without weights.
750 :     //
751 :     // m number of parameters
752 :     // a array of parameters
753 :     // first 1st point number to fit (default =0)
754 :     // last last point number to fit (default=fNpoints-1)
755 :     //
756 :     // based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
757 : brun 11226
758 : couet 13304 const Double_t zero = 0.;
759 :     const Double_t one = 1.;
760 :     const Int_t idim = 20;
761 : brun 11226
762 : couet 13304 Double_t b[400] /* was [20][20] */;
763 :     Int_t i, k, l, ifail, bin;
764 :     Double_t power;
765 :     Double_t da[20], xk, yk;
766 : brun 11226
767 :    
768 : couet 13304 //count the total number of points to fit
769 :     TGraph *g;
770 :     TIter next(fGraphs);
771 :     Double_t *px, *py;
772 :     Int_t n=0;
773 :     Int_t npp;
774 :     while ((g = (TGraph*) next())) {
775 :     px=g->GetX();
776 :     npp=g->GetN();
777 :     for (bin=0; bin<npp; bin++){
778 :     xk=px[bin];
779 :     if (xk < xmin || xk > xmax) continue;
780 :     n++;
781 :     }
782 :     }
783 :     if (m <= 2) {
784 :     LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
785 :     return;
786 :     }
787 :     if (m > idim || m > n) return;
788 :     da[0] = zero;
789 :     for (l = 2; l <= m; ++l) {
790 :     b[l-1] = zero;
791 :     b[m + l*20 - 21] = zero;
792 :     da[l-1] = zero;
793 :     }
794 :     Int_t np = 0;
795 : brun 11226
796 : couet 13304 next.Reset();
797 :     while ((g = (TGraph*) next())) {
798 :     px=g->GetX();
799 :     py=g->GetY();
800 :     npp=g->GetN();
801 :    
802 :     for (k = 0; k <= npp; ++k) {
803 :     xk = px[k];
804 :     if (xk < xmin || xk > xmax) continue;
805 :     np++;
806 :     yk = py[k];
807 :     power = one;
808 :     da[0] += yk;
809 :     for (l = 2; l <= m; ++l) {
810 :     power *= xk;
811 : brun 12643 b[l-1] += power;
812 :     da[l-1] += power*yk;
813 : couet 13304 }
814 :     for (l = 2; l <= m; ++l) {
815 :     power *= xk;
816 :     b[m + l*20 - 21] += power;
817 :     }
818 :     }
819 :     }
820 :     b[0] = Double_t(np);
821 :     for (i = 3; i <= m; ++i) {
822 :     for (k = i; k <= m; ++k) {
823 :     b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
824 :     }
825 :     }
826 :     H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
827 : brun 11226
828 : couet 13304 if (ifail < 0) {
829 :     //a[0] = fY[0];
830 :     py=((TGraph *)fGraphs->First())->GetY();
831 :     a[0]=py[0];
832 :     for (i=1; i<m; ++i) a[i] = 0;
833 :     return;
834 :     }
835 :     for (i=0; i<m; ++i) a[i] = da[i];
836 : brun 11226 }
837 :    
838 : couet 13304
839 : brun 11226 //______________________________________________________________________________
840 :     void TMultiGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
841 :     {
842 : couet 13304 // Least square linear fit without weights.
843 :     //
844 :     // Fit a straight line (a0 + a1*x) to the data in this graph.
845 :     // ndata: number of points to fit
846 :     // first: first point number to fit
847 :     // last: last point to fit O(ndata should be last-first
848 :     // ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
849 :     //
850 :     // extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
851 : brun 11226
852 : couet 13304 Double_t xbar, ybar, x2bar;
853 :     Int_t i;
854 :     Double_t xybar;
855 :     Double_t fn, xk, yk;
856 :     Double_t det;
857 : brun 11226
858 : couet 13304 ifail = -2;
859 :     xbar = ybar = x2bar = xybar = 0;
860 :     Int_t np = 0;
861 :     TGraph *g;
862 :     TIter next(fGraphs);
863 :     Double_t *px, *py;
864 :     Int_t npp;
865 :     while ((g = (TGraph*) next())) {
866 :     px=g->GetX();
867 :     py=g->GetY();
868 :     npp=g->GetN();
869 :     for (i = 0; i < npp; ++i) {
870 :     xk = px[i];
871 :     if (xk < xmin || xk > xmax) continue;
872 :     np++;
873 :     yk = py[i];
874 :     if (ndata < 0) {
875 :     if (yk <= 0) yk = 1e-9;
876 :     yk = TMath::Log(yk);
877 :     }
878 :     xbar += xk;
879 :     ybar += yk;
880 :     x2bar += xk*xk;
881 :     xybar += xk*yk;
882 :     }
883 :     }
884 :     fn = Double_t(np);
885 :     det = fn*x2bar - xbar*xbar;
886 :     ifail = -1;
887 :     if (det <= 0) {
888 :     if (fn > 0) a0 = ybar/fn;
889 :     else a0 = 0;
890 :     a1 = 0;
891 :     return;
892 :     }
893 :     ifail = 0;
894 :     a0 = (x2bar*ybar - xbar*xybar) / det;
895 :     a1 = (fn*xybar - xbar*ybar) / det;
896 :     }
897 : brun 11226
898 :    
899 :     //______________________________________________________________________________
900 : couet 36493 Int_t TMultiGraph::IsInside(Double_t x, Double_t y) const
901 :     {
902 :     // Return 1 if the point (x,y) is inside one of the graphs 0 otherwise.
903 :    
904 :     Int_t in = 0;
905 :     if (!fGraphs) return in;
906 :     TGraph *g;
907 :     TIter next(fGraphs);
908 :     while ((g = (TGraph*) next())) {
909 :     in = g->IsInside(x, y);
910 :     if (in) return in;
911 :     }
912 :     return in;
913 :     }
914 :    
915 :    
916 :     //______________________________________________________________________________
917 : brun 1205 TH1F *TMultiGraph::GetHistogram() const
918 : brun 754 {
919 : couet 13304 // Returns a pointer to the histogram used to draw the axis
920 :     // Takes into account the two following cases.
921 :     // 1- option 'A' was specified in TMultiGraph::Draw. Return fHistogram
922 :     // 2- user had called TPad::DrawFrame. return pointer to hframe histogram
923 : brun 754
924 :     if (fHistogram) return fHistogram;
925 :     if (!gPad) return 0;
926 :     gPad->Modified();
927 :     gPad->Update();
928 :     if (fHistogram) return fHistogram;
929 :     TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
930 :     return h1;
931 :     }
932 :    
933 : couet 13304
934 : brun 754 //______________________________________________________________________________
935 : brun 11226 TF1 *TMultiGraph::GetFunction(const char *name) const
936 :     {
937 : couet 13304 // Return pointer to function with name.
938 :     //
939 :     // Functions such as TGraph::Fit store the fitted function in the list of
940 :     // functions of this graph.
941 : brun 11226
942 :     if (!fFunctions) return 0;
943 :     return (TF1*)fFunctions->FindObject(name);
944 :     }
945 :    
946 : moneta 25487 //______________________________________________________________________________
947 : couet 26029 TList *TMultiGraph::GetListOfFunctions()
948 : moneta 25487 {
949 :     // Return pointer to list of functions
950 :     // if pointer is null create the list
951 : couet 13304
952 : moneta 25487 if (!fFunctions) fFunctions = new TList();
953 : couet 26029 return fFunctions;
954 : moneta 25487 }
955 :    
956 :    
957 : brun 11226 //______________________________________________________________________________
958 : brun 1205 TAxis *TMultiGraph::GetXaxis() const
959 : brun 754 {
960 :     // Get x axis of the graph.
961 :    
962 :     if (!gPad) return 0;
963 : brun 4083 TH1 *h = GetHistogram();
964 :     if (!h) return 0;
965 :     return h->GetXaxis();
966 : brun 754 }
967 :    
968 : couet 13304
969 : brun 754 //______________________________________________________________________________
970 : brun 1205 TAxis *TMultiGraph::GetYaxis() const
971 : brun 754 {
972 :     // Get y axis of the graph.
973 :    
974 :     if (!gPad) return 0;
975 : brun 4083 TH1 *h = GetHistogram();
976 :     if (!h) return 0;
977 :     return h->GetYaxis();
978 : brun 754 }
979 :    
980 : couet 13304
981 : brun 754 //______________________________________________________________________________
982 :     void TMultiGraph::Paint(Option_t *option)
983 :     {
984 : couet 43622 // Paint all the graphs of this multigraph
985 :    
986 : rdm 41541 const TPickerStackGuard pushGuard(this);
987 : brun 754
988 : couet 33513 if (!fGraphs) return;
989 : couet 13304 if (fGraphs->GetSize() == 0) return;
990 : rdm 11750
991 : couet 13304 char *l;
992 :     static char chopt[33];
993 :     Int_t nch = strlen(option);
994 :     Int_t i;
995 :     for (i=0;i<nch;i++) chopt[i] = toupper(option[i]);
996 :     chopt[nch] = 0;
997 : couet 43622
998 : couet 43627 l = (char*)strstr(chopt,"3D");
999 : couet 43622 if (l) {
1000 : couet 43627 l = (char*)strstr(chopt,"L");
1001 : couet 43622 if (l) PaintPolyLine3D(chopt);
1002 :     return;
1003 :     }
1004 :    
1005 : couet 13304 TGraph *g;
1006 : rdm 3742
1007 : couet 43627 l = (char*)strstr(chopt,"A");
1008 : couet 13304 if (l) {
1009 :     *l = ' ';
1010 :     TIter next(fGraphs);
1011 :     Int_t npt = 100;
1012 :     Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
1013 :     rwxmin = gPad->GetUxmin();
1014 :     rwxmax = gPad->GetUxmax();
1015 :     rwymin = gPad->GetUymin();
1016 :     rwymax = gPad->GetUymax();
1017 :     char *xtitle = 0;
1018 :     char *ytitle = 0;
1019 :     Int_t firstx = 0;
1020 :     Int_t lastx = 0;
1021 : couet 40956 Bool_t timedisplay = kFALSE;
1022 :     char *timeformat = 0;
1023 : rdm 11750
1024 : couet 13304 if (fHistogram) {
1025 :     //cleanup in case of a previous unzoom
1026 :     if (fHistogram->GetMinimum() >= fHistogram->GetMaximum()) {
1027 : brun 23506 nch = strlen(fHistogram->GetXaxis()->GetTitle());
1028 : couet 13304 firstx = fHistogram->GetXaxis()->GetFirst();
1029 :     lastx = fHistogram->GetXaxis()->GetLast();
1030 : couet 40956 timedisplay = fHistogram->GetXaxis()->GetTimeDisplay();
1031 : couet 13304 if (nch) {
1032 :     xtitle = new char[nch+1];
1033 : brun 35498 strlcpy(xtitle,fHistogram->GetXaxis()->GetTitle(),nch+1);
1034 : couet 13304 }
1035 :     nch = strlen(fHistogram->GetYaxis()->GetTitle());
1036 :     if (nch) {
1037 :     ytitle = new char[nch+1];
1038 : brun 35498 strlcpy(ytitle,fHistogram->GetYaxis()->GetTitle(),nch+1);
1039 : couet 13304 }
1040 : couet 40956 nch = strlen(fHistogram->GetXaxis()->GetTimeFormat());
1041 :     if (nch) {
1042 :     timeformat = new char[nch+1];
1043 :     strlcpy(timeformat,fHistogram->GetXaxis()->GetTimeFormat(),nch+1);
1044 :     }
1045 : couet 13304 delete fHistogram;
1046 :     fHistogram = 0;
1047 :     }
1048 :     }
1049 :     if (fHistogram) {
1050 :     minimum = fHistogram->GetYaxis()->GetXmin();
1051 :     maximum = fHistogram->GetYaxis()->GetXmax();
1052 :     uxmin = gPad->PadtoX(rwxmin);
1053 :     uxmax = gPad->PadtoX(rwxmax);
1054 :     } else {
1055 : couet 21573 g = (TGraph*) next();
1056 : brun 32499 if (g) g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1057 : couet 13304 while ((g = (TGraph*) next())) {
1058 : couet 20433 Double_t rx1,ry1,rx2,ry2;
1059 :     g->ComputeRange(rx1, ry1, rx2, ry2);
1060 :     if (rx1 < rwxmin) rwxmin = rx1;
1061 :     if (ry1 < rwymin) rwymin = ry1;
1062 :     if (rx2 > rwxmax) rwxmax = rx2;
1063 :     if (ry2 > rwymax) rwymax = ry2;
1064 : couet 13304 if (g->GetN() > npt) npt = g->GetN();
1065 :     }
1066 :     if (rwxmin == rwxmax) rwxmax += 1.;
1067 :     if (rwymin == rwymax) rwymax += 1.;
1068 :     dx = 0.05*(rwxmax-rwxmin);
1069 :     dy = 0.05*(rwymax-rwymin);
1070 :     uxmin = rwxmin - dx;
1071 :     uxmax = rwxmax + dx;
1072 :     if (gPad->GetLogy()) {
1073 :     if (rwymin <= 0) rwymin = 0.001*rwymax;
1074 :     minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
1075 :     maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
1076 :     } else {
1077 :     minimum = rwymin - dy;
1078 :     maximum = rwymax + dy;
1079 :     }
1080 :     if (minimum < 0 && rwymin >= 0) minimum = 0;
1081 :     if (maximum > 0 && rwymax <= 0) maximum = 0;
1082 :     }
1083 : brun 754
1084 : couet 13304 if (fMinimum != -1111) rwymin = minimum = fMinimum;
1085 :     if (fMaximum != -1111) rwymax = maximum = fMaximum;
1086 :     if (uxmin < 0 && rwxmin >= 0) {
1087 :     if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
1088 :     //else uxmin = 0;
1089 :     }
1090 :     if (uxmax > 0 && rwxmax <= 0) {
1091 :     if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
1092 :     //else uxmax = 0;
1093 :     }
1094 :     if (minimum < 0 && rwymin >= 0) {
1095 :     if(gPad->GetLogy()) minimum = 0.9*rwymin;
1096 :     //else minimum = 0;
1097 :     }
1098 :     if (maximum > 0 && rwymax <= 0) {
1099 :     if(gPad->GetLogy()) maximum = 1.1*rwymax;
1100 :     //else maximum = 0;
1101 :     }
1102 :     if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
1103 :     if (uxmin <= 0 && gPad->GetLogx()) {
1104 :     if (uxmax > 1000) uxmin = 1;
1105 :     else uxmin = 0.001*uxmax;
1106 :     }
1107 :     rwymin = minimum;
1108 :     rwymax = maximum;
1109 :     if (fHistogram) {
1110 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1111 :     }
1112 : brun 3575
1113 : couet 13304 // Create a temporary histogram to draw the axis
1114 :     if (!fHistogram) {
1115 :     // the graph is created with at least as many channels as there are points
1116 :     // to permit zooming on the full range
1117 :     rwxmin = uxmin;
1118 :     rwxmax = uxmax;
1119 :     fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1120 :     if (!fHistogram) return;
1121 :     fHistogram->SetMinimum(rwymin);
1122 :     fHistogram->SetBit(TH1::kNoStats);
1123 :     fHistogram->SetMaximum(rwymax);
1124 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1125 :     fHistogram->SetDirectory(0);
1126 :     if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
1127 :     if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
1128 :     if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
1129 : couet 40956 if (timedisplay) {fHistogram->GetXaxis()->SetTimeDisplay(timedisplay);}
1130 :     if (timeformat) {fHistogram->GetXaxis()->SetTimeFormat(timeformat); delete [] timeformat;}
1131 : couet 13304 }
1132 :     fHistogram->Paint("0");
1133 : brun 3575 }
1134 : rdm 3742
1135 : couet 21962 TGraph *gfit = 0;
1136 : brun 754 if (fGraphs) {
1137 : brun 4037 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
1138 : couet 21962 TObject *obj = 0;
1139 : brun 4037
1140 :     while (lnk) {
1141 : rdm 41541
1142 : brun 4037 obj = lnk->GetObject();
1143 : rdm 41541
1144 :     gPad->PushSelectableObject(obj);
1145 :    
1146 :     if (!gPad->PadInHighlightMode() || (gPad->PadInHighlightMode() && obj == gPad->GetSelected())) {
1147 :     if (strlen(lnk->GetOption()))
1148 :     obj->Paint(lnk->GetOption());
1149 :     else
1150 :     obj->Paint(chopt);
1151 :     }
1152 :    
1153 : brun 4037 lnk = (TObjOptLink*)lnk->Next();
1154 :     }
1155 : rdm 41541
1156 : couet 21962 gfit = (TGraph*)obj; // pick one TGraph in the list to paint the fit parameters.
1157 : brun 754 }
1158 : brun 11226
1159 :     TObject *f;
1160 : couet 21962 TF1 *fit = 0;
1161 : brun 11226 if (fFunctions) {
1162 : couet 13304 TIter next(fFunctions);
1163 :     while ((f = (TObject*) next())) {
1164 :     if (f->InheritsFrom(TF1::Class())) {
1165 :     if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
1166 : couet 21962 fit = (TF1*)f;
1167 : couet 13304 } else {
1168 :     f->Paint();
1169 :     }
1170 : brun 11226 }
1171 :     }
1172 : couet 21962
1173 : brun 24084 if (fit) gfit->PaintStats(fit);
1174 : couet 13304 }
1175 : brun 11226
1176 :    
1177 : brun 754 //______________________________________________________________________________
1178 : couet 43622 void TMultiGraph::PaintPolyLine3D(Option_t *option)
1179 :     {
1180 :     // Paint all the graphs of this multigraph as 3D lines
1181 :    
1182 : couet 43629 Int_t i, npt=0;
1183 : couet 43622 char *l;
1184 : couet 43634 Double_t rwxmin=0., rwxmax=0., rwymin=0., rwymax=0.;
1185 : couet 43622 TIter next(fGraphs);
1186 :     TGraph *g;
1187 :    
1188 :     g = (TGraph*) next();
1189 :     if (g) {
1190 :     g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1191 :     npt = g->GetN();
1192 :     }
1193 :    
1194 :     while ((g = (TGraph*) next())) {
1195 :     Double_t rx1,ry1,rx2,ry2;
1196 :     g->ComputeRange(rx1, ry1, rx2, ry2);
1197 :     if (rx1 < rwxmin) rwxmin = rx1;
1198 :     if (ry1 < rwymin) rwymin = ry1;
1199 :     if (rx2 > rwxmax) rwxmax = rx2;
1200 :     if (ry2 > rwymax) rwymax = ry2;
1201 :     if (g->GetN() > npt) npt = g->GetN();
1202 :     }
1203 :    
1204 :     Int_t ndiv = fGraphs->GetSize();
1205 :     TH2F* frame = new TH2F("frame","", ndiv, 0., (Double_t)(ndiv),
1206 :     10, rwxmin, rwxmax);
1207 :    
1208 :     TAxis *Xaxis = frame->GetXaxis();
1209 :     Xaxis->SetNdivisions(-ndiv);
1210 :     next.Reset();
1211 :     for (i=ndiv; i>=1; i--) {
1212 :     g = (TGraph*) next();
1213 :     Xaxis->SetBinLabel(i, g->GetTitle());
1214 :     }
1215 :    
1216 :     frame->SetStats(kFALSE);
1217 :     frame->SetMinimum(rwymin);
1218 :     frame->SetMaximum(rwymax);
1219 :    
1220 : couet 43627 l = (char*)strstr(option,"A");
1221 : couet 43622 if (l) frame->Paint("lego0,fb,bb");
1222 : couet 43627 l = (char*)strstr(option,"BB");
1223 : couet 43622 if (!l) frame->Paint("lego0,fb,a,same");
1224 :    
1225 :     Double_t *x, *y;
1226 :     Double_t xyz1[3], xyz2[3];
1227 :    
1228 :     next.Reset();
1229 :     Int_t j = ndiv;
1230 :     while ((g = (TGraph*) next())) {
1231 :     npt = g->GetN();
1232 :     x = g->GetX();
1233 :     y = g->GetY();
1234 :     gPad->SetLineColor(g->GetLineColor());
1235 :     gPad->SetLineWidth(g->GetLineWidth());
1236 :     gPad->SetLineStyle(g->GetLineStyle());
1237 :     gPad->TAttLine::Modify();
1238 :     for (i=0; i<npt-1; i++) {
1239 :     xyz1[0] = j-0.5;
1240 :     xyz1[1] = x[i];
1241 :     xyz1[2] = y[i];
1242 :     xyz2[0] = j-0.5;
1243 :     xyz2[1] = x[i+1];
1244 :     xyz2[2] = y[i+1];
1245 :     gPad->PaintLine3D(xyz1, xyz2);
1246 :     }
1247 :     j--;
1248 :     }
1249 :    
1250 : couet 43627 l = (char*)strstr(option,"FB");
1251 : couet 43622 if (!l) frame->Paint("lego0,bb,a,same");
1252 :     delete frame;
1253 :     }
1254 :    
1255 :    
1256 :     //______________________________________________________________________________
1257 : brun 1205 void TMultiGraph::Print(Option_t *option) const
1258 : brun 754 {
1259 : couet 13304 // Print the list of graphs
1260 : brun 754
1261 :     TGraph *g;
1262 :     if (fGraphs) {
1263 : couet 13304 TIter next(fGraphs);
1264 :     while ((g = (TGraph*) next())) {
1265 :     g->Print(option);
1266 :     }
1267 : brun 754 }
1268 :     }
1269 :    
1270 : couet 13304
1271 : brun 754 //______________________________________________________________________________
1272 : brun 5552 void TMultiGraph::RecursiveRemove(TObject *obj)
1273 :     {
1274 :     // Recursively remove this object from a list. Typically implemented
1275 :     // by classes that can contain mulitple references to a same object.
1276 :    
1277 :     if (!fGraphs) return;
1278 :     TObject *objr = fGraphs->Remove(obj);
1279 :     if (!objr) return;
1280 :     delete fHistogram; fHistogram = 0;
1281 :     if (gPad) gPad->Modified();
1282 :     }
1283 :    
1284 : couet 13304
1285 : brun 5552 //______________________________________________________________________________
1286 : axel 44507 void TMultiGraph::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1287 : brun 754 {
1288 : couet 13304 // Save primitive as a C++ statement(s) on output stream out
1289 : brun 754
1290 :     char quote = '"';
1291 : axel 44507 out<<" "<<std::endl;
1292 : brun 754 if (gROOT->ClassSaved(TMultiGraph::Class())) {
1293 : couet 13304 out<<" ";
1294 : brun 754 } else {
1295 : couet 13304 out<<" TMultiGraph *";
1296 : brun 754 }
1297 : axel 44507 out<<"multigraph = new TMultiGraph();"<<std::endl;
1298 :     out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
1299 :     out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<std::endl;
1300 : brun 754
1301 :     if (fGraphs) {
1302 : brun 10847 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
1303 :     TObject *g;
1304 :    
1305 :     while (lnk) {
1306 :     g = lnk->GetObject();
1307 : couet 20297 g->SavePrimitive(out, Form("multigraph%s",lnk->GetOption()));
1308 : brun 10847 lnk = (TObjOptLink*)lnk->Next();
1309 : couet 13304 }
1310 : brun 754 }
1311 : couet 36530 const char *l = strstr(option,"th2poly");
1312 :     if (l) {
1313 : axel 44507 out<<" "<<l+7<<"->AddBin(multigraph);"<<std::endl;
1314 : couet 36530 } else {
1315 : axel 44507 out<<" multigraph->Draw(" <<quote<<option<<quote<<");"<<std::endl;
1316 : couet 36530 }
1317 : couet 14017 TAxis *xaxis = GetXaxis();
1318 :     TAxis *yaxis = GetYaxis();
1319 : couet 20433
1320 :     if (xaxis) xaxis->SaveAttributes(out, "multigraph","->GetXaxis()");
1321 :     if (yaxis) yaxis->SaveAttributes(out, "multigraph","->GetYaxis()");
1322 : brun 754 }
1323 :    
1324 : couet 13304
1325 : brun 754 //______________________________________________________________________________
1326 :     void TMultiGraph::SetMaximum(Double_t maximum)
1327 :     {
1328 : couet 13304 // Set multigraph maximum.
1329 :    
1330 : brun 754 fMaximum = maximum;
1331 :     if (fHistogram) fHistogram->SetMaximum(maximum);
1332 :     }
1333 :    
1334 : couet 13304
1335 : brun 754 //______________________________________________________________________________
1336 :     void TMultiGraph::SetMinimum(Double_t minimum)
1337 :     {
1338 : couet 13304 // Set multigraph minimum.
1339 :    
1340 : brun 754 fMinimum = minimum;
1341 :     if (fHistogram) fHistogram->SetMinimum(minimum);
1342 :     }

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9