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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9