Belle II Software development
KlId.h
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#pragma once
10
11#include <framework/datastore/StoreArray.h>
12#include <mdst/dataobjects/MCParticle.h>
13#include <mdst/dataobjects/KLMCluster.h>
14#include <mdst/dataobjects/ECLCluster.h>
15#include <mdst/dataobjects/TrackFitResult.h>
16
17#include <framework/gearbox/Const.h>
18#include <framework/geometry/VectorUtil.h>
19#include <framework/logging/Logger.h>
20#include <tracking/dataobjects/RecoTrack.h>
21#include <genfit/Exception.h>
22#include <utility>
23#include <math.h>
24#include <Math/VectorUtil.h>
25
27namespace Belle2::KlongId {
28
30 int BelleTrackFlag(const Belle2::KLMCluster& cluster, const float angle = 0.26)
31 {
32 const ROOT::Math::XYZVector& pos = cluster.getClusterPosition();
33
35 for (const Belle2::TrackFitResult& track : tracks) {
36 const ROOT::Math::XYZVector& trackPos = track.getPosition();
37
38 double trackAngle = ROOT::Math::VectorUtil::Angle(trackPos, pos);
39 if (trackAngle < angle) {
40 B2DEBUG(20, "BelleFlagTracklAngle::" << trackAngle);
41 return 1;
42 }
43 }
44 return 0;
45 }
46
47
49 int BelleECLFlag(const Belle2::KLMCluster& cluster, const float angle = 0.26)
50 {
51 const ROOT::Math::XYZVector& pos = cluster.getClusterPosition();
53
54 for (const Belle2::ECLCluster& eclcluster : eclclusters) {
55
56 const ROOT::Math::XYZVector& clusterPos = eclcluster.getClusterPosition();
57
58 double clusterAngle = ROOT::Math::VectorUtil::Angle(clusterPos, pos);
59 if (clusterAngle < angle) {
60 B2DEBUG(20, "BelleFlagECLAngle::" << clusterAngle);
61 return 1;
62 }
63 }
64 return 0;
65 }
66
69 {
70 if (part == nullptr) {
71 return 1;
72 } else {
73 return 0;
74 }
75 }
76
79 {
80 unsigned int hirachy_counter = 0;
81 if (mcParticleIsBeamBKG(part)) {
82 return 0;
83 }
84 bool stop = false;
85 while (!stop) {
86 ++hirachy_counter;
87 if (part -> getPDG() == Const::Klong.getPDGCode()) {
88 return hirachy_counter;
89 }
90 if ((part -> getMother() == nullptr)) {
91 stop = true;
92 } else {
93 part = part -> getMother();
94 }
95 }
96 return 0;
97 }
98
100 bool isKLMClusterSignal(const Belle2::KLMCluster& cluster, float mcWeigthCut = 0.66)
101 {
102 const auto mcParticleWeightPair = cluster.getRelatedToWithWeight<Belle2::MCParticle>();
103 Belle2::MCParticle* part = mcParticleWeightPair.first;
104 if (!part) {return false; }
105 float mcWeight = mcParticleWeightPair.second;
106 if (mcParticleIsKlong(part) && (mcWeight > mcWeigthCut)) {
107 return true;
108 } else {
109 return false;
110 }
111 }
112
114 bool isECLClusterSignal(const Belle2::ECLCluster& cluster, float mcWeigthCut = 0.66)
115 {
116 const auto mcParticleWeightPair = cluster.getRelatedToWithWeight<Belle2::MCParticle>();
117 Belle2::MCParticle* part = mcParticleWeightPair.first;
118 if (!part) {return false; }
119 float mcWeight = mcParticleWeightPair.second;
120 if (mcParticleIsKlong(part) && (mcWeight > mcWeigthCut)) {
121 return true;
122 } else {
123 return false;
124 }
125 }
126
129 {
130 bool stop = false;
131 while (!stop) {
132 if (std::abs(part -> getPDG()) == pdg) {
133 return true;
134 }
135 if ((part -> getMother() == nullptr)) {
136 stop = true;
137 } else {
138 part = part -> getMother();
139 }
140 }
141 return false;
142 }
143
144
152 {
153
154 if (mcParticleIsBeamBKG(part)) {
155 return -999;
156 }
157
158 bool stop = false;
159 while (!stop) {
160 if (isMCParticlePDG(part, Const::Klong.getPDGCode())) {
161 return Const::Klong.getPDGCode();
162 }
163 if (isMCParticlePDG(part, Const::Kshort.getPDGCode())) {
164 return Const::Kshort.getPDGCode();
165 }
166 if (isMCParticlePDG(part, Const::kaon.getPDGCode())) {
167 return Const::kaon.getPDGCode();
168 }
169 if (isMCParticlePDG(part, Const::neutron.getPDGCode())) {
170 return Const::neutron.getPDGCode();
171 }
172 if (isMCParticlePDG(part, Const::proton.getPDGCode())) {
173 return Const::proton.getPDGCode();
174 }
175 if (isMCParticlePDG(part, Const::pion.getPDGCode())) {
176 return Const::pion.getPDGCode();
177 }
178 if (isMCParticlePDG(part, Const::pi0.getPDGCode())) {
179 return Const::pi0.getPDGCode();
180 }
181 if (isMCParticlePDG(part, Const::muon.getPDGCode())) {
182 return Const::muon.getPDGCode();
183 }
184 if (isMCParticlePDG(part, Const::electron.getPDGCode())) {
186 }
187 if (isMCParticlePDG(part, Const::photon.getPDGCode())) {
188 return Const::photon.getPDGCode();
189 }
190 if ((part -> getMother() == nullptr)) {
191 stop = true;
192 } else {
193 part = part -> getMother();
194 }
195 }
196 return part ->getPDG();
197 }
198
203 std::pair<Belle2::ECLCluster*, double> findClosestECLCluster(const ROOT::Math::XYZVector& klmClusterPosition,
205 {
206
207 Belle2::ECLCluster* closestECL = nullptr;
208 double closestECLAngleDist = 1e10;
210
211 if (eclClusters.getEntries() > 0) {
212 int index = 0;
213 int indexOfClosestCluster = -1;
214 for (Belle2::ECLCluster& eclcluster : eclClusters) {
215
216 if (eclcluster.hasHypothesis(eclhypothesis)) {
217
218 const ROOT::Math::XYZVector& eclclusterPos = eclcluster.getClusterPosition();
219 double angularDist = ROOT::Math::VectorUtil::Angle(
220 eclclusterPos, klmClusterPosition);
221 if (angularDist < closestECLAngleDist) {
222 closestECLAngleDist = angularDist;
223 indexOfClosestCluster = index;
224 }
225 }
226 ++index;
227 }
228 if (indexOfClosestCluster > -1)
229 closestECL = eclClusters[indexOfClosestCluster];
230 }
231 return std::make_pair(closestECL, closestECLAngleDist);
232 }
233
234
236 std::tuple<const Belle2::KLMCluster*, double, double> findClosestKLMCluster(const ROOT::Math::XYZVector& klmClusterPosition)
237 {
238
240 const Belle2::KLMCluster* closestKLM = nullptr;
241 double closestKLMDist = 1e10;
242 double avInterClusterDist = 0;
243 double nKLMCluster = klmClusters.getEntries();
244
245 if (nKLMCluster > 1) {
246
247 unsigned int index = 0;
248 unsigned int indexOfClosestCluster = 0;
249 for (const Belle2::KLMCluster& nextCluster : klmClusters) {
250
251 const ROOT::Math::XYZVector& nextClusterPos = nextCluster.getClusterPosition();
252 const ROOT::Math::XYZVector& clustDistanceVec = nextClusterPos - klmClusterPosition;
253
254 double nextClusterDist = clustDistanceVec.Mag2();
255 avInterClusterDist = avInterClusterDist + nextClusterDist;
256
257 if ((nextClusterDist < closestKLMDist) and not(nextClusterDist == 0)) {
258 closestKLMDist = nextClusterDist ;
259 indexOfClosestCluster = index;
260 }
261 ++index;
262 }// for next_cluster
263
264 closestKLM = klmClusters[indexOfClosestCluster];
265 avInterClusterDist = avInterClusterDist / (1. * nKLMCluster);
266 }
267
268 return std::make_tuple(closestKLM, closestKLMDist, avInterClusterDist);
269 }
270
271
273 std::tuple<Belle2::RecoTrack*, double, std::unique_ptr<const ROOT::Math::XYZVector> > findClosestTrack(
274 const ROOT::Math::XYZVector& clusterPosition,
275 float cutAngle)
276 {
278 double oldDistance = 10000000;//dont wanna use infty cos that kills tmva...
279 Belle2::RecoTrack* closestTrack = nullptr;
280 ROOT::Math::XYZVector poca = ROOT::Math::XYZVector(0, 0, 0);
281
282
283 for (Belle2::RecoTrack& track : genfitTracks) {
284
285 try {
286 genfit::MeasuredStateOnPlane state;
287 genfit::MeasuredStateOnPlane state_for_cut;
288 state_for_cut = track.getMeasuredStateOnPlaneFromFirstHit();
289
290 // only use tracks that are close cos the extrapolation takes ages
291 double angle = ROOT::Math::VectorUtil::Angle(clusterPosition, state_for_cut.getPos());
292 if (angle < cutAngle) {
293 state = track.getMeasuredStateOnPlaneFromLastHit();
294 state.extrapolateToPoint(XYZToTVector(clusterPosition));
295 const ROOT::Math::XYZVector& trackPos = ROOT::Math::XYZVector(state.getPos());
296
297 const ROOT::Math::XYZVector& distanceVecCluster = clusterPosition - trackPos;
298 double newDistance = distanceVecCluster.Mag2();
299
300 // overwrite old distance
301 if (newDistance < oldDistance) {
302 oldDistance = newDistance;
303 closestTrack = &track;
304 poca = trackPos;
305 }
306 }
307 } catch (genfit::Exception& e) {
308 }// try
309 }// for gftrack
310
311 if (not closestTrack) {
312 return std::make_tuple(nullptr, oldDistance, std::unique_ptr<const ROOT::Math::XYZVector>(nullptr));
313 } else {
314 // actually this is fine because of the datastore
315 // cppcheck-suppress returnDanglingLifetime
316 return std::make_tuple(closestTrack, oldDistance, std::unique_ptr<const ROOT::Math::XYZVector>(new ROOT::Math::XYZVector(poca)));
317 }
318 }
319}//end namespace
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType neutron
neutron particle
Definition: Const.h:675
static const ParticleType pi0
neutral pion particle
Definition: Const.h:674
static const ChargedStable muon
muon particle
Definition: Const.h:660
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ParticleType Klong
K^0_L particle.
Definition: Const.h:678
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ParticleType Kshort
K^0_S particle.
Definition: Const.h:677
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:662
static const ParticleType photon
photon particle
Definition: Const.h:673
static const ChargedStable electron
electron particle
Definition: Const.h:659
ECL cluster data.
Definition: ECLCluster.h:27
EHypothesisBit
The hypothesis bits for this ECLCluster (Connected region (CR) is split using this hypothesis.
Definition: ECLCluster.h:31
@ c_neutralHadron
CR is reconstructed as a neutral hadron (N2)
KLM cluster data.
Definition: KLMCluster.h:28
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Values of the result of a track fit with a given particle hypothesis.
static constexpr auto XYZToTVector
Helper function to convert XYZVector to TVector3.
Definition: VectorUtil.h:24
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
int BelleECLFlag(const Belle2::KLMCluster &cluster, const float angle=0.26)
get Belle stle ECL flag
Definition: KlId.h:49
int BelleTrackFlag(const Belle2::KLMCluster &cluster, const float angle=0.26)
get Belle stle track flag
Definition: KlId.h:30
std::tuple< Belle2::RecoTrack *, double, std::unique_ptr< const ROOT::Math::XYZVector > > findClosestTrack(const ROOT::Math::XYZVector &clusterPosition, float cutAngle)
find nearest genfit track and return it and its distance
Definition: KlId.h:273
int mcParticleIsBeamBKG(const Belle2::MCParticle *part)
return if MCparticle is beambkg
Definition: KlId.h:68
std::tuple< const Belle2::KLMCluster *, double, double > findClosestKLMCluster(const ROOT::Math::XYZVector &klmClusterPosition)
find nearest KLMCluster, tis distance and the av intercluster distance
Definition: KlId.h:236
int isMCParticlePDG(Belle2::MCParticle *part, int pdg)
return if mc particle has a certain pdg in the decay chain
Definition: KlId.h:128
bool isECLClusterSignal(const Belle2::ECLCluster &cluster, float mcWeigthCut=0.66)
checks if a cluster is signal under the mcWeightcondition (mcWeight = energy deposition)
Definition: KlId.h:114
std::pair< Belle2::ECLCluster *, double > findClosestECLCluster(const ROOT::Math::XYZVector &klmClusterPosition, const Belle2::ECLCluster::EHypothesisBit eclhypothesis=Belle2::ECLCluster::EHypothesisBit::c_neutralHadron)
Find the closest ECLCluster with a neutral hadron hypothesis, and return it with its distance.
Definition: KlId.h:203
int getPrimaryPDG(Belle2::MCParticle *part)
return if mc particles primary pdg.
Definition: KlId.h:151
int mcParticleIsKlong(Belle2::MCParticle *part)
return the mc hirachy of the klong 0:not a klong 1:final particle, 2: klong is mother etc
Definition: KlId.h:78