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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9