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

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

Subversion Admin
ViewVC Help
Powered by ViewVC 1.0.9