ATLAS Offline Software
Loading...
Searching...
No Matches
MuonHoughDataNtuple.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
8#include <unordered_set>
9
10using namespace Muon::MuonStationIndex;
13 truthPart{truthPart_}{
14 for (const std::string hitColl : {"truthMdtHits", "truthRpcHits", "truthTgcHits", "truthCscHits"}) {
16 if (!acc.isAvailable(*truthPart)) continue;
17 for (const long long unsigned id : acc(*truthPart)){
18 assocHits.emplace(id);
19 }
20 }
21 }
23 std::unordered_set<Identifier> assocHits{};
24};
25
26MuonHoughDataNtuple::MuonHoughDataNtuple(const std::string& name, ISvcLocator* pSvcLocator) :
27 AthHistogramAlgorithm(name, pSvcLocator) {}
28
30
32{
34 ATH_CHECK(m_truthMuonKey.initialize());
35 ATH_CHECK(m_truthSegmentsKey.initialize());
36 m_truth_tree.addBranch(std::make_unique<MuonVal::EventInfoBranch>(m_truth_tree, 0));
37 m_eta_hit_tree.addBranch(std::make_unique<MuonVal::EventInfoBranch>(m_eta_hit_tree, 0));
38 m_phi_hit_tree.addBranch(std::make_unique<MuonVal::EventInfoBranch>(m_phi_hit_tree, 0));
39
40
41 ATH_CHECK(m_truth_tree.init(this));
42
43 ATH_CHECK(m_eta_hit_tree.init(this));
44 ATH_CHECK(m_phi_hit_tree.init(this));
45 ATH_CHECK(m_idHelperSvc.retrieve());
46 ATH_MSG_DEBUG("MuonHoughDataNtuple succesfully initialised");
47
48 return StatusCode::SUCCESS;
49}
50
52{
53 ATH_CHECK(m_truth_tree.write());
56 return StatusCode::SUCCESS;
57}
58
60{
61 // retrieve containers
62 const EventContext & context = Gaudi::Hive::currentContext();
64
66
67 if(!houghSecVec.isValid()){
68 ATH_MSG_FATAL("Failed to retrieve Hough data per sector vector");
69 return StatusCode::FAILURE;
70 }
71
72 std::vector<HitTruthMatching> truthMuons{};
73
75
76 if(!truthMuonContainer.isValid()){
77 ATH_MSG_FATAL("Failed to retrieve truth muon container");
78 return StatusCode::FAILURE;
79 }
80 for (const xAOD::TruthParticle* truth_mu : *truthMuonContainer) {
81 truthMuons.emplace_back(truth_mu);
82 }
83
84 if(!truthSegContainer.isValid()){
85 ATH_MSG_FATAL("Failed to retrieve truth segment container");
86 return StatusCode::FAILURE;
87 }
88 ATH_MSG_DEBUG("Houghidipuff retrieved all input collections! Start puffing...");
89
90 // filling variables for eta/phi hit with maximum entries for hough
91 for(unsigned int sec_i=0; sec_i<houghSecVec->vec.size(); ++sec_i){ // looping sector
92 // Eta Hit
93 // starting from all layers
94 for (const std::vector<std::shared_ptr<MuonHough::MuonLayerHough::Maximum>>& maxvec :
95 houghSecVec->vec[sec_i].maxVec){
96 // then to all Maxima in the same layer
97 for(const std::shared_ptr<MuonHough::MuonLayerHough::Maximum>& max : maxvec){ // looping maximum collection
98 m_maxHit_sector = sec_i;
99 m_maxHit_z0 = max->pos;
100 m_maxHit_theta = max->theta;
101 m_maxHit_region = toInt(max->refregion);
102 bool firstHit = true; // to store station information
103 // now record all hits in the maximum
104 int hit_index = 0;
105 for (const std::shared_ptr<MuonHough::Hit>& hit : max->hits){
106 m_hit_x.push_back(hit->x);
107 m_hit_ymin.push_back(hit->ymin);
108 m_hit_ymax.push_back(hit->ymax);
109 m_hit_w.push_back(hit->w);
110 // hit debug info
111 const MuonHough::HitDebugInfo* hdi = hit->debugInfo();
112 m_hit_tech.push_back(hdi->type);
113 // more information from prd
114 if(hit->prd){
115 const Trk::PrepRawData* prd = hit->prd;
116 const Identifier measId = prd->identify();
117 // station info
118 if(firstHit){
119 m_maxHit_stationIndex = toInt(m_idHelperSvc->stationIndex(measId));
120 m_maxHit_stationEta = m_idHelperSvc->stationEta(measId);
121 m_maxHit_stationPhi = m_idHelperSvc->stationPhi(measId);
122 firstHit = false;
123 }
124 // local position
125 m_hit_local.push_back(prd->localPosition());
126 // global position -- which is the real one???
127 const Muon::MdtPrepData* DC{dynamic_cast<const Muon::MdtPrepData*>(prd)};
128 if (DC) m_hit_global.push_back(DC->globalPosition());
129 const Muon::MuonCluster* CL{dynamic_cast<const Muon::MuonCluster*>(prd)};
130 if (CL) m_hit_global.push_back(CL->globalPosition());
131 // information from identifiers
132 if (m_idHelperSvc->isMdt(measId)) {
133 m_hit_mdtId.push_back(measId);
134 m_hit_mdtIndex.push_back(hit_index);
135 }
136 else if (m_idHelperSvc->isCsc(measId)) {
137 m_hit_cscId.push_back(measId);
138 m_hit_cscIndex.push_back(hit_index);
139 }
140 else if (m_idHelperSvc->isTgc(measId)) {
141 m_hit_tgcId.push_back(measId);
142 m_hit_tgcIndex.push_back(hit_index);
143 }
144 else if (m_idHelperSvc->isRpc(measId)) {
145 m_hit_rpcId.push_back(measId);
146 m_hit_rpcIndex.push_back(hit_index);
147 }
148 else if (m_idHelperSvc->isMM(measId)) {
149 m_hit_mmId.push_back(measId);
150 m_hit_mmIndex.push_back(hit_index);
151 }
152 else if (m_idHelperSvc->issTgc(measId)) {
153 m_hit_stgcId.push_back(measId);
154 m_hit_stgcIndex.push_back(hit_index);
155 }
156 hit_index++;
157 // truth matching to muons
158 bool matched = std::find_if(truthMuons.begin(), truthMuons.end(),
159 [prd](const HitTruthMatching& hit){
160 return hit.assocHits.count(prd->identify());
161 }) != truthMuons.end();
162
163 m_hit_truthMatched.push_back(matched);
164 }
165 }
166 ATH_CHECK(m_eta_hit_tree.fill(context));
167 }
168 }
169 // PhiHit vec
170 // starting from all layers
171 for (const std::vector<std::shared_ptr<MuonHough::MuonPhiLayerHough::Maximum>>& phimaxvec :
172 houghSecVec->vec[sec_i].phiMaxVec){
173 // then to all Maxima in the same layer
174 for(const std::shared_ptr<MuonHough::MuonPhiLayerHough::Maximum>& max : phimaxvec){ // looping maximum collection
175 m_maxPhiHit_sector = sec_i;
176 m_maxPhiHit_z0 = max->pos;
177 bool firstHit = true; // to store station information
178 // now record all hits in the maximum
179 int phiHit_index = 0;
180 for (const std::shared_ptr<MuonHough::PhiHit>& phi_hit : max->hits){
181 m_phiHit_x.push_back(phi_hit->r);
182 m_phiHit_ymin.push_back(phi_hit->phimin);
183 m_phiHit_ymax.push_back(phi_hit->phimax);
184 m_phiHit_w.push_back(phi_hit->w);
185 // phi hit debug info
186 const MuonHough::HitDebugInfo* hdi = phi_hit->debugInfo();
187 m_phiHit_tech.push_back(hdi->type);
188 // information from prd
189 if(phi_hit->prd){
190 const Trk::PrepRawData* prd = phi_hit->prd;
191 const Identifier measId = prd->identify();
192 // station info
193 if(firstHit){
194 m_maxPhiHit_stationIndex = toInt(m_idHelperSvc->stationIndex(prd->identify()));
195 m_maxPhiHit_stationEta = m_idHelperSvc->stationEta(prd->identify());
196 m_maxPhiHit_stationPhi = m_idHelperSvc->stationPhi(prd->identify());
197 firstHit = false;
198 }
199 // local position
200 m_phiHit_local.push_back(prd->localPosition());
201 // global position -- which is the real one???
202 const Muon::MdtPrepData* DC{dynamic_cast<const Muon::MdtPrepData*>(prd)};
203 if (DC) m_phiHit_global.push_back(DC->globalPosition());
204 const Muon::MuonCluster* CL{dynamic_cast<const Muon::MuonCluster*>(prd)};
205 if (CL) m_phiHit_global.push_back(CL->globalPosition());
206 // information from identifiers
207 if (m_idHelperSvc->isMdt(measId)) {
208 m_phiHit_mdtId.push_back(measId);
209 m_phiHit_mdtIndex.push_back(phiHit_index);
210 }
211 else if (m_idHelperSvc->isCsc(measId)) {
212 m_phiHit_cscId.push_back(measId);
213 m_phiHit_cscIndex.push_back(phiHit_index);
214 }
215 else if (m_idHelperSvc->isTgc(measId)) {
216 m_phiHit_tgcId.push_back(measId);
217 m_phiHit_tgcIndex.push_back(phiHit_index);
218 }
219 else if (m_idHelperSvc->isRpc(measId)) {
220 m_phiHit_rpcId.push_back(measId);
221 m_phiHit_rpcIndex.push_back(phiHit_index);
222 }
223 else if (m_idHelperSvc->isMM(measId)) {
224 m_phiHit_mmId.push_back(measId);
225 m_phiHit_mmIndex.push_back(phiHit_index);
226 }
227 else if (m_idHelperSvc->issTgc(measId)) {
228 m_phiHit_stgcId.push_back(measId);
229 m_phiHit_stgcIndex.push_back(phiHit_index);
230 }
231 phiHit_index++;
232
233 // truth matching to muons
234 bool matched = std::find_if(truthMuons.begin(), truthMuons.end(),
235 [prd](const HitTruthMatching& hit){
236 return hit.assocHits.count(prd->identify());
237 }) != truthMuons.end();
238
239 m_phiHit_truthMatched.push_back(matched);
240 }
241 }
242 ATH_CHECK(m_phi_hit_tree.fill(context));
243 }
244 }
245 }
246
247 // filling truth values
248 for(const xAOD::TruthParticle* truthMu: *truthMuonContainer){
249 m_truth_pdgId = truthMu->pdgId();
251
252 m_truth_pt = truthMu->pt();
253 m_truth_eta = truthMu->eta();
254 m_truth_phi = truthMu->phi();
255 // access the truth segment from truth particle
256 static const SG::AuxElement::ConstAccessor<ElementLink<xAOD::TruthParticleContainer>> truthMuLink("truthParticleLink");
257 // filling truth segment values
258 for(const xAOD::MuonSegment* truthSeg: *truthSegContainer){
259 if(!truthMuLink.isAvailable(*truthSeg)) { continue; } // if segment isn't linked to truth muon, skip
260 ElementLink<xAOD::TruthParticleContainer> truthLink = truthMuLink(*truthSeg);
261 const xAOD::TruthParticle* truthParticle = *truthLink;
262 if(truthParticle!=truthMu) { continue; } // i have no idea if this works
263 m_truth_seg_pos.push_back(truthSeg->x(), truthSeg->y(), truthSeg->z());
264 m_truth_seg_p.push_back(truthSeg->px(), truthSeg->py(), truthSeg->pz());
265
266 m_truth_seg_nPrecisionHits.push_back(truthSeg->nPrecisionHits());
267 m_truth_seg_sector.push_back(truthSeg->sector());
268
269 int nTriggerHits = truthSeg->nPhiLayers() + truthSeg->nTrigEtaLayers();
270 m_truth_seg_nTriggerHits.push_back(nTriggerHits);
271 }
272 ATH_CHECK(m_truth_tree.fill(context));
273 }
274
275
276 return StatusCode::SUCCESS;
277}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
#define max(a, b)
Definition cfImp.cxx:41
AthHistogramAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
MuonVal::ThreeVectorBranch m_hit_global
MuonVal::TgcIdentifierBranch m_phiHit_tgcId
SG::ReadHandleKey< xAOD::MuonSegmentContainer > m_truthSegmentsKey
MuonVal::VectorBranch< float > & m_hit_ymin
MuonVal::VectorBranch< int > & m_hit_tech
MuonVal::CscIdentifierBranch m_phiHit_cscId
MuonVal::TgcIdentifierBranch m_hit_stgcId
MuonVal::VectorBranch< int > & m_truth_seg_nTriggerHits
MuonVal::VectorBranch< float > & m_phiHit_ymin
MuonVal::VectorBranch< int > & m_truth_seg_sector
MuonVal::ScalarBranch< int > & m_maxHit_stationIndex
MuonVal::VectorBranch< int > & m_hit_tgcIndex
MuonVal::MuonTesterTree m_eta_hit_tree
MuonVal::ScalarBranch< float > & m_maxHit_sector
MuonVal::TgcIdentifierBranch m_phiHit_stgcId
MuonVal::MuonTesterTree m_phi_hit_tree
MuonVal::TwoVectorBranch m_phiHit_local
MuonVal::VectorBranch< float > & m_phiHit_ymax
MuonVal::ScalarBranch< float > & m_maxPhiHit_z0
MuonVal::ScalarBranch< int > & m_maxHit_stationEta
MuonVal::CscIdentifierBranch m_hit_cscId
MuonVal::VectorBranch< int > & m_phiHit_mmIndex
MuonVal::RpcIdentifierBranch m_hit_rpcId
MuonVal::VectorBranch< int > & m_hit_stgcIndex
MuonVal::ScalarBranch< int > & m_maxPhiHit_stationEta
MuonVal::ScalarBranch< float > & m_maxPhiHit_sector
MuonVal::ScalarBranch< int > & m_maxHit_region
MuonVal::ScalarBranch< float > & m_maxHit_z0
SG::ReadHandleKey< Muon::HoughDataPerSectorVec > m_houghDataPerSectorVecKey
MuonVal::ThreeVectorBranch m_truth_seg_pos
ServiceHandle< Muon::IMuonIdHelperSvc > m_idHelperSvc
MuonVal::TwoVectorBranch m_hit_local
MuonVal::VectorBranch< int > & m_hit_mmIndex
MuonVal::VectorBranch< float > & m_hit_x
MuonVal::VectorBranch< int > & m_hit_mdtIndex
MuonVal::MdtIdentifierBranch m_phiHit_mdtId
MuonVal::VectorBranch< bool > & m_phiHit_truthMatched
MuonVal::ScalarBranch< float > & m_truth_pt
MuonVal::VectorBranch< int > & m_hit_cscIndex
MuonVal::VectorBranch< float > & m_phiHit_w
MuonVal::ScalarBranch< int > & m_truth_barcode
MuonVal::MdtIdentifierBranch m_hit_mdtId
MuonVal::VectorBranch< float > & m_phiHit_x
MuonVal::VectorBranch< bool > & m_hit_truthMatched
MuonVal::RpcIdentifierBranch m_phiHit_rpcId
MuonVal::VectorBranch< int > & m_hit_rpcIndex
MuonVal::TgcIdentifierBranch m_phiHit_mmId
MuonVal::VectorBranch< float > & m_hit_w
MuonVal::ScalarBranch< float > & m_truth_phi
MuonVal::TgcIdentifierBranch m_hit_mmId
virtual StatusCode finalize() override
MuonVal::VectorBranch< int > & m_phiHit_rpcIndex
MuonVal::VectorBranch< int > & m_phiHit_mdtIndex
MuonVal::ScalarBranch< int > & m_maxPhiHit_stationPhi
MuonHoughDataNtuple(const std::string &name, ISvcLocator *pSvcLocator)
MuonVal::VectorBranch< float > & m_hit_ymax
MuonVal::VectorBranch< int > & m_phiHit_tech
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthMuonKey
MuonVal::VectorBranch< int > & m_phiHit_stgcIndex
MuonVal::VectorBranch< int > & m_truth_seg_nPrecisionHits
MuonVal::ScalarBranch< int > & m_maxPhiHit_stationIndex
MuonVal::TgcIdentifierBranch m_hit_tgcId
MuonVal::ScalarBranch< int > & m_truth_pdgId
MuonVal::VectorBranch< int > & m_phiHit_cscIndex
virtual StatusCode execute() override
MuonVal::ThreeVectorBranch m_phiHit_global
virtual ~MuonHoughDataNtuple() override
virtual StatusCode initialize() override
MuonVal::VectorBranch< int > & m_phiHit_tgcIndex
MuonVal::ThreeVectorBranch m_truth_seg_p
MuonVal::MuonTesterTree m_truth_tree
MuonVal::ScalarBranch< float > & m_maxHit_theta
MuonVal::ScalarBranch< float > & m_truth_eta
MuonVal::ScalarBranch< int > & m_maxHit_stationPhi
struct containing additional debug information on the hits that is not needed for the actual alg but ...
Class to represent measurements from the Monitored Drift Tubes.
Definition MdtPrepData.h:33
virtual const Amg::Vector3D & globalPosition() const
Returns the global position of the CENTER of the drift tube (i.e.
virtual const Amg::Vector3D & globalPosition() const =0
Returns the global position of the measurement (calculated on the fly)
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:569
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const Amg::Vector2D & localPosition() const
return the local position reference
Identifier identify() const
return the identifier
int uniqueID(const T &p)
constexpr int toInt(const EnumType enumVal)
TruthParticle_v1 TruthParticle
Typedef to implementation.
MuonSegment_v1 MuonSegment
Reference the current persistent version:
std::unordered_set< Identifier > assocHits
HitTruthMatching(const xAOD::TruthParticle *truthPart_)
const xAOD::TruthParticle * truthPart