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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9