[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 13304 - (view) (download) (as text)
Original Path: trunk/graf/src/TMultiGraph.cxx

1 : couet 13304 // @(#)root/graf:$Name: $:$Id: TMultiGraph.cxx,v 1.25 2005/09/05 07:25:22 brun Exp $
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 : brun 754
20 : brun 11226
21 : brun 754 #include <ctype.h>
22 :    
23 : brun 11226 extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
24 : brun 754
25 :     ClassImp(TMultiGraph)
26 :    
27 : couet 13304
28 : brun 754 //______________________________________________________________________________
29 :     //
30 :     // A TMultiGraph is a collection of TGraph (or derived) objects
31 :     // Use TMultiGraph::Add to add a new graph to the list.
32 :     // The TMultiGraph owns the objects in the list.
33 :     // Drawing options are the same as for TGraph
34 :     // Example;
35 :     // TGraph *gr1 = new TGraph(...
36 :     // TGraphErrors *gr2 = new TGraphErrors(...
37 :     // TMultiGraph *mg = new TMultiGraph();
38 : brun 4037 // mg->Add(gr1,"lp");
39 :     // mg->Add(gr2,"cp");
40 :     // mg->Draw("a");
41 :     //
42 :     // The drawing option for each TGraph may be specified as an optional
43 :     // second argument of the Add function.
44 :     // If a draw option is specified, it will be used to draw the graph,
45 :     // otherwise the graph will be drawn with the option specified in
46 :     // TMultiGraph::Draw
47 : brun 754
48 : couet 13304
49 : brun 754 //______________________________________________________________________________
50 :     TMultiGraph::TMultiGraph(): TNamed()
51 :     {
52 : couet 13304 // TMultiGraph default constructor
53 : brun 754
54 :     fGraphs = 0;
55 : brun 11226 fFunctions = 0;
56 : brun 754 fHistogram = 0;
57 :     fMaximum = -1111;
58 :     fMinimum = -1111;
59 :     }
60 :    
61 : couet 13304
62 : brun 754 //______________________________________________________________________________
63 :     TMultiGraph::TMultiGraph(const char *name, const char *title)
64 :     : TNamed(name,title)
65 :     {
66 : couet 13304 // constructor with name and title
67 :    
68 : brun 754 fGraphs = 0;
69 : brun 11226 fFunctions = 0;
70 : brun 754 fHistogram = 0;
71 :     fMaximum = -1111;
72 :     fMinimum = -1111;
73 :     }
74 :    
75 : couet 13304
76 : brun 754 //______________________________________________________________________________
77 :     TMultiGraph::~TMultiGraph()
78 :     {
79 : couet 13304 // TMultiGraph destructor
80 : brun 754
81 :     if (!fGraphs) return;
82 : brun 5552 TGraph *g;
83 :     TIter next(fGraphs);
84 :     while ((g = (TGraph*) next())) {
85 : couet 13304 g->ResetBit(kMustCleanup);
86 : brun 5552 }
87 : brun 754 fGraphs->Delete();
88 :     delete fGraphs;
89 :     fGraphs = 0;
90 :     delete fHistogram;
91 :     fHistogram = 0;
92 : brun 11226 if (fFunctions) {
93 :     fFunctions->SetBit(kInvalidObject);
94 :     //special logic to support the case where the same object is
95 :     //added multiple times in fFunctions.
96 :     //This case happens when the same object is added with different
97 :     //drawing modes
98 :     TObject *obj;
99 :     while ((obj = fFunctions->First())) {
100 :     while(fFunctions->Remove(obj));
101 :     delete obj;
102 :     }
103 :     delete fFunctions;
104 :     }
105 : brun 754 }
106 :    
107 : couet 13304
108 : brun 754 //______________________________________________________________________________
109 : brun 4037 void TMultiGraph::Add(TGraph *graph, Option_t *chopt)
110 : brun 754 {
111 :     // add a new graph to the list of graphs
112 : brun 6639 // note that the graph is now owned by the TMultigraph.
113 :     // Deleting the TMultiGraph object will automatically delete the graphs.
114 :     // You should not delete the graphs when the TMultigraph is still active.
115 : rdm 3742
116 : brun 754 if (!fGraphs) fGraphs = new TList();
117 : brun 5552 graph->SetBit(kMustCleanup);
118 : brun 4037 fGraphs->Add(graph,chopt);
119 : brun 754 }
120 :    
121 : couet 13304
122 : brun 754 //______________________________________________________________________________
123 :     void TMultiGraph::Browse(TBrowser *)
124 :     {
125 : couet 13304 // Browse multigraph.
126 :    
127 :     Draw("alp");
128 :     gPad->Update();
129 : brun 754 }
130 :    
131 : couet 13304
132 : brun 754 //______________________________________________________________________________
133 :     Int_t TMultiGraph::DistancetoPrimitive(Int_t px, Int_t py)
134 :     {
135 : couet 13304 // Compute distance from point px,py to each graph
136 : brun 754
137 : couet 13304 // Are we on the axis?
138 : brun 754 const Int_t kMaxDiff = 10;
139 :     Int_t distance = 9999;
140 :     if (fHistogram) {
141 :     distance = fHistogram->DistancetoPrimitive(px,py);
142 :     if (distance <= 0) return distance;
143 :     }
144 :    
145 : couet 13304 // Loop on the list of graphs
146 : brun 754 if (!fGraphs) return distance;
147 :     TGraph *g;
148 :     TIter next(fGraphs);
149 :     while ((g = (TGraph*) next())) {
150 :     Int_t dist = g->DistancetoPrimitive(px,py);
151 : brun 4915 if (dist <= 0) return 0;
152 : brun 754 if (dist < kMaxDiff) {gPad->SetSelected(g); return dist;}
153 :     }
154 :     return distance;
155 :     }
156 :    
157 : couet 13304
158 : brun 754 //______________________________________________________________________________
159 :     void TMultiGraph::Draw(Option_t *option)
160 :     {
161 : couet 13304 // Draw this multigraph with its current attributes.
162 :     //
163 :     // Options to draw a graph are described in TGraph::PainGraph
164 :     //
165 :     // The drawing option for each TGraph may be specified as an optional
166 :     // second argument of the Add function. You can use GetGraphDrawOption
167 :     // to return this option.
168 :     // If a draw option is specified, it will be used to draw the graph,
169 :     // otherwise the graph will be drawn with the option specified in
170 :     // TMultiGraph::Draw. Use GetDrawOption to return the option specified
171 :     // when drawin the TMultiGraph.
172 : brun 754
173 : couet 13304 AppendPad(option);
174 : brun 754 }
175 :    
176 : couet 13304
177 : brun 754 //______________________________________________________________________________
178 : brun 11226 Int_t TMultiGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
179 :     {
180 : couet 13304 // Fit this graph with function with name fname.
181 :     //
182 :     // interface to TF1::Fit(TF1 *f1...
183 : brun 11226
184 :     char *linear;
185 : brun 11231 linear= (char*)strstr(fname, "++");
186 : brun 11226 TF1 *f1=0;
187 :     if (linear)
188 : brun 12643 f1=new TF1(fname, fname, xmin, xmax);
189 : brun 11226 else {
190 :     f1 = (TF1*)gROOT->GetFunction(fname);
191 :     if (!f1) { Printf("Unknown function: %s",fname); return -1; }
192 :     }
193 :    
194 :     return Fit(f1,option,"",xmin,xmax);
195 :     }
196 :    
197 : couet 13304
198 : brun 11226 //______________________________________________________________________________
199 :     Int_t TMultiGraph::Fit(TF1 *f1, Option_t *option, Option_t *, Axis_t rxmin, Axis_t rxmax)
200 :     {
201 : couet 13304 // Fit this multigraph with function f1.
202 :     //
203 :     // In this function all graphs of the multigraph are fitted simultaneously
204 :     //
205 :     // f1 is an already predefined function created by TF1.
206 :     // Predefined functions such as gaus, expo and poln are automatically
207 :     // created by ROOT.
208 :     //
209 :     // The list of fit options is given in parameter option.
210 :     // option = "W" Set all errors to 1
211 :     // = "U" Use a User specified fitting algorithm (via SetFCN)
212 :     // = "Q" Quiet mode (minimum printing)
213 :     // = "V" Verbose mode (default is between Q and V)
214 :     // = "B" Use this option when you want to fix one or more parameters
215 :     // and the fitting function is like "gaus","expo","poln","landau".
216 :     // = "R" Use the Range specified in the function range
217 :     // = "N" Do not store the graphics function, do not draw
218 :     // = "0" Do not plot the result of the fit. By default the fitted function
219 :     // is drawn unless the option"N" above is specified.
220 :     // = "+" Add this new fitted function to the list of fitted functions
221 :     // (by default, any previous function is deleted)
222 :     // = "C" In case of linear fitting, not calculate the chisquare
223 :     // (saves time)
224 :     // = "F" If fitting a polN, switch to minuit fitter
225 :     // = "ROB" In case of linear fitting, compute the LTS regression
226 :     // coefficients (robust(resistant) regression), using
227 :     // the default fraction of good points
228 :     // "ROB=0.x" - compute the LTS regression coefficients, using
229 :     // 0.x as a fraction of good points
230 :     //
231 :     // When the fit is drawn (by default), the parameter goption may be used
232 :     // to specify a list of graphics options. See TGraph::Paint for a complete
233 :     // list of these options.
234 :     //
235 :     // In order to use the Range option, one must first create a function
236 :     // with the expression to be fitted. For example, if your graph
237 :     // has a defined range between -4 and 4 and you want to fit a gaussian
238 :     // only in the interval 1 to 3, you can do:
239 :     // TF1 *f1 = new TF1("f1","gaus",1,3);
240 :     // graph->Fit("f1","R");
241 :     //
242 :     //
243 :     // who is calling this function
244 :     // ============================
245 :     // Note that this function is called when calling TGraphErrors::Fit
246 :     // or TGraphAsymmErrors::Fit ot TGraphBentErrors::Fit
247 :     // see the discussion below on the errors calulation.
248 :     //
249 :     // Setting initial conditions
250 :     // ==========================
251 :     // Parameters must be initialized before invoking the Fit function.
252 :     // The setting of the parameter initial values is automatic for the
253 :     // predefined functions : poln, expo, gaus, landau. One can however disable
254 :     // this automatic computation by specifying the option "B".
255 :     // You can specify boundary limits for some or all parameters via
256 :     // f1->SetParLimits(p_number, parmin, parmax);
257 :     // if parmin>=parmax, the parameter is fixed
258 :     // Note that you are not forced to fix the limits for all parameters.
259 :     // For example, if you fit a function with 6 parameters, you can do:
260 :     // func->SetParameters(0,3.1,1.e-6,0.1,-8,100);
261 :     // func->SetParLimits(4,-10,-4);
262 :     // func->SetParLimits(5, 1,1);
263 :     // With this setup, parameters 0->3 can vary freely
264 :     // Parameter 4 has boundaries [-10,-4] with initial value -8
265 :     // Parameter 5 is fixed to 100.
266 :     //
267 :     // Fit range
268 :     // =========
269 :     // The fit range can be specified in two ways:
270 :     // - specify rxmax > rxmin (default is rxmin=rxmax=0)
271 :     // - specify the option "R". In this case, the function will be taken
272 :     // instead of the full graph range.
273 :     //
274 :     // Changing the fitting function
275 :     // =============================
276 :     // By default the fitting function GraphFitChisquare is used.
277 :     // To specify a User defined fitting function, specify option "U" and
278 :     // call the following functions:
279 :     // TVirtualFitter::Fitter(mygraph)->SetFCN(MyFittingFunction)
280 :     // where MyFittingFunction is of type:
281 :     // extern void MyFittingFunction(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
282 :     //
283 :     // How errors are used in the chisquare function (see TFitter GraphFitChisquare)// Access to the fit results
284 :     // ============================================
285 :     // In case of a TGraphErrors object, ex, the error along x, is projected
286 :     // along the y-direction by calculating the function at the points x-exlow and
287 :     // x+exhigh.
288 :     //
289 :     // The chisquare is computed as the sum of the quantity below at each point:
290 :     //
291 :     // (y - f(x))**2
292 :     // -----------------------------------
293 :     // ey**2 + ((f(x+exhigh) - f(x-exlow))/2)**2
294 :     //
295 :     // where x and y are the point coordinates.
296 :     //
297 :     // In case the function lies below (above) the data point, ey is ey_low (ey_high).
298 :     //
299 :     // thanks to Andy Haas (haas@yahoo.com) for adding the case with TGraphasymmerrors
300 :     // University of Washington
301 :     //
302 :     // a little different approach to approximating the uncertainty in y because of the
303 :     // errors in x, is to make it equal the error in x times the slope of the line.
304 :     // The improvement, compared to the first method (f(x+ exhigh) - f(x-exlow))/2
305 :     // is of (error of x)**2 order. This approach is called "effective variance method".
306 :     // This improvement has been made in version 4.00/08 by Anna Kreshuk.
307 :     //
308 :     // Associated functions
309 :     // ====================
310 :     // One or more object (typically a TF1*) can be added to the list
311 :     // of functions (fFunctions) associated to each graph.
312 :     // When TGraph::Fit is invoked, the fitted function is added to this list.
313 :     // Given a graph gr, one can retrieve an associated function
314 :     // with: TF1 *myfunc = gr->GetFunction("myfunc");
315 :     //
316 :     // If the graph is made persistent, the list of
317 :     // associated functions is also persistent. Given a pointer (see above)
318 :     // to an associated function myfunc, one can retrieve the function/fit
319 :     // parameters with calls such as:
320 :     // Double_t chi2 = myfunc->GetChisquare();
321 :     // Double_t par0 = myfunc->GetParameter(0); //value of 1st parameter
322 :     // Double_t err0 = myfunc->GetParError(0); //error on first parameter
323 :     //
324 :     // Fit Statistics
325 :     // ==============
326 :     // You can change the statistics box to display the fit parameters with
327 :     // the TStyle::SetOptFit(mode) method. This mode has four digits.
328 :     // mode = pcev (default = 0111)
329 :     // v = 1; print name/values of parameters
330 :     // e = 1; print errors (if e=1, v must be 1)
331 :     // c = 1; print Chisquare/Number of degress of freedom
332 :     // p = 1; print Probability
333 :     //
334 :     // For example: gStyle->SetOptFit(1011);
335 :     // prints the fit probability, parameter names/values, and errors.
336 :     // You can change the position of the statistics box with these lines
337 :     // (where g is a pointer to the TGraph):
338 :     //
339 :     // Root > TPaveStats *st = (TPaveStats*)g->GetListOfFunctions()->FindObject("stats")
340 :     // Root > st->SetX1NDC(newx1); //new x start position
341 :     // Root > st->SetX2NDC(newx2); //new x end position
342 : brun 11226
343 :     Int_t fitResult = 0;
344 :     Double_t xmin, xmax, ymin, ymax;
345 :     Int_t i, npar,nvpar,nparx;
346 :     Double_t par, we, al, bl;
347 :     Double_t eplus,eminus,eparab,globcc,amin,edm,errdef,werr;
348 :     Int_t np;
349 :     TF1 *fnew1;
350 : rdm 11750
351 : brun 11226 // Check validity of function
352 :     if (!f1) {
353 :     Error("Fit", "function may not be null pointer");
354 :     return 0;
355 :     }
356 :     if (f1->IsZombie()) {
357 :     Error("Fit", "function is zombie");
358 :     return 0;
359 :     }
360 : rdm 11750
361 : brun 11226 npar = f1->GetNpar();
362 :     if (npar <= 0) {
363 :     Error("Fit", "function %s has illegal number of parameters = %d", f1->GetName(), npar);
364 :     return 0;
365 :     }
366 : rdm 11750
367 : brun 11226 // Check that function has same dimension as graph
368 :     if (f1->GetNdim() > 1) {
369 :     Error("Fit", "function %s is not 1-D", f1->GetName());
370 :     return 0;
371 :     }
372 :    
373 :     TGraph *g;
374 :     TIter next(fGraphs);
375 :    
376 :     Double_t *arglist = new Double_t[100];
377 :     // Decode string choptin and fill fitOption structure
378 :     Foption_t fitOption;
379 : rdm 11750
380 : brun 11226 TString opt = option;
381 :     opt.ToUpper();
382 : brun 12136 Double_t h=0;
383 :     opt.ReplaceAll("ROB", "H");
384 :     //for robust fitting, see if # of good points is defined
385 :     if (opt.Contains("H=0.")) {
386 :     int start = opt.Index("H=0.");
387 :     int numpos = start + strlen("H=0.");
388 :     int numlen = 0;
389 :     int len = opt.Length();
390 :     while( (numpos+numlen<len) && isdigit(opt[numpos+numlen]) ) numlen++;
391 :     TString num = opt(numpos,numlen);
392 :     opt.Remove(start+strlen("H"),strlen("=0.")+numlen);
393 :     h = atof(num.Data());
394 :     h*=TMath::Power(10, -numlen);
395 :     }
396 : brun 11226
397 :     if (opt.Contains("U")) fitOption.User = 1;
398 :     if (opt.Contains("Q")) fitOption.Quiet = 1;
399 :     if (opt.Contains("V")){fitOption.Verbose = 1; fitOption.Quiet = 0;}
400 :     if (opt.Contains("W")) fitOption.W1 = 1;
401 :     if (opt.Contains("E")) fitOption.Errors = 1;
402 :     if (opt.Contains("R")) fitOption.Range = 1;
403 :     if (opt.Contains("N")) fitOption.Nostore = 1;
404 :     if (opt.Contains("0")) fitOption.Nograph = 1;
405 :     if (opt.Contains("+")) fitOption.Plus = 1;
406 :     if (opt.Contains("B")) fitOption.Bound = 1;
407 :     if (opt.Contains("C")) fitOption.Nochisq = 1;
408 : brun 11725 if (opt.Contains("F")) fitOption.Minuit = 1;
409 : brun 12136 if (opt.Contains("H")) fitOption.Robust = 1;
410 : brun 11226
411 :     if (rxmax > rxmin) {
412 :     xmin = rxmin;
413 :     xmax = rxmax;
414 :     } else {
415 :     g=(TGraph *)fGraphs->First();
416 :     if (!g) {
417 : brun 12643 Error("Fit", "No graphs in the multigraph");
418 :     return 0;
419 : brun 11226 }
420 :     Double_t *px, *py;
421 :     np=g->GetN();
422 :     px=g->GetX();
423 :     py=g->GetY();
424 :     xmin=px[0];
425 :     xmax=py[np-1];
426 :     ymin=px[0];
427 :     ymax=py[np-1];
428 :     Double_t err0=g->GetErrorX(0);
429 :     Double_t errn=g->GetErrorX(np-1);
430 :     if (err0 > 0) xmin -= 2*err0;
431 :     if (errn > 0) xmax += 2*errn;
432 :    
433 :     next.Reset();
434 :     while ((g = (TGraph*) next())) {
435 : brun 12643 np=g->GetN();
436 :     px=g->GetX();
437 :     py=g->GetY();
438 :     for (i=0; i<np; i++) {
439 :     if (px[i] < xmin) xmin = px[i];
440 :     if (px[i] > xmax) xmax = px[i];
441 :     if (py[i] < ymin) ymin = py[i];
442 :     if (py[i] > ymax) ymax = py[i];
443 :     }
444 : brun 11226 }
445 :     }
446 :    
447 :     ///////////////
448 :     //set the fitter
449 :     //////////////
450 : brun 11605
451 : brun 11226 Int_t special=f1->GetNumber();
452 :     Bool_t linear = f1->IsLinear();
453 :     if (special==299+npar)
454 :     linear=kTRUE;
455 : brun 11605 if (fitOption.Bound || fitOption.User || fitOption.Errors || fitOption.Minuit)
456 :     linear = kFALSE;
457 : brun 11226
458 :     char l[]="TLinearFitter";
459 :     Int_t strdiff = 0;
460 : brun 12551 Bool_t isSet = kFALSE;
461 : brun 11226 if (TVirtualFitter::GetFitter()){
462 :     //Is a fitter already set? Is it linear?
463 : brun 12551 isSet = kTRUE;
464 : brun 11226 strdiff = strcmp(TVirtualFitter::GetFitter()->IsA()->GetName(), l);
465 :     }
466 :     if (linear){
467 :     TClass *cl = gROOT->GetClass("TLinearFitter");
468 : brun 12551 if (isSet && strdiff!=0) {
469 : brun 12643 delete TVirtualFitter::GetFitter();
470 :     isSet=kFALSE;
471 : brun 11226 }
472 : brun 12551 if (!isSet) {
473 : brun 12643 TVirtualFitter::SetFitter((TVirtualFitter *)cl->New());
474 : brun 11226 }
475 :     } else {
476 : brun 12551 if (isSet && strdiff==0){
477 : brun 12643 delete TVirtualFitter::GetFitter();
478 :     isSet=kFALSE;
479 : brun 11226 }
480 : brun 12551 if (!isSet)
481 : brun 12643 TVirtualFitter::SetFitter(0);
482 : brun 11226 }
483 :    
484 :     TVirtualFitter *grFitter = TVirtualFitter::Fitter(this, f1->GetNpar());
485 :     grFitter->Clear();
486 : rdm 11750
487 : couet 13304 // Get pointer to the function by searching in the list of functions in ROOT
488 : brun 11226 grFitter->SetUserFunc(f1);
489 :     grFitter->SetFitOption(fitOption);
490 :    
491 : couet 13304 // Is a Fit range specified?
492 : brun 11226 if (fitOption.Range) {
493 :     f1->GetRange(xmin, xmax);
494 :     } else {
495 :     f1->SetRange(xmin, xmax);
496 :     }
497 :    
498 : brun 11605 if (linear){
499 : brun 12136 if (fitOption.Robust)
500 :     grFitter->ExecuteCommand("FitMultiGraph", &h, 0);
501 :     else
502 :     grFitter->ExecuteCommand("FitMultiGraph", 0, 0);
503 : brun 11226 } else {
504 :    
505 :     //Int_t special = f1->GetNumber();
506 :     if (fitOption.Bound) special = 0;
507 :     if (special == 100) InitGaus(xmin,xmax);
508 :     else if (special == 400) InitGaus(xmin,xmax);
509 :     else if (special == 200) InitExpo(xmin,xmax);
510 :     else if (special == 299+npar) InitPolynom(xmin,xmax);
511 : rdm 11750
512 : brun 11226 //*-*- Some initialisations
513 :     if (!fitOption.Verbose) {
514 :     arglist[0] = -1;
515 :     grFitter->ExecuteCommand("SET PRINT", arglist,1);
516 :     arglist[0] = 0;
517 :     grFitter->ExecuteCommand("SET NOW", arglist,0);
518 :     }
519 : rdm 11750
520 :     /////////////////////////////////////////////////////////
521 : brun 11226 //*-*- Set error criterion for chisquare
522 :     arglist[0] = TVirtualFitter::GetErrorDef();
523 :     if (!fitOption.User) grFitter->SetFitMethod("MultiGraphFitChisquare");
524 :    
525 :    
526 :     fitResult = grFitter->ExecuteCommand("SET ERR",arglist,1);
527 :     if (fitResult != 0) {
528 : brun 12643 // Abnormal termination, MIGRAD might not have converged on a
529 :     // minimum.
530 :     if (!fitOption.Quiet) {
531 :     Warning("Fit","Abnormal termination of minimization.");
532 :     }
533 :     delete [] arglist;
534 :     return fitResult;
535 : brun 11226 }
536 : rdm 11750
537 : brun 11226 //*-*- Transfer names and initial values of parameters to Minuit
538 :     Int_t nfixed = 0;
539 :     for (i=0;i<npar;i++) {
540 : brun 12643 par = f1->GetParameter(i);
541 :     f1->GetParLimits(i,al,bl);
542 :     if (al*bl != 0 && al >= bl) {
543 :     al = bl = 0;
544 :     arglist[nfixed] = i+1;
545 :     nfixed++;
546 :     }
547 :     we = 0.3*TMath::Abs(par);
548 :     if (we <= TMath::Abs(par)*1e-6) we = 1;
549 :     grFitter->SetParameter(i,f1->GetParName(i),par,we,al,bl);
550 : brun 11226 }
551 :     if(nfixed > 0)grFitter->ExecuteCommand("FIX",arglist,nfixed); // Otto
552 : rdm 11750
553 : brun 11226 //*-*- Reset Print level
554 :     if (!fitOption.Quiet) {
555 : brun 12643 if (fitOption.Verbose) { arglist[0] = 2; grFitter->ExecuteCommand("SET PRINT", arglist,1); }
556 :     else { arglist[0] = 0; grFitter->ExecuteCommand("SET PRINT", arglist,1); }
557 : brun 11226 }
558 :     //*-*- Compute sum of squares of errors in the bin range
559 :     Bool_t hasErrors = kFALSE;
560 :     Double_t ex, ey, sumw2=0;
561 :     next.Reset();
562 :     while ((g = (TGraph*) next())) {
563 : brun 12643 np=g->GetN();
564 :     for (i=0; i<np; i++){
565 :     ex=g->GetErrorX(i);
566 :     ey=g->GetErrorY(i);
567 :     if (ex > 0 || ey > 0) hasErrors=kTRUE;
568 :     sumw2+=ey*ey;
569 :     }
570 : brun 11226 }
571 :    
572 :     //*-*- Perform minimization
573 :    
574 :     arglist[0] = TVirtualFitter::GetMaxIterations();
575 :     arglist[1] = sumw2*TVirtualFitter::GetPrecision();
576 :     grFitter->ExecuteCommand("MIGRAD",arglist,2);
577 :     if (fitOption.Errors) {
578 : brun 12643 grFitter->ExecuteCommand("HESSE",arglist,0);
579 :     grFitter->ExecuteCommand("MINOS",arglist,0);
580 : brun 11226 }
581 : rdm 11750
582 : brun 11226 grFitter->GetStats(amin,edm,errdef,nvpar,nparx);
583 :     f1->SetChisquare(amin);
584 :     Int_t ndf = f1->GetNumberFitPoints()-npar+nfixed;
585 :     f1->SetNDF(ndf);
586 :    
587 :     //*-*- Get return status
588 :     char parName[50];
589 :     for (i=0;i<npar;i++) {
590 : brun 12643 grFitter->GetParameter(i,parName, par,we,al,bl);
591 :     if (!fitOption.Errors) werr = we;
592 :     else {
593 :     grFitter->GetErrors(i,eplus,eminus,eparab,globcc);
594 :     if (eplus > 0 && eminus < 0) werr = 0.5*(eplus-eminus);
595 :     else werr = we;
596 :     }
597 :     if (!hasErrors && ndf > 1) werr *= TMath::Sqrt(amin/(ndf-1));
598 :     f1->SetParameter(i,par);
599 :     f1->SetParError(i,werr);
600 : brun 11226 }
601 :     }
602 : couet 13304 // Print final values of parameters.
603 : brun 11226 if (!fitOption.Quiet) {
604 :     if (fitOption.Errors) grFitter->PrintResults(4,amin);
605 :     else grFitter->PrintResults(3,amin);
606 :     }
607 :     delete [] arglist;
608 :    
609 : couet 13304 // Store fitted function in histogram functions list and draw
610 : brun 11226
611 :     if (!fitOption.Nostore) {
612 :     if (!fFunctions) fFunctions = new TList;
613 :     if (!fitOption.Plus) {
614 :     TIter next2(fFunctions, kIterBackward);
615 :     TObject *obj;
616 :     while ((obj = next2())) {
617 :     if (obj->InheritsFrom(TF1::Class())){
618 : brun 12643 obj = fFunctions->Remove(obj);
619 :     delete obj;
620 :     }
621 : brun 11226 }
622 :     }
623 :     fnew1 = new TF1();
624 :     f1->Copy(*fnew1);
625 :     fFunctions->Add(fnew1);
626 :     fnew1->SetParent(this);
627 :     fnew1->Save(xmin,xmax,0,0,0,0);
628 :     if (fitOption.Nograph) fnew1->SetBit(TF1::kNotDraw);
629 :     fnew1->SetBit(TFormula::kNotGlobal);
630 :    
631 :     if (TestBit(kCanDelete)) return fitResult;
632 :     if (gPad) gPad->Modified();
633 :     }
634 :    
635 :     return fitResult;
636 :     }
637 :    
638 :    
639 :     //______________________________________________________________________________
640 : brun 11584 Option_t *TMultiGraph::GetGraphDrawOption(const TGraph *gr) const
641 :     {
642 : couet 13304 // Return the draw option for the TGraph gr in this TMultiGraph
643 :     // The return option is the one specified when calling TMultiGraph::Add(gr,option).
644 : rdm 11750
645 : brun 11584 if (!fGraphs || !gr) return "";
646 :     TListIter next(fGraphs);
647 :     TObject *obj;
648 :     while ((obj = next())) {
649 : brun 11600 if (obj == (TObject*)gr) return next.GetOption();
650 : brun 11584 }
651 :     return "";
652 :     }
653 :    
654 : couet 13304
655 : brun 11584 //______________________________________________________________________________
656 : brun 11226 void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
657 :     {
658 : couet 13304 // Compute Initial values of parameters for a gaussian.
659 : brun 11226
660 :     Double_t allcha, sumx, sumx2, x, val, rms, mean;
661 :     Int_t bin;
662 :     const Double_t sqrtpi = 2.506628;
663 :    
664 : couet 13304 // Compute mean value and RMS of the graph in the given range
665 : brun 11226 Int_t np = 0;
666 :     allcha = sumx = sumx2 = 0;
667 :     TGraph *g;
668 :     TIter next(fGraphs);
669 :     Double_t *px, *py;
670 :     Int_t npp; //number of points in each graph
671 :     while ((g = (TGraph*) next())) {
672 :     px=g->GetX();
673 :     py=g->GetY();
674 :     npp=g->GetN();
675 :     for (bin=0; bin<npp; bin++){
676 : brun 12643 x=px[bin];
677 :     if (x<xmin || x>xmax) continue;
678 :     np++;
679 :     val=py[bin];
680 :     sumx+=val*x;
681 :     sumx2+=val*x*x;
682 :     allcha+=val;
683 : brun 11226 }
684 :     }
685 :     if (np == 0 || allcha == 0) return;
686 :     mean = sumx/allcha;
687 :     rms = TMath::Sqrt(sumx2/allcha - mean*mean);
688 :    
689 :     Double_t binwidx = TMath::Abs((xmax-xmin)/np);
690 :     if (rms == 0) rms = 1;
691 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
692 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
693 :     f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
694 :     f1->SetParameter(1,mean);
695 :     f1->SetParameter(2,rms);
696 :     f1->SetParLimits(2,0,10*rms);
697 :     }
698 :    
699 : couet 13304
700 : brun 11226 //______________________________________________________________________________
701 :     void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
702 :     {
703 : couet 13304 // Compute Initial values of parameters for an exponential.
704 : brun 11226
705 :     Double_t constant, slope;
706 :     Int_t ifail;
707 :    
708 :     LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
709 :    
710 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
711 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
712 :     f1->SetParameter(0,constant);
713 :     f1->SetParameter(1,slope);
714 :     }
715 :    
716 : couet 13304
717 : brun 11226 //______________________________________________________________________________
718 :     void TMultiGraph::InitPolynom(Double_t xmin, Double_t xmax)
719 :     {
720 : couet 13304 // Compute Initial values of parameters for a polynom.
721 : brun 11226
722 :     Double_t fitpar[25];
723 :    
724 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
725 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
726 :     Int_t npar = f1->GetNpar();
727 :    
728 :     LeastSquareFit(npar, fitpar, xmin, xmax);
729 :    
730 :     for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
731 :     }
732 :    
733 : couet 13304
734 : brun 11226 //______________________________________________________________________________
735 :     void TMultiGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
736 :     {
737 : couet 13304 // Least squares lpolynomial fitting without weights.
738 :     //
739 :     // m number of parameters
740 :     // a array of parameters
741 :     // first 1st point number to fit (default =0)
742 :     // last last point number to fit (default=fNpoints-1)
743 :     //
744 :     // based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
745 : brun 11226
746 : couet 13304 const Double_t zero = 0.;
747 :     const Double_t one = 1.;
748 :     const Int_t idim = 20;
749 : brun 11226
750 : couet 13304 Double_t b[400] /* was [20][20] */;
751 :     Int_t i, k, l, ifail, bin;
752 :     Double_t power;
753 :     Double_t da[20], xk, yk;
754 : brun 11226
755 :    
756 : couet 13304 //count the total number of points to fit
757 :     TGraph *g;
758 :     TIter next(fGraphs);
759 :     Double_t *px, *py;
760 :     Int_t n=0;
761 :     Int_t npp;
762 :     while ((g = (TGraph*) next())) {
763 :     px=g->GetX();
764 :     py=g->GetY();
765 :     npp=g->GetN();
766 :     for (bin=0; bin<npp; bin++){
767 :     xk=px[bin];
768 :     if (xk < xmin || xk > xmax) continue;
769 :     n++;
770 :     }
771 :     }
772 :     if (m <= 2) {
773 :     LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
774 :     return;
775 :     }
776 :     if (m > idim || m > n) return;
777 :     da[0] = zero;
778 :     for (l = 2; l <= m; ++l) {
779 :     b[l-1] = zero;
780 :     b[m + l*20 - 21] = zero;
781 :     da[l-1] = zero;
782 :     }
783 :     Int_t np = 0;
784 : brun 11226
785 : couet 13304 next.Reset();
786 :     while ((g = (TGraph*) next())) {
787 :     px=g->GetX();
788 :     py=g->GetY();
789 :     npp=g->GetN();
790 :    
791 :     for (k = 0; k <= npp; ++k) {
792 :     xk = px[k];
793 :     if (xk < xmin || xk > xmax) continue;
794 :     np++;
795 :     yk = py[k];
796 :     power = one;
797 :     da[0] += yk;
798 :     for (l = 2; l <= m; ++l) {
799 :     power *= xk;
800 : brun 12643 b[l-1] += power;
801 :     da[l-1] += power*yk;
802 : couet 13304 }
803 :     for (l = 2; l <= m; ++l) {
804 :     power *= xk;
805 :     b[m + l*20 - 21] += power;
806 :     }
807 :     }
808 :     }
809 :     b[0] = Double_t(np);
810 :     for (i = 3; i <= m; ++i) {
811 :     for (k = i; k <= m; ++k) {
812 :     b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
813 :     }
814 :     }
815 :     H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
816 : brun 11226
817 : couet 13304 if (ifail < 0) {
818 :     //a[0] = fY[0];
819 :     py=((TGraph *)fGraphs->First())->GetY();
820 :     a[0]=py[0];
821 :     for (i=1; i<m; ++i) a[i] = 0;
822 :     return;
823 :     }
824 :     for (i=0; i<m; ++i) a[i] = da[i];
825 : brun 11226 }
826 :    
827 : couet 13304
828 : brun 11226 //______________________________________________________________________________
829 :     void TMultiGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
830 :     {
831 : couet 13304 // Least square linear fit without weights.
832 :     //
833 :     // Fit a straight line (a0 + a1*x) to the data in this graph.
834 :     // ndata: number of points to fit
835 :     // first: first point number to fit
836 :     // last: last point to fit O(ndata should be last-first
837 :     // ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
838 :     //
839 :     // extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
840 : brun 11226
841 : couet 13304 Double_t xbar, ybar, x2bar;
842 :     Int_t i;
843 :     Double_t xybar;
844 :     Double_t fn, xk, yk;
845 :     Double_t det;
846 : brun 11226
847 : couet 13304 ifail = -2;
848 :     xbar = ybar = x2bar = xybar = 0;
849 :     Int_t np = 0;
850 :     TGraph *g;
851 :     TIter next(fGraphs);
852 :     Double_t *px, *py;
853 :     Int_t npp;
854 :     while ((g = (TGraph*) next())) {
855 :     px=g->GetX();
856 :     py=g->GetY();
857 :     npp=g->GetN();
858 :     for (i = 0; i < npp; ++i) {
859 :     xk = px[i];
860 :     if (xk < xmin || xk > xmax) continue;
861 :     np++;
862 :     yk = py[i];
863 :     if (ndata < 0) {
864 :     if (yk <= 0) yk = 1e-9;
865 :     yk = TMath::Log(yk);
866 :     }
867 :     xbar += xk;
868 :     ybar += yk;
869 :     x2bar += xk*xk;
870 :     xybar += xk*yk;
871 :     }
872 :     }
873 :     fn = Double_t(np);
874 :     det = fn*x2bar - xbar*xbar;
875 :     ifail = -1;
876 :     if (det <= 0) {
877 :     if (fn > 0) a0 = ybar/fn;
878 :     else a0 = 0;
879 :     a1 = 0;
880 :     return;
881 :     }
882 :     ifail = 0;
883 :     a0 = (x2bar*ybar - xbar*xybar) / det;
884 :     a1 = (fn*xybar - xbar*ybar) / det;
885 :     }
886 : brun 11226
887 :    
888 :     //______________________________________________________________________________
889 : brun 1205 TH1F *TMultiGraph::GetHistogram() const
890 : brun 754 {
891 : couet 13304 // Returns a pointer to the histogram used to draw the axis
892 :     // Takes into account the two following cases.
893 :     // 1- option 'A' was specified in TMultiGraph::Draw. Return fHistogram
894 :     // 2- user had called TPad::DrawFrame. return pointer to hframe histogram
895 : brun 754
896 :     if (fHistogram) return fHistogram;
897 :     if (!gPad) return 0;
898 :     gPad->Modified();
899 :     gPad->Update();
900 :     if (fHistogram) return fHistogram;
901 :     TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
902 :     return h1;
903 :     }
904 :    
905 : couet 13304
906 : brun 754 //______________________________________________________________________________
907 : brun 11226 TF1 *TMultiGraph::GetFunction(const char *name) const
908 :     {
909 : couet 13304 // Return pointer to function with name.
910 :     //
911 :     // Functions such as TGraph::Fit store the fitted function in the list of
912 :     // functions of this graph.
913 : brun 11226
914 :     if (!fFunctions) return 0;
915 :     return (TF1*)fFunctions->FindObject(name);
916 :     }
917 :    
918 : couet 13304
919 : brun 11226 //______________________________________________________________________________
920 : brun 1205 TAxis *TMultiGraph::GetXaxis() const
921 : brun 754 {
922 :     // Get x axis of the graph.
923 :    
924 :     if (!gPad) return 0;
925 : brun 4083 TH1 *h = GetHistogram();
926 :     if (!h) return 0;
927 :     return h->GetXaxis();
928 : brun 754 }
929 :    
930 : couet 13304
931 : brun 754 //______________________________________________________________________________
932 : brun 1205 TAxis *TMultiGraph::GetYaxis() const
933 : brun 754 {
934 :     // Get y axis of the graph.
935 :    
936 :     if (!gPad) return 0;
937 : brun 4083 TH1 *h = GetHistogram();
938 :     if (!h) return 0;
939 :     return h->GetYaxis();
940 : brun 754 }
941 :    
942 : couet 13304
943 : brun 754 //______________________________________________________________________________
944 :     void TMultiGraph::Paint(Option_t *option)
945 :     {
946 : couet 13304 // paint all the graphs of this multigraph
947 : brun 754
948 : couet 13304 if (fGraphs->GetSize() == 0) return;
949 : rdm 11750
950 : couet 13304 char *l;
951 :     static char chopt[33];
952 :     Int_t nch = strlen(option);
953 :     Int_t i;
954 :     for (i=0;i<nch;i++) chopt[i] = toupper(option[i]);
955 :     chopt[nch] = 0;
956 :     Double_t *x, *y;
957 :     TGraph *g;
958 : rdm 3742
959 : couet 13304 l = strstr(chopt,"A");
960 :     if (l) {
961 :     *l = ' ';
962 :     TIter next(fGraphs);
963 :     Int_t npt = 100;
964 :     Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
965 :     rwxmin = gPad->GetUxmin();
966 :     rwxmax = gPad->GetUxmax();
967 :     rwymin = gPad->GetUymin();
968 :     rwymax = gPad->GetUymax();
969 :     char *xtitle = 0;
970 :     char *ytitle = 0;
971 :     Int_t firstx = 0;
972 :     Int_t lastx = 0;
973 : rdm 11750
974 : couet 13304 if (fHistogram) {
975 :     //cleanup in case of a previous unzoom
976 :     if (fHistogram->GetMinimum() >= fHistogram->GetMaximum()) {
977 :     Int_t nch = strlen(fHistogram->GetXaxis()->GetTitle());
978 :     firstx = fHistogram->GetXaxis()->GetFirst();
979 :     lastx = fHistogram->GetXaxis()->GetLast();
980 :     if (nch) {
981 :     xtitle = new char[nch+1];
982 :     strcpy(xtitle,fHistogram->GetXaxis()->GetTitle());
983 :     }
984 :     nch = strlen(fHistogram->GetYaxis()->GetTitle());
985 :     if (nch) {
986 :     ytitle = new char[nch+1];
987 :     strcpy(ytitle,fHistogram->GetYaxis()->GetTitle());
988 :     }
989 :     delete fHistogram;
990 :     fHistogram = 0;
991 :     }
992 :     }
993 :     if (fHistogram) {
994 :     minimum = fHistogram->GetYaxis()->GetXmin();
995 :     maximum = fHistogram->GetYaxis()->GetXmax();
996 :     uxmin = gPad->PadtoX(rwxmin);
997 :     uxmax = gPad->PadtoX(rwxmax);
998 :     } else {
999 :     rwxmin = 1e100;
1000 :     rwxmax = -rwxmin;
1001 :     rwymin = rwxmin;
1002 :     rwymax = -rwymin;
1003 :     while ((g = (TGraph*) next())) {
1004 :     Int_t npoints = g->GetN();
1005 :     x = g->GetX();
1006 :     y = g->GetY();
1007 :     for (i=0;i<npoints;i++) {
1008 :     if (x[i] < rwxmin) rwxmin = x[i];
1009 :     if (x[i] > rwxmax) rwxmax = x[i];
1010 :     if (y[i] > rwymax) rwymax = y[i];
1011 :     if (y[i] < rwymin) rwymin = y[i];
1012 :     }
1013 :     g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
1014 :     if (g->GetN() > npt) npt = g->GetN();
1015 :     }
1016 :     if (rwxmin == rwxmax) rwxmax += 1.;
1017 :     if (rwymin == rwymax) rwymax += 1.;
1018 :     dx = 0.05*(rwxmax-rwxmin);
1019 :     dy = 0.05*(rwymax-rwymin);
1020 :     uxmin = rwxmin - dx;
1021 :     uxmax = rwxmax + dx;
1022 :     if (gPad->GetLogy()) {
1023 :     if (rwymin <= 0) rwymin = 0.001*rwymax;
1024 :     minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
1025 :     maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
1026 :     } else {
1027 :     minimum = rwymin - dy;
1028 :     maximum = rwymax + dy;
1029 :     }
1030 :     if (minimum < 0 && rwymin >= 0) minimum = 0;
1031 :     if (maximum > 0 && rwymax <= 0) maximum = 0;
1032 :     }
1033 : brun 754
1034 : couet 13304 if (fMinimum != -1111) rwymin = minimum = fMinimum;
1035 :     if (fMaximum != -1111) rwymax = maximum = fMaximum;
1036 :     if (uxmin < 0 && rwxmin >= 0) {
1037 :     if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
1038 :     //else uxmin = 0;
1039 :     }
1040 :     if (uxmax > 0 && rwxmax <= 0) {
1041 :     if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
1042 :     //else uxmax = 0;
1043 :     }
1044 :     if (minimum < 0 && rwymin >= 0) {
1045 :     if(gPad->GetLogy()) minimum = 0.9*rwymin;
1046 :     //else minimum = 0;
1047 :     }
1048 :     if (maximum > 0 && rwymax <= 0) {
1049 :     if(gPad->GetLogy()) maximum = 1.1*rwymax;
1050 :     //else maximum = 0;
1051 :     }
1052 :     if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
1053 :     if (uxmin <= 0 && gPad->GetLogx()) {
1054 :     if (uxmax > 1000) uxmin = 1;
1055 :     else uxmin = 0.001*uxmax;
1056 :     }
1057 :     rwymin = minimum;
1058 :     rwymax = maximum;
1059 :     if (fHistogram) {
1060 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1061 :     }
1062 : brun 3575
1063 : couet 13304 // Create a temporary histogram to draw the axis
1064 :     if (!fHistogram) {
1065 :     // the graph is created with at least as many channels as there are points
1066 :     // to permit zooming on the full range
1067 :     rwxmin = uxmin;
1068 :     rwxmax = uxmax;
1069 :     fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1070 :     if (!fHistogram) return;
1071 :     fHistogram->SetMinimum(rwymin);
1072 :     fHistogram->SetBit(TH1::kNoStats);
1073 :     fHistogram->SetMaximum(rwymax);
1074 :     fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1075 :     fHistogram->SetDirectory(0);
1076 :     if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
1077 :     if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
1078 :     if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
1079 :     }
1080 :     fHistogram->Paint("0");
1081 : brun 3575 }
1082 : rdm 3742
1083 : brun 754 if (fGraphs) {
1084 : brun 4037 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
1085 :     TObject *obj;
1086 :    
1087 :     while (lnk) {
1088 :     obj = lnk->GetObject();
1089 :     if (strlen(lnk->GetOption())) obj->Paint(lnk->GetOption());
1090 :     else obj->Paint(chopt);
1091 :     lnk = (TObjOptLink*)lnk->Next();
1092 :     }
1093 : brun 754 }
1094 : brun 11226
1095 :     TObject *f;
1096 :     if (fFunctions) {
1097 : couet 13304 TIter next(fFunctions);
1098 :     while ((f = (TObject*) next())) {
1099 :     if (f->InheritsFrom(TF1::Class())) {
1100 :     if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
1101 :     } else {
1102 :     f->Paint();
1103 :     }
1104 : brun 11226 }
1105 :     }
1106 : couet 13304 }
1107 : brun 11226
1108 :    
1109 : brun 754 //______________________________________________________________________________
1110 : brun 1205 void TMultiGraph::Print(Option_t *option) const
1111 : brun 754 {
1112 : couet 13304 // Print the list of graphs
1113 : brun 754
1114 :     TGraph *g;
1115 :     if (fGraphs) {
1116 : couet 13304 TIter next(fGraphs);
1117 :     while ((g = (TGraph*) next())) {
1118 :     g->Print(option);
1119 :     }
1120 : brun 754 }
1121 :     }
1122 :    
1123 : couet 13304
1124 : brun 754 //______________________________________________________________________________
1125 : brun 5552 void TMultiGraph::RecursiveRemove(TObject *obj)
1126 :     {
1127 :     // Recursively remove this object from a list. Typically implemented
1128 :     // by classes that can contain mulitple references to a same object.
1129 :    
1130 :     if (!fGraphs) return;
1131 :     TObject *objr = fGraphs->Remove(obj);
1132 :     if (!objr) return;
1133 :     delete fHistogram; fHistogram = 0;
1134 :     if (gPad) gPad->Modified();
1135 :     }
1136 :    
1137 : couet 13304
1138 : brun 5552 //______________________________________________________________________________
1139 : brun 754 void TMultiGraph::SavePrimitive(ofstream &out, Option_t *option)
1140 :     {
1141 : couet 13304 // Save primitive as a C++ statement(s) on output stream out
1142 : brun 754
1143 :     char quote = '"';
1144 :     out<<" "<<endl;
1145 :     if (gROOT->ClassSaved(TMultiGraph::Class())) {
1146 : couet 13304 out<<" ";
1147 : brun 754 } else {
1148 : couet 13304 out<<" TMultiGraph *";
1149 : brun 754 }
1150 :     out<<"multigraph = new TMultiGraph();"<<endl;
1151 :     out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<endl;
1152 :     out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
1153 :    
1154 :     if (fGraphs) {
1155 : brun 10847 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
1156 :     TObject *g;
1157 :    
1158 :     while (lnk) {
1159 :     g = lnk->GetObject();
1160 :     g->SavePrimitive(out,"multigraph");
1161 :    
1162 :     out<<" multigraph->Add(graph,"<<quote<<lnk->GetOption()<<quote<<");"<<endl;
1163 :    
1164 :     lnk = (TObjOptLink*)lnk->Next();
1165 :    
1166 : couet 13304 }
1167 : brun 754 }
1168 :     out<<" multigraph->Draw("
1169 :     <<quote<<option<<quote<<");"<<endl;
1170 :     }
1171 :    
1172 : couet 13304
1173 : brun 754 //______________________________________________________________________________
1174 :     void TMultiGraph::SetMaximum(Double_t maximum)
1175 :     {
1176 : couet 13304 // Set multigraph maximum.
1177 :    
1178 : brun 754 fMaximum = maximum;
1179 :     if (fHistogram) fHistogram->SetMaximum(maximum);
1180 :     }
1181 :    
1182 : couet 13304
1183 : brun 754 //______________________________________________________________________________
1184 :     void TMultiGraph::SetMinimum(Double_t minimum)
1185 :     {
1186 : couet 13304 // Set multigraph minimum.
1187 :    
1188 : brun 754 fMinimum = minimum;
1189 :     if (fHistogram) fHistogram->SetMinimum(minimum);
1190 :     }

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9