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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9