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