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

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