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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9