[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 26029 - (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 :     #include "TVirtualPad.h"
17 : rdm 3748 #include "Riostream.h"
18 : brun 11226 #include "TVirtualFitter.h"
19 : moneta 25487 #include "TPluginManager.h"
20 : pcanal 14336 #include "TClass.h"
21 : brun 17338 #include "TMath.h"
22 : moneta 25487 #include "TSystem.h"
23 : rdm 22419 #include <stdlib.h>
24 : brun 754
25 :     #include <ctype.h>
26 :    
27 : brun 11226 extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
28 : brun 754
29 :     ClassImp(TMultiGraph)
30 :    
31 : couet 13304
32 : brun 754 //______________________________________________________________________________
33 : couet 20297 /* Begin_Html
34 :     <center><h2>TMultiGraph class</h2></center>
35 :     A TMultiGraph is a collection of TGraph (or derived) objects
36 :     Use <tt>TMultiGraph::Add</tt> to add a new graph to the list.
37 :     The TMultiGraph owns the objects in the list.
38 :     Drawing options are the same as for TGraph.
39 :     <p>
40 :     Example:
41 : couet 20433 <pre>
42 : couet 20297 TGraph *gr1 = new TGraph(...
43 :     TGraphErrors *gr2 = new TGraphErrors(...
44 :     TMultiGraph *mg = new TMultiGraph();
45 :     mg->Add(gr1,"lp");
46 :     mg->Add(gr2,"cp");
47 :     mg->Draw("a");
48 :     </pre>
49 :     The drawing option for each TGraph may be specified as an optional
50 :     second argument of the Add function.
51 :     If a draw option is specified, it will be used to draw the graph,
52 :     otherwise the graph will be drawn with the option specified in
53 :     <tt>TMultiGraph::Draw</tt>.
54 :     End_Html */
55 : brun 754
56 : couet 13304
57 : brun 754 //______________________________________________________________________________
58 :     TMultiGraph::TMultiGraph(): TNamed()
59 :     {
60 : couet 13304 // TMultiGraph default constructor
61 : brun 754
62 :     fGraphs = 0;
63 : brun 11226 fFunctions = 0;
64 : brun 754 fHistogram = 0;
65 :     fMaximum = -1111;
66 :     fMinimum = -1111;
67 :     }
68 :    
69 : couet 13304
70 : brun 754 //______________________________________________________________________________
71 :     TMultiGraph::TMultiGraph(const char *name, const char *title)
72 :     : TNamed(name,title)
73 :     {
74 : couet 13304 // constructor with name and title
75 :    
76 : brun 754 fGraphs = 0;
77 : brun 11226 fFunctions = 0;
78 : brun 754 fHistogram = 0;
79 :     fMaximum = -1111;
80 :     fMinimum = -1111;
81 :     }
82 :    
83 : couet 20433
84 : brun 15134 //______________________________________________________________________________
85 :     TMultiGraph::TMultiGraph(const TMultiGraph& mg) :
86 :     TNamed (mg),
87 :     fGraphs(mg.fGraphs),
88 :     fFunctions(mg.fFunctions),
89 :     fHistogram(mg.fHistogram),
90 :     fMaximum(mg.fMaximum),
91 :     fMinimum(mg.fMinimum)
92 : couet 20433 {
93 : brun 15171 //copy constructor
94 :     }
95 : couet 13304
96 : couet 20433
97 : brun 754 //______________________________________________________________________________
98 : brun 15134 TMultiGraph& TMultiGraph::operator=(const TMultiGraph& mg)
99 :     {
100 : brun 15740 //assignement operator
101 : brun 15171 if(this!=&mg) {
102 :     TNamed::operator=(mg);
103 :     fGraphs=mg.fGraphs;
104 :     fFunctions=mg.fFunctions;
105 :     fHistogram=mg.fHistogram;
106 :     fMaximum=mg.fMaximum;
107 :     fMinimum=mg.fMinimum;
108 : couet 20433 }
109 : brun 15171 return *this;
110 : brun 15134 }
111 :    
112 : couet 20433
113 : brun 15134 //______________________________________________________________________________
114 : brun 754 TMultiGraph::~TMultiGraph()
115 :     {
116 : couet 13304 // TMultiGraph destructor
117 : brun 754
118 :     if (!fGraphs) return;
119 : brun 5552 TGraph *g;
120 :     TIter next(fGraphs);
121 :     while ((g = (TGraph*) next())) {
122 : couet 13304 g->ResetBit(kMustCleanup);
123 : brun 5552 }
124 : brun 754 fGraphs->Delete();
125 :     delete fGraphs;
126 :     fGraphs = 0;
127 :     delete fHistogram;
128 :     fHistogram = 0;
129 : brun 11226 if (fFunctions) {
130 :     fFunctions->SetBit(kInvalidObject);
131 :     //special logic to support the case where the same object is
132 :     //added multiple times in fFunctions.
133 :     //This case happens when the same object is added with different
134 :     //drawing modes
135 :     TObject *obj;
136 :     while ((obj = fFunctions->First())) {
137 : rdm 22546 while(fFunctions->Remove(obj)) { }
138 : brun 11226 delete obj;
139 :     }
140 :     delete fFunctions;
141 :     }
142 : brun 754 }
143 :    
144 : couet 13304
145 : brun 754 //______________________________________________________________________________
146 : brun 4037 void TMultiGraph::Add(TGraph *graph, Option_t *chopt)
147 : brun 754 {
148 :     // add a new graph to the list of graphs
149 : brun 6639 // note that the graph is now owned by the TMultigraph.
150 :     // Deleting the TMultiGraph object will automatically delete the graphs.
151 :     // You should not delete the graphs when the TMultigraph is still active.
152 : rdm 3742
153 : brun 754 if (!fGraphs) fGraphs = new TList();
154 : brun 5552 graph->SetBit(kMustCleanup);
155 : brun 4037 fGraphs->Add(graph,chopt);
156 : brun 754 }
157 :    
158 : couet 13304
159 : brun 754 //______________________________________________________________________________
160 : couet 16436 void TMultiGraph::Add(TMultiGraph *multigraph, Option_t *chopt)
161 :     {
162 :     // add all the graphs in "multigraph" to the list of graphs.
163 :    
164 :     TList *graphlist = multigraph->GetListOfGraphs();
165 :     if (!graphlist) return;
166 :    
167 :     if (!fGraphs) fGraphs = new TList();
168 : couet 20433
169 : couet 16436 TGraph *gr;
170 :     gr = (TGraph*)graphlist->First();
171 :     fGraphs->Add(gr,chopt);
172 : couet 22379 for(Int_t i = 1; i < graphlist->GetSize(); i++){
173 : couet 16436 gr = (TGraph*)graphlist->After(gr);
174 :     fGraphs->Add(gr,chopt);
175 :     }
176 :     }
177 :    
178 :    
179 :     //______________________________________________________________________________
180 : brun 754 void TMultiGraph::Browse(TBrowser *)
181 :     {
182 : couet 13304 // Browse multigraph.
183 :    
184 :     Draw("alp");
185 :     gPad->Update();
186 : brun 754 }
187 :    
188 : couet 13304
189 : brun 754 //______________________________________________________________________________
190 :     Int_t TMultiGraph::DistancetoPrimitive(Int_t px, Int_t py)
191 :     {
192 : couet 13304 // Compute distance from point px,py to each graph
193 : brun 754
194 : couet 13304 // Are we on the axis?
195 : brun 754 const Int_t kMaxDiff = 10;
196 :     Int_t distance = 9999;
197 :     if (fHistogram) {
198 :     distance = fHistogram->DistancetoPrimitive(px,py);
199 :     if (distance <= 0) return distance;
200 :     }
201 :    
202 : couet 13304 // Loop on the list of graphs
203 : brun 754 if (!fGraphs) return distance;
204 :     TGraph *g;
205 :     TIter next(fGraphs);
206 :     while ((g = (TGraph*) next())) {
207 :     Int_t dist = g->DistancetoPrimitive(px,py);
208 : brun 4915 if (dist <= 0) return 0;
209 : brun 754 if (dist < kMaxDiff) {gPad->SetSelected(g); return dist;}
210 :     }
211 :     return distance;
212 :     }
213 :    
214 : couet 13304
215 : brun 754 //______________________________________________________________________________
216 :     void TMultiGraph::Draw(Option_t *option)
217 :     {
218 : couet 13304 // Draw this multigraph with its current attributes.
219 :     //
220 :     // Options to draw a graph are described in TGraph::PainGraph
221 :     //
222 :     // The drawing option for each TGraph may be specified as an optional
223 :     // second argument of the Add function. You can use GetGraphDrawOption
224 :     // to return this option.
225 :     // If a draw option is specified, it will be used to draw the graph,
226 :     // otherwise the graph will be drawn with the option specified in
227 :     // TMultiGraph::Draw. Use GetDrawOption to return the option specified
228 :     // when drawin the TMultiGraph.
229 : brun 754
230 : couet 26029 if (!fGraphs) {
231 :     Error("Draw", "Cannot draw an empty TMultiGraph");
232 :     return;
233 :     }
234 : couet 13304 AppendPad(option);
235 : brun 754 }
236 :    
237 : couet 13304
238 : brun 754 //______________________________________________________________________________
239 : brun 11226 Int_t TMultiGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
240 :     {
241 : couet 13304 // Fit this graph with function with name fname.
242 :     //
243 :     // interface to TF1::Fit(TF1 *f1...
244 : brun 11226
245 :     char *linear;
246 : brun 11231 linear= (char*)strstr(fname, "++");
247 : brun 11226 TF1 *f1=0;
248 :     if (linear)
249 : brun 12643 f1=new TF1(fname, fname, xmin, xmax);
250 : brun 11226 else {
251 :     f1 = (TF1*)gROOT->GetFunction(fname);
252 :     if (!f1) { Printf("Unknown function: %s",fname); return -1; }
253 :     }
254 :    
255 :     return Fit(f1,option,"",xmin,xmax);
256 :     }
257 :    
258 : couet 13304
259 : brun 11226 //______________________________________________________________________________
260 : moneta 25487 Int_t TMultiGraph::Fit(TF1 *f1, Option_t *option, Option_t *goption, Axis_t rxmin, Axis_t rxmax)
261 : brun 11226 {
262 : couet 13304 // Fit this multigraph with function f1.
263 :     //
264 :     // In this function all graphs of the multigraph are fitted simultaneously
265 :     //
266 :     // f1 is an already predefined function created by TF1.
267 :     // Predefined functions such as gaus, expo and poln are automatically
268 :     // created by ROOT.
269 :     //
270 :     // The list of fit options is given in parameter option.
271 :     // option = "W" Set all errors to 1
272 :     // = "U" Use a User specified fitting algorithm (via SetFCN)
273 :     // = "Q" Quiet mode (minimum printing)
274 :     // = "V" Verbose mode (default is between Q and V)
275 :     // = "B" Use this option when you want to fix one or more parameters
276 :     // and the fitting function is like "gaus","expo","poln","landau".
277 :     // = "R" Use the Range specified in the function range
278 :     // = "N" Do not store the graphics function, do not draw
279 :     // = "0" Do not plot the result of the fit. By default the fitted function
280 :     // is drawn unless the option"N" above is specified.
281 :     // = "+" Add this new fitted function to the list of fitted functions
282 :     // (by default, any previous function is deleted)
283 :     // = "C" In case of linear fitting, not calculate the chisquare
284 :     // (saves time)
285 :     // = "F" If fitting a polN, switch to minuit fitter
286 :     // = "ROB" In case of linear fitting, compute the LTS regression
287 : couet 20433 // coefficients (robust(resistant) regression), using
288 : couet 13304 // the default fraction of good points
289 :     // "ROB=0.x" - compute the LTS regression coefficients, using
290 :     // 0.x as a fraction of good points
291 :     //
292 :     // When the fit is drawn (by default), the parameter goption may be used
293 :     // to specify a list of graphics options. See TGraph::Paint for a complete
294 :     // list of these options.
295 :     //
296 :     // In order to use the Range option, one must first create a function
297 :     // with the expression to be fitted. For example, if your graph
298 :     // has a defined range between -4 and 4 and you want to fit a gaussian
299 :     // only in the interval 1 to 3, you can do:
300 :     // TF1 *f1 = new TF1("f1","gaus",1,3);
301 :     // graph->Fit("f1","R");
302 :     //
303 :     //
304 :     // who is calling this function
305 :     // ============================
306 :     // Note that this function is called when calling TGraphErrors::Fit
307 :     // or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
308 :     // see the discussion below on the errors calulation.
309 :     //
310 :     // Setting initial conditions
311 :     // ==========================
312 :     // Parameters must be initialized before invoking the Fit function.
313 :     // The setting of the parameter initial values is automatic for the
314 :     // predefined functions : poln, expo, gaus, landau. One can however disable
315 :     // this automatic computation by specifying the option "B".
316 :     // You can specify boundary limits for some or all parameters via
317 :     // f1->SetParLimits(p_number, parmin, parmax);
318 :     // if parmin>=parmax, the parameter is fixed
319 :     // Note that you are not forced to fix the limits for all parameters.
320 :     // For example, if you fit a function with 6 parameters, you can do:
321 :     // func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
322 :     // func->SetParLimits(4,-10,-4);
323 :     // func->SetParLimits(5, 1,1);
324 :     // With this setup, parameters 0->3 can vary freely
325 :     // Parameter 4 has boundaries [-10,-4] with initial value -8
326 :     // Parameter 5 is fixed to 100.
327 :     //
328 :     // Fit range
329 :     // =========
330 :     // The fit range can be specified in two ways:
331 :     // - specify rxmax > rxmin (default is rxmin=rxmax=0)
332 :     // - specify the option "R". In this case, the function will be taken
333 :     // instead of the full graph range.
334 :     //
335 :     // Changing the fitting function
336 :     // =============================
337 :     // By default the fitting function GraphFitChisquare is used.
338 :     // To specify a User defined fitting function, specify option "U" and
339 :     // call the following functions:
340 :     // TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
341 :     // where MyFittingFunction is of type:
342 :     // extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
343 :     //
344 :     // How errors are used in the chisquare function (see TFitter GraphFitChisquare)// Access to the fit results
345 :     // ============================================
346 :     // In case of a TGraphErrors object, ex, the error along x, is projected
347 :     // along the y-direction by calculating the function at the points x-exlow and
348 :     // x+exhigh.
349 :     //
350 :     // The chisquare is computed as the sum of the quantity below at each point:
351 :     //
352 :     // (y - f(x))**2
353 :     // -----------------------------------
354 :     // ey**2 + ((f(x+exhigh) - f(x-exlow))/2)**2
355 :     //
356 :     // where x and y are the point coordinates.
357 :     //
358 :     // In case the function lies below (above) the data point, ey is ey_low (ey_high).
359 :     //
360 :     // thanks to Andy Haas (haas@yahoo.com) for adding the case with TGraphasymmerrors
361 :     // University of Washington
362 :     //
363 :     // a little different approach to approximating the uncertainty in y because of the
364 :     // errors in x, is to make it equal the error in x times the slope of the line.
365 :     // The improvement, compared to the first method (f(x+ exhigh) - f(x-exlow))/2
366 :     // is of (error of x)**2 order. This approach is called "effective variance method".
367 :     // This improvement has been made in version 4.00/08 by Anna Kreshuk.
368 :     //
369 :     // Associated functions
370 :     // ====================
371 :     // One or more object (typically a TF1*) can be added to the list
372 :     // of functions (fFunctions) associated to each graph.
373 :     // When TGraph::Fit is invoked, the fitted function is added to this list.
374 :     // Given a graph gr, one can retrieve an associated function
375 :     // with: TF1 *myfunc = gr->GetFunction("myfunc");
376 :     //
377 :     // If the graph is made persistent, the list of
378 :     // associated functions is also persistent. Given a pointer (see above)
379 :     // to an associated function myfunc, one can retrieve the function/fit
380 :     // parameters with calls such as:
381 :     // Double_t chi2 = myfunc->GetChisquare();
382 :     // Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
383 :     // Double_t err0 = myfunc->GetParError(0); //error on first parameter
384 :     //
385 :     // Fit Statistics
386 :     // ==============
387 :     // You can change the statistics box to display the fit parameters with
388 :     // the TStyle::SetOptFit(mode) method. This mode has four digits.
389 :     // mode = pcev (default = 0111)
390 :     // v = 1; print name/values of parameters
391 :     // e = 1; print errors (if e=1, v must be 1)
392 :     // c = 1; print Chisquare/Number of degress of freedom
393 :     // p = 1; print Probability
394 :     //
395 :     // For example: gStyle->SetOptFit(1011);
396 :     // prints the fit probability, parameter names/values, and errors.
397 :     // You can change the position of the statistics box with these lines
398 :     // (where g is a pointer to the TGraph):
399 :     //
400 :     // Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
401 :     // Root > st->SetX1NDC(newx1); //new x start position
402 :     // Root > st->SetX2NDC(newx2); //new x end position
403 : brun 11226
404 : moneta 25487 return DoFit(f1,option,goption,rxmin,rxmax); // implemented in HFitImpl.cxx
405 : rdm 11750
406 : moneta 25487 }
407 : rdm 11750
408 : moneta 25487 //______________________________________________________________________________
409 :     void TMultiGraph::FitPanel()
410 :     {
411 :     // -*-*-*-*-*Display a panel with all histogram fit options*-*-*-*-*-*
412 :     // ==============================================
413 :     //
414 :     // See class TFitPanel for example
415 : rdm 11750
416 : moneta 25487 if (!gPad)
417 :     gROOT->MakeDefCanvas();
418 : brun 11226
419 : moneta 25487 if (!gPad) {
420 :     Error("FitPanel", "Unable to create a default canvas");
421 :     return;
422 : brun 12136 }
423 : brun 11226
424 : moneta 25487 // use plugin manager to create instance of TFitEditor
425 :     TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
426 :     if (handler && handler->LoadPlugin() != -1) {
427 :     if (handler->ExecPlugin(2, gPad, this) == 0)
428 :     Error("FitPanel", "Unable to crate the FitPanel");
429 : brun 11226 }
430 : couet 26029 else
431 : moneta 25487 Error("FitPanel", "Unable to find the FitPanel plug-in");
432 : brun 11226 }
433 :    
434 :     //______________________________________________________________________________
435 : brun 11584 Option_t *TMultiGraph::GetGraphDrawOption(const TGraph *gr) const
436 :     {
437 : couet 13304 // Return the draw option for the TGraph gr in this TMultiGraph
438 :     // The return option is the one specified when calling TMultiGraph::Add(gr,option).
439 : rdm 11750
440 : brun 11584 if (!fGraphs || !gr) return "";
441 :     TListIter next(fGraphs);
442 :     TObject *obj;
443 :     while ((obj = next())) {
444 : brun 11600 if (obj == (TObject*)gr) return next.GetOption();
445 : brun 11584 }
446 :     return "";
447 :     }
448 :    
449 : couet 13304
450 : brun 11584 //______________________________________________________________________________
451 : brun 11226 void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
452 :     {
453 : couet 13304 // Compute Initial values of parameters for a gaussian.
454 : brun 11226
455 :     Double_t allcha, sumx, sumx2, x, val, rms, mean;
456 :     Int_t bin;
457 :     const Double_t sqrtpi = 2.506628;
458 :    
459 : couet 13304 // Compute mean value and RMS of the graph in the given range
460 : brun 11226 Int_t np = 0;
461 :     allcha = sumx = sumx2 = 0;
462 :     TGraph *g;
463 :     TIter next(fGraphs);
464 :     Double_t *px, *py;
465 :     Int_t npp; //number of points in each graph
466 :     while ((g = (TGraph*) next())) {
467 :     px=g->GetX();
468 :     py=g->GetY();
469 :     npp=g->GetN();
470 :     for (bin=0; bin<npp; bin++){
471 : brun 12643 x=px[bin];
472 :     if (x<xmin || x>xmax) continue;
473 :     np++;
474 :     val=py[bin];
475 :     sumx+=val*x;
476 :     sumx2+=val*x*x;
477 :     allcha+=val;
478 : brun 11226 }
479 :     }
480 :     if (np == 0 || allcha == 0) return;
481 :     mean = sumx/allcha;
482 :     rms = TMath::Sqrt(sumx2/allcha - mean*mean);
483 :    
484 :     Double_t binwidx = TMath::Abs((xmax-xmin)/np);
485 :     if (rms == 0) rms = 1;
486 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
487 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
488 :     f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
489 :     f1->SetParameter(1,mean);
490 :     f1->SetParameter(2,rms);
491 :     f1->SetParLimits(2,0,10*rms);
492 :     }
493 :    
494 : couet 13304
495 : brun 11226 //______________________________________________________________________________
496 :     void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
497 :     {
498 : couet 13304 // Compute Initial values of parameters for an exponential.
499 : brun 11226
500 :     Double_t constant, slope;
501 :     Int_t ifail;
502 :    
503 :     LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
504 :    
505 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
506 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
507 :     f1->SetParameter(0,constant);
508 :     f1->SetParameter(1,slope);
509 :     }
510 :    
511 : couet 13304
512 : brun 11226 //______________________________________________________________________________
513 :     void TMultiGraph::InitPolynom(Double_t xmin, Double_t xmax)
514 :     {
515 : couet 13304 // Compute Initial values of parameters for a polynom.
516 : brun 11226
517 :     Double_t fitpar[25];
518 :    
519 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
520 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
521 :     Int_t npar = f1->GetNpar();
522 :    
523 :     LeastSquareFit(npar, fitpar, xmin, xmax);
524 :    
525 :     for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
526 :     }
527 :    
528 : couet 13304
529 : brun 11226 //______________________________________________________________________________
530 :     void TMultiGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
531 :     {
532 : couet 13304 // Least squares lpolynomial fitting without weights.
533 :     //
534 :     // m number of parameters
535 :     // a array of parameters
536 :     // first 1st point number to fit (default =0)
537 :     // last last point number to fit (default=fNpoints-1)
538 :     //
539 :     // based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
540 : brun 11226
541 : couet 13304 const Double_t zero = 0.;
542 :     const Double_t one = 1.;
543 :     const Int_t idim = 20;
544 : brun 11226
545 : couet 13304 Double_t b[400] /* was [20][20] */;
546 :     Int_t i, k, l, ifail, bin;
547 :     Double_t power;
548 :     Double_t da[20], xk, yk;
549 : brun 11226
550 :    
551 : couet 13304 //count the total number of points to fit
552 :     TGraph *g;
553 :     TIter next(fGraphs);
554 :     Double_t *px, *py;
555 :     Int_t n=0;
556 :     Int_t npp;
557 :     while ((g = (TGraph*) next())) {
558 :     px=g->GetX();
559 :     py=g->GetY();
560 :     npp=g->GetN();
561 :     for (bin=0; bin<npp; bin++){
562 :     xk=px[bin];
563 :     if (xk < xmin || xk > xmax) continue;
564 :     n++;
565 :     }
566 :     }
567 :     if (m <= 2) {
568 :     LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
569 :     return;
570 :     }
571 :     if (m > idim || m > n) return;
572 :     da[0] = zero;
573 :     for (l = 2; l <= m; ++l) {
574 :     b[l-1] = zero;
575 :     b[m + l*20 - 21] = zero;
576 :     da[l-1] = zero;
577 :     }
578 :     Int_t np = 0;
579 : brun 11226
580 : couet 13304 next.Reset();
581 :     while ((g = (TGraph*) next())) {
582 :     px=g->GetX();
583 :     py=g->GetY();
584 :     npp=g->GetN();
585 :    
586 :     for (k = 0; k <= npp; ++k) {
587 :     xk = px[k];
588 :     if (xk < xmin || xk > xmax) continue;
589 :     np++;
590 :     yk = py[k];
591 :     power = one;
592 :     da[0] += yk;
593 :     for (l = 2; l <= m; ++l) {
594 :     power *= xk;
595 : brun 12643 b[l-1] += power;
596 :     da[l-1] += power*yk;
597 : couet 13304 }
598 :     for (l = 2; l <= m; ++l) {
599 :     power *= xk;
600 :     b[m + l*20 - 21] += power;
601 :     }
602 :     }
603 :     }
604 :     b[0] = Double_t(np);
605 :     for (i = 3; i <= m; ++i) {
606 :     for (k = i; k <= m; ++k) {
607 :     b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
608 :     }
609 :     }
610 :     H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
611 : brun 11226
612 : couet 13304 if (ifail < 0) {
613 :     //a[0] = fY[0];
614 :     py=((TGraph *)fGraphs->First())->GetY();
615 :     a[0]=py[0];
616 :     for (i=1; i<m; ++i) a[i] = 0;
617 :     return;
618 :     }
619 :     for (i=0; i<m; ++i) a[i] = da[i];
620 : brun 11226 }
621 :    
622 : couet 13304
623 : brun 11226 //______________________________________________________________________________
624 :     void TMultiGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
625 :     {
626 : couet 13304 // Least square linear fit without weights.
627 :     //
628 :     // Fit a straight line (a0 + a1*x) to the data in this graph.
629 :     // ndata: number of points to fit
630 :     // first: first point number to fit
631 :     // last: last point to fit O(ndata should be last-first
632 :     // ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
633 :     //
634 :     // extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
635 : brun 11226
636 : couet 13304 Double_t xbar, ybar, x2bar;
637 :     Int_t i;
638 :     Double_t xybar;
639 :     Double_t fn, xk, yk;
640 :     Double_t det;
641 : brun 11226
642 : couet 13304 ifail = -2;
643 :     xbar = ybar = x2bar = xybar = 0;
644 :     Int_t np = 0;
645 :     TGraph *g;
646 :     TIter next(fGraphs);
647 :     Double_t *px, *py;
648 :     Int_t npp;
649 :     while ((g = (TGraph*) next())) {
650 :     px=g->GetX();
651 :     py=g->GetY();
652 :     npp=g->GetN();
653 :     for (i = 0; i < npp; ++i) {
654 :     xk = px[i];
655 :     if (xk < xmin || xk > xmax) continue;
656 :     np++;
657 :     yk = py[i];
658 :     if (ndata < 0) {
659 :     if (yk <= 0) yk = 1e-9;
660 :     yk = TMath::Log(yk);
661 :     }
662 :     xbar += xk;
663 :     ybar += yk;
664 :     x2bar += xk*xk;
665 :     xybar += xk*yk;
666 :     }
667 :     }
668 :     fn = Double_t(np);
669 :     det = fn*x2bar - xbar*xbar;
670 :     ifail = -1;
671 :     if (det <= 0) {
672 :     if (fn > 0) a0 = ybar/fn;
673 :     else a0 = 0;
674 :     a1 = 0;
675 :     return;
676 :     }
677 :     ifail = 0;
678 :     a0 = (x2bar*ybar - xbar*xybar) / det;
679 :     a1 = (fn*xybar - xbar*ybar) / det;
680 :     }
681 : brun 11226
682 :    
683 :     //______________________________________________________________________________
684 : brun 1205 TH1F *TMultiGraph::GetHistogram() const
685 : brun 754 {
686 : couet 13304 // Returns a pointer to the histogram used to draw the axis
687 :     // Takes into account the two following cases.
688 :     // 1- option 'A' was specified in TMultiGraph::Draw. Return fHistogram
689 :     // 2- user had called TPad::DrawFrame. return pointer to hframe histogram
690 : brun 754
691 :     if (fHistogram) return fHistogram;
692 :     if (!gPad) return 0;
693 :     gPad->Modified();
694 :     gPad->Update();
695 :     if (fHistogram) return fHistogram;
696 :     TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
697 :     return h1;
698 :     }
699 :    
700 : couet 13304
701 : brun 754 //______________________________________________________________________________
702 : brun 11226 TF1 *TMultiGraph::GetFunction(const char *name) const
703 :     {
704 : couet 13304 // Return pointer to function with name.
705 :     //
706 :     // Functions such as TGraph::Fit store the fitted function in the list of
707 :     // functions of this graph.
708 : brun 11226
709 :     if (!fFunctions) return 0;
710 :     return (TF1*)fFunctions->FindObject(name);
711 :     }
712 :    
713 : moneta 25487 //______________________________________________________________________________
714 : couet 26029 TList *TMultiGraph::GetListOfFunctions()
715 : moneta 25487 {
716 :     // Return pointer to list of functions
717 :     // if pointer is null create the list
718 : couet 13304
719 : moneta 25487 if (!fFunctions) fFunctions = new TList();
720 : couet 26029 return fFunctions;
721 : moneta 25487 }
722 :    
723 :    
724 : brun 11226 //______________________________________________________________________________
725 : brun 1205 TAxis *TMultiGraph::GetXaxis() const
726 : brun 754 {
727 :     // Get x axis of the graph.
728 :    
729 :     if (!gPad) return 0;
730 : brun 4083 TH1 *h = GetHistogram();
731 :     if (!h) return 0;
732 :     return h->GetXaxis();
733 : brun 754 }
734 :    
735 : couet 13304
736 : brun 754 //______________________________________________________________________________
737 : brun 1205 TAxis *TMultiGraph::GetYaxis() const
738 : brun 754 {
739 :     // Get y axis of the graph.
740 :    
741 :     if (!gPad) return 0;
742 : brun 4083 TH1 *h = GetHistogram();
743 :     if (!h) return 0;
744 :     return h->GetYaxis();
745 : brun 754 }
746 :    
747 : couet 13304
748 : brun 754 //______________________________________________________________________________
749 :     void TMultiGraph::Paint(Option_t *option)
750 :     {
751 : couet 13304 // paint all the graphs of this multigraph
752 : brun 754
753 : couet 13304 if (fGraphs->GetSize() == 0) return;
754 : rdm 11750
755 : couet 13304 char *l;
756 :     static char chopt[33];
757 :     Int_t nch = strlen(option);
758 :     Int_t i;
759 :     for (i=0;i<nch;i++) chopt[i] = toupper(option[i]);
760 :     chopt[nch] = 0;
761 :     TGraph *g;
762 : rdm 3742
763 : couet 13304 l = strstr(chopt,"A");
764 :     if (l) {
765 :     *l = ' ';
766 :     TIter next(fGraphs);
767 :     Int_t npt = 100;
768 :     Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
769 :     rwxmin = gPad->GetUxmin();
770 :     rwxmax = gPad->GetUxmax();
771 :     rwymin = gPad->GetUymin();
772 :     rwymax = gPad->GetUymax();
773 :     char *xtitle = 0;
774 :     char *ytitle = 0;
775 :     Int_t firstx = 0;
776 :     Int_t lastx = 0;
777 : rdm 11750
778 : couet 13304 if (fHistogram) {
779 :     //cleanup in case of a previous unzoom
780 :     if (fHistogram->GetMinimum() >= fHistogram->GetMaximum()) {
781 : brun 23506 nch = strlen(fHistogram->GetXaxis()->GetTitle());
782 : couet 13304 firstx = fHistogram->GetXaxis()->GetFirst();
783 :     lastx = fHistogram->GetXaxis()->GetLast();
784 :     if (nch) {
785 :     xtitle = new char[nch+1];
786 :     strcpy(xtitle,fHistogram->GetXaxis()->GetTitle());
787 :     }
788 :     nch = strlen(fHistogram->GetYaxis()->GetTitle());
789 :     if (nch) {
790 :     ytitle = new char[nch+1];
791 :     strcpy(ytitle,fHistogram->GetYaxis()->GetTitle());
792 :     }
793 :     delete fHistogram;
794 :     fHistogram = 0;
795 :     }
796 :     }
797 :     if (fHistogram) {
798 :     minimum = fHistogram->GetYaxis()->GetXmin();
799 :     maximum = fHistogram->GetYaxis()->GetXmax();
800 :     uxmin = gPad->PadtoX(rwxmin);
801 :     uxmax = gPad->PadtoX(rwxmax);
802 :     } else {
803 : couet 21573 g = (TGraph*) next();
804 :     g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
805 : couet 13304 while ((g = (TGraph*) next())) {
806 : couet 20433 Double_t rx1,ry1,rx2,ry2;
807 :     g->ComputeRange(rx1, ry1, rx2, ry2);
808 :     if (rx1 < rwxmin) rwxmin = rx1;
809 :     if (ry1 < rwymin) rwymin = ry1;
810 :     if (rx2 > rwxmax) rwxmax = rx2;
811 :     if (ry2 > rwymax) rwymax = ry2;
812 : couet 13304 if (g->GetN() > npt) npt = g->GetN();
813 :     }
814 :     if (rwxmin == rwxmax) rwxmax += 1.;
815 :     if (rwymin == rwymax) rwymax += 1.;
816 :     dx = 0.05*(rwxmax-rwxmin);
817 :     dy = 0.05*(rwymax-rwymin);
818 :     uxmin = rwxmin - dx;
819 :     uxmax = rwxmax + dx;
820 :     if (gPad->GetLogy()) {
821 :     if (rwymin <= 0) rwymin = 0.001*rwymax;
822 :     minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
823 :     maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
824 :     } else {
825 :     minimum = rwymin - dy;
826 :     maximum = rwymax + dy;
827 :     }
828 :     if (minimum < 0 && rwymin >= 0) minimum = 0;
829 :     if (maximum > 0 && rwymax <= 0) maximum = 0;
830 :     }
831 : brun 754
832 : couet 13304 if (fMinimum != -1111) rwymin = minimum = fMinimum;
833 :     if (fMaximum != -1111) rwymax = maximum = fMaximum;
834 :     if (uxmin < 0 && rwxmin >= 0) {
835 :     if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
836 :     //else uxmin = 0;
837 :     }
838 :     if (uxmax > 0 && rwxmax <= 0) {
839 :     if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
840 :     //else uxmax = 0;
841 :     }
842 :     if (minimum < 0 && rwymin >= 0) {
843 :     if(gPad->GetLogy()) minimum = 0.9*rwymin;
844 :     //else minimum = 0;
845 :     }
846 :     if (maximum > 0 && rwymax <= 0) {
847 :     if(gPad->GetLogy()) maximum = 1.1*rwymax;
848 :     //else maximum = 0;
849 :     }
850 :     if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
851 :     if (uxmin <= 0 && gPad->GetLogx()) {
852 :     if (uxmax > 1000) uxmin = 1;
853 :     else uxmin = 0.001*uxmax;
854 :     }
855 :     rwymin = minimum;
856 :     rwymax = maximum;
857 :     if (fHistogram) {
858 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
859 :     }
860 : brun 3575
861 : couet 13304 // Create a temporary histogram to draw the axis
862 :     if (!fHistogram) {
863 :     // the graph is created with at least as many channels as there are points
864 :     // to permit zooming on the full range
865 :     rwxmin = uxmin;
866 :     rwxmax = uxmax;
867 :     fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
868 :     if (!fHistogram) return;
869 :     fHistogram->SetMinimum(rwymin);
870 :     fHistogram->SetBit(TH1::kNoStats);
871 :     fHistogram->SetMaximum(rwymax);
872 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
873 :     fHistogram->SetDirectory(0);
874 :     if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
875 :     if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
876 :     if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
877 :     }
878 :     fHistogram->Paint("0");
879 : brun 3575 }
880 : rdm 3742
881 : couet 21962 TGraph *gfit = 0;
882 : brun 754 if (fGraphs) {
883 : brun 4037 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
884 : couet 21962 TObject *obj = 0;
885 : brun 4037
886 :     while (lnk) {
887 :     obj = lnk->GetObject();
888 :     if (strlen(lnk->GetOption())) obj->Paint(lnk->GetOption());
889 :     else obj->Paint(chopt);
890 :     lnk = (TObjOptLink*)lnk->Next();
891 :     }
892 : couet 21962 gfit = (TGraph*)obj; // pick one TGraph in the list to paint the fit parameters.
893 : brun 754 }
894 : brun 11226
895 :     TObject *f;
896 : couet 21962 TF1 *fit = 0;
897 : brun 11226 if (fFunctions) {
898 : couet 13304 TIter next(fFunctions);
899 :     while ((f = (TObject*) next())) {
900 :     if (f->InheritsFrom(TF1::Class())) {
901 :     if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
902 : couet 21962 fit = (TF1*)f;
903 : couet 13304 } else {
904 :     f->Paint();
905 :     }
906 : brun 11226 }
907 :     }
908 : couet 21962
909 : brun 24084 if (fit) gfit->PaintStats(fit);
910 : couet 13304 }
911 : brun 11226
912 :    
913 : brun 754 //______________________________________________________________________________
914 : brun 1205 void TMultiGraph::Print(Option_t *option) const
915 : brun 754 {
916 : couet 13304 // Print the list of graphs
917 : brun 754
918 :     TGraph *g;
919 :     if (fGraphs) {
920 : couet 13304 TIter next(fGraphs);
921 :     while ((g = (TGraph*) next())) {
922 :     g->Print(option);
923 :     }
924 : brun 754 }
925 :     }
926 :    
927 : couet 13304
928 : brun 754 //______________________________________________________________________________
929 : brun 5552 void TMultiGraph::RecursiveRemove(TObject *obj)
930 :     {
931 :     // Recursively remove this object from a list. Typically implemented
932 :     // by classes that can contain mulitple references to a same object.
933 :    
934 :     if (!fGraphs) return;
935 :     TObject *objr = fGraphs->Remove(obj);
936 :     if (!objr) return;
937 :     delete fHistogram; fHistogram = 0;
938 :     if (gPad) gPad->Modified();
939 :     }
940 :    
941 : couet 13304
942 : brun 5552 //______________________________________________________________________________
943 : brun 15672 void TMultiGraph::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
944 : brun 754 {
945 : couet 13304 // Save primitive as a C++ statement(s) on output stream out
946 : brun 754
947 :     char quote = '"';
948 :     out<<" "<<endl;
949 :     if (gROOT->ClassSaved(TMultiGraph::Class())) {
950 : couet 13304 out<<" ";
951 : brun 754 } else {
952 : couet 13304 out<<" TMultiGraph *";
953 : brun 754 }
954 :     out<<"multigraph = new TMultiGraph();"<<endl;
955 :     out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<endl;
956 :     out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
957 :    
958 :     if (fGraphs) {
959 : brun 10847 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
960 :     TObject *g;
961 :    
962 :     while (lnk) {
963 :     g = lnk->GetObject();
964 : couet 20297 g->SavePrimitive(out, Form("multigraph%s",lnk->GetOption()));
965 : brun 10847 lnk = (TObjOptLink*)lnk->Next();
966 : couet 13304 }
967 : brun 754 }
968 : couet 20297 out<<" multigraph->Draw(" <<quote<<option<<quote<<");"<<endl;
969 : couet 14017
970 :     TAxis *xaxis = GetXaxis();
971 :     TAxis *yaxis = GetYaxis();
972 : couet 20433
973 :     if (xaxis) xaxis->SaveAttributes(out, "multigraph","->GetXaxis()");
974 :     if (yaxis) yaxis->SaveAttributes(out, "multigraph","->GetYaxis()");
975 : brun 754 }
976 :    
977 : couet 13304
978 : brun 754 //______________________________________________________________________________
979 :     void TMultiGraph::SetMaximum(Double_t maximum)
980 :     {
981 : couet 13304 // Set multigraph maximum.
982 :    
983 : brun 754 fMaximum = maximum;
984 :     if (fHistogram) fHistogram->SetMaximum(maximum);
985 :     }
986 :    
987 : couet 13304
988 : brun 754 //______________________________________________________________________________
989 :     void TMultiGraph::SetMinimum(Double_t minimum)
990 :     {
991 : couet 13304 // Set multigraph minimum.
992 :    
993 : brun 754 fMinimum = minimum;
994 :     if (fHistogram) fHistogram->SetMinimum(minimum);
995 :     }

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9