Belle II Software development
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
24using namespace std;
25using namespace Belle2;
26using namespace Belle2::KlongId;
27
28REG_MODULE(KlongValidation);
29
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)",
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 addParam("contact", m_contact, "Email address of contact for validation plots.", m_contact);
40}
41
43{
44}
45
47{
48 // require existence of necessary datastore obj
49 m_klmClusters.isRequired();
51
52 // initialize root tree to write stuff into
53 m_f = new TFile(m_outputName.c_str(), "recreate");
55 m_Phi_Pass = new TH1F("Phi Passed", "Passed #Phi;#Phi [rad];Count", 32, -3.2, 3.2);
56 m_Phi_all = new TH1F("Phi All", "All #Phi;#Phi [rad];Count", 32, -3.2, 3.2);
57 m_effPhi = new TH1F("Phi Efficiency", "Efficiency #Phi;#Phi [rad];Efficiency", 32, -3.2, 3.2);
58
59 m_Theta_Pass = new TH1F("Theta Passed", "Passed #Theta;#Theta [rad];Count", 32, 0, 3.2);
60 m_Theta_all = new TH1F("Theta All", "All #Theta;#Theta [rad];Count", 32, 0, 3.2);
61 m_effTheta = new TH1F("Theta Efficiency", "Efficiency #Theta;#Theta [rad];Efficiency", 32, 0, 3.2);
62
63 m_Mom_Pass = new TH1F("Momentum Efficiency", "Efficiency Momentum;Momentum [GeV];Efficiency", 25, 0, 5);
64 m_Mom_all = new TH1F("Momentum Efficiency normalise", "Efficiency Momentum;Momentum [GeV];Count", 25, 0, 5);
65 m_Mom_all_plot = new TH1F("Momentum Efficiency all, obtained from clusters", "Efficiency Momentum;Momentum [GeV];Count", 25, 0,
66 5);
67 m_effMom = new TH1F("Momentum Efficiency obtained from cluster", "Efficiency Momentum;Momentum [GeV];Efficiency", 25, 0, 5);
68
69 m_fakePhi_Pass = new TH1F("Phi Fake Passed", "Fake Passed #Phi;#Phi [rad];Count", 32, -3.2, 3.2);
70 m_fakePhi = new TH1F("Phi Fake Rate", "Fake Rate #Phi;#Phi [rad];Fake Rate", 32, -3.2, 3.2);
71
72 m_fakeTheta_Pass = new TH1F("Theta Fake Passed", "Fake Passed #Theta;#Theta [rad];Count", 32, 0, 3.2);
73 m_fakeTheta = new TH1F("Theta Fake Rate", "Fake Rate #Theta;#Theta [rad];Fake Rate", 32, 0, 3.2);
74
75 m_fakeMom_Pass = new TH1F("Momentum Fake Passed", "Momentum Fake Passed;Momentum [GeV];Count", 25, 0, 5);
76 m_fakeMom = new TH1F("Momentum Fake Rate", "Momentum Fake Rate;Momentum [GeV];Fake Rate", 25, 0, 5);
77
78 m_time = new TH1F("KLM Cluster Time", "Cluster Timing;Cluster time [ns];Count", 100, -30, 70);
79 m_trackSep = new TH1F("KLM trackSeperation Distance", "KLM trackSeperation Distance;Distance [mm];Count", 100, 0, 4000);
80 m_nLayer = new TH1F("KLM N-Layer", "N-layer;N-layer;count", 20, 0, 20);
81
82 m_bkgPhi = new TH1F("Phi Beam BKG", "BeamBKG #Phi;#Phi [rad];Background count", 32, -3.2, 3.2);
83 m_bkgTheta = new TH1F("Theta Beam BKG", "BeamBKG #Theta;#Theta [rad];Background count", 32, 0, 3.2);
84 m_bkgMom = new TH1F("Momentum Beam BKG", "BeamBKG #Momentum;#Momentum [GeV];Background count", 25, 0, 5);
85
86 m_innermostLayer = new TH1F("Innermost layer", "Innermost layer;Innermost layer; Count", 20, 0, 20);
87 m_trackFlag = new TH1F("Track flag", "TrackFlag;Flag; Count", 2, 0, 1);
88 m_ECLFlag = new TH1F("ECL flag", "ECLFlag;Flag; Count", 2, 0, 1);
89
90 m_klidFake = new TH1F("Klid Fake", "klid Fake;klid; Count", m_xbins.size() - 1, m_xbins.data());
91 m_klidTrue = new TH1F("Klid True", "Klid True;klid; Count", m_xbins.size() - 1, m_xbins.data());
92 m_klidAll = new TH1F("Klid all", "klid all clusters;klid; Count", m_xbins.size() - 1, m_xbins.data());
93}
94
96{
97}
98
100{
101}
102
104{
105 // the performance of the classifier depends on the cut on the classifier
106 // output. This is generally arbitrary as it depends on what the user wants.
107 // Finding K_L or rejecting?
108 // Therefore an arbitrary cut is chosen so that something happens and we can judge
109 // if the behavior changes.
110 for (const KLMCluster& cluster : m_klmClusters) {
111
112 MCParticle* mcpart = NULL;
113
114 auto klidObj = cluster.getRelatedTo<KlId>();
115
116 if (m_KlIdCut >= 0) {
117 if (klidObj) {
118 if (klidObj->getKlId() > m_KlIdCut) {
119 m_reconstructedAsKl = true;
120 } else {
121 m_reconstructedAsKl = false;
122 }
123 }
124 } else {
125 if (!cluster.getAssociatedTrackFlag()) {
126 m_reconstructedAsKl = true;
127 }
128 }
129 mcpart = cluster.getRelatedTo<MCParticle>();
130 // find mc truth
131 // go thru all mcparts mothers up to highest particle and check if its a Klong
132 m_isKl = 0;
133 if (mcpart == nullptr) {
134 m_isKl = 0; // this is the case for beambkg
135 m_isBeamBKG = 1;
136 } else {
137 m_isBeamBKG = 0;
138 while (!(mcpart -> getMother() == nullptr)) {
139 if (mcpart -> getPDG() == Const::Klong.getPDGCode()) {
140 m_isKl = 1;
141 break;
142 }
143 mcpart = mcpart -> getMother();
144 }// while
145 }
146 // second param is cut on
147 m_isKl = isKLMClusterSignal(cluster, 0);
148
149 double time = cluster.getTime();
150 double nLayer = cluster.getLayers();
151 double trackSep = 1.0e20;
152 auto trackSeperations = cluster.getRelationsTo<TrackClusterSeparation>();
153 for (auto trackSeperation : trackSeperations) {
154 double dist = trackSeperation.getDistance();
155 if (dist < trackSep) {
156 trackSep = dist;
157 }
158 }
159 // only fill if cluster is correctly reconstructed Kl
160 if (m_reconstructedAsKl && m_isKl) {m_passed = true;}
161 else {m_passed = false;}
162
163 // used for fake rate
164 if (m_reconstructedAsKl && (!m_isKl)) {m_faked = true;}
165 else {m_faked = false;}
166
167 m_phi = cluster.getMomentum().Phi();
168 m_theta = cluster.getMomentum().Theta();
169 m_momentum = std::abs(cluster.getMomentumMag());
170
171 if (m_passed) {
172 m_Phi_Pass -> Fill(m_phi);
173 m_Theta_Pass -> Fill(m_theta);
174 m_Mom_Pass -> Fill(m_momentum);
175 }
176
177 if (m_faked) {
178 m_fakePhi_Pass -> Fill(m_phi);
179 m_fakeTheta_Pass -> Fill(m_theta);
180 m_fakeMom_Pass -> Fill(m_momentum);
181 }
182
183 if (m_isBeamBKG) {
184 m_bkgPhi-> Fill(m_phi);
185 m_bkgTheta-> Fill(m_theta);
186 m_bkgMom-> Fill(m_momentum);
187 }
188
189 m_klidAll->SetMinimum(0.);
190 m_Mom_all_plot->SetMinimum(0.);
191 m_effPhi->SetMinimum(0.);
192 m_effTheta->SetMinimum(0.);
193 m_effMom->SetMinimum(0.);
194 m_fakePhi->SetMinimum(0.);
195 m_fakeTheta->SetMinimum(0.);
196 m_fakeMom->SetMinimum(0.);
197 m_time->SetMinimum(0.);
198 m_trackSep->SetMinimum(0.);
199 m_nLayer->SetMinimum(0.);
200 m_bkgPhi->SetMinimum(0.);
201 m_bkgTheta->SetMinimum(0.);
202 m_bkgMom->SetMinimum(0.);
203 m_innermostLayer->SetMinimum(0.);
204 m_trackFlag->SetMinimum(0.);
205 m_ECLFlag->SetMinimum(0.);
206
207 //fill all to normalise later
208 m_Phi_all -> Fill(m_phi);
209 m_Theta_all -> Fill(m_theta);
210 m_Mom_all -> Fill(m_momentum);
211 m_Mom_all_plot -> Fill(m_momentum);
212
213 m_time ->Fill(time);
214 m_trackSep->Fill(trackSep);
215 m_nLayer->Fill(nLayer);
216 m_innermostLayer->Fill(cluster.getInnermostLayer());
217 m_trackFlag->Fill(cluster.getAssociatedTrackFlag());
218 m_ECLFlag->Fill(cluster.getAssociatedEclClusterFlag());
219
220 if (klidObj) {
221 if (!m_isKl) {
222 m_klidFake ->Fill(klidObj->getKlId());
223 }
224 if (m_isKl) {
225 m_klidTrue ->Fill(klidObj->getKlId());
226 }
227 m_klidAll ->Fill(klidObj->getKlId());
228 }
229 }// for klm clusters
230
231}
232
234{
235 // write TEff to root file ,
236 m_f -> cd();
237 // TH1F is not compatible with the validation server
238 m_effPhi->Divide(m_Phi_Pass, m_Phi_all);
240 m_effMom->Divide(m_Mom_Pass, m_Mom_all);
241
245
246 std::vector<double> efficiency(m_xbins.size());
247 std::vector<double> purity(m_xbins.size());
248 std::vector<double> backrej(m_xbins.size());
249 double NallKlong = m_klidTrue->GetEntries();
250 double NallFakes = m_klidFake->GetEntries();
251 for (unsigned int i = 0; i < m_xbins.size() - 1; ++i) {
252 double NtruePositive = m_klidTrue->Integral(i, m_xbins.size() - 1);
253 double NfalsePositive = m_klidFake->Integral(i, m_xbins.size() - 1);
254 if ((NtruePositive + NfalsePositive) <= 0) {
255 purity[i] = 0;
256 } else {
257 purity[i] = NtruePositive / (NtruePositive + NfalsePositive);
258 }
259 if (NallFakes) {
260 backrej[i] = 1 - NfalsePositive / NallFakes;
261 } else {
262 backrej[i] = 0;
263 }
264 if (NallKlong <= 0) {
265 B2WARNING("No Klongs found in ROC calculation.");
266 } else {
267 efficiency[i] = NtruePositive / NallKlong;
268 }
269 }
270
271 m_ROC = new TGraph(m_xbins.size(), efficiency.data(), purity.data());
272 m_backRej = new TGraph(m_xbins.size(), efficiency.data(), backrej.data());
273
274 m_ROC->SetMarkerStyle(3);
275 m_backRej->SetMarkerStyle(3);
276 m_ROC->SetLineColorAlpha(kBlue, 0.0);
277 m_backRej->SetLineColorAlpha(kBlue, 0.0);
278 m_ROC->GetXaxis()->SetTitle("Efficiency");
279 m_backRej->GetXaxis()->SetTitle("Efficiency");
280 m_ROC->GetYaxis()->SetTitle("Purity");
281 m_backRej->GetYaxis()->SetTitle("Background rejection");
282 m_ROC-> SetTitle("Klong Purity vs Efficiency");
283 m_backRej-> SetTitle("Klong background rejection vs Efficiency");
284 m_ROC-> GetListOfFunctions() -> Add(new TNamed("Description",
285 "Purity vs Efficiency each point represents a cut on the klong ID."));
286 m_ROC -> GetListOfFunctions() -> Add(new TNamed("Check", "Should be as high as possible"));
287 m_backRej -> GetListOfFunctions() -> Add(new TNamed("Check", "Should be as high as possible"));
288 m_ROC -> GetListOfFunctions() -> Add(new TNamed("Contact", m_contact));
289 m_backRej -> GetListOfFunctions() -> Add(new TNamed("Contact", m_contact));
290
291 // tuple: pointer to the plot, name of the plot, true for shifter plots
292 std::vector<std::tuple<TH1F*, std::string, std::string, bool>> histograms;
293 std::string defaultCheck{"Nightly result should not differ significantly from the reference"};
294 histograms.push_back(make_tuple(m_klidAll, "KlId distribution", defaultCheck, true));
295 histograms.push_back(make_tuple(m_Mom_all_plot, "All Momenta generated", defaultCheck, false));
296 histograms.push_back(make_tuple(m_effPhi, "KlId efficiency in Phi", defaultCheck, false));
297 histograms.push_back(make_tuple(m_effTheta, "KlId efficiency in Theta",
298 defaultCheck + "; efficiency must be 0 below 0.3 rad and above 2.7 rad", true));
299 histograms.push_back(make_tuple(m_effMom, "KlId efficiency in Momentum", defaultCheck, true));
300 histograms.push_back(make_tuple(m_fakePhi, "KlId fake rate in Phi", defaultCheck, false));
301 histograms.push_back(make_tuple(m_fakeTheta, "KlId fake rate in Theta",
302 defaultCheck + "; fake rate must be 0 below 0.3 rad and above 2.7 rad", true));
303 histograms.push_back(make_tuple(m_fakeMom, "KlId fake rate in Momentum", defaultCheck, true));
304 histograms.push_back(make_tuple(m_time, "KLMClusters time", defaultCheck + "; no secondary peaks", true));
305 histograms.push_back(make_tuple(m_trackSep, "Distance of KLMCluster to closest Track", defaultCheck, true));
306 histograms.push_back(make_tuple(m_nLayer, "KLMClusters nLayer", defaultCheck, false));
307 histograms.push_back(make_tuple(m_innermostLayer, "KLMClusters innermostLayer", defaultCheck, false));
308 histograms.push_back(make_tuple(m_bkgPhi, "Beam bkg in Phi", defaultCheck, false));
309 histograms.push_back(make_tuple(m_bkgTheta, "Beam bkg in Theta",
310 defaultCheck + "; background must be 0 below 0.3 rad and above 2.7 rad", false));
311 histograms.push_back(make_tuple(m_bkgMom, "Beam bkg in Momentum", defaultCheck, false));
312 histograms.push_back(make_tuple(m_trackFlag, "Track flag", defaultCheck, false));
313 histograms.push_back(make_tuple(m_ECLFlag, "ECLCluster flag", defaultCheck, false));
314
315 for (auto hist : histograms) {
316 std::get<0>(hist) -> SetTitle(std::get<1>(hist).c_str());
317 std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("Description", std::get<1>(hist)));
318 std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("Check", std::get<2>(hist).c_str()));
319 std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("Contact", m_contact));
320 if (std::get<3>(hist))
321 std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("MetaOptions", "shifter,pvalue-warn=0.99,pvalue-error=0.02"));
322 else
323 std::get<0>(hist) -> GetListOfFunctions() -> Add(new TNamed("MetaOptions", "pvalue-warn=0.99,pvalue-error=0.1"));
324 std::get<0>(hist) -> Write();
325 }
326
327 //this guy is a tgraph
328 m_ROC -> Write();
329 m_backRej -> Write();
330 m_f -> Close();
331}
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType Klong
K^0_L particle.
Definition: Const.h:678
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
std::string m_contact
contact email address
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.
STL namespace.