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

Last change on this file since 8785 was 8022, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 17.0 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: Acceptance 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// If you are using filtercuts which gives you only two discrete values
59// of the hadronness (0.25 and 0.75) you may want to use MHHadronness(2).
60// If you request Q05() from such a MHHadronness instance you will get
61// Acc_g/sqrt(Acc_h)
62//
63////////////////////////////////////////////////////////////////////////////
64#include "MHHadronness.h"
65
66#include <TPad.h>
67#include <TGraph.h>
68#include <TStyle.h>
69#include <TCanvas.h>
70#include <TMarker.h>
71
72#include "MParList.h"
73#include "MBinning.h"
74#include "MHMatrix.h"
75#include "MParameters.h"
76
77#include "MLog.h"
78#include "MLogManip.h"
79
80#include "MMcEvt.hxx"
81
82ClassImp(MHHadronness);
83
84using namespace std;
85
86// --------------------------------------------------------------------------
87//
88// Setup histograms, nbins is the number of bins used for the evaluation.
89// The default is 100 bins.
90//
91MHHadronness::MHHadronness(Int_t nbins, const char *name, const char *title)
92 : fMatrix(NULL)
93{
94 //
95 // set the name and title of this object
96 //
97 fName = name ? name : "MHHadronness";
98 fTitle = title ? title : "Gamma/Hadron Separation Quality Histograms";
99
100 fGraph = new TGraph;
101 fGraph->SetTitle("Acceptance Gammas vs. Hadrons");
102 fGraph->SetMarkerStyle(kFullDotSmall);
103
104 fGhness = new TH1D("Ghness", "Acceptance vs. Hadronness (Gammas)", nbins, 0, 1.0001);
105 fPhness = new TH1D("Phness", "Acceptance vs. Hadronness (Hadrons)", nbins, 0, 1.0001);
106 fGhness->SetXTitle("Hadronness");
107 fPhness->SetXTitle("Hadronness");
108 fGhness->SetYTitle("Acceptance");
109 fPhness->SetYTitle("Acceptance");
110 fPhness->SetLineColor(kRed);
111 fGhness->SetDirectory(NULL);
112 fPhness->SetDirectory(NULL);
113
114 fIntGhness = new TH1D("AccGammas", "Integral Acceptance vs. Hadronness (Gammas)", nbins, 0, 1);
115 fIntPhness = new TH1D("AccHadrons", "Integral Acceptance vs. Hadronness (Hadrons)", nbins, 0, 1);
116 fIntGhness->SetXTitle("Hadronness");
117 fIntPhness->SetXTitle("Hadronness");
118 fIntGhness->SetYTitle("Acceptance");
119 fIntPhness->SetYTitle("Acceptance");
120 fIntGhness->SetMaximum(1.1);
121 fIntPhness->SetMaximum(1.1);
122 fIntGhness->SetDirectory(NULL);
123 fIntPhness->SetDirectory(NULL);
124 fIntPhness->SetLineColor(kRed);
125 fIntGhness->SetBit(TH1::kNoStats);
126 fIntPhness->SetBit(TH1::kNoStats);
127
128 fQfac = new TGraph;
129 fQfac->SetTitle(" Naive Quality factor ");
130 fQfac->SetMarkerStyle(kFullDotSmall);
131}
132
133// --------------------------------------------------------------------------
134//
135// Delete the histograms.
136//
137MHHadronness::~MHHadronness()
138{
139 delete fGhness;
140 delete fIntGhness;
141 delete fPhness;
142 delete fIntPhness;
143 delete fQfac;
144 delete fGraph;
145}
146
147// --------------------------------------------------------------------------
148//
149// Setup Filling of the histograms. It needs:
150// MMcEvt and MHadronness
151//
152Bool_t MHHadronness::SetupFill(const MParList *plist)
153{
154 if (!fMatrix)
155 {
156 fMcEvt = (MMcEvt*)plist->FindObject(AddSerialNumber("MMcEvt"));
157 if (!fMcEvt)
158 {
159 *fLog << err << dbginf << AddSerialNumber("MMcEvt");
160 *fLog << " not found... aborting." << endl;
161 return kFALSE;
162 }
163 }
164
165 fHadronness = (MParameterD*)plist->FindObject("MHadronness");
166
167// fGhness->Reset();
168// fPhness->Reset();
169
170 /*
171 MBinning* bins = (MBinning*)plist->FindObject("BinningHadronness");
172 if (!bins)
173 {
174 *fLog << err << dbginf << "BinningHadronness [MBinning] not found... aborting." << endl;
175 return kFALSE;
176 }
177
178 SetBinning(&fHist, binsalpha, binsenergy, binstheta);
179
180 fHist.Sumw2();
181 */
182
183 return kTRUE;
184}
185
186// --------------------------------------------------------------------------
187//
188// Fill the Hadronness from a MHadronness container into the corresponding
189// histogram dependant on the particle id.
190//
191// Every particle Id different than kGAMMA is considered a hadron.
192//
193// If you call Fill with a pointer (eg. from MFillH) this container is
194// used as a hadronness (Warning: No type check is done!) otherwise
195// the default is used ("MHadronness")
196//
197// Sometimes a distance is calculated as NaN (not a number). Such events
198// are skipped at the moment.
199//
200Bool_t MHHadronness::Fill(const MParContainer *par, const Stat_t w)
201{
202 // Preliminary Workaround: FIXME!
203 if (!par && !fHadronness)
204 {
205 *fLog << err << "MHHadronness::Fill: No MHadronness container specified!" << endl;
206 return kFALSE;
207 }
208
209 const MParameterD &had = par ? *(MParameterD*)par : *fHadronness;
210
211 const Double_t h = TMath::Min(TMath::Max(had.GetVal(), 0.), 1.);
212
213 if (!TMath::Finite(h))
214 return kCONTINUE;
215
216 const Int_t particleid = fMatrix ? (Int_t)(*fMatrix)[fMap] : fMcEvt->GetPartId();
217
218 if (particleid==MMcEvt::kGAMMA)
219 fGhness->Fill(h, w);
220 else
221 fPhness->Fill(h, w);
222
223 return kTRUE;
224}
225
226// --------------------------------------------------------------------------
227//
228// Returns the quality factor at gamma acceptance 0.5.
229//
230// If the histogram containes only two bins we return the
231// naive quality: ag/sqrt(ah)
232// with ag the acceptance of gammas and ah the acceptance of hadrons.
233//
234// You can use this (nbins=2) in case of some kind of filter cuts giving
235// only a result: gamma yes/no (means discrete values of hadronness 0.25
236// or 0.75)
237//
238// FIXME: In the later case weights cannot be used!
239//
240Float_t MHHadronness::GetQ05() const
241{
242 if (fGhness->GetNbinsX()==2)
243 {
244 // acceptance of all gamma-like gammas (h<0.5)
245 const Double_t ig = fGhness->GetBinContent(1);
246
247 // acceptance of all gamma-like hadrons (h<0.5)
248 const Double_t ip = fPhness->GetBinContent(1);
249
250 if (ip==0)
251 return 0; // FIXME!
252
253 // naive quality factor
254 const Double_t q = ig / sqrt(ip);
255
256 *fLog << all << ip << "/" << ig << ": " << q << endl;
257
258 return q;
259 }
260
261 const Int_t n = fQfac->GetN();
262
263 Double_t val1x=0;
264 Double_t val2x=1;
265
266 Double_t val1y=0;
267 Double_t val2y=0;
268
269 for (Int_t i=1; i<=n; i++)
270 {
271 Double_t x, y;
272
273 fQfac->GetPoint(i, x, y);
274
275 if (x<0.5 && x>val1x)
276 {
277 val1x = x;
278 val1y = y;
279 }
280
281 if (x>0.5 && x<val2x)
282 {
283 val2x = x;
284 val2y = y;
285 }
286 }
287
288 //*fLog << dbg << val1x << "/" << val1y << " " << val2x << "/" << val2y << endl;
289
290 return val2x-val1x == 0 ? 0 : val1y - (val2y-val1y)/(val2x-val1x) * (val1x-0.5);
291}
292
293// --------------------------------------------------------------------------
294//
295// Finalize the histograms:
296// - integrate the hadroness histograms --> acceptance
297// - fill the Minimum Distance histogram (formular see class description)
298// - fill the Quality histogram (formular see class description)
299//
300void MHHadronness::CalcGraph(Double_t sumg, Double_t sump)
301{
302 Int_t n = fGhness->GetNbinsX();
303
304 fGraph->Set(n);
305 fQfac->Set(n);
306
307 // Calculate acceptances
308 Float_t max=0;
309
310 for (Int_t i=1; i<=n; i++)
311 {
312 const Stat_t ip = fPhness->Integral(1, i)/sump;
313 const Stat_t ig = fGhness->Integral(1, i)/sumg;
314
315 fIntPhness->SetBinContent(i, ip);
316 fIntGhness->SetBinContent(i, ig);
317
318 fGraph->SetPoint(i, ip, ig);
319
320 if (ip<=0)
321 continue;
322
323 const Double_t val = ig/sqrt(ip);
324 fQfac->SetPoint(i, ig, val);
325
326 if (val>max)
327 max = val;
328 }
329
330 fQfac->SetMaximum(max*1.05);
331}
332
333Bool_t MHHadronness::Finalize()
334{
335 const Stat_t sumg = fGhness->Integral();
336 const Stat_t sump = fPhness->Integral();
337
338 *fLog << inf << "Sum Hadronness: gammas=" << sumg << " hadrons=" << sump << endl;
339
340 // Normalize photon distribution
341 if (sumg>0)
342 fGhness->Scale(1./sumg);
343 else
344 *fLog << warn << "Cannot calculate hadronness for 'gammas'." << endl;
345
346 // Normalize hadron distribution
347 if (sump>0)
348 fPhness->Scale(1./sump);
349 else
350 *fLog << warn << "Cannot calculate hadronness for 'hadrons'." << endl;
351
352 CalcGraph(1, 1);
353
354 return kTRUE;
355}
356
357void MHHadronness::Paint(Option_t *opt)
358{
359 Stat_t sumg = fGhness->Integral();
360 Stat_t sump = fPhness->Integral();
361
362 // Normalize photon distribution
363 if (sumg<=0)
364 sumg=1;
365
366 // Normalize hadron distribution
367 if (sump<=0)
368 sump=1;
369
370 CalcGraph(sumg, sump);
371}
372
373// --------------------------------------------------------------------------
374//
375// Search the corresponding points for the given hadron acceptance (acchad)
376// and interpolate the tow points (linear)
377//
378Double_t MHHadronness::GetGammaAcceptance(Double_t acchad) const
379{
380 const Int_t n = fGraph->GetN();
381 const Double_t *x = fGraph->GetX();
382 const Double_t *y = fGraph->GetY();
383
384 Int_t i = 0;
385 while (i<n && x[i]<acchad)
386 i++;
387
388 if (i==0 || i==n)
389 return 0;
390
391 if (i==n-1)
392 i--;
393
394 const Double_t x1 = x[i-1];
395 const Double_t y1 = y[i-1];
396
397 const Double_t x2 = x[i];
398 const Double_t y2 = y[i];
399
400 return (y2-y1)/(x2-x1) * (acchad-x2) + y2;
401}
402
403// --------------------------------------------------------------------------
404//
405// Search the corresponding points for the given gamma acceptance (accgam)
406// and interpolate the tow points (linear)
407//
408Double_t MHHadronness::GetHadronAcceptance(Double_t accgam) const
409{
410 const Int_t n = fGraph->GetN();
411 const Double_t *x = fGraph->GetX();
412 const Double_t *y = fGraph->GetY();
413
414 Int_t i = 0;
415 while (i<n && y[i]<accgam)
416 i++;
417
418 if (i==0 || i==n)
419 return 0;
420
421 if (i==n-1)
422 i--;
423
424 const Double_t x1 = y[i-1];
425 const Double_t y1 = x[i-1];
426
427 const Double_t x2 = y[i];
428 const Double_t y2 = x[i];
429
430 return (y2-y1)/(x2-x1) * (accgam-x2) + y2;
431}
432
433// --------------------------------------------------------------------------
434//
435// Search the hadronness corresponding to a given hadron acceptance.
436//
437Double_t MHHadronness::GetHadronness(Double_t acchad) const
438{
439 for (int i=1; i<fIntPhness->GetNbinsX()+1; i++)
440 if (fIntPhness->GetBinContent(i)>acchad)
441 return fIntPhness->GetBinLowEdge(i);
442
443 return -1;
444}
445
446// --------------------------------------------------------------------------
447//
448// Print the corresponding Gammas Acceptance for a hadron acceptance of
449// 10%, 20%, 30%, 40% and 50%. Also the minimum distance to the optimum
450// acceptance and the corresponding acceptances and hadroness value is
451// printed, together with the maximum Q-factor.
452//
453void MHHadronness::Print(Option_t *) const
454{
455 *fLog << all;
456 *fLog << underline << GetDescriptor() << endl;
457
458 if (fGraph->GetN()==0)
459 {
460 *fLog << " <No Entries>" << endl;
461 return;
462 }
463
464 *fLog << "Used " << fGhness->GetEntries() << " Gammas and " << fPhness->GetEntries() << " Hadrons." << endl;
465 *fLog << "acc(hadron) acc(gamma) acc(g)/acc(h) h" << endl <<endl;
466
467 *fLog << " 0.005 " << Form("%6.3f", GetGammaAcceptance(0.005)) << " " << Form("%6.3f", GetGammaAcceptance(0.005)/0.005) << " " << GetHadronness(0.005) << endl;
468 *fLog << " 0.02 " << Form("%6.3f", GetGammaAcceptance(0.02)) << " " << Form("%6.3f", GetGammaAcceptance(0.02)/0.02) << " " << GetHadronness(0.02) << endl;
469 *fLog << " 0.05 " << Form("%6.3f", GetGammaAcceptance(0.05)) << " " << Form("%6.3f", GetGammaAcceptance(0.05)/0.05) << " " << GetHadronness(0.05) << endl;
470 *fLog << " 0.1 " << Form("%6.3f", GetGammaAcceptance(0.1 )) << " " << Form("%6.3f", GetGammaAcceptance(0.1)/0.1) << " " << GetHadronness(0.1) << endl;
471 *fLog << " 0.2 " << Form("%6.3f", GetGammaAcceptance(0.2 )) << " " << Form("%6.3f", GetGammaAcceptance(0.2)/0.2) << " " << GetHadronness(0.2) << endl;
472 *fLog << " 0.3 " << Form("%6.3f", GetGammaAcceptance(0.3 )) << " " << Form("%6.3f", GetGammaAcceptance(0.3)/0.3) << " " << GetHadronness(0.3) << endl;
473 *fLog << Form("%6.3f", GetHadronAcceptance(0.1)) << " 0.1 " << endl;
474 *fLog << Form("%6.3f", GetHadronAcceptance(0.2)) << " 0.2 " << endl;
475 *fLog << Form("%6.3f", GetHadronAcceptance(0.3)) << " 0.3 " << endl;
476 *fLog << Form("%6.3f", GetHadronAcceptance(0.4)) << " 0.4 " << endl;
477 *fLog << Form("%6.3f", GetHadronAcceptance(0.5)) << " 0.5 " << endl;
478 *fLog << Form("%6.3f", GetHadronAcceptance(0.6)) << " 0.6 " << endl;
479 *fLog << Form("%6.3f", GetHadronAcceptance(0.7)) << " 0.7 " << endl;
480 *fLog << Form("%6.3f", GetHadronAcceptance(0.8)) << " 0.8 " << endl;
481 *fLog << endl;
482
483 *fLog << "Q-Factor @ Acc Gammas=0.5: Q(0.5)=" << Form("%.1f", GetQ05()) << endl;
484 *fLog << " Acc Hadrons = " << Form("%5.1f", GetHadronAcceptance(0.5)*100) << "%" << endl;
485 *fLog << endl;
486}
487
488// --------------------------------------------------------------------------
489//
490// Draw all histograms. (For the Meaning see class description)
491//
492void MHHadronness::Draw(Option_t *)
493{
494 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas("Hadronness", fTitle);
495 pad->SetBorderMode(0);
496
497 AppendPad("");
498
499 pad->Divide(2, 2);
500
501 TH1 *h;
502
503 pad->cd(1);
504 gPad->SetBorderMode(0);
505 //gStyle->SetOptStat(10);
506 MH::DrawSame(*fGhness, *fPhness, "Hadronness"); // Displ both stat boxes
507
508 pad->cd(2);
509 gPad->SetBorderMode(0);
510 gPad->SetGridx();
511 gPad->SetGridy();
512 fIntGhness->Draw();
513 fIntPhness->Draw("same");
514
515 pad->cd(3);
516 gPad->SetBorderMode(0);
517 gPad->SetGridx();
518 gPad->SetGridy();
519 fQfac->Draw("A*");
520 gPad->Modified();
521 gPad->Update();
522 if ((h=fQfac->GetHistogram()))
523 {
524 h->GetXaxis()->SetRangeUser(0, 1);
525 h->SetXTitle("Acceptance Gammas");
526 h->SetYTitle("Quality");
527 fQfac->Draw("P");
528 gPad->Modified();
529 gPad->Update();
530 }
531
532 pad->cd(4);
533 gPad->SetBorderMode(0);
534 gPad->SetGridx();
535 gPad->SetGridy();
536 fGraph->Draw("AC");
537 gPad->Modified();
538 gPad->Update();
539 if ((h=fGraph->GetHistogram()))
540 {
541 h->GetXaxis()->SetRangeUser(0, 1);
542 h->SetXTitle("Acceptance Hadrons");
543 h->SetYTitle("Acceptance Gammas");
544 fGraph->SetMaximum(1);
545 fGraph->Draw("P");
546 gPad->Modified();
547 gPad->Update();
548 }
549}
550
551// --------------------------------------------------------------------------
552//
553// You can use this function if you want to use a MHMatrix instead of
554// MMcEvt. This function adds all necessary columns to the
555// given matrix. Afterward you should fill the matrix with the corresponding
556// data (eg from a file by using MHMatrix::Fill). If you now loop
557// through the matrix (eg using MMatrixLoop) MHHadronness::Fill
558// will take the values from the matrix instead of the containers.
559//
560void MHHadronness::InitMapping(MHMatrix *mat)
561{
562 if (fMatrix)
563 return;
564
565 fMatrix = mat;
566
567 TString str = AddSerialNumber("MMcEvt");
568 str += ".fPartId";
569
570 fMap = fMatrix->AddColumn(str);
571}
572
573void MHHadronness::StopMapping()
574{
575 fMatrix = NULL;
576}
Note: See TracBrowser for help on using the repository browser.