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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9