source: trunk/MagicSoft/Mars/mhist/MHHadronness.cc@ 2194

Last change on this file since 2194 was 2180, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 13.8 KB
Line 
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-2003
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MHHadronness
29//
30// This is histogram is a way to evaluate the quality of a gamma/hadron
31// seperation method. It is filled from a MHadronness 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 (Ag) vs. the hadroness
46// - red: acceptance of non gammas (Ah) vs. the hadroness
47// * Bottom Left Corner:
48// Naive quality factor: Ag/sqrt(Ah)
49// * Bottom Right Corner:
50// - black: Acceprtance Gammas vs. Acceptance Hadrons
51// - blue cross: minimum of distance to (0, 1)
52//
53// As a default MHHadronness searches for the container "MHadronness".
54// This can be overwritten by a different pointer specified in the
55// Fill function (No type check is done!) Use a different name in
56// MFillH.
57//
58////////////////////////////////////////////////////////////////////////////
59#include "MHHadronness.h"
60
61#include <TPad.h>
62#include <TGraph.h>
63#include <TStyle.h>
64#include <TCanvas.h>
65#include <TMarker.h>
66
67#include "MLog.h"
68#include "MLogManip.h"
69
70#include "MParList.h"
71#include "MBinning.h"
72#include "MHadronness.h"
73
74#include "MMcEvt.hxx"
75
76ClassImp(MHHadronness);
77
78using namespace std;
79
80// --------------------------------------------------------------------------
81//
82// Setup histograms, nbins is the number of bins used for the evaluation.
83// The default is 100 bins.
84//
85MHHadronness::MHHadronness(Int_t nbins, const char *name, const char *title)
86{
87 //
88 // set the name and title of this object
89 //
90 fName = name ? name : "MHHadronness";
91 fTitle = title ? title : "Gamma/Hadron Separation Quality Histograms";
92
93 fGraph = new TGraph;
94 fGraph->SetTitle("Acceptance Gammas vs. Hadrons");
95 fGraph->SetMarkerStyle(kFullDotSmall);
96
97 fGhness = new TH1D("Ghness", "H. Gammas", nbins, 0, 1);
98 fPhness = new TH1D("Phness", "H. Hadrons", nbins, 0, 1);
99 fGhness->SetXTitle("Hadronness");
100 fPhness->SetXTitle("Hadronness");
101 fGhness->SetYTitle("Counts");
102 fPhness->SetYTitle("Counts");
103 fPhness->SetLineColor(kRed);
104 fGhness->SetDirectory(NULL);
105 fPhness->SetDirectory(NULL);
106
107 fIntGhness = new TH1D("AccGammas", "Acceptance", nbins, 0, 1);
108 fIntPhness = new TH1D("AccHadrons", "Acceptance", nbins, 0, 1);
109 fIntGhness->SetXTitle("Hadronness");
110 fIntPhness->SetXTitle("Hadronness");
111 fIntGhness->SetYTitle("Acceptance");
112 fIntPhness->SetYTitle("Acceptance");
113 fIntGhness->SetMaximum(1.1);
114 fIntPhness->SetMaximum(1.1);
115 fIntGhness->SetDirectory(NULL);
116 fIntPhness->SetDirectory(NULL);
117 fIntPhness->SetLineColor(kRed);
118 fIntGhness->SetBit(TH1::kNoStats);
119 fIntPhness->SetBit(TH1::kNoStats);
120
121 fQfac = new TGraph;
122 fQfac->SetTitle(" Naive Quality factor ");
123 fQfac->SetMarkerStyle(kFullDotSmall);
124}
125
126// --------------------------------------------------------------------------
127//
128// Delete the histograms.
129//
130MHHadronness::~MHHadronness()
131{
132 delete fGhness;
133 delete fIntGhness;
134 delete fPhness;
135 delete fIntPhness;
136 delete fQfac;
137 delete fGraph;
138}
139
140// --------------------------------------------------------------------------
141//
142// Setup Filling of the histograms. It needs:
143// MMcEvt and MHadronness
144//
145Bool_t MHHadronness::SetupFill(const MParList *plist)
146{
147 fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
148 if (!fMcEvt)
149 {
150 *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
151 return kFALSE;
152 }
153
154 fHadronness = (MHadronness*)plist->FindObject("MHadronness");
155
156 /*
157 MBinning* bins = (MBinning*)plist->FindObject("BinningHadronness");
158 if (!bins)
159 {
160 *fLog << err << dbginf << "BinningHadronness [MBinning] not found... aborting." << endl;
161 return kFALSE;
162 }
163
164 SetBinning(&fHist, binsalpha, binsenergy, binstheta);
165
166 fHist.Sumw2();
167 */
168
169 return kTRUE;
170}
171
172// --------------------------------------------------------------------------
173//
174// Fill the Hadronness from a MHadronness container into the corresponding
175// histogram dependant on the particle id.
176//
177// Every particle Id different than kGAMMA is considered a hadron.
178//
179// If you call Fill with a pointer (eg. from MFillH) this container is
180// used as a hadronness (Warning: No type check is done!) otherwise
181// the default is used ("MHadronness")
182//
183// Sometimes a distance is calculated as NaN (not a number). Such events
184// are skipped at the moment.
185//
186Bool_t MHHadronness::Fill(const MParContainer *par, const Stat_t w)
187{
188 // Preliminary Workaround: FIXME!
189 if (!par && !fHadronness)
190 {
191 *fLog << err << "MHHadronness::Fill: No MHadronness container specified!" << endl;
192 return kFALSE;
193 }
194
195 const MHadronness &had = par ? *(MHadronness*)par : *fHadronness;
196
197 const Double_t h = had.GetHadronness();
198
199 if (TMath::IsNaN(h))
200 return kCONTINUE;
201
202 if (fMcEvt->GetPartId()==kGAMMA)
203 fGhness->Fill(h, w);
204 else
205 fPhness->Fill(h, w);
206
207 return kTRUE;
208}
209
210Float_t MHHadronness::GetQ05() const
211{
212 Int_t n = fQfac->GetN();
213
214 Double_t val1x=0;
215 Double_t val2x=1;
216
217 Double_t val1y=0;
218 Double_t val2y=0;
219
220 for (Int_t i=1; i<=n; i++)
221 {
222 Double_t x, y;
223
224 fQfac->GetPoint(i, x, y);
225
226 if (x<0.5 && x>val1x)
227 {
228 val1x = x;
229 val1y = y;
230 }
231
232 if (x>0.5 && x<val2x)
233 {
234 val2x = x;
235 val2y = y;
236 }
237 }
238
239 return val2x-val1x == 0 ? 0 : val1y - (val2y-val1y)/(val2x-val1x) * (val1x-0.5);
240}
241
242// --------------------------------------------------------------------------
243//
244// Finalize the histograms:
245// - integrate the hadroness histograms --> acceptance
246// - fill the Minimum Distance histogram (formular see class description)
247// - fill the Quality histogram (formular see class description)
248//
249void MHHadronness::CalcGraph(Double_t sumg, Double_t sump)
250{
251 Int_t n = fGhness->GetNbinsX();
252
253 fGraph->Set(n);
254 fQfac->Set(n);
255
256 // Calculate acceptances
257 Float_t max=0;
258
259 for (Int_t i=1; i<=n; i++)
260 {
261 const Stat_t ip = fPhness->Integral(1, i)/sump;
262 const Stat_t ig = fGhness->Integral(1, i)/sumg;
263
264 fIntPhness->SetBinContent(i, ip);
265 fIntGhness->SetBinContent(i, ig);
266
267 fGraph->SetPoint(i, ip, ig);
268
269 if (ip<=0)
270 continue;
271
272 const Double_t val = ig/sqrt(ip);
273 fQfac->SetPoint(i, ig, val);
274
275 if (val>max)
276 max = val;
277 }
278
279 fQfac->SetMaximum(max*1.05);
280}
281
282Bool_t MHHadronness::Finalize()
283{
284 const Stat_t sumg = fGhness->Integral();
285 const Stat_t sump = fPhness->Integral();
286
287 *fLog << inf << "Sum Hadronness: gammas=" << sumg << " hadrons=" << sump << endl;
288
289 // Normalize photon distribution
290 if (sumg>0)
291 fGhness->Scale(1./sumg);
292 else
293 *fLog << warn << "Cannot calculate hadronness for 'gammas'." << endl;
294
295 // Normalize hadron distribution
296 if (sump>0)
297 fPhness->Scale(1./sump);
298 else
299 *fLog << warn << "Cannot calculate hadronness for 'hadrons'." << endl;
300
301 CalcGraph(1, 1);
302
303 return kTRUE;
304}
305
306void MHHadronness::Paint(Option_t *opt)
307{
308 Stat_t sumg = fGhness->Integral();
309 Stat_t sump = fPhness->Integral();
310
311 // Normalize photon distribution
312 if (sumg<=0)
313 sumg=1;
314
315 // Normalize hadron distribution
316 if (sump<=0)
317 sump=1;
318
319 CalcGraph(sumg, sump);
320}
321
322// --------------------------------------------------------------------------
323//
324// Search the corresponding points for the given hadron acceptance (acchad)
325// and interpolate the tow points (linear)
326//
327Double_t MHHadronness::GetGammaAcceptance(Double_t acchad) const
328{
329 const Int_t n = fGraph->GetN();
330 const Double_t *x = fGraph->GetX();
331 const Double_t *y = fGraph->GetY();
332
333 Int_t i = 0;
334 while (i<n && x[i]<acchad)
335 i++;
336
337 if (i==0 || i==n)
338 return 0;
339
340 if (i==n-1)
341 i--;
342
343 const Double_t x1 = x[i-1];
344 const Double_t y1 = y[i-1];
345
346 const Double_t x2 = x[i];
347 const Double_t y2 = y[i];
348
349 return (y2-y1)/(x2-x1) * (acchad-x2) + y2;
350}
351
352// --------------------------------------------------------------------------
353//
354// Search the corresponding points for the given gamma acceptance (accgam)
355// and interpolate the tow points (linear)
356//
357Double_t MHHadronness::GetHadronAcceptance(Double_t accgam) const
358{
359 const Int_t n = fGraph->GetN();
360 const Double_t *x = fGraph->GetX();
361 const Double_t *y = fGraph->GetY();
362
363 Int_t i = 0;
364 while (i<n && y[i]<accgam)
365 i++;
366
367 if (i==0 || i==n)
368 return 0;
369
370 if (i==n-1)
371 i--;
372
373 const Double_t x1 = y[i-1];
374 const Double_t y1 = x[i-1];
375
376 const Double_t x2 = y[i];
377 const Double_t y2 = x[i];
378
379 return (y2-y1)/(x2-x1) * (accgam-x2) + y2;
380}
381
382// --------------------------------------------------------------------------
383//
384// Print the corresponding Gammas Acceptance for a hadron acceptance of
385// 10%, 20%, 30%, 40% and 50%. Also the minimum distance to the optimum
386// acceptance and the corresponding acceptances and hadroness value is
387// printed, together with the maximum Q-factor.
388//
389void MHHadronness::Print(Option_t *) const
390{
391 *fLog << all;
392 *fLog << "Hadronness histograms:" << endl;
393 *fLog << "---------------------" << endl;
394
395 if (fGraph->GetN()==0)
396 {
397 *fLog << " <No Entries>" << endl;
398 return;
399 }
400
401 *fLog << "Used " << fGhness->GetEntries() << " Gammas and " << fPhness->GetEntries() << " Hadrons." << endl;
402 *fLog << "acc(hadron) acc(gamma) acc(g)/acc(h)" << endl <<endl;
403
404 *fLog << " 0.005 " << Form("%6.3f", GetGammaAcceptance(0.005)) << " " << Form("%6.3f", GetGammaAcceptance(0.005)/0.005) << endl;
405 *fLog << " 0.02 " << Form("%6.3f", GetGammaAcceptance(0.02)) << " " << Form("%6.3f", GetGammaAcceptance(0.02)/0.02) << endl;
406 *fLog << " 0.05 " << Form("%6.3f", GetGammaAcceptance(0.05)) << " " << Form("%6.3f", GetGammaAcceptance(0.05)/0.05) << endl;
407 *fLog << " 0.1 " << Form("%6.3f", GetGammaAcceptance(0.1 )) << " " << Form("%6.3f", GetGammaAcceptance(0.1)/0.1) << endl;
408 *fLog << " 0.2 " << Form("%6.3f", GetGammaAcceptance(0.2 )) << " " << Form("%6.3f", GetGammaAcceptance(0.2)/0.2) << endl;
409 *fLog << " 0.3 " << Form("%6.3f", GetGammaAcceptance(0.3 )) << " " << Form("%6.3f", GetGammaAcceptance(0.3)/0.3) << endl;
410 *fLog << Form("%6.3f", GetHadronAcceptance(0.1)) << " 0.1 " << endl;
411 *fLog << Form("%6.3f", GetHadronAcceptance(0.2)) << " 0.2 " << endl;
412 *fLog << Form("%6.3f", GetHadronAcceptance(0.3)) << " 0.3 " << endl;
413 *fLog << Form("%6.3f", GetHadronAcceptance(0.4)) << " 0.4 " << endl;
414 *fLog << Form("%6.3f", GetHadronAcceptance(0.5)) << " 0.5 " << endl;
415 *fLog << Form("%6.3f", GetHadronAcceptance(0.6)) << " 0.6 " << endl;
416 *fLog << Form("%6.3f", GetHadronAcceptance(0.7)) << " 0.7 " << endl;
417 *fLog << Form("%6.3f", GetHadronAcceptance(0.8)) << " 0.8 " << endl;
418 *fLog << endl;
419
420 *fLog << "Q-Factor @ Acc Gammas=0.5: Q(0.5)=" << Form("%.1f", GetQ05()) << endl;
421 *fLog << " Acc Hadrons = " << Form("%5.1f", GetHadronAcceptance(0.5)*100) << "%" << endl;
422 *fLog << endl;
423}
424
425// --------------------------------------------------------------------------
426//
427// Draw all histograms. (For the Meaning see class description)
428//
429void MHHadronness::Draw(Option_t *)
430{
431 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas("Hadronness", fTitle);
432 pad->SetBorderMode(0);
433
434 AppendPad("");
435
436 pad->Divide(2, 2);
437
438 TH1 *h;
439
440 pad->cd(1);
441 gPad->SetBorderMode(0);
442 //gStyle->SetOptStat(10);
443 MH::Draw(*fGhness, *fPhness, "Hadronness"); // Displ both stat boxes
444
445 pad->cd(2);
446 gPad->SetBorderMode(0);
447 fIntGhness->Draw();
448 fIntPhness->Draw("same");
449
450 pad->cd(3);
451 gPad->SetBorderMode(0);
452 fQfac->Draw("A*");
453 gPad->Modified();
454 gPad->Update();
455 if ((h=fQfac->GetHistogram()))
456 {
457 h->GetXaxis()->SetRangeUser(0, 1);
458 h->SetXTitle("Acceptance Gammas");
459 h->SetYTitle("Quality");
460 fQfac->Draw("P");
461 gPad->Modified();
462 gPad->Update();
463 }
464
465 pad->cd(4);
466 gPad->SetBorderMode(0);
467 fGraph->Draw("AC");
468 gPad->Modified();
469 gPad->Update();
470 if ((h=fGraph->GetHistogram()))
471 {
472 h->GetXaxis()->SetRangeUser(0, 1);
473 h->SetXTitle("Acceptance Hadrons");
474 h->SetYTitle("Acceptance Gammas");
475 fGraph->SetMaximum(1);
476 fGraph->Draw("P");
477 gPad->Modified();
478 gPad->Update();
479 }
480}
Note: See TracBrowser for help on using the repository browser.