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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9