ATLAS Offline Software
Loading...
Searching...
No Matches
HGTD::ClusterTruthTool Class Reference

#include <ClusterTruthTool.h>

Inheritance diagram for HGTD::ClusterTruthTool:
Collaboration diagram for HGTD::ClusterTruthTool:

Public Member Functions

 ClusterTruthTool (const std::string &, const std::string &, const IInterface *)
virtual HGTD::ClusterTruthInfo classifyCluster (const HGTD_Cluster *cluster, const xAOD::TruthParticle *tp, const InDetSimDataCollection *sim_data, const HepMC::GenEvent *hard_scatter_evnt=nullptr) const override final
 The InDetSimDataCollection is a map, connecting Identifiers from the RDOs to InDetSimData objects, that keep a vector of HepMcParticleLink, float pairs, where the float keeps the energy of the given energy deposit.

Detailed Description

Definition at line 25 of file ClusterTruthTool.h.

Constructor & Destructor Documentation

◆ ClusterTruthTool()

HGTD::ClusterTruthTool::ClusterTruthTool ( const std::string & t,
const std::string & n,
const IInterface * p )

Definition at line 23 of file ClusterTruthTool.cxx.

26 : base_class(t, n, p) {}

Member Function Documentation

◆ classifyCluster()

HGTD::ClusterTruthInfo HGTD::ClusterTruthTool::classifyCluster ( const HGTD_Cluster * cluster,
const xAOD::TruthParticle * tp,
const InDetSimDataCollection * sim_data,
const HepMC::GenEvent * hard_scatter_evnt = nullptr ) const
finaloverridevirtual

The InDetSimDataCollection is a map, connecting Identifiers from the RDOs to InDetSimData objects, that keep a vector of HepMcParticleLink, float pairs, where the float keeps the energy of the given energy deposit.

For HGTD, we store the time of the deposit instead, and can then select the first deposit, since only this one would contribute to the time read out by the ASIC.

Parameters
[in]clusterHit in HGTD, for which we want the truth information.
[in]tpTruth particle that is potentially matched to the cluster.
[in]sim_dataSDO collection storing the truth links of the charged diodes.
[in]hard_scatter_evntIf given, a cluster can be categorised as coming from the HS event, even if the direct match with the truth particle fails.
Returns
Struct combining the relevant truth information.

Definition at line 28 of file ClusterTruthTool.cxx.

