| 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): A. Moralejo 3/2003 <mailto:moralejo@pd.infn.it>
|
|---|
| 19 | !
|
|---|
| 20 | ! Copyright: MAGIC Software Development, 2000-2003
|
|---|
| 21 | !
|
|---|
| 22 | !
|
|---|
| 23 | \* ======================================================================== */
|
|---|
| 24 |
|
|---|
| 25 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 26 | // //
|
|---|
| 27 | // MHMcCT1CollectionArea //
|
|---|
| 28 | // //
|
|---|
| 29 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 30 |
|
|---|
| 31 | #include "MHMcCT1CollectionArea.h"
|
|---|
| 32 |
|
|---|
| 33 | #include <TH2.h>
|
|---|
| 34 | #include <TCanvas.h>
|
|---|
| 35 |
|
|---|
| 36 | #include "MMcEvt.hxx"
|
|---|
| 37 | #include "MH.h"
|
|---|
| 38 | #include "MBinning.h"
|
|---|
| 39 | #include "MParList.h"
|
|---|
| 40 | #include "MLog.h"
|
|---|
| 41 | #include "MLogManip.h"
|
|---|
| 42 |
|
|---|
| 43 | ClassImp(MHMcCT1CollectionArea);
|
|---|
| 44 |
|
|---|
| 45 | using namespace std;
|
|---|
| 46 |
|
|---|
| 47 | // --------------------------------------------------------------------------
|
|---|
| 48 | //
|
|---|
| 49 | // Creates the three necessary histograms:
|
|---|
| 50 | // - selected showers (input)
|
|---|
| 51 | // - all showers (input)
|
|---|
| 52 | // - collection area (result)
|
|---|
| 53 | //
|
|---|
| 54 | MHMcCT1CollectionArea::MHMcCT1CollectionArea(const char *name, const char *title)
|
|---|
| 55 | {
|
|---|
| 56 | //
|
|---|
| 57 | // nbins, minEnergy, maxEnergy defaults:
|
|---|
| 58 | // we set the energy range from 100 Gev to 30000 GeV (in log, 3.5 orders
|
|---|
| 59 | // of magnitude) and for each order we take 10 subdivisions --> 35 xbins
|
|---|
| 60 | // we set the theta range from 12.5 to 48 deg, with 6 bins (the latter
|
|---|
| 61 | // choice has been done to make the bin centers as close as possible to
|
|---|
| 62 | // the actual zenith angles in the CT1 MC sample).
|
|---|
| 63 | //
|
|---|
| 64 |
|
|---|
| 65 | fName = name ? name : "MHMcCT1CollectionArea";
|
|---|
| 66 | fTitle = title ? title : "Collection Area vs. log10 Energy";
|
|---|
| 67 |
|
|---|
| 68 | fHistAll = new TH2D;
|
|---|
| 69 | fHistSel = new TH2D;
|
|---|
| 70 | fHistCol = new TH2D;
|
|---|
| 71 |
|
|---|
| 72 | fHistCol->SetName(fName);
|
|---|
| 73 | fHistAll->SetName("AllEvents");
|
|---|
| 74 | fHistSel->SetName("SelectedEvents");
|
|---|
| 75 |
|
|---|
| 76 | fHistCol->SetTitle(fTitle);
|
|---|
| 77 | fHistAll->SetTitle("All showers - Theta vs log10 Energy distribution");
|
|---|
| 78 | fHistSel->SetTitle("Selected showers - Theta vs log10 Energy distribution");
|
|---|
| 79 |
|
|---|
| 80 | fHistAll->SetDirectory(NULL);
|
|---|
| 81 | fHistSel->SetDirectory(NULL);
|
|---|
| 82 | fHistCol->SetDirectory(NULL);
|
|---|
| 83 |
|
|---|
| 84 | fHistAll->SetXTitle("log10 E [GeV]");
|
|---|
| 85 | fHistAll->SetYTitle("theta [deg]");
|
|---|
| 86 | fHistAll->SetZTitle("N");
|
|---|
| 87 |
|
|---|
| 88 | fHistSel->SetXTitle("log10 E [GeV]");
|
|---|
| 89 | fHistSel->SetYTitle("theta [deg]");
|
|---|
| 90 | fHistSel->SetZTitle("N");
|
|---|
| 91 |
|
|---|
| 92 | fHistCol->SetXTitle("log10 E [GeV]");
|
|---|
| 93 | fHistCol->SetYTitle("theta [deg]");
|
|---|
| 94 | fHistCol->SetZTitle("A [m^{2}]");
|
|---|
| 95 | }
|
|---|
| 96 |
|
|---|
| 97 | // --------------------------------------------------------------------------
|
|---|
| 98 | //
|
|---|
| 99 | // Delete the three histograms
|
|---|
| 100 | //
|
|---|
| 101 | MHMcCT1CollectionArea::~MHMcCT1CollectionArea()
|
|---|
| 102 | {
|
|---|
| 103 | delete fHistAll;
|
|---|
| 104 | delete fHistSel;
|
|---|
| 105 | delete fHistCol;
|
|---|
| 106 | }
|
|---|
| 107 |
|
|---|
| 108 | // --------------------------------------------------------------------------
|
|---|
| 109 | //
|
|---|
| 110 | // Set the binnings and prepare the filling of the histograms
|
|---|
| 111 | //
|
|---|
| 112 | Bool_t MHMcCT1CollectionArea::SetupFill(const MParList *plist)
|
|---|
| 113 | {
|
|---|
| 114 | const MBinning* binsenergy = (MBinning*)plist->FindObject("BinningE");
|
|---|
| 115 | const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
|
|---|
| 116 |
|
|---|
| 117 | if (!binsenergy || !binstheta)
|
|---|
| 118 | {
|
|---|
| 119 | *fLog << err << dbginf << "At least one MBinning not found... aborting.";
|
|---|
| 120 | *fLog << endl;
|
|---|
| 121 | return kFALSE;
|
|---|
| 122 | }
|
|---|
| 123 |
|
|---|
| 124 | SetBinning(fHistAll, binsenergy, binstheta);
|
|---|
| 125 | SetBinning(fHistSel, binsenergy, binstheta);
|
|---|
| 126 | SetBinning(fHistCol, binsenergy, binstheta);
|
|---|
| 127 |
|
|---|
| 128 | fHistAll->Sumw2();
|
|---|
| 129 | fHistSel->Sumw2();
|
|---|
| 130 | fHistCol->Sumw2();
|
|---|
| 131 |
|
|---|
| 132 | return kTRUE;
|
|---|
| 133 | }
|
|---|
| 134 |
|
|---|
| 135 |
|
|---|
| 136 | // --------------------------------------------------------------------------
|
|---|
| 137 | //
|
|---|
| 138 | // Fill data into the histogram which contains the selected showers
|
|---|
| 139 | //
|
|---|
| 140 | Bool_t MHMcCT1CollectionArea::Fill(const MParContainer *par, const Stat_t w)
|
|---|
| 141 | {
|
|---|
| 142 | MMcEvt &mcevt = *(MMcEvt*)par;
|
|---|
| 143 |
|
|---|
| 144 | fHistSel->Fill(log10(mcevt.GetEnergy()), kRad2Deg*mcevt.GetTelescopeTheta(), w);
|
|---|
| 145 | return kTRUE;
|
|---|
| 146 | }
|
|---|
| 147 |
|
|---|
| 148 | // --------------------------------------------------------------------------
|
|---|
| 149 | //
|
|---|
| 150 | // Draw the histogram with all showers
|
|---|
| 151 | //
|
|---|
| 152 | void MHMcCT1CollectionArea::DrawAll(Option_t* option)
|
|---|
| 153 | {
|
|---|
| 154 | if (!gPad)
|
|---|
| 155 | MH::MakeDefCanvas(fHistAll);
|
|---|
| 156 |
|
|---|
| 157 | fHistAll->Draw(option);
|
|---|
| 158 |
|
|---|
| 159 | gPad->Modified();
|
|---|
| 160 | gPad->Update();
|
|---|
| 161 | }
|
|---|
| 162 |
|
|---|
| 163 | // --------------------------------------------------------------------------
|
|---|
| 164 | //
|
|---|
| 165 | // Draw the histogram with the selected showers only.
|
|---|
| 166 | //
|
|---|
| 167 | void MHMcCT1CollectionArea::DrawSel(Option_t* option)
|
|---|
| 168 | {
|
|---|
| 169 | if (!gPad)
|
|---|
| 170 | MH::MakeDefCanvas(fHistSel);
|
|---|
| 171 |
|
|---|
| 172 | fHistSel->Draw(option);
|
|---|
| 173 |
|
|---|
| 174 | gPad->Modified();
|
|---|
| 175 | gPad->Update();
|
|---|
| 176 | }
|
|---|
| 177 |
|
|---|
| 178 | // --------------------------------------------------------------------------
|
|---|
| 179 | //
|
|---|
| 180 | // Creates a new canvas and draws the histogram into it.
|
|---|
| 181 | // Be careful: The histogram belongs to this object and won't get deleted
|
|---|
| 182 | // together with the canvas.
|
|---|
| 183 | //
|
|---|
| 184 | TObject *MHMcCT1CollectionArea::DrawClone(Option_t* option) const
|
|---|
| 185 | {
|
|---|
| 186 | TCanvas &c = *MakeDefCanvas("CollArea", "Collection area plots", 600, 600);
|
|---|
| 187 | c.Divide(2,2);
|
|---|
| 188 |
|
|---|
| 189 | //
|
|---|
| 190 | // This is necessary to get the expected behaviour of DrawClone
|
|---|
| 191 | //
|
|---|
| 192 | gROOT->SetSelectedPad(NULL);
|
|---|
| 193 |
|
|---|
| 194 | c.cd(1);
|
|---|
| 195 | fHistCol->SetDirectory(NULL);
|
|---|
| 196 | fHistCol->DrawCopy(option);
|
|---|
| 197 |
|
|---|
| 198 | c.cd(2);
|
|---|
| 199 | fHistSel->SetDirectory(NULL);
|
|---|
| 200 | fHistSel->DrawCopy(option);
|
|---|
| 201 |
|
|---|
| 202 | c.cd(3);
|
|---|
| 203 | fHistAll->SetDirectory(NULL);
|
|---|
| 204 | fHistAll->DrawCopy(option);
|
|---|
| 205 |
|
|---|
| 206 |
|
|---|
| 207 | c.Modified();
|
|---|
| 208 | c.Update();
|
|---|
| 209 |
|
|---|
| 210 | return &c;
|
|---|
| 211 | }
|
|---|
| 212 |
|
|---|
| 213 | void MHMcCT1CollectionArea::Draw(Option_t* option)
|
|---|
| 214 | {
|
|---|
| 215 | if (!gPad)
|
|---|
| 216 | MH::MakeDefCanvas(fHistCol);
|
|---|
| 217 |
|
|---|
| 218 | fHistCol->Draw(option);
|
|---|
| 219 |
|
|---|
| 220 | gPad->Modified();
|
|---|
| 221 | gPad->Update();
|
|---|
| 222 | }
|
|---|
| 223 |
|
|---|
| 224 | //
|
|---|
| 225 | // Calculate the Efficiency (collection area) for the CT1 MC sample
|
|---|
| 226 | // and set the 'ReadyToSave' flag
|
|---|
| 227 | //
|
|---|
| 228 | void MHMcCT1CollectionArea::CalcEfficiency()
|
|---|
| 229 | {
|
|---|
| 230 | //
|
|---|
| 231 | // Here we estimate the total number of showers in each energy bin
|
|---|
| 232 | // from the known the energy range and spectral index of the generated
|
|---|
| 233 | // showers. This procedure is intended for the CT1 MC files. The total
|
|---|
| 234 | // number of generated events, collection area, spectral index etc will be
|
|---|
| 235 | // set here by hand, so make sure that the MC sample you are using is the
|
|---|
| 236 | // right one (check all these quantities in your files and compare with
|
|---|
| 237 | // what is written below. In some theta bins, there are two different
|
|---|
| 238 | // productions, with different energy limits but with the same spectral
|
|---|
| 239 | // slope. We account for this when calculating the original number of
|
|---|
| 240 | // events in each energy bin.
|
|---|
| 241 | //
|
|---|
| 242 |
|
|---|
| 243 | for (Int_t thetabin = 1; thetabin <= fHistAll->GetNbinsY(); thetabin++)
|
|---|
| 244 | {
|
|---|
| 245 | // This theta is not exactly the one of the MC events, just about
|
|---|
| 246 | // the same:
|
|---|
| 247 | Float_t theta = fHistAll->GetYaxis()->GetBinCenter(thetabin);
|
|---|
| 248 |
|
|---|
| 249 | Float_t emin1, emax1, emin2, emax2;
|
|---|
| 250 | Float_t index, expo, k1, k2;
|
|---|
| 251 | Float_t numevts1, numevts2;
|
|---|
| 252 | Float_t r1, r2; // Impact parameter range (on ground).
|
|---|
| 253 |
|
|---|
| 254 | emin1 = 0; emax1 = 0; emin2 = 0; emax2 = 0;
|
|---|
| 255 | expo = 0.; k1 = 0.; k2 = 0.; r1 = 0.; r2 = 0.;
|
|---|
| 256 | numevts1 = 0; numevts2 = 0;
|
|---|
| 257 |
|
|---|
| 258 | if (theta > 14 && theta < 16) // 15 deg
|
|---|
| 259 | {
|
|---|
| 260 | r1 = 0.;
|
|---|
| 261 | r2 = 250.; //meters
|
|---|
| 262 | emin1 = 300.;
|
|---|
| 263 | emax1 = 400.; // Energies in GeV.
|
|---|
| 264 | emin2 = 400.;
|
|---|
| 265 | emax2 = 30000.;
|
|---|
| 266 | numevts1 = 4000.;
|
|---|
| 267 | numevts2 = 25740.;
|
|---|
| 268 | }
|
|---|
| 269 | else if (theta > 20 && theta < 21) // 20.5 deg
|
|---|
| 270 | {
|
|---|
| 271 | r1 = 0.;
|
|---|
| 272 | r2 = 263.; //meters
|
|---|
| 273 | emin1 = 300.;
|
|---|
| 274 | emax1 = 400.; // Energies in GeV.
|
|---|
| 275 | emin2 = 400.;
|
|---|
| 276 | emax2 = 30000.;
|
|---|
| 277 | numevts1 = 6611.;
|
|---|
| 278 | numevts2 = 24448.;
|
|---|
| 279 | }
|
|---|
| 280 | else if (theta > 26 && theta < 27) // 26.5 degrees
|
|---|
| 281 | {
|
|---|
| 282 | r1 = 0.;
|
|---|
| 283 | r2 = 290.; //meters
|
|---|
| 284 | emin1 = 300.;
|
|---|
| 285 | emax1 = 400.; // Energies in GeV.
|
|---|
| 286 | emax2 = emax1; emin2 = 400.;
|
|---|
| 287 | emax2 = 30000.;
|
|---|
| 288 | numevts1 = 4000.;
|
|---|
| 289 | numevts2 = 26316.;
|
|---|
| 290 | }
|
|---|
| 291 | else if (theta > 32 && theta < 33) // 32.5 degrees
|
|---|
| 292 | {
|
|---|
| 293 | r1 = 0.;
|
|---|
| 294 | r2 = 350.; //meters
|
|---|
| 295 | emin1 = 300.;
|
|---|
| 296 | emax1 = 30000.; // Energies in GeV.
|
|---|
| 297 | emax2 = emax1;
|
|---|
| 298 | numevts1 = 33646.;
|
|---|
| 299 | }
|
|---|
| 300 | else if (theta > 38 && theta < 39) // 38.75 degrees
|
|---|
| 301 | {
|
|---|
| 302 | r1 = 0.;
|
|---|
| 303 | r2 = 380.; //meters
|
|---|
| 304 | emin1 = 300.;
|
|---|
| 305 | emax1 = 30000.; // Energies in GeV.
|
|---|
| 306 | emax2 = emax1;
|
|---|
| 307 | numevts1 = 38415.;
|
|---|
| 308 | }
|
|---|
| 309 | else if (theta > 44 && theta < 46) // 45 degrees
|
|---|
| 310 | {
|
|---|
| 311 | r1 = 0.;
|
|---|
| 312 | r2 = 565.; //meters
|
|---|
| 313 | emin1 = 300.;
|
|---|
| 314 | emax1 = 50000.; // Energies in GeV.
|
|---|
| 315 | emax2 = emax1;
|
|---|
| 316 | numevts1 = 30197.;
|
|---|
| 317 | }
|
|---|
| 318 |
|
|---|
| 319 | index = 1.5; // Differential spectral Index.
|
|---|
| 320 | expo = 1.-index;
|
|---|
| 321 | k1 = numevts1 / (pow(emax1,expo) - pow(emin1,expo));
|
|---|
| 322 | k2 = numevts2 / (pow(emax2,expo) - pow(emin2,expo));
|
|---|
| 323 |
|
|---|
| 324 | for (Int_t i=1; i <= fHistAll->GetNbinsX(); i++)
|
|---|
| 325 | {
|
|---|
| 326 | const Float_t e1 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i));
|
|---|
| 327 | const Float_t e2 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i+1));
|
|---|
| 328 |
|
|---|
| 329 | if (e1 < emin1 || e2 > emax2)
|
|---|
| 330 | continue;
|
|---|
| 331 |
|
|---|
| 332 | Float_t events;
|
|---|
| 333 |
|
|---|
| 334 | if (e2 <= emax1)
|
|---|
| 335 | events = k1 * (pow(e2, expo) - pow(e1, expo));
|
|---|
| 336 | else if (e1 >= emin2)
|
|---|
| 337 | events = k2 * (pow(e2, expo) - pow(e1, expo));
|
|---|
| 338 | else
|
|---|
| 339 | events =
|
|---|
| 340 | k1 * (pow(emax1, expo) - pow(e1, expo))+
|
|---|
| 341 | k2 * (pow(e2, expo) - pow(emin2, expo));
|
|---|
| 342 |
|
|---|
| 343 | fHistAll->SetBinContent(i, thetabin, events);
|
|---|
| 344 | fHistAll->SetBinError(i, thetabin, sqrt(events));
|
|---|
| 345 | }
|
|---|
| 346 |
|
|---|
| 347 | // -----------------------------------------------------------
|
|---|
| 348 |
|
|---|
| 349 | const Float_t dr = TMath::Pi() * (r2*r2 - r1*r1);
|
|---|
| 350 |
|
|---|
| 351 | for (Int_t ix = 1; ix <= fHistAll->GetNbinsX(); ix++)
|
|---|
| 352 | {
|
|---|
| 353 | const Float_t Na = fHistAll->GetBinContent(ix,thetabin);
|
|---|
| 354 |
|
|---|
| 355 | if (Na <= 0)
|
|---|
| 356 | {
|
|---|
| 357 | //
|
|---|
| 358 | // If energy is large, this case means that no or very few events
|
|---|
| 359 | // were generated at this energy bin. In this case we assign it
|
|---|
| 360 | // the effective area of the bin below it in energy. If energy is
|
|---|
| 361 | // below 1E4, it means that no events triggered -> eff area = 0
|
|---|
| 362 | //
|
|---|
| 363 |
|
|---|
| 364 | if (fHistSel->GetXaxis()->GetBinLowEdge(ix) > 4.)
|
|---|
| 365 | {
|
|---|
| 366 | fHistCol->SetBinContent(ix, thetabin, fHistCol->GetBinContent(ix-1, thetabin));
|
|---|
| 367 | fHistCol->SetBinError(ix, thetabin, fHistCol->GetBinError(ix-1, thetabin));
|
|---|
| 368 | }
|
|---|
| 369 | continue;
|
|---|
| 370 | }
|
|---|
| 371 |
|
|---|
| 372 | const Float_t Ns = fHistSel->GetBinContent(ix,thetabin);
|
|---|
| 373 |
|
|---|
| 374 | // Since Na is an estimate of the total number of showers generated
|
|---|
| 375 | // in the energy bin, it may happen that Ns (triggered showers) is
|
|---|
| 376 | // larger than Na. In that case, the bin is skipped:
|
|---|
| 377 |
|
|---|
| 378 | if (Na < Ns)
|
|---|
| 379 | continue;
|
|---|
| 380 |
|
|---|
| 381 | const Double_t eff = Ns/Na;
|
|---|
| 382 | const Double_t efferr = sqrt((1.-eff)*Ns)/Na;
|
|---|
| 383 |
|
|---|
| 384 |
|
|---|
| 385 | const Float_t area = dr * cos(theta*TMath::Pi()/180.);
|
|---|
| 386 |
|
|---|
| 387 | fHistCol->SetBinContent(ix, thetabin, eff*area);
|
|---|
| 388 | fHistCol->SetBinError(ix, thetabin, efferr*area);
|
|---|
| 389 |
|
|---|
| 390 | }
|
|---|
| 391 | }
|
|---|
| 392 |
|
|---|
| 393 | SetReadyToSave();
|
|---|
| 394 | }
|
|---|