| 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, 5/2002 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 19 | ! Author(s): Abelardo Moralejo <mailto:moralejo@pd.infn.it>
|
|---|
| 20 | !
|
|---|
| 21 | ! Copyright: MAGIC Software Development, 2000-2002
|
|---|
| 22 | !
|
|---|
| 23 | !
|
|---|
| 24 | \* ======================================================================== */
|
|---|
| 25 |
|
|---|
| 26 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 27 | //
|
|---|
| 28 | // MHHadroness
|
|---|
| 29 | //
|
|---|
| 30 | // This is histogram is a way to evaluate the quality of a gamma/hadron
|
|---|
| 31 | // seperation method. It is filled from a MHadroness container, which
|
|---|
| 32 | // stores a hadroness for the current event. The Value must be in the
|
|---|
| 33 | // range [0,1]. To fill the histograms correctly the information
|
|---|
| 34 | // whether it is a gamma or hadron (not a gamma) must be available from
|
|---|
| 35 | // a MMcEvt container.
|
|---|
| 36 | //
|
|---|
| 37 | // In the constructor you can change the number of used bns for the
|
|---|
| 38 | // evaluation.
|
|---|
| 39 | //
|
|---|
| 40 | // The meaning of the histograms (Draw, DrawClone) are the following:
|
|---|
| 41 | // * Upper Left Corner:
|
|---|
| 42 | // - black: histogram of all hadronesses for gammas
|
|---|
| 43 | // - red: histogram of all hadronesses for non gammas
|
|---|
| 44 | // * Upper Right Corner:
|
|---|
| 45 | // - black: acceptance of gammas vs. the hadroness
|
|---|
| 46 | // - red: acceptance of non gammas vs. the hadroness
|
|---|
| 47 | // - blue: 2D distance of (acceptance_hadrons, acceptances_gammas)
|
|---|
| 48 | // to optimum (0, 1)
|
|---|
| 49 | // 1-sqrt(Ag*Ag + (1-Ah)*(1-Ah))
|
|---|
| 50 | // * Bottom Left Corner:
|
|---|
| 51 | // Naive quality factor: Ag/sqrt(Ah)
|
|---|
| 52 | // * Bottom Right Corner:
|
|---|
| 53 | // - black: Acceprtance Gammas vs. Acceptance Hadrons
|
|---|
| 54 | // - blue cross: minimum of distance to (0, 1)
|
|---|
| 55 | //
|
|---|
| 56 | ////////////////////////////////////////////////////////////////////////////
|
|---|
| 57 | #include "MHHadroness.h"
|
|---|
| 58 |
|
|---|
| 59 | #include <TPad.h>
|
|---|
| 60 | #include <TGraph.h>
|
|---|
| 61 | #include <TStyle.h>
|
|---|
| 62 | #include <TCanvas.h>
|
|---|
| 63 | #include <TMarker.h>
|
|---|
| 64 |
|
|---|
| 65 | #include "MLog.h"
|
|---|
| 66 | #include "MLogManip.h"
|
|---|
| 67 |
|
|---|
| 68 | #include "MParList.h"
|
|---|
| 69 | #include "MBinning.h"
|
|---|
| 70 | #include "MHadroness.h"
|
|---|
| 71 |
|
|---|
| 72 | #include "MMcEvt.hxx"
|
|---|
| 73 |
|
|---|
| 74 | ClassImp(MHHadroness);
|
|---|
| 75 |
|
|---|
| 76 | // --------------------------------------------------------------------------
|
|---|
| 77 | //
|
|---|
| 78 | // Setup histograms, nbins is the number of bins used for the evaluation.
|
|---|
| 79 | // The default is 100 bins.
|
|---|
| 80 | //
|
|---|
| 81 | MHHadroness::MHHadroness(Int_t nbins, const char *name, const char *title)
|
|---|
| 82 | {
|
|---|
| 83 | //
|
|---|
| 84 | // set the name and title of this object
|
|---|
| 85 | //
|
|---|
| 86 | fName = name ? name : "MHHadroness";
|
|---|
| 87 | fTitle = title ? title : "Gamma/Hadron Separation Quality Histograms";
|
|---|
| 88 |
|
|---|
| 89 | fGraph = new TGraph;
|
|---|
| 90 | fGraph->SetTitle("Acceptance Gammas vs. Hadrons");
|
|---|
| 91 | fGraph->SetMaximum(1);
|
|---|
| 92 |
|
|---|
| 93 | fGhness = new TH1D("Ghness", "Hadronness", nbins, 0, 1);
|
|---|
| 94 | fPhness = new TH1D("Phness", "Hadronness", nbins, 0, 1);
|
|---|
| 95 | fGhness->SetXTitle("Hadroness");
|
|---|
| 96 | fPhness->SetXTitle("Hadroness");
|
|---|
| 97 | fGhness->SetYTitle("Counts");
|
|---|
| 98 | fPhness->SetYTitle("Counts");
|
|---|
| 99 |
|
|---|
| 100 | fIntGhness = new TH1D("AccGammas", "Acceptance", nbins, 0, 1);
|
|---|
| 101 | fIntPhness = new TH1D("AccHadrons", "Acceptance", nbins, 0, 1);
|
|---|
| 102 | fIntGhness->SetXTitle("Hadroness");
|
|---|
| 103 | fIntPhness->SetXTitle("Hadroness");
|
|---|
| 104 | fIntGhness->SetYTitle("Acceptance");
|
|---|
| 105 | fIntPhness->SetYTitle("Acceptance");
|
|---|
| 106 |
|
|---|
| 107 | fQfac = new TH1D("Qfac", "Naive Quality factor", nbins, 0, 1);
|
|---|
| 108 | fQfac->SetXTitle("Hadroness");
|
|---|
| 109 | fQfac->SetYTitle("Quality");
|
|---|
| 110 |
|
|---|
| 111 | fMinDist = new TH1D("MinDist", "Minimum Dist to (0, 1)", nbins, 0, 1);
|
|---|
| 112 | fMinDist->SetXTitle("Hadroness");
|
|---|
| 113 | fMinDist->SetYTitle("Distance");
|
|---|
| 114 |
|
|---|
| 115 | fQfac->SetDirectory(NULL);
|
|---|
| 116 | fGhness->SetDirectory(NULL);
|
|---|
| 117 | fPhness->SetDirectory(NULL);
|
|---|
| 118 | fMinDist->SetDirectory(NULL);
|
|---|
| 119 | fIntGhness->SetDirectory(NULL);
|
|---|
| 120 | fIntPhness->SetDirectory(NULL);
|
|---|
| 121 | }
|
|---|
| 122 |
|
|---|
| 123 | // --------------------------------------------------------------------------
|
|---|
| 124 | //
|
|---|
| 125 | // Delete the histograms.
|
|---|
| 126 | //
|
|---|
| 127 | MHHadroness::~MHHadroness()
|
|---|
| 128 | {
|
|---|
| 129 | delete fGhness;
|
|---|
| 130 | delete fIntGhness;
|
|---|
| 131 | delete fPhness;
|
|---|
| 132 | delete fIntPhness;
|
|---|
| 133 | delete fQfac;
|
|---|
| 134 | delete fMinDist;
|
|---|
| 135 | delete fGraph;
|
|---|
| 136 | }
|
|---|
| 137 |
|
|---|
| 138 | // --------------------------------------------------------------------------
|
|---|
| 139 | //
|
|---|
| 140 | // Setup Filling of the histograms. It needs:
|
|---|
| 141 | // MMcEvt and MHadroness
|
|---|
| 142 | //
|
|---|
| 143 | Bool_t MHHadroness::SetupFill(const MParList *plist)
|
|---|
| 144 | {
|
|---|
| 145 | fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
|
|---|
| 146 | if (!fMcEvt)
|
|---|
| 147 | {
|
|---|
| 148 | *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
|
|---|
| 149 | return kFALSE;
|
|---|
| 150 | }
|
|---|
| 151 |
|
|---|
| 152 | fHadroness = (MHadroness*)plist->FindObject("MHadroness");
|
|---|
| 153 | if (!fHadroness)
|
|---|
| 154 | {
|
|---|
| 155 | *fLog << err << dbginf << "MHadroness not found... aborting." << endl;
|
|---|
| 156 | return kFALSE;
|
|---|
| 157 | }
|
|---|
| 158 |
|
|---|
| 159 | /*
|
|---|
| 160 | MBinning* bins = (MBinning*)plist->FindObject("BinningHadroness");
|
|---|
| 161 | if (!bins)
|
|---|
| 162 | {
|
|---|
| 163 | *fLog << err << dbginf << "BinningHadroness [MBinning] not found... aborting." << endl;
|
|---|
| 164 | return kFALSE;
|
|---|
| 165 | }
|
|---|
| 166 |
|
|---|
| 167 | SetBinning(&fHist, binsalpha, binsenergy, binstheta);
|
|---|
| 168 |
|
|---|
| 169 | fHist.Sumw2();
|
|---|
| 170 | */
|
|---|
| 171 |
|
|---|
| 172 | return kTRUE;
|
|---|
| 173 | }
|
|---|
| 174 |
|
|---|
| 175 | // --------------------------------------------------------------------------
|
|---|
| 176 | //
|
|---|
| 177 | // Fill the Hadroness from a MHadroness container into the corresponding
|
|---|
| 178 | // histogram dependant on the particle id.
|
|---|
| 179 | //
|
|---|
| 180 | Bool_t MHHadroness::Fill(const MParContainer *par)
|
|---|
| 181 | {
|
|---|
| 182 | if (fMcEvt->GetPartId()==kGAMMA)
|
|---|
| 183 | fGhness->Fill(fHadroness->GetHadroness());
|
|---|
| 184 | else
|
|---|
| 185 | fPhness->Fill(fHadroness->GetHadroness());
|
|---|
| 186 |
|
|---|
| 187 | return kTRUE;
|
|---|
| 188 | }
|
|---|
| 189 |
|
|---|
| 190 | // --------------------------------------------------------------------------
|
|---|
| 191 | //
|
|---|
| 192 | // Finalize the histograms:
|
|---|
| 193 | // - integrate the hadroness histograms --> acceptance
|
|---|
| 194 | // - fill the Minimum Distance histogram (formular see class description)
|
|---|
| 195 | // - fill the Quality histogram (formular see class description)
|
|---|
| 196 | //
|
|---|
| 197 | Bool_t MHHadroness::Finalize()
|
|---|
| 198 | {
|
|---|
| 199 | Int_t n = fGhness->GetNbinsX();
|
|---|
| 200 |
|
|---|
| 201 | fGraph->Set(n);
|
|---|
| 202 |
|
|---|
| 203 | const Stat_t sumg = fGhness->Integral(1, n+1);
|
|---|
| 204 | const Stat_t sump = fPhness->Integral(1, n+1);
|
|---|
| 205 |
|
|---|
| 206 | for (Int_t i=1; i<=n; i++)
|
|---|
| 207 | {
|
|---|
| 208 | const Stat_t ip = fPhness->Integral(1, i)/sump;
|
|---|
| 209 | const Stat_t ig = fGhness->Integral(1, i)/sumg;
|
|---|
| 210 |
|
|---|
| 211 | fIntPhness->SetBinContent(i, ip);
|
|---|
| 212 | fIntGhness->SetBinContent(i, ig);
|
|---|
| 213 |
|
|---|
| 214 | fGraph->SetPoint(i, ip, ig);
|
|---|
| 215 |
|
|---|
| 216 | if (ip<=0)
|
|---|
| 217 | continue;
|
|---|
| 218 |
|
|---|
| 219 | fMinDist->SetBinContent(i, 1-sqrt(ip*ip + (ig-1)*(ig-1)));
|
|---|
| 220 | fQfac->SetBinContent(i, ig/sqrt(ip));
|
|---|
| 221 | }
|
|---|
| 222 |
|
|---|
| 223 | return kTRUE;
|
|---|
| 224 | }
|
|---|
| 225 |
|
|---|
| 226 | // --------------------------------------------------------------------------
|
|---|
| 227 | //
|
|---|
| 228 | // Search the corresponding points for the given hadron acceptance (acchad)
|
|---|
| 229 | // and interpolate the to pointsd (linear)
|
|---|
| 230 | //
|
|---|
| 231 | Double_t MHHadroness::GetGammaAcceptance(Double_t acchad) const
|
|---|
| 232 | {
|
|---|
| 233 | const Int_t n = fGraph->GetN();
|
|---|
| 234 | const Double_t *x = fGraph->GetX();
|
|---|
| 235 | const Double_t *y = fGraph->GetY();
|
|---|
| 236 |
|
|---|
| 237 | Int_t i = 0;
|
|---|
| 238 | while (i<n && x[i]<acchad)
|
|---|
| 239 | i++;
|
|---|
| 240 |
|
|---|
| 241 | if (i==0 || i==n)
|
|---|
| 242 | return 0;
|
|---|
| 243 |
|
|---|
| 244 | if (i==n-1)
|
|---|
| 245 | i--;
|
|---|
| 246 |
|
|---|
| 247 | const Double_t x1 = x[i-1];
|
|---|
| 248 | const Double_t y1 = y[i-1];
|
|---|
| 249 |
|
|---|
| 250 | const Double_t x2 = x[i];
|
|---|
| 251 | const Double_t y2 = y[i];
|
|---|
| 252 |
|
|---|
| 253 | return (y2-y1)/(x2-x1) * (acchad-x2) + y2;
|
|---|
| 254 | }
|
|---|
| 255 |
|
|---|
| 256 | // --------------------------------------------------------------------------
|
|---|
| 257 | //
|
|---|
| 258 | // Print the corresponding Gammas Acceptance for a hadron acceptance of
|
|---|
| 259 | // 10%, 20%, 30%, 40% and 50%. Also the minimum distance to the optimum
|
|---|
| 260 | // acceptance and the corresponding acceptances and hadroness value is
|
|---|
| 261 | // printed, together with the maximum Q-factor.
|
|---|
| 262 | //
|
|---|
| 263 | void MHHadroness::Print(Option_t *) const
|
|---|
| 264 | {
|
|---|
| 265 | *fLog << all;
|
|---|
| 266 | *fLog << "Hadroness histograms:" << endl;
|
|---|
| 267 | *fLog << "---------------------" << endl;
|
|---|
| 268 | *fLog << "Acc Gammas at 1% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.01)*100) << "%" << endl;
|
|---|
| 269 | *fLog << "Acc Gammas at 2% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.02)*100) << "%" << endl;
|
|---|
| 270 | *fLog << "Acc Gammas at 5% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.05)*100) << "%" << endl;
|
|---|
| 271 | *fLog << "Acc Gammas at 10% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.1)*100) << "%" << endl;
|
|---|
| 272 | *fLog << "Acc Gammas at 20% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.2)*100) << "%" << endl;
|
|---|
| 273 | *fLog << "Acc Gammas at 30% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.3)*100) << "%" << endl;
|
|---|
| 274 | *fLog << "Acc Gammas at 40% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.4)*100) << "%" << endl;
|
|---|
| 275 | *fLog << "Acc Gammas at 50% Hadron-acc: " << Form("%3.0f", GetGammaAcceptance(0.5)*100) << "%" << endl;
|
|---|
| 276 | const Int_t h = fMinDist->GetMaximumBin();
|
|---|
| 277 | *fLog << "Minimum Distance to (0, 1): " << Form("%.2f", 1-fMinDist->GetMaximum()) << " @ H=" << fMinDist->GetBinCenter(h) << endl;
|
|---|
| 278 | *fLog << " Acc Gammas = " << Form("%3.0f", fIntGhness->GetBinContent(h)*100) << "%, ";
|
|---|
| 279 | *fLog << "Acc Hadrons = " << Form("%3.0f", fIntPhness->GetBinContent(h)*100) << "%" << endl;
|
|---|
| 280 | const Int_t q = GetQfac()->GetMaximumBin();
|
|---|
| 281 | *fLog << "Maximum Q-Factor: " << GetQfac()->GetMaximum() << " @ H=";
|
|---|
| 282 | *fLog << GetQfac()->GetBinCenter(q) << endl;;
|
|---|
| 283 | *fLog << " Acc Gammas = " << Form("%3.0f", fIntGhness->GetBinContent(q)*100) << "%, ";
|
|---|
| 284 | *fLog << "Acc Hadrons = " << Form("%3.0f", fIntPhness->GetBinContent(q)*100) << "%" << endl;
|
|---|
| 285 | *fLog << endl;
|
|---|
| 286 | }
|
|---|
| 287 |
|
|---|
| 288 | // --------------------------------------------------------------------------
|
|---|
| 289 | //
|
|---|
| 290 | // Draw clone of all histograms. (For the Meaning see class description)
|
|---|
| 291 | //
|
|---|
| 292 | TObject *MHHadroness::DrawClone(Option_t *opt) const
|
|---|
| 293 | {
|
|---|
| 294 | TCanvas &c = *MakeDefCanvas("Hadroness", fTitle);
|
|---|
| 295 | c.Divide(2, 2);
|
|---|
| 296 |
|
|---|
| 297 | gROOT->SetSelectedPad(NULL);
|
|---|
| 298 |
|
|---|
| 299 | c.cd(1);
|
|---|
| 300 | gStyle->SetOptStat(10);
|
|---|
| 301 | Getghness()->DrawCopy();
|
|---|
| 302 | Getphness()->SetLineColor(kRed);
|
|---|
| 303 | Getphness()->DrawCopy("same");
|
|---|
| 304 | c.cd(4);
|
|---|
| 305 | TGraph &g = (TGraph&)*fGraph->DrawClone("AC");
|
|---|
| 306 | g.SetBit(kCanDelete);
|
|---|
| 307 | gPad->Modified();
|
|---|
| 308 | gPad->Update();
|
|---|
| 309 | if (g.GetHistogram())
|
|---|
| 310 | {
|
|---|
| 311 | g.GetXaxis()->SetRangeUser(0, 1);
|
|---|
| 312 | g.GetXaxis()->SetTitle("Acceptance Hadrons");
|
|---|
| 313 | g.GetYaxis()->SetTitle("Acceptance Gammas");
|
|---|
| 314 | g.SetMarkerStyle(kFullDotSmall);
|
|---|
| 315 | g.Draw("P");
|
|---|
| 316 |
|
|---|
| 317 | gPad->Modified();
|
|---|
| 318 | gPad->Update();
|
|---|
| 319 | }
|
|---|
| 320 | const Int_t h = fMinDist->GetMaximumBin();
|
|---|
| 321 | TMarker *m = new TMarker(fIntPhness->GetBinContent(h),
|
|---|
| 322 | fIntGhness->GetBinContent(h), kStar);
|
|---|
| 323 | m->SetMarkerColor(kBlue);
|
|---|
| 324 | m->SetBit(kCanDelete);
|
|---|
| 325 | m->Draw();
|
|---|
| 326 | c.cd(2);
|
|---|
| 327 | gStyle->SetOptStat(0);
|
|---|
| 328 | Getighness()->DrawCopy();
|
|---|
| 329 | Getiphness()->SetLineColor(kRed);
|
|---|
| 330 | Getiphness()->DrawCopy("same");
|
|---|
| 331 | fMinDist->SetLineColor(kBlue);
|
|---|
| 332 | fMinDist->DrawCopy("same");
|
|---|
| 333 | c.cd(3);
|
|---|
| 334 | GetQfac()->DrawCopy();
|
|---|
| 335 |
|
|---|
| 336 | return &c;
|
|---|
| 337 | }
|
|---|
| 338 |
|
|---|
| 339 | // --------------------------------------------------------------------------
|
|---|
| 340 | //
|
|---|
| 341 | // Draw all histograms. (For the Meaning see class description)
|
|---|
| 342 | //
|
|---|
| 343 | void MHHadroness::Draw(Option_t *)
|
|---|
| 344 | {
|
|---|
| 345 | if (!gPad)
|
|---|
| 346 | MakeDefCanvas("Hadroness", fTitle);
|
|---|
| 347 |
|
|---|
| 348 | gPad->Divide(2, 2);
|
|---|
| 349 |
|
|---|
| 350 | gPad->cd(1);
|
|---|
| 351 | gStyle->SetOptStat(10);
|
|---|
| 352 | Getghness()->Draw();
|
|---|
| 353 | Getphness()->SetLineColor(kRed);
|
|---|
| 354 | Getphness()->Draw("same");
|
|---|
| 355 | gPad->cd(4);
|
|---|
| 356 | fGraph->Draw("AC");
|
|---|
| 357 | gPad->Modified();
|
|---|
| 358 | gPad->Update();
|
|---|
| 359 | if (fGraph->GetHistogram())
|
|---|
| 360 | {
|
|---|
| 361 | fGraph->GetXaxis()->SetRangeUser(0, 1);
|
|---|
| 362 | fGraph->GetXaxis()->SetTitle("Acceptance Hadrons");
|
|---|
| 363 | fGraph->GetYaxis()->SetTitle("Acceptance Gammas");
|
|---|
| 364 | fGraph->SetMarkerStyle(kFullDotSmall);
|
|---|
| 365 | fGraph->Draw("P");
|
|---|
| 366 | gPad->Modified();
|
|---|
| 367 | gPad->Update();
|
|---|
| 368 | }
|
|---|
| 369 | const Int_t h = fMinDist->GetMaximumBin();
|
|---|
| 370 | TMarker *m = new TMarker(fIntPhness->GetBinContent(h),
|
|---|
| 371 | fIntGhness->GetBinContent(h), kStar);
|
|---|
| 372 | m->SetMarkerColor(kBlue);
|
|---|
| 373 | m->SetBit(kCanDelete);
|
|---|
| 374 | m->Draw();
|
|---|
| 375 | gPad->cd(2);
|
|---|
| 376 | gStyle->SetOptStat(0);
|
|---|
| 377 | Getighness()->Draw();
|
|---|
| 378 | Getiphness()->SetLineColor(kRed);
|
|---|
| 379 | Getiphness()->Draw("same");
|
|---|
| 380 | fMinDist->SetLineColor(kBlue);
|
|---|
| 381 | fMinDist->Draw("same");
|
|---|
| 382 | gPad->cd(3);
|
|---|
| 383 | GetQfac()->Draw();
|
|---|
| 384 | }
|
|---|