31 {
32
33 if (not sim_data or not tp) {
34 return {HGTD::ClusterTruthOrigin::UNIDENTIFIED, false, false};
35 }
36
37 ATH_MSG_DEBUG("[ClusterTruthTool] " << sim_data->size() << " SDO elements");
38 ATH_MSG_DEBUG("[ClusterTruthTool] all SimData IDs: ");
39 for (const auto& elem : *sim_data) {
40 ATH_MSG_DEBUG(elem.first);
41 }
42
43 const std::vector<Identifier>& rdo_id_list = cluster->rdoList();
44 // keep record of the cluster origins and if they are shadowed or not
45 std::vector<std::pair<HGTD::ClusterTruthOrigin, bool>> shadowed_origins;
46
47 for (const auto& rdo_id : rdo_id_list) {
48
49 ATH_MSG_DEBUG("[ClusterTruthTool] looking for ID:" << rdo_id);
50
51 auto pos = sim_data->find(rdo_id);
52 // the InDetSimData contains a std::pair<HepMcParticleLink, float>, where
53 // the second entry in the pair holds the time of the SiChargedDiode
54 if (pos == sim_data->end()) {
56 "[HGTD::ClusterTruthTool::classifyCluster] ID not found in SDO "
57 "map, going to next ID");
58 continue;
59 }
60 // collect deposits, sorted with first deposit at start of map bu default
61 std::map<float, HGTD::ClusterTruthOrigin> sorted_deposits;
62
64 // the following is taken from 20.20 as is
65 ATH_MSG_DEBUG("[ClusterTruthTool] going through "
66 << pos->second.getdeposits().size() << " deposits");
67 for (const auto& deposit : pos->second.getdeposits()) {
68 const HepMcParticleLink& particle_link = deposit.first;
69
70 // check for identity with original particle
71 HepMC::ConstGenParticlePtr gen_part = particle_link.cptr();
72
73 // HepMcParticleLinks with uniqueIDs equal to
74 // HepMC::UNDEFINED_ID (used for detector noise, delta rays and
75 // random energy deposits) or secondaries generated by GEANT4.
76 if (HepMC::uniqueID(particle_link) == HepMC::UNDEFINED_ID ||
77 (gen_part && HepMC::is_simulation_particle(gen_part))) {
78 sorted_deposits.emplace(deposit.second, HGTD::ClusterTruthOrigin::SECONDARY);
79 continue;
80 }
81 if (gen_part) {
82 TLorentzVector l4(gen_part->momentum().px(), gen_part->momentum().py(),
83 gen_part->momentum().pz(), gen_part->momentum().e());
84 // flag if this deposit came from the tested truth particle
85 if (HepMC::is_same_particle(gen_part,tp) && tp->p4().DeltaR(l4) < 0.05) {
86 sorted_deposits.emplace(deposit.second,
88 // otherwise if given, the parent event can be checked
89 } else if (hard_scatter_evnt and
90 gen_part->parent_event() == hard_scatter_evnt) {
91 sorted_deposits.emplace(deposit.second,
93 // otherwise, a particle that was generated but doesn't come from the
94 // hard scatter is considered to be originating from a pileup
95 // interaction
96 } else {
97 sorted_deposits.emplace(deposit.second,
99 }
100 } else {
101 // if there is no gen particle, we can guess based on the event index
102 if (particle_link.eventIndex() == 0) {
103 sorted_deposits.emplace(deposit.second,
105 } else {
106 // If the gen was not kept, we assume it is cut away from truth record
107 // and is pileup
108 sorted_deposits.emplace(deposit.second,
110 }
111 }
112 } // END lOOP over the deposits
113
114 ATH_MSG_DEBUG("[ClusterTruthTool] " << sorted_deposits.size()
115 << " sorted_deposits");
116
117 // check for shadowing, which means that a deposit was left by the truth
118 // particle, but it was not the first deposit and is thus not used for the
119 // time measurement -> I will have an incorrect time
120 bool is_shadowed = false;
122
123 std::map<float, HGTD::ClusterTruthOrigin>::iterator elem =
124 sorted_deposits.begin();
125
126 for (; elem != sorted_deposits.end(); ++elem) {
127 if (elem == sorted_deposits.begin()) {
128 current_origin = elem->second;
129 } else {
130 // if one of the later deposits originates from the truth particle, then
131 // the hit I want to find was shadowed by something else, and I will
132 // reconstruct an incorrect time
133 if (current_origin != HGTD::ClusterTruthOrigin::TRUTH_PARTICLE &&
135 is_shadowed = true;
136 }
137 }
138 } // END LOOP over the time sorted deposits
139 shadowed_origins.emplace_back(current_origin, is_shadowed);
140 } // END LOOP over RDO identifiers
141
142 HGTD::ClusterTruthInfo result;
143
144 if (shadowed_origins.empty()) {
145 ATH_MSG_DEBUG("did not manage to understand any RDOs...");
147 result.is_shadowed = false;
148 // A cluster is considered to be merged if more than one particle deposited
149 // energy in a given pad.
150 result.is_merged = false;
151 } else {
152 result.is_merged = false;
153 result.origin = shadowed_origins.at(0).first;
154 result.is_shadowed = shadowed_origins.at(0).second;
155 for (size_t i = 1; i < shadowed_origins.size(); ++i) {
156 if (shadowed_origins.at(i).first != result.origin) {
157 result.is_merged = true;
158 } else {
159 result.is_shadowed |= shadowed_origins.at(i).second;
160 }
161 }
162 }
163 return result;
164}
#define ATH_MSG_DEBUG(x)
const std::vector< Identifier > & rdoList() const
return the List of rdo identifiers (pointers)
int uniqueID(const T &p)
constexpr int UNDEFINED_ID
bool is_same_particle(const T1 &p1, const T2 &p2)
Method to establish if two particles in the GenEvent actually represent the same particle.
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38

The documentation for this class was generated from the following files: