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

1 : rdm 11750 // @(#)root/graf:$Name: $:$Id: TMultiGraph.cxx,v 1.21 2005/05/02 21:45:05 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 : rdm 11750 //
196 : brun 11226 // 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 : rdm 11750 // is of (error of x)**2 order. This approach is called "effective variance method".
292 : brun 11226 // 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 : rdm 11750
339 : brun 11226 // 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 : rdm 11750
349 : brun 11226 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 : rdm 11750
355 : brun 11226 // 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 : rdm 11750
368 : brun 11226 TString opt = option;
369 :     opt.ToUpper();
370 :    
371 :     if (opt.Contains("U")) fitOption.User = 1;
372 :     if (opt.Contains("Q")) fitOption.Quiet = 1;
373 :     if (opt.Contains("V")){fitOption.Verbose = 1; fitOption.Quiet = 0;}
374 :     if (opt.Contains("W")) fitOption.W1 = 1;
375 :     if (opt.Contains("E")) fitOption.Errors = 1;
376 :     if (opt.Contains("R")) fitOption.Range = 1;
377 :     if (opt.Contains("N")) fitOption.Nostore = 1;
378 :     if (opt.Contains("0")) fitOption.Nograph = 1;
379 :     if (opt.Contains("+")) fitOption.Plus = 1;
380 :     if (opt.Contains("B")) fitOption.Bound = 1;
381 :     if (opt.Contains("C")) fitOption.Nochisq = 1;
382 : brun 11725 if (opt.Contains("F")) fitOption.Minuit = 1;
383 : brun 11226
384 :     if (rxmax > rxmin) {
385 :     xmin = rxmin;
386 :     xmax = rxmax;
387 :     } else {
388 :     g=(TGraph *)fGraphs->First();
389 :     if (!g) {
390 :     Error("Fit", "No graphs in the multigraph");
391 :     return 0;
392 :     }
393 :     Double_t *px, *py;
394 :     np=g->GetN();
395 :     px=g->GetX();
396 :     py=g->GetY();
397 :     xmin=px[0];
398 :     xmax=py[np-1];
399 :     ymin=px[0];
400 :     ymax=py[np-1];
401 :     Double_t err0=g->GetErrorX(0);
402 :     Double_t errn=g->GetErrorX(np-1);
403 :     if (err0 > 0) xmin -= 2*err0;
404 :     if (errn > 0) xmax += 2*errn;
405 :    
406 :     next.Reset();
407 :     while ((g = (TGraph*) next())) {
408 :     np=g->GetN();
409 :     px=g->GetX();
410 :     py=g->GetY();
411 :     for (i=0; i<np; i++) {
412 :     if (px[i] < xmin) xmin = px[i];
413 :     if (px[i] > xmax) xmax = px[i];
414 :     if (py[i] < ymin) ymin = py[i];
415 :     if (py[i] > ymax) ymax = py[i];
416 :     }
417 :     }
418 :     }
419 :    
420 :     ///////////////
421 :     //set the fitter
422 :     //////////////
423 : brun 11605
424 : brun 11226 Int_t special=f1->GetNumber();
425 :     Bool_t linear = f1->IsLinear();
426 :     if (special==299+npar)
427 :     linear=kTRUE;
428 : brun 11605 if (fitOption.Bound || fitOption.User || fitOption.Errors || fitOption.Minuit)
429 :     linear = kFALSE;
430 : brun 11226
431 :     char l[]="TLinearFitter";
432 :     Int_t strdiff = 0;
433 :     Bool_t IsSet = kFALSE;
434 :     if (TVirtualFitter::GetFitter()){
435 :     //Is a fitter already set? Is it linear?
436 :     IsSet = kTRUE;
437 :     strdiff = strcmp(TVirtualFitter::GetFitter()->IsA()->GetName(), l);
438 :     }
439 :     if (linear){
440 :     TClass *cl = gROOT->GetClass("TLinearFitter");
441 :     if (IsSet && strdiff!=0) {
442 :     delete TVirtualFitter::GetFitter();
443 :     IsSet=kFALSE;
444 :     }
445 :     if (!IsSet) {
446 :     TVirtualFitter::SetFitter((TVirtualFitter *)cl->New());
447 :     }
448 :     } else {
449 :     if (IsSet && strdiff==0){
450 :     delete TVirtualFitter::GetFitter();
451 :     IsSet=kFALSE;
452 :     }
453 : rdm 11750 if (!IsSet)
454 :     TVirtualFitter::SetFitter(0);
455 : brun 11226 }
456 :    
457 :     TVirtualFitter *grFitter = TVirtualFitter::Fitter(this, f1->GetNpar());
458 :     grFitter->Clear();
459 : rdm 11750
460 : brun 11226 //*-*- Get pointer to the function by searching in the list of functions in ROOT
461 :     grFitter->SetUserFunc(f1);
462 :     grFitter->SetFitOption(fitOption);
463 :    
464 :     //*-*- Is a Fit range specified?
465 :     if (fitOption.Range) {
466 :     f1->GetRange(xmin, xmax);
467 :     } else {
468 :     f1->SetRange(xmin, xmax);
469 :     }
470 :    
471 : brun 11605 if (linear){
472 : brun 11226 grFitter->ExecuteCommand("FitMultiGraph", 0, 0);
473 : rdm 11750
474 : brun 11226 } else {
475 :    
476 :     //Int_t special = f1->GetNumber();
477 :     if (fitOption.Bound) special = 0;
478 :     if (special == 100) InitGaus(xmin,xmax);
479 :     else if (special == 400) InitGaus(xmin,xmax);
480 :     else if (special == 200) InitExpo(xmin,xmax);
481 :     else if (special == 299+npar) InitPolynom(xmin,xmax);
482 : rdm 11750
483 : brun 11226 //*-*- Some initialisations
484 :     if (!fitOption.Verbose) {
485 :     arglist[0] = -1;
486 :     grFitter->ExecuteCommand("SET PRINT", arglist,1);
487 :     arglist[0] = 0;
488 :     grFitter->ExecuteCommand("SET NOW", arglist,0);
489 :     }
490 : rdm 11750
491 :     /////////////////////////////////////////////////////////
492 : brun 11226 //*-*- Set error criterion for chisquare
493 :     arglist[0] = TVirtualFitter::GetErrorDef();
494 :     if (!fitOption.User) grFitter->SetFitMethod("MultiGraphFitChisquare");
495 :    
496 :    
497 :     fitResult = grFitter->ExecuteCommand("SET ERR",arglist,1);
498 :     if (fitResult != 0) {
499 :     // Abnormal termination, MIGRAD might not have converged on a
500 :     // minimum.
501 :     if (!fitOption.Quiet) {
502 :     Warning("Fit","Abnormal termination of minimization.");
503 :     }
504 :     delete [] arglist;
505 :     return fitResult;
506 :     }
507 : rdm 11750
508 : brun 11226 //*-*- Transfer names and initial values of parameters to Minuit
509 :     Int_t nfixed = 0;
510 :     for (i=0;i<npar;i++) {
511 :     par = f1->GetParameter(i);
512 :     f1->GetParLimits(i,al,bl);
513 :     if (al*bl != 0 && al >= bl) {
514 :     al = bl = 0;
515 :     arglist[nfixed] = i+1;
516 :     nfixed++;
517 :     }
518 :     we = 0.3*TMath::Abs(par);
519 :     if (we <= TMath::Abs(par)*1e-6) we = 1;
520 :     grFitter->SetParameter(i,f1->GetParName(i),par,we,al,bl);
521 :     }
522 :     if(nfixed > 0)grFitter->ExecuteCommand("FIX",arglist,nfixed); // Otto
523 : rdm 11750
524 : brun 11226 //*-*- Reset Print level
525 :     if (!fitOption.Quiet) {
526 :     if (fitOption.Verbose) { arglist[0] = 2; grFitter->ExecuteCommand("SET PRINT", arglist,1); }
527 :     else { arglist[0] = 0; grFitter->ExecuteCommand("SET PRINT", arglist,1); }
528 :     }
529 :     //*-*- Compute sum of squares of errors in the bin range
530 :     Bool_t hasErrors = kFALSE;
531 :     Double_t ex, ey, sumw2=0;
532 :     next.Reset();
533 :     while ((g = (TGraph*) next())) {
534 :     np=g->GetN();
535 :     for (i=0; i<np; i++){
536 :     ex=g->GetErrorX(i);
537 :     ey=g->GetErrorY(i);
538 :     if (ex > 0 || ey > 0) hasErrors=kTRUE;
539 :     sumw2+=ey*ey;
540 :     }
541 :     }
542 :    
543 :     //*-*- Perform minimization
544 :    
545 :     arglist[0] = TVirtualFitter::GetMaxIterations();
546 :     arglist[1] = sumw2*TVirtualFitter::GetPrecision();
547 :     grFitter->ExecuteCommand("MIGRAD",arglist,2);
548 :     if (fitOption.Errors) {
549 :     grFitter->ExecuteCommand("HESSE",arglist,0);
550 :     grFitter->ExecuteCommand("MINOS",arglist,0);
551 :     }
552 : rdm 11750
553 : brun 11226 grFitter->GetStats(amin,edm,errdef,nvpar,nparx);
554 :     f1->SetChisquare(amin);
555 :     Int_t ndf = f1->GetNumberFitPoints()-npar+nfixed;
556 :     f1->SetNDF(ndf);
557 :    
558 :     //*-*- Get return status
559 :     char parName[50];
560 :     for (i=0;i<npar;i++) {
561 :     grFitter->GetParameter(i,parName, par,we,al,bl);
562 :     if (!fitOption.Errors) werr = we;
563 :     else {
564 :     grFitter->GetErrors(i,eplus,eminus,eparab,globcc);
565 :     if (eplus > 0 && eminus < 0) werr = 0.5*(eplus-eminus);
566 :     else werr = we;
567 :     }
568 :     if (!hasErrors && ndf > 1) werr *= TMath::Sqrt(amin/(ndf-1));
569 :     f1->SetParameter(i,par);
570 :     f1->SetParError(i,werr);
571 :     }
572 :     }
573 :     //*-*- Print final values of parameters.
574 :     if (!fitOption.Quiet) {
575 :     if (fitOption.Errors) grFitter->PrintResults(4,amin);
576 :     else grFitter->PrintResults(3,amin);
577 :     }
578 :     delete [] arglist;
579 :    
580 :     //*-*- Store fitted function in histogram functions list and draw
581 :    
582 :     if (!fitOption.Nostore) {
583 :     if (!fFunctions) fFunctions = new TList;
584 :     if (!fitOption.Plus) {
585 :     TIter next2(fFunctions, kIterBackward);
586 :     TObject *obj;
587 :     while ((obj = next2())) {
588 :     if (obj->InheritsFrom(TF1::Class())){
589 : rdm 11750 obj = fFunctions->Remove(obj);
590 : brun 11226 delete obj;
591 :     }
592 :     }
593 :     }
594 :     fnew1 = new TF1();
595 :     f1->Copy(*fnew1);
596 :     fFunctions->Add(fnew1);
597 :     fnew1->SetParent(this);
598 :     fnew1->Save(xmin,xmax,0,0,0,0);
599 :     if (fitOption.Nograph) fnew1->SetBit(TF1::kNotDraw);
600 :     fnew1->SetBit(TFormula::kNotGlobal);
601 :    
602 :     if (TestBit(kCanDelete)) return fitResult;
603 :     if (gPad) gPad->Modified();
604 :     }
605 :    
606 :    
607 :     return fitResult;
608 :    
609 :     }
610 :    
611 :    
612 :     //______________________________________________________________________________
613 : brun 11584 Option_t *TMultiGraph::GetGraphDrawOption(const TGraph *gr) const
614 :     {
615 :     // Return the draw option for the TGraph gr in this TMultiGraph
616 :     // The return option is the one specified when calling TMultiGraph::Add(gr,option).
617 : rdm 11750
618 : brun 11584 if (!fGraphs || !gr) return "";
619 :     TListIter next(fGraphs);
620 :     TObject *obj;
621 :     while ((obj = next())) {
622 : brun 11600 if (obj == (TObject*)gr) return next.GetOption();
623 : brun 11584 }
624 :     return "";
625 :     }
626 :    
627 :     //______________________________________________________________________________
628 : brun 11226 void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
629 :     {
630 :     //*-*-*-*-*-*Compute Initial values of parameters for a gaussian*-*-*-*-*-*-*
631 :     //*-* ===================================================
632 :    
633 :     Double_t allcha, sumx, sumx2, x, val, rms, mean;
634 :     Int_t bin;
635 :     const Double_t sqrtpi = 2.506628;
636 :    
637 :     //*-*- Compute mean value and RMS of the graph in the given range
638 :     Int_t np = 0;
639 :     allcha = sumx = sumx2 = 0;
640 :     TGraph *g;
641 :     TIter next(fGraphs);
642 :     Double_t *px, *py;
643 :     Int_t npp; //number of points in each graph
644 :     while ((g = (TGraph*) next())) {
645 :     px=g->GetX();
646 :     py=g->GetY();
647 :     npp=g->GetN();
648 :     for (bin=0; bin<npp; bin++){
649 :     x=px[bin];
650 :     if (x<xmin || x>xmax) continue;
651 :     np++;
652 :     val=py[bin];
653 :     sumx+=val*x;
654 :     sumx2+=val*x*x;
655 :     allcha+=val;
656 :     }
657 :     }
658 :     if (np == 0 || allcha == 0) return;
659 :     mean = sumx/allcha;
660 :     rms = TMath::Sqrt(sumx2/allcha - mean*mean);
661 :    
662 :     Double_t binwidx = TMath::Abs((xmax-xmin)/np);
663 :     if (rms == 0) rms = 1;
664 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
665 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
666 :     f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
667 :     f1->SetParameter(1,mean);
668 :     f1->SetParameter(2,rms);
669 :     f1->SetParLimits(2,0,10*rms);
670 :     }
671 :    
672 :     //______________________________________________________________________________
673 :     void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
674 :     {
675 :     //*-*-*-*-*-*Compute Initial values of parameters for an exponential*-*-*-*-*
676 :     //*-* =======================================================
677 :    
678 :     Double_t constant, slope;
679 :     Int_t ifail;
680 :    
681 :     LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
682 :    
683 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
684 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
685 :     f1->SetParameter(0,constant);
686 :     f1->SetParameter(1,slope);
687 :    
688 :     }
689 :    
690 :     //______________________________________________________________________________
691 :     void TMultiGraph::InitPolynom(Double_t xmin, Double_t xmax)
692 :     {
693 :     //*-*-*-*-*-*Compute Initial values of parameters for a polynom*-*-*-*-*-*-*
694 :     //*-* ===================================================
695 :    
696 :     Double_t fitpar[25];
697 :    
698 :     TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
699 :     TF1 *f1 = (TF1*)grFitter->GetUserFunc();
700 :     Int_t npar = f1->GetNpar();
701 :    
702 :     LeastSquareFit(npar, fitpar, xmin, xmax);
703 :    
704 :     for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
705 : brun 11605
706 : brun 11226 }
707 :    
708 :     //______________________________________________________________________________
709 :     void TMultiGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
710 :     {
711 :     //*-*-*-*-*-*-*-*Least squares lpolynomial fitting without weights*-*-*-*-*-*-*
712 :     //*-* =================================================
713 :     //
714 :     // m number of parameters
715 :     // a array of parameters
716 :     // first 1st point number to fit (default =0)
717 :     // last last point number to fit (default=fNpoints-1)
718 :     //
719 :     // based on CERNLIB routine LSQ: Translated to C++ by Rene Brun
720 :     //
721 :     //
722 :     const Double_t zero = 0.;
723 :     const Double_t one = 1.;
724 :     const Int_t idim = 20;
725 :    
726 :     Double_t b[400] /* was [20][20] */;
727 :     Int_t i, k, l, ifail, bin;
728 :     Double_t power;
729 :     Double_t da[20], xk, yk;
730 :    
731 :    
732 :     //count the total number of points to fit
733 :     TGraph *g;
734 :     TIter next(fGraphs);
735 :     Double_t *px, *py;
736 : rdm 11750 Int_t n=0;
737 : brun 11226 Int_t npp;
738 :     while ((g = (TGraph*) next())) {
739 :     px=g->GetX();
740 :     py=g->GetY();
741 :     npp=g->GetN();
742 :     for (bin=0; bin<npp; bin++){
743 :     xk=px[bin];
744 :     if (xk < xmin || xk > xmax) continue;
745 :     n++;
746 :     }
747 :     }
748 :     if (m <= 2) {
749 :     LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
750 :     return;
751 :     }
752 :     if (m > idim || m > n) return;
753 :     da[0] = zero;
754 :     for (l = 2; l <= m; ++l) {
755 :     b[l-1] = zero;
756 :     b[m + l*20 - 21] = zero;
757 :     da[l-1] = zero;
758 :     }
759 :     Int_t np = 0;
760 :    
761 :     next.Reset();
762 :     while ((g = (TGraph*) next())) {
763 :     px=g->GetX();
764 :     py=g->GetY();
765 :     npp=g->GetN();
766 :    
767 :     for (k = 0; k <= npp; ++k) {
768 :     xk = px[k];
769 :     if (xk < xmin || xk > xmax) continue;
770 :     np++;
771 :     yk = py[k];
772 :     power = one;
773 :     da[0] += yk;
774 :     for (l = 2; l <= m; ++l) {
775 :     power *= xk;
776 :     b[l-1] += power;
777 :     da[l-1] += power*yk;
778 :     }
779 :     for (l = 2; l <= m; ++l) {
780 :     power *= xk;
781 :     b[m + l*20 - 21] += power;
782 :     }
783 :     }
784 :     }
785 :     b[0] = Double_t(np);
786 :     for (i = 3; i <= m; ++i) {
787 :     for (k = i; k <= m; ++k) {
788 :     b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
789 :     }
790 :     }
791 :     H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
792 :    
793 :     if (ifail < 0) {
794 :     //a[0] = fY[0];
795 :     py=((TGraph *)fGraphs->First())->GetY();
796 :     a[0]=py[0];
797 :     for (i=1; i<m; ++i) a[i] = 0;
798 :     return;
799 :     }
800 :     for (i=0; i<m; ++i) a[i] = da[i];
801 :    
802 :     }
803 :    
804 :     //______________________________________________________________________________
805 :     void TMultiGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
806 :     {
807 :     //*-*-*-*-*-*-*-*-*-*Least square linear fit without weights*-*-*-*-*-*-*-*-*
808 :     //*-* =======================================
809 :     //
810 :     // Fit a straight line (a0 + a1*x) to the data in this graph.
811 :     // ndata: number of points to fit
812 :     // first: first point number to fit
813 :     // last: last point to fit O(ndata should be last-first
814 :     // ifail: return parameter indicating the status of the fit (ifail=0, fit is OK)
815 :     //
816 :     // extracted from CERNLIB LLSQ: Translated to C++ by Rene Brun
817 :     //
818 :    
819 :     Double_t xbar, ybar, x2bar;
820 :     Int_t i;
821 :     Double_t xybar;
822 :     Double_t fn, xk, yk;
823 :     Double_t det;
824 :    
825 :     ifail = -2;
826 :     xbar = ybar = x2bar = xybar = 0;
827 :     Int_t np = 0;
828 :     TGraph *g;
829 :     TIter next(fGraphs);
830 :     Double_t *px, *py;
831 : rdm 11750 Int_t npp;
832 : brun 11226 while ((g = (TGraph*) next())) {
833 :     px=g->GetX();
834 :     py=g->GetY();
835 :     npp=g->GetN();
836 :     for (i = 0; i < npp; ++i) {
837 :     xk = px[i];
838 :     if (xk < xmin || xk > xmax) continue;
839 :     np++;
840 :     yk = py[i];
841 :     if (ndata < 0) {
842 :     if (yk <= 0) yk = 1e-9;
843 :     yk = TMath::Log(yk);
844 :     }
845 :     xbar += xk;
846 :     ybar += yk;
847 :     x2bar += xk*xk;
848 :     xybar += xk*yk;
849 :     }
850 :     }
851 :     fn = Double_t(np);
852 :     det = fn*x2bar - xbar*xbar;
853 :     ifail = -1;
854 :     if (det <= 0) {
855 :     if (fn > 0) a0 = ybar/fn;
856 :     else a0 = 0;
857 :     a1 = 0;
858 :     return;
859 :     }
860 :     ifail = 0;
861 :     a0 = (x2bar*ybar - xbar*xybar) / det;
862 :     a1 = (fn*xybar - xbar*ybar) / det;
863 :    
864 :    
865 :    
866 :     }
867 :    
868 :     //______________________________________________________________________________
869 : brun 1205 TH1F *TMultiGraph::GetHistogram() const
870 : brun 754 {
871 :     // Returns a pointer to the histogram used to draw the axis
872 :     // Takes into account the two following cases.
873 :     // 1- option 'A' was specified in TMultiGraph::Draw. Return fHistogram
874 :     // 2- user had called TPad::DrawFrame. return pointer to hframe histogram
875 :    
876 :     if (fHistogram) return fHistogram;
877 :     if (!gPad) return 0;
878 :     gPad->Modified();
879 :     gPad->Update();
880 :     if (fHistogram) return fHistogram;
881 :     TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
882 :     return h1;
883 :     }
884 :    
885 :     //______________________________________________________________________________
886 : brun 11226 TF1 *TMultiGraph::GetFunction(const char *name) const
887 :     {
888 :     //*-*-*-*-*Return pointer to function with name*-*-*-*-*-*-*-*-*-*-*-*-*
889 :     //*-* ===================================
890 :     //
891 :     // Functions such as TGraph::Fit store the fitted function in the list of
892 :     // functions of this graph.
893 :    
894 :     if (!fFunctions) return 0;
895 :     return (TF1*)fFunctions->FindObject(name);
896 :     }
897 :    
898 :     //______________________________________________________________________________
899 : brun 1205 TAxis *TMultiGraph::GetXaxis() const
900 : brun 754 {
901 :     // Get x axis of the graph.
902 :    
903 :     if (!gPad) return 0;
904 : brun 4083 TH1 *h = GetHistogram();
905 :     if (!h) return 0;
906 :     return h->GetXaxis();
907 : brun 754 }
908 :    
909 :     //______________________________________________________________________________
910 : brun 1205 TAxis *TMultiGraph::GetYaxis() const
911 : brun 754 {
912 :     // Get y axis of the graph.
913 :    
914 :     if (!gPad) return 0;
915 : brun 4083 TH1 *h = GetHistogram();
916 :     if (!h) return 0;
917 :     return h->GetYaxis();
918 : brun 754 }
919 :    
920 :     //______________________________________________________________________________
921 :     void TMultiGraph::Paint(Option_t *option)
922 :     {
923 :     // paint all the graphs of this multigraph
924 :    
925 : brun 5552 if (fGraphs->GetSize() == 0) return;
926 : rdm 11750
927 : brun 754 char *l;
928 :     static char chopt[33];
929 :     Int_t nch = strlen(option);
930 : brun 1500 Int_t i;
931 :     for (i=0;i<nch;i++) chopt[i] = toupper(option[i]);
932 : brun 754 chopt[nch] = 0;
933 : brun 1500 Double_t *x, *y;
934 : brun 3575 TGraph *g;
935 : rdm 3742
936 : brun 754 l = strstr(chopt,"A");
937 :     if (l) {
938 :     *l = ' ';
939 :     TIter next(fGraphs);
940 :     Int_t npt = 100;
941 :     Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
942 : brun 3575 rwxmin = gPad->GetUxmin();
943 :     rwxmax = gPad->GetUxmax();
944 :     rwymin = gPad->GetUymin();
945 :     rwymax = gPad->GetUymax();
946 : brun 3799 char *xtitle = 0;
947 :     char *ytitle = 0;
948 :     Int_t firstx = 0;
949 :     Int_t lastx = 0;
950 : rdm 11750
951 : brun 754 if (fHistogram) {
952 : brun 3575 //cleanup in case of a previous unzoom
953 :     if (fHistogram->GetMinimum() >= fHistogram->GetMaximum()) {
954 : brun 3799 Int_t nch = strlen(fHistogram->GetXaxis()->GetTitle());
955 :     firstx = fHistogram->GetXaxis()->GetFirst();
956 :     lastx = fHistogram->GetXaxis()->GetLast();
957 :     if (nch) {
958 :     xtitle = new char[nch+1];
959 :     strcpy(xtitle,fHistogram->GetXaxis()->GetTitle());
960 :     }
961 :     nch = strlen(fHistogram->GetYaxis()->GetTitle());
962 :     if (nch) {
963 :     ytitle = new char[nch+1];
964 :     strcpy(ytitle,fHistogram->GetYaxis()->GetTitle());
965 :     }
966 : brun 3575 delete fHistogram;
967 :     fHistogram = 0;
968 :     }
969 :     }
970 :     if (fHistogram) {
971 :     minimum = fHistogram->GetYaxis()->GetXmin();
972 :     maximum = fHistogram->GetYaxis()->GetXmax();
973 :     uxmin = gPad->PadtoX(rwxmin);
974 :     uxmax = gPad->PadtoX(rwxmax);
975 : brun 754 } else {
976 :     rwxmin = 1e100;
977 :     rwxmax = -rwxmin;
978 :     rwymin = rwxmin;
979 :     rwymax = -rwymin;
980 :     while ((g = (TGraph*) next())) {
981 : brun 1500 Int_t npoints = g->GetN();
982 :     x = g->GetX();
983 :     y = g->GetY();
984 :     for (i=0;i<npoints;i++) {
985 :     if (x[i] < rwxmin) rwxmin = x[i];
986 :     if (x[i] > rwxmax) rwxmax = x[i];
987 : brun 9245 if (y[i] > rwymax) rwymax = y[i];
988 : brun 1500 if (y[i] < rwymin) rwymin = y[i];
989 :     }
990 : brun 754 g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
991 :     if (g->GetN() > npt) npt = g->GetN();
992 :     }
993 :     if (rwxmin == rwxmax) rwxmax += 1.;
994 :     if (rwymin == rwymax) rwymax += 1.;
995 : brun 3575 dx = 0.05*(rwxmax-rwxmin);
996 :     dy = 0.05*(rwymax-rwymin);
997 : brun 754 uxmin = rwxmin - dx;
998 :     uxmax = rwxmax + dx;
999 : brun 9245 if (gPad->GetLogy()) {
1000 :     if (rwymin <= 0) rwymin = 0.001*rwymax;
1001 :     minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
1002 :     maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
1003 :     } else {
1004 :     minimum = rwymin - dy;
1005 :     maximum = rwymax + dy;
1006 :     }
1007 : brun 9010 if (minimum < 0 && rwymin >= 0) minimum = 0;
1008 :     if (maximum > 0 && rwymax <= 0) maximum = 0;
1009 : brun 754 }
1010 :    
1011 :     if (fMinimum != -1111) rwymin = minimum = fMinimum;
1012 :     if (fMaximum != -1111) rwymax = maximum = fMaximum;
1013 :     if (uxmin < 0 && rwxmin >= 0) {
1014 :     if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
1015 : brun 3575 //else uxmin = 0;
1016 : brun 754 }
1017 :     if (uxmax > 0 && rwxmax <= 0) {
1018 :     if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
1019 : brun 3575 //else uxmax = 0;
1020 : brun 754 }
1021 :     if (minimum < 0 && rwymin >= 0) {
1022 :     if(gPad->GetLogy()) minimum = 0.9*rwymin;
1023 : brun 3575 //else minimum = 0;
1024 : brun 754 }
1025 :     if (maximum > 0 && rwymax <= 0) {
1026 :     if(gPad->GetLogy()) maximum = 1.1*rwymax;
1027 : brun 3575 //else maximum = 0;
1028 : brun 754 }
1029 :     if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
1030 :     if (uxmin <= 0 && gPad->GetLogx()) {
1031 :     if (uxmax > 1000) uxmin = 1;
1032 :     else uxmin = 0.001*uxmax;
1033 :     }
1034 :     rwymin = minimum;
1035 :     rwymax = maximum;
1036 :     if (fHistogram) {
1037 : brun 3575 fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1038 :     }
1039 :    
1040 :     //*-*- Create a temporary histogram to draw the axis
1041 :     if (!fHistogram) {
1042 :     // the graph is created with at least as many channels as there are points
1043 :     // to permit zooming on the full range
1044 :     rwxmin = uxmin;
1045 :     rwxmax = uxmax;
1046 :     fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
1047 :     if (!fHistogram) return;
1048 : brun 754 fHistogram->SetMinimum(rwymin);
1049 : brun 3575 fHistogram->SetBit(TH1::kNoStats);
1050 : brun 754 fHistogram->SetMaximum(rwymax);
1051 : brun 3575 fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
1052 :     fHistogram->SetDirectory(0);
1053 : brun 3799 if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
1054 :     if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
1055 :     if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
1056 : brun 754 }
1057 : brun 4037 fHistogram->Paint("0");
1058 : brun 3575 }
1059 : rdm 3742
1060 : brun 754 if (fGraphs) {
1061 : brun 4037 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
1062 :     TObject *obj;
1063 :    
1064 :     while (lnk) {
1065 :     obj = lnk->GetObject();
1066 :     if (strlen(lnk->GetOption())) obj->Paint(lnk->GetOption());
1067 :     else obj->Paint(chopt);
1068 :     lnk = (TObjOptLink*)lnk->Next();
1069 :     }
1070 : brun 754 }
1071 : brun 11226
1072 :     TObject *f;
1073 :     if (fFunctions) {
1074 :     TIter next(fFunctions);
1075 :     while ((f = (TObject*) next())) {
1076 :     if (f->InheritsFrom(TF1::Class())) {
1077 :     if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
1078 :     } else {
1079 :     f->Paint();
1080 :     }
1081 :     }
1082 :     }
1083 :    
1084 :    
1085 : brun 754 }
1086 :    
1087 :     //______________________________________________________________________________
1088 : brun 1205 void TMultiGraph::Print(Option_t *option) const
1089 : brun 754 {
1090 :     // Print the list of graphs
1091 :    
1092 :     TGraph *g;
1093 :     if (fGraphs) {
1094 :     TIter next(fGraphs);
1095 :     while ((g = (TGraph*) next())) {
1096 :     g->Print(option);
1097 :     }
1098 :     }
1099 :     }
1100 :    
1101 :     //______________________________________________________________________________
1102 : brun 5552 void TMultiGraph::RecursiveRemove(TObject *obj)
1103 :     {
1104 :     // Recursively remove this object from a list. Typically implemented
1105 :     // by classes that can contain mulitple references to a same object.
1106 :    
1107 :     if (!fGraphs) return;
1108 :     TObject *objr = fGraphs->Remove(obj);
1109 :     if (!objr) return;
1110 :     delete fHistogram; fHistogram = 0;
1111 :     if (gPad) gPad->Modified();
1112 :     }
1113 :    
1114 :     //______________________________________________________________________________
1115 : brun 754 void TMultiGraph::SavePrimitive(ofstream &out, Option_t *option)
1116 :     {
1117 :     // Save primitive as a C++ statement(s) on output stream out
1118 :    
1119 :     char quote = '"';
1120 :     out<<" "<<endl;
1121 :     if (gROOT->ClassSaved(TMultiGraph::Class())) {
1122 :     out<<" ";
1123 :     } else {
1124 :     out<<" TMultiGraph *";
1125 :     }
1126 :     out<<"multigraph = new TMultiGraph();"<<endl;
1127 :     out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<endl;
1128 :     out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
1129 :    
1130 :     if (fGraphs) {
1131 : brun 10847 TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
1132 :     TObject *g;
1133 :    
1134 :     while (lnk) {
1135 :     g = lnk->GetObject();
1136 :     g->SavePrimitive(out,"multigraph");
1137 :    
1138 :     out<<" multigraph->Add(graph,"<<quote<<lnk->GetOption()<<quote<<");"<<endl;
1139 :    
1140 :     lnk = (TObjOptLink*)lnk->Next();
1141 :    
1142 : brun 754 }
1143 :     }
1144 :     out<<" multigraph->Draw("
1145 :     <<quote<<option<<quote<<");"<<endl;
1146 :     }
1147 :    
1148 :     //______________________________________________________________________________
1149 :     void TMultiGraph::SetMaximum(Double_t maximum)
1150 :     {
1151 :     fMaximum = maximum;
1152 :     if (fHistogram) fHistogram->SetMaximum(maximum);
1153 :     }
1154 :    
1155 :     //______________________________________________________________________________
1156 :     void TMultiGraph::SetMinimum(Double_t minimum)
1157 :     {
1158 :     fMinimum = minimum;
1159 :     if (fHistogram) fHistogram->SetMinimum(minimum);
1160 :     }

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9