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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9