[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 33513 - (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 : moneta 31207 #include "HFitInterface.h"
26 :     #include "Fit/DataRange.h"
27 :     #include "Math/MinimizerOptions.h"
28 :    
29 : brun 754 #include <ctype.h>
30 :    
31 : brun 11226 extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
32 : brun 754
33 :     ClassImp(TMultiGraph)
34 :    
35 : couet 13304
36 : brun 754 //______________________________________________________________________________
37 : couet 20297 /* Begin_Html
38 :     <center><h2>TMultiGraph class</h2></center>
39 :     A TMultiGraph is a collection of TGraph (or derived) objects
40 :     Use <tt>TMultiGraph::Add</tt> to add a new graph to the list.
41 :     The TMultiGraph owns the objects in the list.
42 :     Drawing options are the same as for TGraph.
43 :     <p>
44 :     Example:
45 : couet 20433 <pre>
46 : couet 20297 TGraph *gr1 = new TGraph(...
47 :     TGraphErrors *gr2 = new TGraphErrors(...
48 :     TMultiGraph *mg = new TMultiGraph();
49 :     mg->Add(gr1,"lp");
50 :     mg->Add(gr2,"cp");
51 :     mg->Draw("a");
52 :     </pre>
53 :     The drawing option for each TGraph may be specified as an optional
54 :     second argument of the Add function.
55 :     If a draw option is specified, it will be used to draw the graph,
56 :     otherwise the graph will be drawn with the option specified in
57 :     <tt>TMultiGraph::Draw</tt>.
58 :     End_Html */
59 : brun 754
60 : couet 13304
61 : brun 754 //______________________________________________________________________________
62 :     TMultiGraph::TMultiGraph(): TNamed()
63 :     {
64 : couet 13304 // TMultiGraph default constructor
65 : brun 754
66 :     fGraphs = 0;
67 : brun 11226 fFunctions = 0;
68 : brun 754 fHistogram = 0;
69 :     fMaximum = -1111;
70 :     fMinimum = -1111;
71 :     }
72 :    
73 : couet 13304
74 : brun 754 //______________________________________________________________________________
75 :     TMultiGraph::TMultiGraph(const char *name, const char *title)
76 :     : TNamed(name,title)
77 :     {
78 : couet 13304 // constructor with name and title
79 :    
80 : brun 754 fGraphs = 0;
81 : brun 11226 fFunctions = 0;
82 : brun 754 fHistogram = 0;
83 :     fMaximum = -1111;
84 :     fMinimum = -1111;
85 :     }
86 :    
87 : couet 20433
88 : brun 15134 //______________________________________________________________________________
89 :     TMultiGraph::TMultiGraph(const TMultiGraph& mg) :
90 :     TNamed (mg),
91 :     fGraphs(mg.fGraphs),
92 :     fFunctions(mg.fFunctions),
93 :     fHistogram(mg.fHistogram),
94 :     fMaximum(mg.fMaximum),
95 :     fMinimum(mg.fMinimum)
96 : couet 20433 {
97 : brun 15171 //copy constructor
98 :     }
99 : couet 13304
100 : couet 20433
101 : brun 754 //______________________________________________________________________________
102 : brun 15134 TMultiGraph& TMultiGraph::operator=(const TMultiGraph& mg)
103 :     {
104 : brun 15740 //assignement operator
105 : brun 15171 if(this!=&mg) {
106 :     TNamed::operator=(mg);
107 :     fGraphs=mg.fGraphs;
108 :     fFunctions=mg.fFunctions;
109 :     fHistogram=mg.fHistogram;
110 :     fMaximum=mg.fMaximum;
111 :     fMinimum=mg.fMinimum;
112 : couet 20433 }
113 : brun 15171 return *this;
114 : brun 15134 }
115 :    
116 : couet 20433
117 : brun 15134 //______________________________________________________________________________
118 : brun 754 TMultiGraph::~TMultiGraph()
119 :     {
120 : couet 13304 // TMultiGraph destructor
121 : brun 754
122 :     if (!fGraphs) return;
123 : brun 5552 TGraph *g;
124 :     TIter next(fGraphs);
125 :     while ((g = (TGraph*) next())) {
126 : couet 13304 g->ResetBit(kMustCleanup);
127 : brun 5552 }
128 : brun 754 fGraphs->Delete();
129 :     delete fGraphs;
130 :     fGraphs = 0;
131 :     delete fHistogram;
132 :     fHistogram = 0;
133 : brun 11226 if (fFunctions) {
134 :     fFunctions->SetBit(kInvalidObject);
135 :     //special logic to support the case where the same object is
136 :     //added multiple times in fFunctions.
137 :     //This case happens when the same object is added with different
138 :     //drawing modes
139 :     TObject *obj;
140 :     while ((obj = fFunctions->First())) {
141 : rdm 22546 while(fFunctions->Remove(obj)) { }
142 : brun 11226 delete obj;
143 :     }
144 :     delete fFunctions;
145 :     }
146 : brun 754 }
147 :    
148 : couet 13304
149 : brun 754 //______________________________________________________________________________
150 : brun 4037 void TMultiGraph::Add(TGraph *graph, Option_t *chopt)
151 : brun 754 {
152 :     // add a new graph to the list of graphs
153 : brun 6639 // note that the graph is now owned by the TMultigraph.
154 :     // Deleting the TMultiGraph object will automatically delete the graphs.
155 :     // You should not delete the graphs when the TMultigraph is still active.
156 : rdm 3742
157 : brun 754 if (!fGraphs) fGraphs = new TList();
158 : brun 5552 graph->SetBit(kMustCleanup);
159 : brun 4037 fGraphs->Add(graph,chopt);
160 : brun 754 }
161 :    
162 : couet 13304
163 : brun 754 //______________________________________________________________________________
164 : couet 16436 void TMultiGraph::Add(TMultiGraph *multigraph, Option_t *chopt)
165 :     {
166 :     // add all the graphs in "multigraph" to the list of graphs.
167 :    
168 :     TList *graphlist = multigraph->GetListOfGraphs();
169 :     if (!graphlist) return;
170 :    
171 :     if (!fGraphs) fGraphs = new TList();
172 : couet 20433
173 : couet 16436 TGraph *gr;
174 :     gr = (TGraph*)graphlist->First();
175 :     fGraphs->Add(gr,chopt);
176 : couet 22379 for(Int_t i = 1; i < graphlist->GetSize(); i++){
177 : couet 16436 gr = (TGraph*)graphlist->After(gr);
178 :     fGraphs->Add(gr,chopt);
179 :     }
180 :     }
181 :    
182 :    
183 :     //______________________________________________________________________________
184 : brun 754 void TMultiGraph::Browse(TBrowser *)
185 :     {
186 : couet 13304 // Browse multigraph.
187 :    
188 :     Draw("alp");
189 :     gPad->Update();
190 : brun 754 }
191 :    
192 : couet 13304
193 : brun 754 //______________________________________________________________________________
194 :     Int_t TMultiGraph::DistancetoPrimitive(Int_t px, Int_t py)
195 :     {
196 : couet 13304 // Compute distance from point px,py to each graph
197 : brun 754
198 : couet 13304 // Are we on the axis?
199 : brun 754 const Int_t kMaxDiff = 10;
200 :     Int_t distance = 9999;
201 :     if (fHistogram) {
202 :     distance = fHistogram->DistancetoPrimitive(px,py);
203 :     if (distance <= 0) return distance;
204 :     }
205 :    
206 : couet 13304 // Loop on the list of graphs
207 : brun 754 if (!fGraphs) return distance;
208 :     TGraph *g;
209 :     TIter next(fGraphs);
210 :     while ((g = (TGraph*) next())) {
211 :     Int_t dist = g->DistancetoPrimitive(px,py);
212 : brun 4915 if (dist <= 0) return 0;
213 : brun 754 if (dist < kMaxDiff) {gPad->SetSelected(g); return dist;}
214 :     }
215 :     return distance;
216 :     }
217 :    
218 : couet 13304
219 : brun 754 //______________________________________________________________________________
220 :     void TMultiGraph::Draw(Option_t *option)
221 :     {
222 : couet 13304 // Draw this multigraph with its current attributes.
223 :     //
224 :     // Options to draw a graph are described in TGraph::PainGraph
225 :     //
226 :     // The drawing option for each TGraph may be specified as an optional
227 :     // second argument of the Add function. You can use GetGraphDrawOption
228 :     // to return this option.
229 :     // If a draw option is specified, it will be used to draw the graph,
230 :     // otherwise the graph will be drawn with the option specified in
231 :     // TMultiGraph::Draw. Use GetDrawOption to return the option specified
232 :     // when drawin the TMultiGraph.
233 : brun 754
234 : couet 13304 AppendPad(option);
235 : brun 754 }
236 :    
237 : couet 13304
238 : brun 754 //______________________________________________________________________________
239 : moneta 31491 TFitResultPtr TMultiGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
240 : brun 11226 {
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 31491 TFitResultPtr 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 31207 // internal multigraph fitting methods
405 :     Foption_t fitOption;
406 :     ROOT::Fit::FitOptionsMake(option,fitOption);
407 : rdm 11750
408 : moneta 31207 // create range and minimizer options with default values
409 :     ROOT::Fit::DataRange range(rxmin,rxmax);
410 :     ROOT::Math::MinimizerOptions minOption;
411 :     return ROOT::Fit::FitObject(this, f1 , fitOption , minOption, goption, range);
412 :    
413 : moneta 25487 }
414 : rdm 11750
415 : moneta 25487 //______________________________________________________________________________
416 :     void TMultiGraph::FitPanel()
417 :     {
418 :     // -*-*-*-*-*Display a panel with all histogram fit options*-*-*-*-*-*
419 :     // ==============================================
420 :     //
421 :     // See class TFitPanel for example
422 : rdm 11750
423 : moneta 25487 if (!gPad)
424 :     gROOT->MakeDefCanvas();
425 : brun 11226
426 : moneta 25487 if (!gPad) {
427 :     Error("FitPanel", "Unable to create a default canvas");
428 :     return;
429 : brun 12136 }
430 : brun 11226
431 : moneta 25487 // use plugin manager to create instance of TFitEditor
432 :     TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
433 :     if (handler && handler->LoadPlugin() != -1) {
434 :     if (handler->ExecPlugin(2, gPad, this) == 0)
435 :     Error("FitPanel", "Unable to crate the FitPanel");
436 : brun 11226 }
437 : couet 26029 else
438 : moneta 25487 Error("FitPanel", "Unable to find the FitPanel plug-in");
439 : brun 11226 }
440 :    
441 :     //______________________________________________________________________________
442 : brun 11584 Option_t *TMultiGraph::GetGraphDrawOption(const TGraph *gr) const
443 :     {
444 : couet 13304 // Return the draw option for the TGraph gr in this TMultiGraph
445 :     // The return option is the one specified when calling TMultiGraph::Add(gr,option).
446 : rdm 11750
447 : brun 11584 if (!fGraphs || !gr) return "";
448 :     TListIter next(fGraphs);
449 :     TObject *obj;
450 :     while ((obj = next())) {
451 : brun 11600 if (obj == (TObject*)gr) return next.GetOption();
452 : brun 11584 }
453 :     return "";
454 :     }
455 :    
456 : couet 13304
457 : brun 11584 //______________________________________________________________________________
458 : brun 11226 void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
459 :     {
460 : couet 13304 // Compute Initial values of parameters for a gaussian.
461 : brun 11226
462 :     Double_t allcha, sumx, sumx2, x, val, rms, mean;
463 :     Int_t bin;
464 :     const Double_t sqrtpi = 2.506628;
465 :    
466 : couet 13304 // Compute mean value and RMS of the graph in the given range
467 : brun 11226 Int_t np = 0;
468 :     allcha = sumx = sumx2 = 0;
469 :     TGraph *g;
470 :     TIter next(fGraphs);
471 :     Double_t *px, *py;
472 :     Int_t npp; //number of points in each graph
473 :     while ((g = (TGraph*) next())) {
474 :     px=g->GetX();
475 :     py=g->GetY();
476 :     npp=g->GetN();
477 :     for (bin=0; bin<npp; bin++){
478 : brun 12643 x=px[bin];
479 :     if (x<xmin || x>xmax) continue;
480 :     np++;
481 :     val=py[bin];
482 :     sumx+=val*x;
483 :     sumx2+=val*x*x;
484 :     allcha+=val;
485 : brun 11226 }
486 :     }
487 :     if (np == 0 || allcha == 0) return;
488 :     mean = sumx/allcha;
489 :     rms = TMath::Sqrt(sumx2/allcha - mean*mean);
490 :    
491 :     Double_t binwidx = TMath::Abs((xmax-xmin)/np);
492 :     if (rms == 0) rms = 1;
493 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
494 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
495 :     f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
496 :     f1->SetParameter(1,mean);
497 :     f1->SetParameter(2,rms);
498 :     f1->SetParLimits(2,0,10*rms);
499 :     }
500 :    
501 : couet 13304
502 : brun 11226 //______________________________________________________________________________
503 :     void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
504 :     {
505 : couet 13304 // Compute Initial values of parameters for an exponential.
506 : brun 11226
507 :     Double_t constant, slope;
508 :     Int_t ifail;
509 :    
510 :     LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
511 :    
512 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
513 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
514 :     f1->SetParameter(0,constant);
515 :     f1->SetParameter(1,slope);
516 :     }
517 :    
518 : couet 13304
519 : brun 11226 //______________________________________________________________________________
520 :     void TMultiGraph::InitPolynom(Double_t xmin, Double_t xmax)
521 :     {
522 : couet 13304 // Compute Initial values of parameters for a polynom.
523 : brun 11226
524 :     Double_t fitpar[25];
525 :    
526 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
527 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
528 :     Int_t npar = f1->GetNpar();
529 :    
530 :     LeastSquareFit(npar, fitpar, xmin, xmax);
531 :    
532 :     for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
533 :     }
534 :    
535 : couet 13304
536 : brun 11226 //______________________________________________________________________________
537 :     void TMultiGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
538 :     {
539 : couet 13304 // Least squares lpolynomial fitting without weights.
540 :     //
541 :     // m number of parameters
542 :     // a array of parameters
543 :     // first 1st point number to fit (default =0)
544 :     // last last point number to fit (default=fNpoints-1)
545 :     //
546 :     // based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
547 : brun 11226
548 : couet 13304 const Double_t zero = 0.;
549 :     const Double_t one = 1.;
550 :     const Int_t idim = 20;
551 : brun 11226
552 : couet 13304 Double_t b[400] /* was [20][20] */;
553 :     Int_t i, k, l, ifail, bin;
554 :     Double_t power;
555 :     Double_t da[20], xk, yk;
556 : brun 11226
557 :    
558 : couet 13304 //count the total number of points to fit
559 :     TGraph *g;
560 :     TIter next(fGraphs);
561 :     Double_t *px, *py;
562 :     Int_t n=0;
563 :     Int_t npp;
564 :     while ((g = (TGraph*) next())) {
565 :     px=g->GetX();
566 :     py=g->GetY();
567 :     npp=g->GetN();
568 :     for (bin=0; bin<npp; bin++){
569 :     xk=px[bin];
570 :     if (xk < xmin || xk > xmax) continue;
571 :     n++;
572 :     }
573 :     }
574 :     if (m <= 2) {
575 :     LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
576 :     return;
577 :     }
578 :     if (m > idim || m > n) return;
579 :     da[0] = zero;
580 :     for (l = 2; l <= m; ++l) {
581 :     b[l-1] = zero;
582 :     b[m + l*20 - 21] = zero;
583 :     da[l-1] = zero;
584 :     }
585 :     Int_t np = 0;
586 : brun 11226
587 : couet 13304 next.Reset();
588 :     while ((g = (TGraph*) next())) {
589 :     px=g->GetX();
590 :     py=g->GetY();
591 :     npp=g->GetN();
592 :    
593 :     for (k = 0; k <= npp; ++k) {
594 :     xk = px[k];
595 :     if (xk < xmin || xk > xmax) continue;
596 :     np++;
597 :     yk = py[k];
598 :     power = one;
599 :     da[0] += yk;
600 :     for (l = 2; l <= m; ++l) {
601 :     power *= xk;
602 : brun 12643 b[l-1] += power;
603 :     da[l-1] += power*yk;
604 : couet 13304 }
605 :     for (l = 2; l <= m; ++l) {
606 :     power *= xk;
607 :     b[m + l*20 - 21] += power;
608 :     }
609 :     }
610 :     }
611 :     b[0] = Double_t(np);
612 :     for (i = 3; i <= m; ++i) {
613 :     for (k = i; k <= m; ++k) {
614 :     b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
615 :     }
616 :     }
617 :     H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
618 : brun 11226
619 : couet 13304 if (ifail < 0) {
620 :     //a[0] = fY[0];
621 :     py=((TGraph *)fGraphs->First())->GetY();
622 :     a[0]=py[0];
623 :     for (i=1; i<m; ++i) a[i] = 0;
624 :     return;
625 :     }
626 :     for (i=0; i<m; ++i) a[i] = da[i];
627 : brun 11226 }
628 :    
629 : couet 13304
630 : brun 11226 //______________________________________________________________________________
631 :     void TMultiGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
632 :     {
633 : couet 13304 // Least square linear fit without weights.
634 :     //
635 :     // Fit a straight line (a0 + a1*x) to the data in this graph.
636 :     // ndata: number of points to fit
637 :     // first: first point number to fit
638 :     // last: last point to fit O(ndata should be last-first
639 :     // ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
640 :     //
641 :     // extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
642 : brun 11226
643 : couet 13304 Double_t xbar, ybar, x2bar;
644 :     Int_t i;
645 :     Double_t xybar;
646 :     Double_t fn, xk, yk;
647 :     Double_t det;
648 : brun 11226
649 : couet 13304 ifail = -2;
650 :     xbar = ybar = x2bar = xybar = 0;
651 :     Int_t np = 0;
652 :     TGraph *g;
653 :     TIter next(fGraphs);
654 :     Double_t *px, *py;
655 :     Int_t npp;
656 :     while ((g = (TGraph*) next())) {
657 :     px=g->GetX();
658 :     py=g->GetY();
659 :     npp=g->GetN();
660 :     for (i = 0; i < npp; ++i) {
661 :     xk = px[i];
662 :     if (xk < xmin || xk > xmax) continue;
663 :     np++;
664 :     yk = py[i];
665 :     if (ndata < 0) {
666 :     if (yk <= 0) yk = 1e-9;
667 :     yk = TMath::Log(yk);
668 :     }
669 :     xbar += xk;
670 :     ybar += yk;
671 :     x2bar += xk*xk;
672 :     xybar += xk*yk;
673 :     }
674 :     }
675 :     fn = Double_t(np);
676 :     det = fn*x2bar - xbar*xbar;
677 :     ifail = -1;
678 :     if (det <= 0) {
679 :     if (fn > 0) a0 = ybar/fn;
680 :     else a0 = 0;
681 :     a1 = 0;
682 :     return;
683 :     }
684 :     ifail = 0;
685 :     a0 = (x2bar*ybar - xbar*xybar) / det;
686 :     a1 = (fn*xybar - xbar*ybar) / det;
687 :     }
688 : brun 11226
689 :    
690 :     //______________________________________________________________________________
691 : brun 1205 TH1F *TMultiGraph::GetHistogram() const
692 : brun 754 {
693 : couet 13304 // Returns a pointer to the histogram used to draw the axis
694 :     // Takes into account the two following cases.
695 :     // 1- option 'A' was specified in TMultiGraph::Draw. Return fHistogram
696 :     // 2- user had called TPad::DrawFrame. return pointer to hframe histogram
697 : brun 754
698 :     if (fHistogram) return fHistogram;
699 :     if (!gPad) return 0;
700 :     gPad->Modified();
701 :     gPad->Update();
702 :     if (fHistogram) return fHistogram;
703 :     TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
704 :     return h1;
705 :     }
706 :    
707 : couet 13304
708 : brun 754 //______________________________________________________________________________
709 : brun 11226 TF1 *TMultiGraph::GetFunction(const char *name) const
710 :     {
711 : couet 13304 // Return pointer to function with name.
712 :     //
713 :     // Functions such as TGraph::Fit store the fitted function in the list of
714 :     // functions of this graph.
715 : brun 11226
716 :     if (!fFunctions) return 0;
717 :     return (TF1*)fFunctions->FindObject(name);
718 :     }
719 :    
720 : moneta 25487 //______________________________________________________________________________
721 : couet 26029 TList *TMultiGraph::GetListOfFunctions()
722 : moneta 25487 {
723 :     // Return pointer to list of functions
724 :     // if pointer is null create the list
725 : couet 13304
726 : moneta 25487 if (!fFunctions) fFunctions = new TList();
727 : couet 26029 return fFunctions;
728 : moneta 25487 }
729 :    
730 :    
731 : brun 11226 //______________________________________________________________________________
732 : brun 1205 TAxis *TMultiGraph::GetXaxis() const
733 : brun 754 {
734 :     // Get x axis of the graph.
735 :    
736 :     if (!gPad) return 0;
737 : brun 4083 TH1 *h = GetHistogram();
738 :     if (!h) return 0;
739 :     return h->GetXaxis();
740 : brun 754 }
741 :    
742 : couet 13304
743 : brun 754 //______________________________________________________________________________
744 : brun 1205 TAxis *TMultiGraph::GetYaxis() const
745 : brun 754 {
746 :     // Get y axis of the graph.
747 :    
748 :     if (!gPad) return 0;
749 : brun 4083 TH1 *h = GetHistogram();
750 :     if (!h) return 0;
751 :     return h->GetYaxis();
752 : brun 754 }
753 :    
754 : couet 13304
755 : brun 754 //______________________________________________________________________________
756 :     void TMultiGraph::Paint(Option_t *option)
757 :     {
758 : couet 13304 // paint all the graphs of this multigraph
759 : brun 754
760 : couet 33513 if (!fGraphs) return;
761 : couet 13304 if (fGraphs->GetSize() == 0) return;
762 : rdm 11750
763 : couet 13304 char *l;
764 :     static char chopt[33];
765 :     Int_t nch = strlen(option);
766 :     Int_t i;
767 :     for (i=0;i<nch;i++) chopt[i] = toupper(option[i]);
768 :     chopt[nch] = 0;
769 :     TGraph *g;
770 : rdm 3742
771 : couet 13304 l = strstr(chopt,"A");
772 :     if (l) {
773 :     *l = ' ';
774 :     TIter next(fGraphs);
775 :     Int_t npt = 100;
776 :     Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
777 :     rwxmin = gPad->GetUxmin();
778 :     rwxmax = gPad->GetUxmax();
779 :     rwymin = gPad->GetUymin();
780 :     rwymax = gPad->GetUymax();
781 :     char *xtitle = 0;
782 :     char *ytitle = 0;
783 :     Int_t firstx = 0;
784 :     Int_t lastx = 0;
785 : rdm 11750
786 : couet 13304 if (fHistogram) {
787 :     //cleanup in case of a previous unzoom
788 :     if (fHistogram->GetMinimum() >= fHistogram->GetMaximum()) {
789 : brun 23506 nch = strlen(fHistogram->GetXaxis()->GetTitle());
790 : couet 13304 firstx = fHistogram->GetXaxis()->GetFirst();
791 :     lastx = fHistogram->GetXaxis()->GetLast();
792 :     if (nch) {
793 :     xtitle = new char[nch+1];
794 :     strcpy(xtitle,fHistogram->GetXaxis()->GetTitle());
795 :     }
796 :     nch = strlen(fHistogram->GetYaxis()->GetTitle());
797 :     if (nch) {
798 :     ytitle = new char[nch+1];
799 :     strcpy(ytitle,fHistogram->GetYaxis()->GetTitle());
800 :     }
801 :     delete fHistogram;
802 :     fHistogram = 0;
803 :     }
804 :     }
805 :     if (fHistogram) {
806 :     minimum = fHistogram->GetYaxis()->GetXmin();
807 :     maximum = fHistogram->GetYaxis()->GetXmax();
808 :     uxmin = gPad->PadtoX(rwxmin);
809 :     uxmax = gPad->PadtoX(rwxmax);
810 :     } else {
811 : couet 21573 g = (TGraph*) next();
812 : brun 32499 if (g) g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
813 : couet 13304 while ((g = (TGraph*) next())) {
814 : couet 20433 Double_t rx1,ry1,rx2,ry2;
815 :     g->ComputeRange(rx1, ry1, rx2, ry2);
816 :     if (rx1 < rwxmin) rwxmin = rx1;
817 :     if (ry1 < rwymin) rwymin = ry1;
818 :     if (rx2 > rwxmax) rwxmax = rx2;
819 :     if (ry2 > rwymax) rwymax = ry2;
820 : couet 13304 if (g->GetN() > npt) npt = g->GetN();
821 :     }
822 :     if (rwxmin == rwxmax) rwxmax += 1.;
823 :     if (rwymin == rwymax) rwymax += 1.;
824 :     dx = 0.05*(rwxmax-rwxmin);
825 :     dy = 0.05*(rwymax-rwymin);
826 :     uxmin = rwxmin - dx;
827 :     uxmax = rwxmax + dx;
828 :     if (gPad->GetLogy()) {
829 :     if (rwymin <= 0) rwymin = 0.001*rwymax;
830 :     minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
831 :     maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
832 :     } else {
833 :     minimum = rwymin - dy;
834 :     maximum = rwymax + dy;
835 :     }
836 :     if (minimum < 0 && rwymin >= 0) minimum = 0;
837 :     if (maximum > 0 && rwymax <= 0) maximum = 0;
838 :     }
839 : brun 754
840 : couet 13304 if (fMinimum != -1111) rwymin = minimum = fMinimum;
841 :     if (fMaximum != -1111) rwymax = maximum = fMaximum;
842 :     if (uxmin < 0 && rwxmin >= 0) {
843 :     if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
844 :     //else uxmin = 0;
845 :     }
846 :     if (uxmax > 0 && rwxmax <= 0) {
847 :     if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
848 :     //else uxmax = 0;
849 :     }
850 :     if (minimum < 0 && rwymin >= 0) {
851 :     if(gPad->GetLogy()) minimum = 0.9*rwymin;
852 :     //else minimum = 0;
853 :     }
854 :     if (maximum > 0 && rwymax <= 0) {
855 :     if(gPad->GetLogy()) maximum = 1.1*rwymax;
856 :     //else maximum = 0;
857 :     }
858 :     if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
859 :     if (uxmin <= 0 && gPad->GetLogx()) {
860 :     if (uxmax > 1000) uxmin = 1;
861 :     else uxmin = 0.001*uxmax;
862 :     }
863 :     rwymin = minimum;
864 :     rwymax = maximum;
865 :     if (fHistogram) {
866 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
867 :     }
868 : brun 3575
869 : couet 13304 // Create a temporary histogram to draw the axis
870 :     if (!fHistogram) {
871 :     // the graph is created with at least as many channels as there are points
872 :     // to permit zooming on the full range
873 :     rwxmin = uxmin;
874 :     rwxmax = uxmax;
875 :     fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
876 :     if (!fHistogram) return;
877 :     fHistogram->SetMinimum(rwymin);
878 :     fHistogram->SetBit(TH1::kNoStats);
879 :     fHistogram->SetMaximum(rwymax);
880 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
881 :     fHistogram->SetDirectory(0);
882 :     if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
883 :     if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
884 :     if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
885 :     }
886 :     fHistogram->Paint("0");
887 : brun 3575 }
888 : rdm 3742
889 : couet 21962 TGraph *gfit = 0;
890 : brun 754 if (fGraphs) {
891 : brun 4037 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
892 : couet 21962 TObject *obj = 0;
893 : brun 4037
894 :     while (lnk) {
895 :     obj = lnk->GetObject();
896 :     if (strlen(lnk->GetOption())) obj->Paint(lnk->GetOption());
897 :     else obj->Paint(chopt);
898 :     lnk = (TObjOptLink*)lnk->Next();
899 :     }
900 : couet 21962 gfit = (TGraph*)obj; // pick one TGraph in the list to paint the fit parameters.
901 : brun 754 }
902 : brun 11226
903 :     TObject *f;
904 : couet 21962 TF1 *fit = 0;
905 : brun 11226 if (fFunctions) {
906 : couet 13304 TIter next(fFunctions);
907 :     while ((f = (TObject*) next())) {
908 :     if (f->InheritsFrom(TF1::Class())) {
909 :     if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
910 : couet 21962 fit = (TF1*)f;
911 : couet 13304 } else {
912 :     f->Paint();
913 :     }
914 : brun 11226 }
915 :     }
916 : couet 21962
917 : brun 24084 if (fit) gfit->PaintStats(fit);
918 : couet 13304 }
919 : brun 11226
920 :    
921 : brun 754 //______________________________________________________________________________
922 : brun 1205 void TMultiGraph::Print(Option_t *option) const
923 : brun 754 {
924 : couet 13304 // Print the list of graphs
925 : brun 754
926 :     TGraph *g;
927 :     if (fGraphs) {
928 : couet 13304 TIter next(fGraphs);
929 :     while ((g = (TGraph*) next())) {
930 :     g->Print(option);
931 :     }
932 : brun 754 }
933 :     }
934 :    
935 : couet 13304
936 : brun 754 //______________________________________________________________________________
937 : brun 5552 void TMultiGraph::RecursiveRemove(TObject *obj)
938 :     {
939 :     // Recursively remove this object from a list. Typically implemented
940 :     // by classes that can contain mulitple references to a same object.
941 :    
942 :     if (!fGraphs) return;
943 :     TObject *objr = fGraphs->Remove(obj);
944 :     if (!objr) return;
945 :     delete fHistogram; fHistogram = 0;
946 :     if (gPad) gPad->Modified();
947 :     }
948 :    
949 : couet 13304
950 : brun 5552 //______________________________________________________________________________
951 : brun 15672 void TMultiGraph::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
952 : brun 754 {
953 : couet 13304 // Save primitive as a C++ statement(s) on output stream out
954 : brun 754
955 :     char quote = '"';
956 :     out<<" "<<endl;
957 :     if (gROOT->ClassSaved(TMultiGraph::Class())) {
958 : couet 13304 out<<" ";
959 : brun 754 } else {
960 : couet 13304 out<<" TMultiGraph *";
961 : brun 754 }
962 :     out<<"multigraph = new TMultiGraph();"<<endl;
963 :     out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<endl;
964 :     out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
965 :    
966 :     if (fGraphs) {
967 : brun 10847 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
968 :     TObject *g;
969 :    
970 :     while (lnk) {
971 :     g = lnk->GetObject();
972 : couet 20297 g->SavePrimitive(out, Form("multigraph%s",lnk->GetOption()));
973 : brun 10847 lnk = (TObjOptLink*)lnk->Next();
974 : couet 13304 }
975 : brun 754 }
976 : couet 20297 out<<" multigraph->Draw(" <<quote<<option<<quote<<");"<<endl;
977 : couet 14017
978 :     TAxis *xaxis = GetXaxis();
979 :     TAxis *yaxis = GetYaxis();
980 : couet 20433
981 :     if (xaxis) xaxis->SaveAttributes(out, "multigraph","->GetXaxis()");
982 :     if (yaxis) yaxis->SaveAttributes(out, "multigraph","->GetYaxis()");
983 : brun 754 }
984 :    
985 : couet 13304
986 : brun 754 //______________________________________________________________________________
987 :     void TMultiGraph::SetMaximum(Double_t maximum)
988 :     {
989 : couet 13304 // Set multigraph maximum.
990 :    
991 : brun 754 fMaximum = maximum;
992 :     if (fHistogram) fHistogram->SetMaximum(maximum);
993 :     }
994 :    
995 : couet 13304
996 : brun 754 //______________________________________________________________________________
997 :     void TMultiGraph::SetMinimum(Double_t minimum)
998 :     {
999 : couet 13304 // Set multigraph minimum.
1000 :    
1001 : brun 754 fMinimum = minimum;
1002 :     if (fHistogram) fHistogram->SetMinimum(minimum);
1003 :     }

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9