ATLAS Offline Software
Loading...
Searching...
No Matches
HGTD_RDOAnalysis.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "HGTD_RDOAnalysis.h"
9
11
12#include "TTree.h"
13#include "TString.h"
14
15#include <algorithm>
16#include <math.h>
17#include <functional>
18#include <iostream>
19
20
21
23 ATH_MSG_DEBUG( "Initializing HGTD_RDOAnalysis" );
24
25 // This will check that the properties were initialized
26 // properly by job configuration.
27 ATH_CHECK( m_inputKey.initialize() );
28 ATH_CHECK( m_inputTruthKey.initialize() );
30
31 // Grab HGTDID helper
32 ATH_CHECK(detStore()->retrieve(m_HGTD_Manager, m_HGTD_Name.value()));
33 ATH_CHECK(detStore()->retrieve(m_HGTD_ID, m_HGTDID_Name.value()));
34
35 m_tree = new TTree(m_ntupleName.value().c_str(), "HGTD_RDOAnalysis");
36 ATH_CHECK(histSvc()->regTree(m_ntuplePath.value() + m_ntupleName.value(), m_tree));
37 m_tree->Branch("m_rdo_module_layer", &m_rdo_module_layer);
38 m_tree->Branch("m_rdo_module_x", &m_rdo_module_x);
39 m_tree->Branch("m_rdo_module_y", &m_rdo_module_y);
40 m_tree->Branch("m_rdo_module_z", &m_rdo_module_z);
41 m_tree->Branch("m_rdo_module_ID", &m_rdo_module_ID);
42 m_tree->Branch("m_rdo_hit_x", &m_rdo_hit_x);
43 m_tree->Branch("m_rdo_hit_y", &m_rdo_hit_y);
44 m_tree->Branch("m_rdo_hit_z", &m_rdo_hit_z);
45 m_tree->Branch("m_rdo_hit_toa", &m_rdo_hit_toa);
46 m_tree->Branch("m_rdo_hit_sdo_toa", &m_rdo_hit_sdo_toa);
47 m_tree->Branch("m_rdo_hit_sdo_truth_category", &m_rdo_hit_sdo_truth_category);
48 /*
49 // HISTOGRAMS
50
52 m_h_rdoID = new TH1F("h_rdoID", "rdoID", 100, 0, 10e17);
53 m_h_rdoID->StatOverflows();
54 ATH_CHECK(histSvc()->regHist(m_histPath + m_h_rdoID->GetName(), m_h_rdoID));
55 */
56
57 return StatusCode::SUCCESS;
58}
59
61 ATH_MSG_DEBUG(" In HGTD_RDOAnalysis::execute()" );
62 m_rdo_module_layer.clear();
63 m_rdo_module_x.clear();
64 m_rdo_module_y.clear();
65 m_rdo_module_z.clear();
66 m_rdo_module_ID.clear();
67 m_rdo_hit_x.clear();
68 m_rdo_hit_y.clear();
69 m_rdo_hit_z.clear();
70 m_rdo_hit_toa.clear();
71 m_rdo_hit_sdo_toa.clear();
73
74 const EventContext& ctx{Gaudi::Hive::currentContext()};
75 // Raw HGTD Data
76 const HGTD_RDO_Container* p_RDO_cont{nullptr};
77 ATH_CHECK(SG::get(p_RDO_cont, m_inputKey, ctx));
78
79 //Adding SimMap and McEvent here for added truthMatching checks
80 const InDetSimDataCollection* simDataMap{nullptr};
81 ATH_CHECK(SG::get(simDataMap, m_inputTruthKey, ctx));
82
83 const McEventCollection* mcEventCollection{nullptr};
84 ATH_CHECK(SG::get(mcEventCollection, m_inputMcEventCollectionKey, ctx));
85
86 bool doTruthMatching = true;
87 const HepMC::GenEvent* hardScatterEvent(nullptr);
88
89 if (mcEventCollection->size()==0){
90 ATH_MSG_WARNING("Failed to retrieve a nonzero sized truth event collection, disabling truthMatching");
91 doTruthMatching = false;
92 }
93 if(doTruthMatching) hardScatterEvent = mcEventCollection->at(0);
94
95 // loop over RDO container
96 for ( HGTD_RDO_Container::const_iterator rdoCont_itr = p_RDO_cont->begin(); rdoCont_itr != p_RDO_cont->end(); ++rdoCont_itr ) {
97
98 const HGTD_RDO_Collection* p_RDO_coll(*rdoCont_itr);
99 const Identifier rdoIDColl((*rdoCont_itr)->identify());
100
101 const InDetDD::HGTD_DetectorElement *elementColl = m_HGTD_Manager->getDetectorElement(rdoIDColl);
102 InDetDD::SiLocalPosition localPosColl = elementColl->rawLocalPositionOfCell(rdoIDColl);
103 Amg::Vector3D globalPosColl = elementColl->globalPosition(localPosColl);
104
105
106 for ( HGTD_RDO_Collection::const_iterator rdo_itr= p_RDO_coll->begin(); rdo_itr != p_RDO_coll->end(); ++rdo_itr ) {
107
108 const Identifier rdoID((*rdo_itr)->identify());
109
110 const InDetDD::HGTD_DetectorElement *rdo_element = m_HGTD_Manager->getDetectorElement(rdoID);
111 InDetDD::SiLocalPosition localPos_hit = rdo_element->rawLocalPositionOfCell(rdoID);
112 Amg::Vector3D globalPos_hit = rdo_element->globalPosition(localPos_hit);
113
114 m_rdo_hit_x.push_back(globalPos_hit[Amg::x]);
115 m_rdo_hit_y.push_back(globalPos_hit[Amg::y]);
116 m_rdo_hit_z.push_back(globalPos_hit[Amg::z]);
117
118 // Move global position storage for RDO module from outer to inner loop
119 m_rdo_module_x.push_back(globalPosColl[Amg::x]);
120 m_rdo_module_y.push_back(globalPosColl[Amg::y]);
121 m_rdo_module_z.push_back(globalPosColl[Amg::z]);
122
123 m_rdo_module_layer.push_back(m_HGTD_ID->layer(rdoIDColl));
124 m_rdo_module_ID.push_back(rdoIDColl.get_compact());
125
126 // Get the Time of Arrival which include the tof and the digitisation smearing of 25 ps
127 float rdo_toa = (*rdo_itr)->getTOA();
128 m_rdo_hit_toa.push_back(rdo_toa);
129
130 // For truth matching studies we need to get the truth nature of each deposit.
131 // For each hit, there can be up to 7 deposit from different nature accessed through the SDO map.
132 if(doTruthMatching){
133 if(simDataMap) {
134 InDetSimDataCollection::const_iterator iter = (*simDataMap).find((*rdo_itr)->identify());
135 std::vector<SdoInfo> sdo_info;
136 if ( iter != (*simDataMap).end() ) {
137 const InDetSimData& sdo = iter->second;
138 const std::vector< InDetSimData::Deposit >& deposits = sdo.getdeposits();
139 for( std::vector< InDetSimData::Deposit >::const_iterator nextdeposit = deposits.begin() ; nextdeposit!=deposits.end(); ++nextdeposit) {
140 SdoInfo sdoi;
141 // The times of the deposit also include the smearing of digi of 25 ps
142 sdoi.time = nextdeposit->second;
143 const HepMcParticleLink& particleLink = nextdeposit->first;
144 if(particleLink.isValid()){
145 HepMC::ConstGenParticlePtr genPart(particleLink.cptr());
146 if(isHSGoodParticle(genPart,hardScatterEvent)) sdoi.truth = 1; // signal
147 else sdoi.truth = 2; // pile-up
148 }
149 sdo_info.push_back(sdoi);
150 }
151 } else {
152 // Secondaries are not saved in the SDO map, hence saving the TOA for them
153 SdoInfo sdoi;
154 sdoi.time = rdo_toa;
155 sdoi.truth = 3; // secondary
156 sdo_info.push_back(sdoi);
157 }
158
159 // Only the particle that creates the first deposit is considered to be the correct particle
160 std::sort(sdo_info.begin(), sdo_info.end());
161 if (sdo_info.size() > 0) {
162 m_rdo_hit_sdo_toa.push_back(sdo_info.at(0).time);
163 m_rdo_hit_sdo_truth_category.push_back(sdo_info.at(0).truth);
164 }
165 }
166 else {
167 ATH_MSG_ERROR("No HGTD SDO map found!");
168 }
169 }
170 }
171 }
172
173
174 m_tree->Fill();
175
176 return StatusCode::SUCCESS;
177}
178
179bool HGTD_RDOAnalysis::isHSGoodParticle(HepMC::ConstGenParticlePtr particlePtr,const HepMC::GenEvent* hardScatterGenEvent, float min_pt_cut) {
180 bool decision = false;
181 if( MC::isGenStable(particlePtr) and
182 MC::isCharged(particlePtr) and
183 particlePtr->momentum().perp() >= min_pt_cut and
184 particlePtr->momentum().eta() < 4.0 and
185 particlePtr->parent_event() == hardScatterGenEvent)
186 decision = true;
187
188 return decision;
189}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
ATLAS-specific HepMC functions.
Handle class for reading from StoreGate.
const ServiceHandle< StoreGateSvc > & detStore() const
const ServiceHandle< ITHistSvc > & histSvc() const
The standard THistSvc (for writing histograms and TTrees and more to a root file) Returns (kind of) a...
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
const T * at(size_type n) const
Access an element, as an rvalue.
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
std::vector< float > m_rdo_module_x
bool isHSGoodParticle(HepMC::ConstGenParticlePtr particlePtr, const HepMC::GenEvent *hardScatterEvent, float min_pt_cut=1000.)
std::vector< unsigned long long > m_rdo_module_ID
std::vector< float > m_rdo_hit_sdo_toa
virtual StatusCode initialize() override final
std::vector< float > m_rdo_hit_z
Gaudi::Property< std::string > m_ntupleName
const HGTD_ID * m_HGTD_ID
std::vector< float > m_rdo_hit_x
Gaudi::Property< std::string > m_HGTDID_Name
SG::ReadHandleKey< HGTD_RDO_Container > m_inputKey
SG::ReadHandleKey< InDetSimDataCollection > m_inputTruthKey
std::vector< int > m_rdo_hit_sdo_truth_category
std::vector< float > m_rdo_hit_toa
Gaudi::Property< std::string > m_HGTD_Name
std::vector< float > m_rdo_hit_y
std::vector< float > m_rdo_module_y
std::vector< int > m_rdo_module_layer
SG::ReadHandleKey< McEventCollection > m_inputMcEventCollectionKey
Gaudi::Property< std::string > m_ntuplePath
std::vector< float > m_rdo_module_z
virtual StatusCode execute() override final
const HGTD_DetectorManager * m_HGTD_Manager
const_iterator end() const
return const_iterator for end of container
const_iterator begin() const
return const_iterator for first entry
value_type get_compact() const
Get the compact id.
Class to hold geometrical description of an HGTD detector element.
Class to represent a position in the natural frame of a silicon sensor, for Pixel and SCT For Pixel: ...
HepGeom::Point3D< double > globalPosition(const HepGeom::Point3D< double > &localPos) const
transform a reconstruction local position into a global position (inline):
Amg::Vector2D rawLocalPositionOfCell(const SiCellId &cellId) const
Returns position (center) of cell.
const std::vector< Deposit > & getdeposits() const
This defines the McEventCollection, which is really just an ObjectVector of McEvent objectsFile: Gene...
Eigen::Matrix< double, 3, 1 > Vector3D
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
bool isCharged(const T &p)
bool isGenStable(const T &p)
Determine if the particle is stable at the generator (not det-sim) level,.
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.