| 1 | /* ======================================================================== *\ | 
|---|
| 2 | ! | 
|---|
| 3 | ! * | 
|---|
| 4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction | 
|---|
| 5 | ! * Software. It is distributed to you in the hope that it can be a useful | 
|---|
| 6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. | 
|---|
| 7 | ! * It is distributed WITHOUT ANY WARRANTY. | 
|---|
| 8 | ! * | 
|---|
| 9 | ! * Permission to use, copy, modify and distribute this software and its | 
|---|
| 10 | ! * documentation for any purpose is hereby granted without fee, | 
|---|
| 11 | ! * provided that the above copyright notice appear in all copies and | 
|---|
| 12 | ! * that both that copyright notice and this permission notice appear | 
|---|
| 13 | ! * in supporting documentation. It is provided "as is" without express | 
|---|
| 14 | ! * or implied warranty. | 
|---|
| 15 | ! * | 
|---|
| 16 | ! | 
|---|
| 17 | ! | 
|---|
| 18 | !   Author(s): Wolfgang Wittek 5/2002 <mailto:wittek@mppmu.mpg.de> | 
|---|
| 19 | ! | 
|---|
| 20 | !   Copyright: MAGIC Software Development, 2000-2002 | 
|---|
| 21 | ! | 
|---|
| 22 | ! | 
|---|
| 23 | \* ======================================================================== */ | 
|---|
| 24 |  | 
|---|
| 25 | ////////////////////////////////////////////////////////////////////////////// | 
|---|
| 26 | //                                                                          // | 
|---|
| 27 | //  MHEffOnTime                                                             // | 
|---|
| 28 | //                                                                          // | 
|---|
| 29 | //  calculates the effective on time for each bin in the variable Var       // | 
|---|
| 30 | //             (Var may be time or Theta)                                   // | 
|---|
| 31 | //                                                                          // | 
|---|
| 32 | //  Important : The sample of events used for this should be as close       // | 
|---|
| 33 | //              to the sample of recorded events as possible.               // | 
|---|
| 34 | //              This means that NO additional cuts (be it quality cuts or   // | 
|---|
| 35 | //              gamma/hadron cuts) should be applied (see MAGIC-TDAS 02-02).// | 
|---|
| 36 | //                                                                          // | 
|---|
| 37 | ////////////////////////////////////////////////////////////////////////////// | 
|---|
| 38 |  | 
|---|
| 39 | #include "MHEffOnTime.h" | 
|---|
| 40 |  | 
|---|
| 41 | #include <TStyle.h> | 
|---|
| 42 |  | 
|---|
| 43 | #include <TMinuit.h> | 
|---|
| 44 | #include <TFitter.h> | 
|---|
| 45 |  | 
|---|
| 46 | #include <TF1.h> | 
|---|
| 47 | #include <TH2.h> | 
|---|
| 48 | #include <TCanvas.h> | 
|---|
| 49 |  | 
|---|
| 50 | #include "MBinning.h" | 
|---|
| 51 | #include "MParList.h" | 
|---|
| 52 |  | 
|---|
| 53 | #include "MLog.h" | 
|---|
| 54 | #include "MLogManip.h" | 
|---|
| 55 |  | 
|---|
| 56 | #include "MHTimeDiffTheta.h" | 
|---|
| 57 |  | 
|---|
| 58 | ClassImp(MHEffOnTime); | 
|---|
| 59 |  | 
|---|
| 60 |  | 
|---|
| 61 | // -------------------------------------------------------------------------- | 
|---|
| 62 | // | 
|---|
| 63 | // Default Constructor. It sets the name of the variable ("Time" or "Theta") | 
|---|
| 64 | //                      and the units of the variable ("[s]" or "[\\circ])") | 
|---|
| 65 | // | 
|---|
| 66 | MHEffOnTime::MHEffOnTime(const char *varname, const char *unit) | 
|---|
| 67 | : fHEffOn(), fHProb(), fHLambda(), fHRdead() | 
|---|
| 68 | { | 
|---|
| 69 | fVarname  = varname; | 
|---|
| 70 | fUnit     = unit; | 
|---|
| 71 |  | 
|---|
| 72 | TString strg = fVarname + " " + fUnit; | 
|---|
| 73 |  | 
|---|
| 74 | // | 
|---|
| 75 | //   set the name and title of this object | 
|---|
| 76 | // | 
|---|
| 77 | fName  = TString("MHEffOnTime-")+fVarname; | 
|---|
| 78 | fTitle = "1-D histogram of Eff On Time"; | 
|---|
| 79 |  | 
|---|
| 80 | // effective on time versus Var | 
|---|
| 81 | fHEffOn.SetName("EffOn"); | 
|---|
| 82 | fHEffOn.SetTitle(TString("T_{on, eff} vs. ")+fVarname); | 
|---|
| 83 |  | 
|---|
| 84 | fHEffOn.SetDirectory(NULL); | 
|---|
| 85 |  | 
|---|
| 86 | fHEffOn.SetXTitle(strg); | 
|---|
| 87 | fHEffOn.SetYTitle("T_{on, eff} [s]"); | 
|---|
| 88 |  | 
|---|
| 89 | // chi2-probability versus Var | 
|---|
| 90 | fHProb.SetName("Chi2-prob"); | 
|---|
| 91 | TString strg3("\\chi^{2}-prob of OnTimeFit vs. "); | 
|---|
| 92 | strg3 += fVarname; | 
|---|
| 93 | fHProb.SetTitle(strg3); | 
|---|
| 94 |  | 
|---|
| 95 | fHProb.SetDirectory(NULL); | 
|---|
| 96 |  | 
|---|
| 97 | fHProb.SetXTitle(strg); | 
|---|
| 98 | fHProb.SetYTitle("\\chi^{2}-probability"); | 
|---|
| 99 |  | 
|---|
| 100 | // lambda versus Var | 
|---|
| 101 | fHLambda.SetName("lambda"); | 
|---|
| 102 | fHLambda.SetTitle(TString("\\lambda from OnTimeFit vs. ")+fVarname); | 
|---|
| 103 |  | 
|---|
| 104 | fHLambda.SetDirectory(NULL); | 
|---|
| 105 |  | 
|---|
| 106 | fHLambda.SetXTitle(strg); | 
|---|
| 107 | fHLambda.SetYTitle("\\lambda [Hz]"); | 
|---|
| 108 |  | 
|---|
| 109 | // Rdead versus Var | 
|---|
| 110 | // Rdead is the fraction from real time lost by the dead time | 
|---|
| 111 | fHRdead.SetName("Rdead"); | 
|---|
| 112 |  | 
|---|
| 113 | fHRdead.SetTitle(TString("Rdead of OnTimeFit vs. ")+fVarname); | 
|---|
| 114 |  | 
|---|
| 115 | fHRdead.SetDirectory(NULL); | 
|---|
| 116 |  | 
|---|
| 117 | fHRdead.SetXTitle(strg); | 
|---|
| 118 | fHRdead.SetYTitle("Rdead"); | 
|---|
| 119 |  | 
|---|
| 120 | } | 
|---|
| 121 |  | 
|---|
| 122 | // -------------------------------------------------------------------- | 
|---|
| 123 | // | 
|---|
| 124 | //  determine range (yq[0], yq[1]) of time differences | 
|---|
| 125 | //  where fit should be performed; | 
|---|
| 126 | //  require a fraction >=min of all entries to lie below yq[0] | 
|---|
| 127 | //      and a fraction <=max of all entries to lie below yq[1]; | 
|---|
| 128 | //  within the range (yq[0], yq[1]) there must be no empty bin; | 
|---|
| 129 | // | 
|---|
| 130 | void MHEffOnTime::GetQuantiles(const TH1 &h, Double_t min, Double_t max, Double_t yq[2]) const | 
|---|
| 131 | { | 
|---|
| 132 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,02,07) | 
|---|
| 133 | // WOrkaround for missing const qualifier of TH1::Integral | 
|---|
| 134 | TH1 &h1 = (TH1&)h; | 
|---|
| 135 |  | 
|---|
| 136 | // choose pedestrian approach as long as GetQuantiles is | 
|---|
| 137 | // not available | 
|---|
| 138 | const Int_t    jbins = h1.GetNbinsX(); | 
|---|
| 139 | const Double_t Nm    = h1.Integral(); | 
|---|
| 140 |  | 
|---|
| 141 | const Double_t xq[2] = { 0.15*Nm, 0.98*Nm }; | 
|---|
| 142 |  | 
|---|
| 143 | yq[0] = yq[1] = h1.GetBinLowEdge(jbins+1); | 
|---|
| 144 |  | 
|---|
| 145 | for (int j=1; j<=jbins; j++) | 
|---|
| 146 | if (h1.Integral(2, j) >= xq[0]) | 
|---|
| 147 | { | 
|---|
| 148 | yq[0] = h1.GetBinLowEdge(j); | 
|---|
| 149 | break; | 
|---|
| 150 | } | 
|---|
| 151 |  | 
|---|
| 152 | for (int j=1; j<=jbins; j++) | 
|---|
| 153 | if (h1.Integral(1, j) >= xq[1] || h1.GetBinContent(j)==0) | 
|---|
| 154 | { | 
|---|
| 155 | yq[1] = h1.GetBinLowEdge(j); | 
|---|
| 156 | break; | 
|---|
| 157 | } | 
|---|
| 158 | #else | 
|---|
| 159 | // GetQuantiles doesn't seem to be available in root 3.01/06 | 
|---|
| 160 | Double_t xq[2] = { min, max }; | 
|---|
| 161 | ((TH1&)h).GetQuantiles(2, yq, xq); | 
|---|
| 162 | #endif | 
|---|
| 163 | } | 
|---|
| 164 |  | 
|---|
| 165 | void MHEffOnTime::DrawBin(TH1 &h, Int_t i) const | 
|---|
| 166 | { | 
|---|
| 167 | TString strg1 = fVarname+"-bin #"; | 
|---|
| 168 | strg1 += i; | 
|---|
| 169 |  | 
|---|
| 170 | new TCanvas(strg1, strg1); | 
|---|
| 171 |  | 
|---|
| 172 | gPad->SetLogy(); | 
|---|
| 173 | gStyle->SetOptFit(1011); | 
|---|
| 174 |  | 
|---|
| 175 | TString name="Bin_"; | 
|---|
| 176 | name += i; | 
|---|
| 177 |  | 
|---|
| 178 | h.SetName(name); | 
|---|
| 179 | h.SetXTitle("\\Delta t [s]"); | 
|---|
| 180 | h.SetYTitle("Counts"); | 
|---|
| 181 | h.DrawCopy(); | 
|---|
| 182 | } | 
|---|
| 183 |  | 
|---|
| 184 | Bool_t MHEffOnTime::CalcResults(const TF1 &func, Double_t Nm, Int_t i) | 
|---|
| 185 | { | 
|---|
| 186 | const Double_t lambda = func.GetParameter(0); | 
|---|
| 187 | const Double_t N0     = func.GetParameter(1); | 
|---|
| 188 | const Double_t prob   = func.GetProb(); | 
|---|
| 189 | const Int_t    NDF    = func.GetNDF(); | 
|---|
| 190 |  | 
|---|
| 191 | Double_t xmin, xmax; | 
|---|
| 192 | ((TF1&)func).GetRange(xmin, xmax); | 
|---|
| 193 |  | 
|---|
| 194 | *fLog << inf; | 
|---|
| 195 | *fLog << "Fitted bin #" << i << " from " << xmin << " to " << xmax; | 
|---|
| 196 | *fLog << ",  got: lambda=" << lambda << "Hz N0=" << N0 << endl; | 
|---|
| 197 |  | 
|---|
| 198 | if (prob<=0.01) | 
|---|
| 199 | { | 
|---|
| 200 | *fLog << warn << "WARN - Fit bin#" << i << " gives:"; | 
|---|
| 201 | *fLog << " Chi^2-Probab(" << prob << ")<0.01"; | 
|---|
| 202 | *fLog << " NoFitPts=" << func.GetNumberFitPoints(); | 
|---|
| 203 | *fLog << " Chi^2=" << func.GetChisquare(); | 
|---|
| 204 | } | 
|---|
| 205 |  | 
|---|
| 206 | // was fit successful ? | 
|---|
| 207 | if (NDF<=0 || /*prob<=0.001 ||*/ lambda<=0 || N0<=0) | 
|---|
| 208 | { | 
|---|
| 209 | *fLog << warn << dbginf << "Fit failed bin #" << i << ": "; | 
|---|
| 210 | if (NDF<=0) | 
|---|
| 211 | *fLog << " NDF(Number of Degrees of Freedom)=0"; | 
|---|
| 212 | if (lambda<=0) | 
|---|
| 213 | *fLog << " Parameter#0(lambda)=0"; | 
|---|
| 214 | if (N0<=0) | 
|---|
| 215 | *fLog << " Parameter#1(N0)=0"; | 
|---|
| 216 | *fLog << endl; | 
|---|
| 217 |  | 
|---|
| 218 | return kFALSE; | 
|---|
| 219 | } | 
|---|
| 220 |  | 
|---|
| 221 | // | 
|---|
| 222 | // -------------- start error calculation ---------------- | 
|---|
| 223 | // | 
|---|
| 224 | Double_t emat[2][2]; | 
|---|
| 225 | gMinuit->mnemat(&emat[0][0], 2); | 
|---|
| 226 |  | 
|---|
| 227 | // | 
|---|
| 228 | // Rdead : fraction of real time lost by the dead time | 
|---|
| 229 | // kappa = number of observed events (Nm) divided by | 
|---|
| 230 | //         the number of genuine events (N0) | 
|---|
| 231 | // Teff  : effective on-time | 
|---|
| 232 | // | 
|---|
| 233 | const Double_t Teff  = Nm/lambda; | 
|---|
| 234 | const Double_t kappa = Nm/N0; | 
|---|
| 235 | const Double_t Rdead = 1.0 - kappa; | 
|---|
| 236 |  | 
|---|
| 237 | const Double_t dldl   = emat[0][0]; | 
|---|
| 238 | const Double_t dN0dN0 = emat[1][1]; | 
|---|
| 239 |  | 
|---|
| 240 | const Double_t dTeff = Teff * sqrt(dldl/(lambda*lambda) + 1.0/Nm); | 
|---|
| 241 | const Double_t dl = sqrt(dldl); | 
|---|
| 242 | const Double_t dRdead = kappa * sqrt(dN0dN0/(N0*N0) + 1.0/Nm); | 
|---|
| 243 | // | 
|---|
| 244 | // -------------- end error calculation ---------------- | 
|---|
| 245 | // | 
|---|
| 246 |  | 
|---|
| 247 | // the effective on time is Nm/lambda | 
|---|
| 248 | fHEffOn.SetBinContent(i,  Teff); | 
|---|
| 249 | fHEffOn.SetBinError  (i, dTeff); | 
|---|
| 250 |  | 
|---|
| 251 | // plot chi2-probability of fit | 
|---|
| 252 | fHProb.SetBinContent(i, prob); | 
|---|
| 253 |  | 
|---|
| 254 | // lambda from fit | 
|---|
| 255 | fHLambda.SetBinContent(i, lambda); | 
|---|
| 256 | fHLambda.SetBinError  (i,     dl); | 
|---|
| 257 |  | 
|---|
| 258 | // Rdead from fit | 
|---|
| 259 | fHRdead.SetBinContent(i, Rdead); | 
|---|
| 260 | fHRdead.SetBinError  (i,dRdead); | 
|---|
| 261 |  | 
|---|
| 262 | return kTRUE; | 
|---|
| 263 | } | 
|---|
| 264 |  | 
|---|
| 265 | void MHEffOnTime::ResetBin(Int_t i) | 
|---|
| 266 | { | 
|---|
| 267 | fHEffOn.SetBinContent (i, 1.e-20); | 
|---|
| 268 | fHProb.SetBinContent  (i, 1.e-20); | 
|---|
| 269 | fHLambda.SetBinContent(i, 1.e-20); | 
|---|
| 270 | fHRdead.SetBinContent (i, 1.e-20); | 
|---|
| 271 | } | 
|---|
| 272 |  | 
|---|
| 273 | // ----------------------------------------------------------------------- | 
|---|
| 274 | // | 
|---|
| 275 | // Calculate the effective on time by fitting the distribution of | 
|---|
| 276 | // time differences | 
|---|
| 277 | // | 
|---|
| 278 | void MHEffOnTime::Calc(const TH2D *hist, const Bool_t draw) | 
|---|
| 279 | { | 
|---|
| 280 | // nbins = number of Var bins | 
|---|
| 281 | const Int_t nbins = hist->GetNbinsY(); | 
|---|
| 282 |  | 
|---|
| 283 | for (int i=1; i<=nbins; i++) | 
|---|
| 284 | { | 
|---|
| 285 | TH1D &h = *((TH2D*)hist)->ProjectionX(TString("Calc-")+fVarname, | 
|---|
| 286 | i, i, "E"); | 
|---|
| 287 | if (draw) | 
|---|
| 288 | DrawBin(h, i); | 
|---|
| 289 |  | 
|---|
| 290 | ResetBin(i); | 
|---|
| 291 |  | 
|---|
| 292 | // Nmdel = Nm * binwidth,  with Nm = number of observed events | 
|---|
| 293 | const Double_t Nm    = h.Integral(); | 
|---|
| 294 | const Double_t Nmdel = h.Integral("width"); | 
|---|
| 295 | if (Nm <= 0) | 
|---|
| 296 | { | 
|---|
| 297 | *fLog << warn << dbginf << "Nm<0 for bin #" << i << endl; | 
|---|
| 298 | delete &h; | 
|---|
| 299 | continue; | 
|---|
| 300 | } | 
|---|
| 301 |  | 
|---|
| 302 | Double_t yq[2]; | 
|---|
| 303 | GetQuantiles(h, 0.15, 0.95, yq); | 
|---|
| 304 |  | 
|---|
| 305 | // | 
|---|
| 306 | // Setup Poisson function for the fit: | 
|---|
| 307 | // lambda [Hz], N0 = ideal no of evts, del = bin width of dt | 
|---|
| 308 | // | 
|---|
| 309 | TF1 func("Poisson", " [1]*[2] * [0] * exp(-[0] *x)", yq[0], yq[1]); | 
|---|
| 310 | func.SetParNames("lambda", "N0", "del"); | 
|---|
| 311 |  | 
|---|
| 312 | func.SetParameter(0, 1); | 
|---|
| 313 | func.SetParameter(1, Nm); | 
|---|
| 314 | func.FixParameter(2, Nmdel/Nm); | 
|---|
| 315 |  | 
|---|
| 316 | // For me (TB) it seems that setting parameter limits isn't necessary | 
|---|
| 317 |  | 
|---|
| 318 | // For a description of the options see TH1::Fit | 
|---|
| 319 | h.Fit("Poisson", "0IRQ"); | 
|---|
| 320 |  | 
|---|
| 321 | // Calc and fill results of fit into the histograms | 
|---|
| 322 | const Bool_t rc = CalcResults(func, Nm, i); | 
|---|
| 323 |  | 
|---|
| 324 | // draw the distribution of time differences if requested | 
|---|
| 325 | if (rc && draw) | 
|---|
| 326 | func.DrawCopy("same"); | 
|---|
| 327 |  | 
|---|
| 328 | delete &h; | 
|---|
| 329 | } | 
|---|
| 330 | } | 
|---|
| 331 |  | 
|---|
| 332 | // ------------------------------------------------------------------------- | 
|---|
| 333 | // | 
|---|
| 334 | // Set the binnings and prepare the filling of the histograms | 
|---|
| 335 | // | 
|---|
| 336 | Bool_t MHEffOnTime::SetupFill(const MParList *plist) | 
|---|
| 337 | { | 
|---|
| 338 | TString bn = "Binning"+fVarname; | 
|---|
| 339 |  | 
|---|
| 340 | const MBinning* binsVar = (MBinning*)plist->FindObject(bn); | 
|---|
| 341 |  | 
|---|
| 342 | if (!binsVar) | 
|---|
| 343 | { | 
|---|
| 344 | *fLog << err << dbginf << bn << " [MBinning] not found... aborting." << endl; | 
|---|
| 345 | return kFALSE; | 
|---|
| 346 | } | 
|---|
| 347 |  | 
|---|
| 348 | SetBinning(&fHEffOn,  binsVar); | 
|---|
| 349 | SetBinning(&fHProb,   binsVar); | 
|---|
| 350 | SetBinning(&fHLambda, binsVar); | 
|---|
| 351 | SetBinning(&fHRdead,  binsVar); | 
|---|
| 352 |  | 
|---|
| 353 | fHEffOn.Sumw2(); | 
|---|
| 354 | fHProb.Sumw2(); | 
|---|
| 355 | fHLambda.Sumw2(); | 
|---|
| 356 | fHRdead.Sumw2(); | 
|---|
| 357 |  | 
|---|
| 358 | return kTRUE; | 
|---|
| 359 | } | 
|---|
| 360 |  | 
|---|
| 361 | // ------------------------------------------------------------------------- | 
|---|
| 362 | // | 
|---|
| 363 | // Draw the histogram | 
|---|
| 364 | // | 
|---|
| 365 | void MHEffOnTime::Draw(Option_t *opt) | 
|---|
| 366 | { | 
|---|
| 367 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this); | 
|---|
| 368 | pad->SetBorderMode(0); | 
|---|
| 369 |  | 
|---|
| 370 | pad->Divide(2,2); | 
|---|
| 371 |  | 
|---|
| 372 | pad->cd(1); | 
|---|
| 373 | gPad->SetBorderMode(0); | 
|---|
| 374 | fHEffOn.Draw(opt); | 
|---|
| 375 |  | 
|---|
| 376 | pad->cd(2); | 
|---|
| 377 | gPad->SetBorderMode(0); | 
|---|
| 378 | fHProb.Draw(opt); | 
|---|
| 379 |  | 
|---|
| 380 | pad->cd(3); | 
|---|
| 381 | gPad->SetBorderMode(0); | 
|---|
| 382 | fHLambda.Draw(opt); | 
|---|
| 383 |  | 
|---|
| 384 | pad->cd(4); | 
|---|
| 385 | gPad->SetBorderMode(0); | 
|---|
| 386 | fHRdead.Draw(opt); | 
|---|
| 387 |  | 
|---|
| 388 | pad->Modified(); | 
|---|
| 389 | pad->Update(); | 
|---|
| 390 | } | 
|---|
| 391 |  | 
|---|
| 392 | void MHEffOnTime::Calc(const MHTimeDiffTheta &hist, const Bool_t Draw) | 
|---|
| 393 | { | 
|---|
| 394 | Calc(hist.GetHist(), Draw); | 
|---|
| 395 | } | 
|---|
| 396 |  | 
|---|