| 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): Thomas Bretz 12/2000 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 19 | ! Author(s): Harald Kornmayer 1/2001
|
|---|
| 20 | !
|
|---|
| 21 | ! Copyright: MAGIC Software Development, 2000-2009
|
|---|
| 22 | !
|
|---|
| 23 | !
|
|---|
| 24 | \* ======================================================================== */
|
|---|
| 25 |
|
|---|
| 26 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 27 | //
|
|---|
| 28 | // MHCollectionArea
|
|---|
| 29 | //
|
|---|
| 30 | // Class Version 2:
|
|---|
| 31 | // ----------------
|
|---|
| 32 | // + added //! to fMcEvt which was missing before
|
|---|
| 33 | // - removed obsolete data member fEnergy
|
|---|
| 34 | //
|
|---|
| 35 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 36 | #include "MHCollectionArea.h"
|
|---|
| 37 |
|
|---|
| 38 | #include <TMath.h>
|
|---|
| 39 |
|
|---|
| 40 | #include <TLatex.h>
|
|---|
| 41 | #include <TCanvas.h>
|
|---|
| 42 | #include <TPaveStats.h>
|
|---|
| 43 |
|
|---|
| 44 | #include "MLog.h"
|
|---|
| 45 | #include "MLogManip.h"
|
|---|
| 46 |
|
|---|
| 47 | #include "MString.h"
|
|---|
| 48 | #include "MBinning.h"
|
|---|
| 49 |
|
|---|
| 50 | #include "MMcEvt.hxx"
|
|---|
| 51 | #include "MMcRunHeader.hxx"
|
|---|
| 52 | #include "MMcCorsikaRunHeader.h"
|
|---|
| 53 |
|
|---|
| 54 | #include "MParList.h"
|
|---|
| 55 | #include "MParameters.h"
|
|---|
| 56 |
|
|---|
| 57 | ClassImp(MHCollectionArea);
|
|---|
| 58 |
|
|---|
| 59 | using namespace std;
|
|---|
| 60 |
|
|---|
| 61 | // --------------------------------------------------------------------------
|
|---|
| 62 | //
|
|---|
| 63 | // Creates the three necessary histograms:
|
|---|
| 64 | // - selected showers (input)
|
|---|
| 65 | // - all showers (input)
|
|---|
| 66 | // - collection area (result)
|
|---|
| 67 | //
|
|---|
| 68 | MHCollectionArea::MHCollectionArea(const char *name, const char *title)
|
|---|
| 69 | : fMcEvt(0), fHeader(0), fMcAreaRadius(-1), fIsExtern(kFALSE)
|
|---|
| 70 | {
|
|---|
| 71 | // initialize the histogram for the distribution r vs E
|
|---|
| 72 | //
|
|---|
| 73 | // we set the energy range from 2 Gev to 20000 GeV (in log 4 orders
|
|---|
| 74 | // of magnitude) and for each order we take 25 subdivision --> 100 xbins
|
|---|
| 75 | //
|
|---|
| 76 | // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
|
|---|
| 77 | //
|
|---|
| 78 | fName = name ? name : "MHCollectionArea";
|
|---|
| 79 | fTitle = title ? title : "Collection Area vs. Energy/Theta";
|
|---|
| 80 |
|
|---|
| 81 | fHistSel.SetName("SelEvts");
|
|---|
| 82 | fHistSel.SetTitle("Number of Events after cuts");
|
|---|
| 83 | fHistSel.SetXTitle("\\Theta [deg]");
|
|---|
| 84 | fHistSel.SetYTitle("E [GeV]");
|
|---|
| 85 | fHistSel.SetDirectory(NULL);
|
|---|
| 86 | fHistSel.UseCurrentStyle();
|
|---|
| 87 | fHistSel.SetLineColor(kBlue);
|
|---|
| 88 |
|
|---|
| 89 | fHistAll.SetName("AllEvts");
|
|---|
| 90 | fHistAll.SetTitle("Number of events produced");
|
|---|
| 91 | fHistAll.SetXTitle("\\Theta [deg]");
|
|---|
| 92 | fHistAll.SetYTitle("E_{mc} [GeV]");
|
|---|
| 93 | fHistAll.SetDirectory(NULL);
|
|---|
| 94 | fHistAll.UseCurrentStyle();
|
|---|
| 95 |
|
|---|
| 96 | fHEnergy.SetName("CollEnergy");
|
|---|
| 97 | fHEnergy.SetTitle("Collection Area");
|
|---|
| 98 | fHEnergy.SetXTitle("E [GeV]");
|
|---|
| 99 | fHEnergy.SetYTitle("A_{eff} [m^{2}]");
|
|---|
| 100 | fHEnergy.SetDirectory(NULL);
|
|---|
| 101 | fHEnergy.UseCurrentStyle();
|
|---|
| 102 |
|
|---|
| 103 | MBinning binsa, binse, binst;
|
|---|
| 104 | binse.SetEdgesLog(21, 6.3, 100000);
|
|---|
| 105 | binst.SetEdgesASin(67, -0.005, 0.665);
|
|---|
| 106 |
|
|---|
| 107 | binse.Apply(fHEnergy);
|
|---|
| 108 |
|
|---|
| 109 | MH::SetBinning(&fHistSel, &binst, &binse);
|
|---|
| 110 | MH::SetBinning(&fHistAll, &binst, &binse);
|
|---|
| 111 |
|
|---|
| 112 | // For some unknown reasons this must be called after
|
|---|
| 113 | // the binning has been initialized at least once
|
|---|
| 114 | fHistSel.Sumw2();
|
|---|
| 115 | fHistAll.Sumw2();
|
|---|
| 116 | fHEnergy.Sumw2();
|
|---|
| 117 | }
|
|---|
| 118 |
|
|---|
| 119 | // --------------------------------------------------------------------------
|
|---|
| 120 | //
|
|---|
| 121 | // Return the Area defined by fMcAreaRadius
|
|---|
| 122 | //
|
|---|
| 123 | Double_t MHCollectionArea::GetCollectionAreaAbs() const
|
|---|
| 124 | {
|
|---|
| 125 | return TMath::Pi()*fMcAreaRadius*fMcAreaRadius;
|
|---|
| 126 | }
|
|---|
| 127 |
|
|---|
| 128 | // --------------------------------------------------------------------------
|
|---|
| 129 | //
|
|---|
| 130 | // Calculate the Efficiency (collection area) and set the 'ReadyToSave'
|
|---|
| 131 | // flag
|
|---|
| 132 | //
|
|---|
| 133 | void MHCollectionArea::CalcEfficiency()
|
|---|
| 134 | {
|
|---|
| 135 | TH1D *hsel = fHistSel.ProjectionY("Spy", -1, -1, "E");;
|
|---|
| 136 | TH1D *hall = fHistAll.ProjectionY("Apy", -1, -1, "E");
|
|---|
| 137 |
|
|---|
| 138 | //
|
|---|
| 139 | // Impact parameter range.
|
|---|
| 140 | //
|
|---|
| 141 | const Float_t totalarea = GetCollectionAreaAbs();//TMath::Pi() * (r2*r2 - r1*r1);
|
|---|
| 142 |
|
|---|
| 143 | // "b" option: calculate binomial errors
|
|---|
| 144 | // Do not use totalarea inside the binomial error calculation:
|
|---|
| 145 | // it is not a weight.
|
|---|
| 146 | fHEnergy.Divide(hsel, hall, 1, 1, "b");
|
|---|
| 147 | #if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
|
|---|
| 148 | MH::SetBinomialErrors(fHEnergy, *hsel, *hall);
|
|---|
| 149 | #endif
|
|---|
| 150 | if (fMcAreaRadius>0)
|
|---|
| 151 | fHEnergy.Scale(totalarea);
|
|---|
| 152 |
|
|---|
| 153 | delete hsel;
|
|---|
| 154 | delete hall;
|
|---|
| 155 | }
|
|---|
| 156 |
|
|---|
| 157 |
|
|---|
| 158 | Bool_t MHCollectionArea::SetupFill(const MParList *pl)
|
|---|
| 159 | {
|
|---|
| 160 | fHistSel.Reset();
|
|---|
| 161 | if (!fIsExtern)
|
|---|
| 162 | fHistAll.Reset();
|
|---|
| 163 |
|
|---|
| 164 | fHeader = (MMcRunHeader*)pl->FindObject("MMcRunHeader");
|
|---|
| 165 | if (!fHeader)
|
|---|
| 166 | {
|
|---|
| 167 | *fLog << err << "MMcRunHeader not found... abort." << endl;
|
|---|
| 168 | return kFALSE;
|
|---|
| 169 | }
|
|---|
| 170 |
|
|---|
| 171 | fMcEvt = (MMcEvt*)pl->FindObject("MMcEvt");
|
|---|
| 172 | if (!fMcEvt)
|
|---|
| 173 | {
|
|---|
| 174 | *fLog << err << "MMcEvt not found... abort." << endl;
|
|---|
| 175 | return kFALSE;
|
|---|
| 176 | }
|
|---|
| 177 | /*
|
|---|
| 178 | fEnergy = (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
|
|---|
| 179 | if (!fEnergy)
|
|---|
| 180 | {
|
|---|
| 181 | *fLog << err << "MEnergyEst [MParameterD] not found... abort." << endl;
|
|---|
| 182 | return kFALSE;
|
|---|
| 183 | }
|
|---|
| 184 | */
|
|---|
| 185 | MBinning binst, binse;
|
|---|
| 186 | binst.SetEdges(fHistAll, 'x');
|
|---|
| 187 | binse.SetEdges(fHistAll, 'y');
|
|---|
| 188 |
|
|---|
| 189 | //if (!fIsExtern)
|
|---|
| 190 | {
|
|---|
| 191 | MBinning *bins = (MBinning*)pl->FindObject("BinningTheta", "MBinning");
|
|---|
| 192 | if (bins)
|
|---|
| 193 | binst.SetEdges(*bins);
|
|---|
| 194 |
|
|---|
| 195 | bins = (MBinning*)pl->FindObject("BinningEnergyEst", "MBinning");
|
|---|
| 196 | if (bins)
|
|---|
| 197 | binse.SetEdges(*bins);
|
|---|
| 198 | }
|
|---|
| 199 |
|
|---|
| 200 | binse.Apply(fHEnergy);
|
|---|
| 201 |
|
|---|
| 202 | MH::SetBinning(&fHistSel, &binst, &binse);
|
|---|
| 203 | MH::SetBinning(&fHistAll, &binst, &binse);
|
|---|
| 204 |
|
|---|
| 205 | fMcAreaRadius = -1;
|
|---|
| 206 | fCorsikaVersion = 0;
|
|---|
| 207 |
|
|---|
| 208 | return kTRUE;
|
|---|
| 209 | }
|
|---|
| 210 |
|
|---|
| 211 | void MHCollectionArea::GetImpactMax()
|
|---|
| 212 | {
|
|---|
| 213 | if (fHeader->GetImpactMax()<=fMcAreaRadius*100 || fHeader->GetImpactMax()<0)
|
|---|
| 214 | return;
|
|---|
| 215 |
|
|---|
| 216 | fMcAreaRadius = 0.01*fHeader->GetImpactMax(); // cm->m
|
|---|
| 217 | *fLog << inf << "Maximum simulated impact: " << fMcAreaRadius << "m" << endl;
|
|---|
| 218 | }
|
|---|
| 219 |
|
|---|
| 220 | Bool_t MHCollectionArea::ReInit(MParList *plist)
|
|---|
| 221 | {
|
|---|
| 222 | GetImpactMax();
|
|---|
| 223 |
|
|---|
| 224 | if (fCorsikaVersion!=0 && fCorsikaVersion!=fHeader->GetCorsikaVersion())
|
|---|
| 225 | {
|
|---|
| 226 | *fLog << warn;
|
|---|
| 227 | *fLog << "Warning - Read files have different Corsika versions..." << endl;
|
|---|
| 228 | *fLog << " Last file=" << fCorsikaVersion << " New file=" << fHeader->GetCorsikaVersion() << endl;
|
|---|
| 229 | }
|
|---|
| 230 | fCorsikaVersion = fHeader->GetCorsikaVersion();
|
|---|
| 231 |
|
|---|
| 232 | if (fIsExtern)
|
|---|
| 233 | return kTRUE;
|
|---|
| 234 |
|
|---|
| 235 | fTotalNumSimulatedShowers += fHeader->GetNumSimulatedShowers();
|
|---|
| 236 | *fLog << inf << "Total Number of Simulated showers: " << fTotalNumSimulatedShowers << endl;
|
|---|
| 237 |
|
|---|
| 238 | fAllEvtsTriggered |= fHeader->GetAllEvtsTriggered();
|
|---|
| 239 | *fLog << inf << "Only triggered events avail: " << (fAllEvtsTriggered?"yes":"no") << endl;
|
|---|
| 240 |
|
|---|
| 241 | MMcCorsikaRunHeader *crh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
|
|---|
| 242 | if (!crh)
|
|---|
| 243 | {
|
|---|
| 244 | *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
|
|---|
| 245 | return kFALSE;
|
|---|
| 246 | }
|
|---|
| 247 |
|
|---|
| 248 | //
|
|---|
| 249 | // Calculate approximately the original number of events in each
|
|---|
| 250 | // energy bin:
|
|---|
| 251 | //
|
|---|
| 252 | const Float_t emin = crh->GetELowLim();
|
|---|
| 253 | const Float_t emax = crh->GetEUppLim();
|
|---|
| 254 | const Float_t expo = 1 + crh->GetSlopeSpec();
|
|---|
| 255 | const Float_t k = fHeader->GetNumSimulatedShowers() /
|
|---|
| 256 | (pow(emax,expo) - pow(emin,expo));
|
|---|
| 257 |
|
|---|
| 258 | const Int_t nbiny = fHistAll.GetNbinsY();
|
|---|
| 259 |
|
|---|
| 260 | TAxis &axe = *fHistAll.GetYaxis();
|
|---|
| 261 | for (Int_t i = 1; i <= nbiny; i++)
|
|---|
| 262 | {
|
|---|
| 263 | const Float_t e1 = axe.GetBinLowEdge(i);
|
|---|
| 264 | const Float_t e2 = axe.GetBinLowEdge(i+1);
|
|---|
| 265 |
|
|---|
| 266 | if (e1 < emin || e2 > emax)
|
|---|
| 267 | continue;
|
|---|
| 268 |
|
|---|
| 269 | const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
|
|---|
| 270 | //
|
|---|
| 271 | // We fill the i-th energy bin, with the total number of events
|
|---|
| 272 | // Second argument of Fill would be impact parameter of each
|
|---|
| 273 | // event, but we don't really need it for the collection area,
|
|---|
| 274 | // so we just put a dummy value (1.)
|
|---|
| 275 | //
|
|---|
| 276 |
|
|---|
| 277 | const Float_t energy = (e1+e2)/2.;
|
|---|
| 278 | for (int j=0; j<TMath::Nint(events); j++)
|
|---|
| 279 | fHistAll.Fill(0., energy);
|
|---|
| 280 | // you have MMcRunHeader.fShowerThetaMin and MMcRunHeader.fShowerThetaMax
|
|---|
| 281 | }
|
|---|
| 282 |
|
|---|
| 283 | return kTRUE;
|
|---|
| 284 | }
|
|---|
| 285 |
|
|---|
| 286 | void MHCollectionArea::Paint(Option_t *option)
|
|---|
| 287 | {
|
|---|
| 288 | if (TString(option)=="paint4" && fMcAreaRadius>0)
|
|---|
| 289 | {
|
|---|
| 290 | const TString txt = MString::Format("r_{max}=%.0fm --> A_{max}=%.0fm^{2}",
|
|---|
| 291 | fMcAreaRadius, GetCollectionAreaAbs());
|
|---|
| 292 |
|
|---|
| 293 | TLatex text(0.31, 0.95, txt);
|
|---|
| 294 | text.SetBit(TLatex::kTextNDC);
|
|---|
| 295 | text.SetTextSize(0.04);
|
|---|
| 296 | text.Paint();
|
|---|
| 297 | return;
|
|---|
| 298 | }
|
|---|
| 299 |
|
|---|
| 300 | TVirtualPad *pad;
|
|---|
| 301 |
|
|---|
| 302 | TPaveStats *st=0;
|
|---|
| 303 | for (int x=0; x<4; x++)
|
|---|
| 304 | {
|
|---|
| 305 | pad=gPad->GetPad(x+1);
|
|---|
| 306 | if (!pad || !(st = (TPaveStats*)pad->GetPrimitive("stats")))
|
|---|
| 307 | continue;
|
|---|
| 308 |
|
|---|
| 309 | if (st->GetOptStat()==11)
|
|---|
| 310 | continue;
|
|---|
| 311 |
|
|---|
| 312 | const Double_t y1 = st->GetY1NDC();
|
|---|
| 313 | const Double_t y2 = st->GetY2NDC();
|
|---|
| 314 | const Double_t x1 = st->GetX1NDC();
|
|---|
| 315 | const Double_t x2 = st->GetX2NDC();
|
|---|
| 316 |
|
|---|
| 317 | st->SetY1NDC((y2-y1)/3+y1);
|
|---|
| 318 | st->SetX1NDC((x2-x1)/3+x1);
|
|---|
| 319 | st->SetOptStat(11);
|
|---|
| 320 | }
|
|---|
| 321 |
|
|---|
| 322 | pad = gPad;
|
|---|
| 323 |
|
|---|
| 324 | TH1 *h1=0, *h2=0;
|
|---|
| 325 |
|
|---|
| 326 | pad->cd(1);
|
|---|
| 327 | if (gPad->FindObject("ProjSelX"))
|
|---|
| 328 | fHistSel.ProjectionX("ProjSelX", -1, -1, "E");
|
|---|
| 329 |
|
|---|
| 330 | pad->cd(2);
|
|---|
| 331 | if (gPad->FindObject("ProjAllY"))
|
|---|
| 332 | h1=fHistAll.ProjectionY("ProjAllY", -1, -1, "E");
|
|---|
| 333 | if (gPad->FindObject("ProjSelY"))
|
|---|
| 334 | h2=fHistSel.ProjectionY("ProjSelY", -1, -1, "E");
|
|---|
| 335 |
|
|---|
| 336 | if (h1 && h1->GetMaximum()>0)
|
|---|
| 337 | {
|
|---|
| 338 | gPad->SetLogx();
|
|---|
| 339 | gPad->SetLogy();
|
|---|
| 340 | }
|
|---|
| 341 |
|
|---|
| 342 | pad->cd(3);
|
|---|
| 343 | TH1 *h=dynamic_cast<TH1*>(gPad->FindObject("Efficiency"));
|
|---|
| 344 | if (h1 && h2 && h)
|
|---|
| 345 | {
|
|---|
| 346 | h->Divide(h2, h1, 1, 1, "b");
|
|---|
| 347 | #if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
|
|---|
| 348 | MH::SetBinomialErrors(*h, *h2, *h1);
|
|---|
| 349 | #endif
|
|---|
| 350 | h->SetMinimum(0);
|
|---|
| 351 | }
|
|---|
| 352 |
|
|---|
| 353 | pad->cd(4);
|
|---|
| 354 | CalcEfficiency();
|
|---|
| 355 | if (fHEnergy.GetMaximum()>0)
|
|---|
| 356 | {
|
|---|
| 357 | gPad->SetLogx();
|
|---|
| 358 | gPad->SetLogy();
|
|---|
| 359 | }
|
|---|
| 360 | }
|
|---|
| 361 |
|
|---|
| 362 | void MHCollectionArea::Draw(Option_t *option)
|
|---|
| 363 | {
|
|---|
| 364 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
|---|
| 365 |
|
|---|
| 366 | // Do the projection before painting the histograms into
|
|---|
| 367 | // the individual pads
|
|---|
| 368 | AppendPad();
|
|---|
| 369 |
|
|---|
| 370 | pad->SetBorderMode(0);
|
|---|
| 371 | pad->Divide(2,2);
|
|---|
| 372 |
|
|---|
| 373 | TH1 *h=0, *h1=0, *h2=0;
|
|---|
| 374 |
|
|---|
| 375 | if (fHistSel.GetNbinsX()>1)
|
|---|
| 376 | {
|
|---|
| 377 | pad->cd(1);
|
|---|
| 378 | gPad->SetBorderMode(0);
|
|---|
| 379 | gPad->SetGridx();
|
|---|
| 380 | gPad->SetGridy();
|
|---|
| 381 | /*
|
|---|
| 382 | h = fHistAll.ProjectionX("ProjAllX", -1, -1, "E");
|
|---|
| 383 | h->SetXTitle("\\Theta [\\circ]");
|
|---|
| 384 | h->SetDirectory(NULL);
|
|---|
| 385 | h->SetLineColor(kGreen);
|
|---|
| 386 | h->SetBit(kCanDelete);
|
|---|
| 387 | h->Draw();
|
|---|
| 388 | */
|
|---|
| 389 | h = fHistSel.ProjectionX("ProjSelX", -1, -1, "E");
|
|---|
| 390 | h->SetXTitle("\\Theta [\\circ]");
|
|---|
| 391 | h->SetDirectory(NULL);
|
|---|
| 392 | h->SetLineColor(kRed);
|
|---|
| 393 | h->SetBit(kCanDelete);
|
|---|
| 394 | h->Draw("hist"/*"same"*/);
|
|---|
| 395 | }
|
|---|
| 396 | else
|
|---|
| 397 | delete pad->GetPad(1);
|
|---|
| 398 |
|
|---|
| 399 | if (fHistSel.GetNbinsY()>1)
|
|---|
| 400 | {
|
|---|
| 401 | pad->cd(2);
|
|---|
| 402 | gPad->SetBorderMode(0);
|
|---|
| 403 | gPad->SetGridx();
|
|---|
| 404 | gPad->SetGridy();
|
|---|
| 405 |
|
|---|
| 406 | h1 = fHistAll.ProjectionY("ProjAllY", -1, -1, "E");
|
|---|
| 407 | h1->SetDirectory(NULL);
|
|---|
| 408 | h1->SetLineColor(kGreen);
|
|---|
| 409 | h1->SetXTitle("E [GeV]");
|
|---|
| 410 | h1->SetBit(kCanDelete);
|
|---|
| 411 | h1->Draw();
|
|---|
| 412 |
|
|---|
| 413 | h2 = fHistSel.ProjectionY("ProjSelY", -1, -1, "E");
|
|---|
| 414 | h2->SetDirectory(NULL);
|
|---|
| 415 | h2->SetLineColor(kRed);
|
|---|
| 416 | h2->SetBit(kCanDelete);
|
|---|
| 417 | h2->Draw("same");
|
|---|
| 418 | }
|
|---|
| 419 | else
|
|---|
| 420 | delete pad->GetPad(2);
|
|---|
| 421 |
|
|---|
| 422 | if (h1 && h2)
|
|---|
| 423 | {
|
|---|
| 424 | pad->cd(3);
|
|---|
| 425 | gPad->SetBorderMode(0);
|
|---|
| 426 | gPad->SetGridx();
|
|---|
| 427 | gPad->SetGridy();
|
|---|
| 428 | gPad->SetLogx();
|
|---|
| 429 | h = h2->DrawCopy();
|
|---|
| 430 | h->Divide(h2, h1, 1, 1, "b");
|
|---|
| 431 | #if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
|
|---|
| 432 | MH::SetBinomialErrors(*h, *h2, *h1);
|
|---|
| 433 | #endif
|
|---|
| 434 | h->SetNameTitle("Efficiency", "Efficiency");
|
|---|
| 435 | h->SetDirectory(NULL);
|
|---|
| 436 | //AppendPad("paint3");
|
|---|
| 437 | }
|
|---|
| 438 | else
|
|---|
| 439 | delete pad->GetPad(4);
|
|---|
| 440 |
|
|---|
| 441 | if (fHEnergy.GetNbinsX()>1)
|
|---|
| 442 | {
|
|---|
| 443 | pad->cd(4);
|
|---|
| 444 | gPad->SetBorderMode(0);
|
|---|
| 445 | gPad->SetGridx();
|
|---|
| 446 | gPad->SetGridy();
|
|---|
| 447 | fHEnergy.Draw();
|
|---|
| 448 | AppendPad("paint4");
|
|---|
| 449 | }
|
|---|
| 450 | else
|
|---|
| 451 | delete pad->GetPad(4);
|
|---|
| 452 | }
|
|---|
| 453 |
|
|---|
| 454 | Int_t MHCollectionArea::Fill(const MParContainer *par, const Stat_t weight)
|
|---|
| 455 | {
|
|---|
| 456 | // This is not perfect because it selects the maximum impact only
|
|---|
| 457 | // from the selected events. Hoever, it will get overwritten
|
|---|
| 458 | // in finalize anyway.
|
|---|
| 459 | const Double_t impact = fMcEvt->GetImpact()*0.01; // cm->m
|
|---|
| 460 | if (impact>fMcAreaRadius)
|
|---|
| 461 | fMcAreaRadius = impact;
|
|---|
| 462 |
|
|---|
| 463 | const Double_t energy = fMcEvt->GetEnergy();
|
|---|
| 464 | const Double_t theta = fMcEvt->GetTelescopeTheta()*TMath::RadToDeg();
|
|---|
| 465 |
|
|---|
| 466 | fHistSel.Fill(theta, energy, weight);
|
|---|
| 467 |
|
|---|
| 468 | return kTRUE;
|
|---|
| 469 | }
|
|---|
| 470 |
|
|---|
| 471 | Bool_t MHCollectionArea::Finalize()
|
|---|
| 472 | {
|
|---|
| 473 | GetImpactMax();
|
|---|
| 474 |
|
|---|
| 475 | *fLog << all << "Maximum simulated impact found: " << fMcAreaRadius << "m" << endl;
|
|---|
| 476 |
|
|---|
| 477 | CalcEfficiency();
|
|---|
| 478 |
|
|---|
| 479 | return kTRUE;
|
|---|
| 480 | }
|
|---|