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

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