Belle II Software  release-08-01-10
KlongValidationModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <reconstruction/modules/KlId/KlongValidation/KlongValidationModule.h>
10 #include <mdst/dataobjects/KlId.h>
11 
12 #include <framework/datastore/StoreArray.h>
13 #include <framework/logging/Logger.h>
14 #include <framework/gearbox/Const.h>
15 
16 // here's where the functions are hidden
17 #include <reconstruction/modules/KlId/KLMExpert/KlId.h>
18 #include <tracking/dataobjects/TrackClusterSeparation.h>
19 
20 #include <TFile.h>
21 #include <TGraph.h>
22 #include <cstring>
23 
24 using namespace std;
25 using namespace Belle2;
26 using namespace Belle2::KlongId;
27 
28 REG_MODULE(KlongValidation);
29 
30 KlongValidationModule::KlongValidationModule(): Module()
31 {
32  setDescription("Module used by the validation server to generate root files for the K0L validation. It calculates fake rates and efficiencies.");
33  addParam("outputName", m_outputName,
34  "Name of the output root file (must end with .root)",
35  m_outputName);
36  addParam("KlIdCut", m_KlIdCut,
37  "If cut < 0, then only the neutral KLMClusters will be used. Otherwise, only the KLMClusters that satisfies the selection will be used.",
38  m_KlIdCut);
39 }
40 
42 {
43 }
44 
46 {
47  // require existence of necessary datastore obj
48  m_klmClusters.isRequired();
50 
51  // initialize root tree to write stuff into
52  m_f = new TFile(m_outputName.c_str(), "recreate");
54  m_Phi_Pass = new TH1F("Phi Passed", "Passed #Phi;#Phi [rad];Count", 32, -3.2, 3.2);
55  m_Phi_all = new TH1F("Phi All", "All #Phi;#Phi [rad];Count", 32, -3.2, 3.2);
56  m_effPhi = new TH1F("Phi Efficiency", "Efficiency #Phi;#Phi [rad];Efficiency", 32, -3.2, 3.2);
57 
58  m_Theta_Pass = new TH1F("Theta Passed", "Passed #Theta;#Theta [rad];Count", 32, 0, 3.2);
59  m_Theta_all = new TH1F("Theta All", "All #Theta;#Theta [rad];Count", 32, 0, 3.2);
60  m_effTheta = new TH1F("Theta Efficiency", "Efficiency #Theta;#Theta [rad];Efficiency", 32, 0, 3.2);
61 
62  m_Mom_Pass = new TH1F("Momentum Efficiency", "Efficiency Momentum;Momentum [GeV];Efficiency", 25, 0, 5);
63  m_Mom_all = new TH1F("Momentum Efficiency normalise", "Efficiency Momentum;Momentum [GeV];Count", 25, 0, 5);
64  m_Mom_all_plot = new TH1F("Momentum Efficiency all, obtained from clusters", "Efficiency Momentum;Momentum [GeV];Count", 25, 0,
65  5);
66  m_effMom = new TH1F("Momentum Efficiency obtained from cluster", "Efficiency Momentum;Momentum [GeV];Efficiency", 25, 0, 5);
67 
68  m_fakePhi_Pass = new TH1F("Phi Fake Passsed", "Fake Passed #Phi;#Phi [rad];Count", 32, -3.2, 3.2);
69  m_fakePhi = new TH1F("Phi Fake Rate", "Fake Rate #Phi;#Phi [rad];Fake Rate", 32, -3.2, 3.2);
70 
71  m_fakeTheta_Pass = new TH1F("Theta Fake Passed", "Fake Passed #Theta;#Theta [rad];Count", 32, 0, 3.2);
72  m_fakeTheta = new TH1F("Theta Fake Rate", "Fake Rate #Theta;#Theta [rad];Fake Rate", 32, 0, 3.2);
73 
74  m_fakeMom_Pass = new TH1F("Momentum Fake Passed", "Momentum Fake Passed;Momentum [GeV];Count", 25, 0, 5);
75  m_fakeMom = new TH1F("Momentum Fake Rate", "Momentum Fake Rate;Momentum [GeV];Fake Rate", 25, 0, 5);
76 
77  m_time = new TH1F("KLM Cluster Time", "Cluster Timing;Cluster time [ns];Count", 200, -30, 70);
78  m_trackSep = new TH1F("KLM trackSeperation Distance", "KLM trackSeperation Distance;Distance [mm];Count", 100, 0, 4000);
79  m_nLayer = new TH1F("KLM N-Layer", "N-layer;N-layer;count", 20, 0, 20);
80 
81  m_bkgPhi = new TH1F("Phi Beam BKG", "BeamBKG #Phi;#Phi [rad];Background count", 32, -3.2, 3.2);
82  m_bkgTheta = new TH1F("Theta Beam BKG", "BeamBKG #Theta;#Theta [rad];Background count", 32, 0, 3.2);
83  m_bkgMom = new TH1F("Momentum Beam BKG", "BeamBKG #Momentum;#Momentum [GeV];Background count", 25, 0, 5);
84 
85  m_innermostLayer = new TH1F("Innermost layer", "Innermost layer;Innermost layer; Count", 20, 0, 20);
86  m_trackFlag = new TH1F("Track flag", "TrackFlag;Flag; Count", 2, 0, 1);
87  m_ECLFlag = new TH1F("ECL flag", "ECLFlag;Flag; Count", 2, 0, 1);
88 
89  m_klidFake = new TH1F("Klid Fake", "klid Fake;klid; Count", m_xbins.size() - 1, m_xbins.data());
90  m_klidTrue = new TH1F("Klid True", "Klid True;klid; Count", m_xbins.size() - 1, m_xbins.data());
91  m_klidAll = new TH1F("Klid all", "klid all clusters;klid; Count", m_xbins.size() - 1, m_xbins.data());
92 }
93 
95 {
96 }
97 
99 {
100 }
101 
103 {
104  // the performance of the classifier depends on the cut on the classifier
105  // output. This is generally arbitrary as it depends on what the user wants.
106  // Finding K_L or rejecting?
107  // Therefore an arbitrary cut is chosen so that something happens and we can judge
108  // if the behavior changes.
109  for (const KLMCluster& cluster : m_klmClusters) {
110 
111  MCParticle* mcpart = NULL;
112 
113  auto klidObj = cluster.getRelatedTo<KlId>();
114 
115  if (m_KlIdCut >= 0) {
116  if (klidObj) {
117  if (klidObj->getKlId() > m_KlIdCut) {
118  m_reconstructedAsKl = true;
119  } else {
120  m_reconstructedAsKl = false;
121  }
122  }
123  } else {
124  if (!cluster.getAssociatedTrackFlag()) {
125  m_reconstructedAsKl = true;
126  }
127  }
128  mcpart = cluster.getRelatedTo<MCParticle>();
129  // find mc truth
130  // go thru all mcparts mothers up to highest particle and check if its a Klong
131  m_isKl = 0;
132  if (mcpart == nullptr) {
133  m_isKl = 0; // this is the case for beambkg
134  m_isBeamBKG = 1;
135  } else {
136  m_isBeamBKG = 0;
137  while (!(mcpart -> getMother() == nullptr)) {
138  if (mcpart -> getPDG() == Const::Klong.getPDGCode()) {
139  m_isKl = 1;
140  break;
141  }
142  mcpart = mcpart -> getMother();
143  }// while
144  }
145  // second param is cut on
146  m_isKl = isKLMClusterSignal(cluster, 0);
147 
148  double time = cluster.getTime();
149  double nLayer = cluster.getLayers();
150  double trackSep = 1.0e20;
151  auto trackSeperations = cluster.getRelationsTo<TrackClusterSeparation>();
152  for (auto trackSeperation : trackSeperations) {
153  double dist = trackSeperation.getDistance();
154  if (dist < trackSep) {
155  trackSep = dist;
156  }
157  }
158  // only fill if cluster is correctly reconstructed Kl
159  if (m_reconstructedAsKl && m_isKl) {m_passed = true;}
160  else {m_passed = false;}
161 
162  // used for fake rate
163  if (m_reconstructedAsKl && (!m_isKl)) {m_faked = true;}
164  else {m_faked = false;}
165 
166  m_phi = cluster.getMomentum().Phi();
167  m_theta = cluster.getMomentum().Theta();
168  m_momentum = std::abs(cluster.getMomentumMag());
169 
170  if (m_passed) {
171  m_Phi_Pass -> Fill(m_phi);
172  m_Theta_Pass -> Fill(m_theta);
173  m_Mom_Pass -> Fill(m_momentum);
174  }
175 
176  if (m_faked) {
177  m_fakePhi_Pass -> Fill(m_phi);
178  m_fakeTheta_Pass -> Fill(m_theta);
179  m_fakeMom_Pass -> Fill(m_momentum);
180  }
181 
182  if (m_isBeamBKG) {
183  m_bkgPhi-> Fill(m_phi);
184  m_bkgTheta-> Fill(m_theta);
185  m_bkgMom-> Fill(m_momentum);
186  }
187 
188  m_klidAll->SetMinimum(0.);
189  m_Mom_all_plot->SetMinimum(0.);
190  m_effPhi->SetMinimum(0.);
191  m_effTheta->SetMinimum(0.);
192  m_effMom->SetMinimum(0.);
193  m_fakePhi->SetMinimum(0.);
194  m_fakeTheta->SetMinimum(0.);
195  m_fakeMom->SetMinimum(0.);
196  m_time->SetMinimum(0.);
197  m_trackSep->SetMinimum(0.);
198  m_nLayer->SetMinimum(0.);
199  m_bkgPhi->SetMinimum(0.);
200  m_bkgTheta->SetMinimum(0.);
201  m_bkgMom->SetMinimum(0.);
202  m_innermostLayer->SetMinimum(0.);
203  m_trackFlag->SetMinimum(0.);
204  m_ECLFlag->SetMinimum(0.);
205 
206  //fil all to normalise later
207  m_Phi_all -> Fill(m_phi);
208  m_Theta_all -> Fill(m_theta);
209  m_Mom_all -> Fill(m_momentum);
210  m_Mom_all_plot -> Fill(m_momentum);
211 
212  m_time ->Fill(time);
213  m_trackSep->Fill(trackSep);
214  m_nLayer->Fill(nLayer);
215  m_innermostLayer->Fill(cluster.getInnermostLayer());
216  m_trackFlag->Fill(cluster.getAssociatedTrackFlag());
217  m_ECLFlag->Fill(cluster.getAssociatedEclClusterFlag());
218 
219  if (klidObj) {
220  if (!m_isKl) {
221  m_klidFake ->Fill(klidObj->getKlId());
222  }
223  if (m_isKl) {
224  m_klidTrue ->Fill(klidObj->getKlId());
225  }
226  m_klidAll ->Fill(klidObj->getKlId());
227  }
228  }// for klm clusters
229 
230 }
231 
233 {
234  // write TEff to root file ,
235  m_f -> cd();
236  // TH1F is not compatible with the validation server
237  m_effPhi->Divide(m_Phi_Pass, m_Phi_all);
239  m_effMom->Divide(m_Mom_Pass, m_Mom_all);
240 
244 
245  std::vector<double> efficiency(m_xbins.size());
246  std::vector<double> purity(m_xbins.size());
247  std::vector<double> backrej(m_xbins.size());
248  double NallKlong = m_klidTrue->GetEntries();
249  double NallFakes = m_klidFake->GetEntries();
250  for (unsigned int i = 0; i < m_xbins.size() - 1; ++i) {
251  double NtruePositive = m_klidTrue->Integral(i, m_xbins.size() - 1);
252  double NfalsePositive = m_klidFake->Integral(i, m_xbins.size() - 1);
253  if ((NtruePositive + NfalsePositive) <= 0) {
254  purity[i] = 0;
255  } else {
256  purity[i] = NtruePositive / (NtruePositive + NfalsePositive);
257  }
258  if (NallFakes) {
259  backrej[i] = 1 - NfalsePositive / NallFakes;
260  } else {
261  backrej[i] = 0;
262  }
263  if (NallKlong <= 0) {
264  B2WARNING("No Klongs found in ROC calculation.");
265  } else {
266  efficiency[i] = NtruePositive / NallKlong;
267  }
268  }
269 
270  m_ROC = new TGraph(m_xbins.size(), efficiency.data(), purity.data());
271  m_backRej = new TGraph(m_xbins.size(), efficiency.data(), backrej.data());
272 
273  m_ROC->SetMarkerStyle(3);
274  m_backRej->SetMarkerStyle(3);
275  m_ROC->SetLineColorAlpha(kBlue, 0.0);
276  m_backRej->SetLineColorAlpha(kBlue, 0.0);
277  m_ROC->GetXaxis()->SetTitle("Efficiency");
278  m_backRej->GetXaxis()->SetTitle("Efficiency");
279  m_ROC->GetYaxis()->SetTitle("Purity");
280  m_backRej->GetYaxis()->SetTitle("Background rejection");
281  m_ROC-> SetTitle("Klong Purity vs Efficiency");
282  m_backRej-> SetTitle("Klong background rejection vs Efficiency");
283  m_ROC-> GetListOfFunctions() -> Add(new TNamed("Description",
284  "Purity vs Efficiency each point represents a cut on the klong ID."));
285  m_ROC -> GetListOfFunctions() -> Add(new TNamed("Check", "Should be as high as possible"));
286  m_backRej -> GetListOfFunctions() -> Add(new TNamed("Check", "Should be as high as possible"));
287  m_ROC -> GetListOfFunctions() -> Add(new TNamed("Contact", "fnc@lnf.infn.it"));
288  m_backRej -> GetListOfFunctions() -> Add(new TNamed("Contact", "fnc@lnf.infn.it"));
289 
290  // tuple: pointer to the plot, name of the plot, true for shifter plots
291  std::vector<std::tuple<TH1F*, std::string, std::string, bool>> histograms;
292  std::string defaultCheck{"Nightly result should not differ significantly from the reference"};
293  histograms.push_back(make_tuple(m_klidAll, "KlId distribution", defaultCheck, true));
294  histograms.push_back(make_tuple(m_Mom_all_plot, "All Momenta generated", defaultCheck, false));
295  histograms.push_back(make_tuple(m_effPhi, "KlId efficiency in Phi", defaultCheck, false));
296  histograms.push_back(make_tuple(m_effTheta, "KlId efficiency in Theta",
297  defaultCheck + "; efficiency must be 0 below 0.3 rad and above 2.7 rad", true));
298  histograms.push_back(make_tuple(m_effMom, "KlId efficiency in Momentum", defaultCheck, true));
299  histograms.push_back(make_tuple(m_fakePhi, "KlId fake rate in Phi", defaultCheck, false));
300  histograms.push_back(make_tuple(m_fakeTheta, "KlId fake rate in Theta",
301  defaultCheck + "; fake rate must be 0 below 0.3 rad and above 2.7 rad", true));
302  histograms.push_back(make_tuple(m_fakeMom, "KlId fake rate in Momentum", defaultCheck, true));
303  histograms.push_back(make_tuple(m_time, "KLMClusters time", defaultCheck + "; no secondary peaks", true));
304  histograms.push_back(make_tuple(m_trackSep, "Distance of KLMCluster to closest Track", defaultCheck, true));
305  histograms.push_back(make_tuple(m_nLayer, "KLMClusters nLayer", defaultCheck, false));
306  histograms.push_back(make_tuple(m_innermostLayer, "KLMClusters innermostLayer", defaultCheck, false));
307  histograms.push_back(make_tuple(m_bkgPhi, "Beam bkg in Phi", defaultCheck, false));
308  histograms.push_back(make_tuple(m_bkgTheta, "Beam bkg in Theta",
309  defaultCheck + "; background must be 0 below 0.3 rad and above 2.7 rad", false));
310  histograms.push_back(make_tuple(m_bkgMom, "Beam bkg in Momentum", defaultCheck, false));
311  histograms.push_back(make_tuple(m_trackFlag, "Track flag", defaultCheck, false));
312  histograms.push_back(make_tuple(m_ECLFlag, "ECLCluster flag", defaultCheck, false));
313 
314  for (auto hist : histograms) {
315  std::get<0>(hist) -> SetTitle(std::get<1>(hist).c_str());
316  std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("Description", std::get<1>(hist)));
317  std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("Check", std::get<2>(hist).c_str()));
318  std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("Contact", "fnc@lnf.infn.it"));
319  if (std::get<3>(hist))
320  std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("MetaOptions", "shifter,pvalue-warn=0.99,pvalue-error=0.02"));
321  else
322  std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("MetaOptions", "pvalue-warn=0.99,pvalue-error=0.1"));
323  std::get<0>(hist) -> Write();
324  }
325 
326  //this guy is a tgraph
327  m_ROC -> Write();
328  m_backRej -> Write();
329  m_f -> Close();
330 }
int getPDGCode() const
PDG code.
Definition: Const.h:464
static const ParticleType Klong
K^0_L particle.
Definition: Const.h:669
KLM cluster data.
Definition: KLMCluster.h:28
Klong identifcation (KlId) datastore object to store results from KlId calculations.
Definition: KlId.h:23
TH1F * m_Mom_Pass
momentum efficiency
TH1F * m_effTheta
efficiency in angle to z
bool m_passed
did cluster pass selection of algorythm?
TH1F * m_trackSep
track separation distance plot
TH1F * m_Mom_all_plot
momentum efficiency
TH1F * m_Theta_Pass
efficiency in angle to z
TH1F * m_fakeMom_Pass
fake momentum plot
StoreArray< KLMCluster > m_klmClusters
storearrays
virtual void initialize() override
initialize
TGraph * m_backRej
background rejection
TH1F * m_fakePhi_Pass
fake phi, angle in x-y
virtual void event() override
process event
const std::vector< double > m_xbins
bins used for the ROC plots
virtual void endRun() override
end run
virtual void terminate() override
terminate
TH1F * m_fakeTheta_Pass
fake theta, angle to z
TH1F * m_Phi_Pass
efficiency in x-y plane
TH1F * m_Phi_all
efficiency in x-y plane
virtual void beginRun() override
beginn run
std::string m_outputName
output file name
TH1F * m_time
cluster timing plot
float m_KlIdCut
of > 0 use Klid else use trackflag as reconstruction criterion
virtual ~KlongValidationModule()
Destructor
StoreArray< MCParticle > m_mcParticles
storearrays
bool m_faked
cluster wrongly reconstructed as K0L?
bool m_reconstructedAsKl
cluster reconstructed as K0L?
TH1F * m_fakeTheta
fake theta, angle to z
TH1F * m_effMom
momentum efficiency
TH1F * m_effPhi
efficiency in x-y plane
TH1F * m_Mom_all
momentum efficiency
TH1F * m_fakeMom
fake momentum plot
TH1F * m_fakePhi
fake phi, angle in x-y
TH1F * m_Theta_all
efficiency in angle to z
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Store one Track-KLMCluster separation as a ROOT object.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Helper functions for all klid modules to improve readability of the code.
Definition: KlId.h:27
bool isKLMClusterSignal(const Belle2::KLMCluster &cluster, float mcWeigthCut=0.66)
checks if a cluster is signal under the mcWeightcondition (mcWeight = energy deposition)
Definition: KlId.h:100
Abstract base class for different kinds of events.