ATLAS Offline Software
Loading...
Searching...
No Matches
SimHitTreeCreator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
5// class header include
6#include "SimHitTreeCreator.h"
7
8// Athena includes
11
12// HepMC includes
14
15// ROOT includes
16#include "TTree.h"
17
18ISF::SimHitTreeCreator::SimHitTreeCreator( const std::string& name, ISvcLocator* pSvcLocator )
19 : ::AthAlgorithm( name, pSvcLocator )
20 , m_thistSvc("THistSvc",name)
21 , m_validationStream("ISFSimHit")
22 , m_t_simHits(nullptr)
23 , m_pileup(-1)
24 , m_type(-1)
25 , m_id(-1)
26 , m_mother(-1)
27 , m_barcode(-1)
28 , m_time(-1.)
29 , m_drift(-1.)
30 , m_edeposit(-1.)
31 , m_momentum(-1.)
32 , m_theta(-1.)
33 , m_phi(-1.)
34 , m_eta(-1.)
35 , m_bcmHits("BCMHits")
36 , m_blmHits("BLMHits")
37 , m_pixHits("PixelHits")
38 , m_sctHits("SCT_Hits")
39 , m_trtHits("TRTUncompressedHits")
40 , m_pixPileupHits("PileupPixelHits")
41 , m_sctPileupHits("PileupSCT_Hits")
42 , m_trtPileupHits("PileupTRTUncompressedHits")
43 , m_mdtHits("MDT_Hits")
44 , m_rpcHits("RPC_Hits")
45 , m_tgcHits("TGC_Hits")
46 , m_cscHits("CSC_Hits")
47{
48 declareProperty("ValidationStreamName",
49 m_validationStream = "ISFSimHit",
50 "Name of the output stream" );
51 declareProperty("THistService",
53 "The THistSvc" );
54
55 declareProperty("BCM_HitCollection", m_bcmHits );
56 declareProperty("BLM_HitCollection", m_blmHits );
57 declareProperty("PixelHitCollection", m_pixHits );
58 declareProperty("SCT_HitCollection", m_sctHits );
59 declareProperty("TRT_HitCollection", m_trtHits );
60 declareProperty("PileupPixelHitCollection", m_pixPileupHits );
61 declareProperty("PileupSCT_HitCollection", m_sctPileupHits );
62 declareProperty("PileupTRT_HitCollection", m_trtPileupHits );
63 declareProperty("MDT_HitCollection", m_mdtHits );
64 declareProperty("RPC_HitCollection", m_rpcHits );
65 declareProperty("TGC_HitCollection", m_tgcHits );
66 declareProperty("CSC_HitCollection", m_cscHits );
67
68}
69
70// Destructor
74
77 // retrieve the histogram service
78 ATH_CHECK( m_thistSvc.retrieve() );
79
81
82 ATH_CHECK( m_bcmHits.initialize() );
83 ATH_CHECK( m_blmHits.initialize() );
84 ATH_CHECK( m_pixHits.initialize() );
85 ATH_CHECK( m_sctHits.initialize() );
86 ATH_CHECK( m_trtHits.initialize() );
87 ATH_CHECK( m_pixPileupHits.initialize() );
88 ATH_CHECK( m_sctPileupHits.initialize() );
89 ATH_CHECK( m_trtPileupHits.initialize() );
90 ATH_CHECK( m_mdtHits.initialize() );
91 ATH_CHECK( m_rpcHits.initialize() );
92 ATH_CHECK( m_tgcHits.initialize() );
93 ATH_CHECK( m_cscHits.initialize() );
94 return StatusCode::SUCCESS;
95}
96
99 return StatusCode::SUCCESS;
100}
101
104 // Create the prefix of histogram names for the THistSvc
105 //std::string prefix = "/" + m_validationStream + "/";
106 const char *treeName="simhits";
107 const std::string prefix = "/" + m_validationStream + "/"+ treeName;
108 m_t_simHits = new TTree( treeName, treeName );
109 m_t_simHits->Branch("pileup" , &m_pileup , "pileup/I" );
110 m_t_simHits->Branch("type" , &m_type , "type/I" );
111 m_t_simHits->Branch("id" , &m_id , "id/I" );
112 m_t_simHits->Branch("mother" , &m_mother , "mother/I" );
113 m_t_simHits->Branch("barcode" , &m_barcode , "barcode/I" );
114 m_t_simHits->Branch("momentum" , &m_momentum , "momentum/F" );
115 m_t_simHits->Branch("time" , &m_time , "time/F" );
116 m_t_simHits->Branch("drift" , &m_drift , "dist/F" );
117 m_t_simHits->Branch("edeposit" , &m_edeposit , "edeposit/F" );
118 m_t_simHits->Branch("theta" , &m_theta , "theta/F" );
119 m_t_simHits->Branch("phi" , &m_phi , "phi/F" );
120 m_t_simHits->Branch("eta" , &m_eta , "eta/F" );
121
122 // register the Tree to the THistSvc and return it's StatusCodes
123 ATH_CHECK( m_thistSvc->regTree( prefix, m_t_simHits) );
124
125 return StatusCode::SUCCESS;
126}
127
130{
131 // loop over collections
132
134 if (mdtHits.isValid() && mdtHits->size()) {
135 MDTSimHitCollection::const_iterator ih=mdtHits->begin();
136 while ( ih!=mdtHits->end()) {
137 m_type = 1;
138 m_id = (*ih).MDTid();
139 m_mother = (*ih).particleEncoding();
140 m_time = (*ih).globalTime();
141 m_drift = (*ih).driftRadius();
142 m_edeposit = (*ih).energyDeposit();
143 m_barcode = (*ih).truthBarcode();
144 HepMcParticleLink HMPL = (*ih).particleLink();
146
147 ++ih;
148 while (ih!=mdtHits->end() && m_id==(*ih).MDTid() && m_barcode==(*ih).truthBarcode() ) {
149 // merge energy deposits and move on
150 m_edeposit += (*ih).energyDeposit();
151 ++ih;
152 }
153 m_t_simHits->Fill();
154 }
155 }
156
158 if (rpcHits.isValid() && rpcHits->size()) {
159 RPCSimHitCollection::const_iterator ih=rpcHits->begin();
160 while (ih!=rpcHits->end()) {
161 m_type = 2;
162 m_id = (*ih).RPCid();
163 m_mother = (*ih).particleEncoding();
164 m_time = (*ih).globalTime();
165 m_drift = 0.;
166 m_edeposit = (*ih).energyDeposit();
167 m_barcode = (*ih).truthBarcode();
168 HepMcParticleLink HMPL = (*ih).particleLink();
170
171 ++ih;
172 while (ih!=rpcHits->end() && m_id==(*ih).RPCid() && m_barcode==(*ih).truthBarcode() ) {
173 // merge energy deposits and move on
174 m_edeposit += (*ih).energyDeposit();
175 ++ih;
176 }
177 m_t_simHits->Fill();
178 }
179 }
180
182 if (tgcHits.isValid() && tgcHits->size()) {
183 TGCSimHitCollection::const_iterator ih=tgcHits->begin();
184 while ( ih!=tgcHits->end()) {
185 m_type = 3;
186 m_id = (*ih).TGCid();
187 m_mother = (*ih).particleEncoding();
188 m_time = (*ih).globalTime();
189 m_drift = 0.;
190 m_edeposit = (*ih).energyDeposit();
191 m_barcode = (*ih).truthBarcode();
192 HepMcParticleLink HMPL = (*ih).particleLink();
194
195 ++ih;
196 while (ih!=tgcHits->end() && m_id==(*ih).TGCid() && m_barcode==(*ih).truthBarcode() ) {
197 // merge energy deposits and move on
198 m_edeposit += (*ih).energyDeposit();
199 ++ih;
200 }
201 m_t_simHits->Fill();
202 }
203 }
204
206 if (cscHits.isValid() && cscHits->size()) {
207 CSCSimHitCollection::const_iterator ih=cscHits->begin();
208 while ( ih!=cscHits->end()) {
209 m_type = 4;
210 m_id = (*ih).CSCid();
211 m_mother = (*ih).particleID();
212 m_time = (*ih).globalTime();
213 m_drift = 0.;
214 m_edeposit = (*ih).energyDeposit();
215 m_barcode = (*ih).truthBarcode();
216 HepMcParticleLink HMPL = (*ih).particleLink();
218
219 ++ih;
220 while (ih!=cscHits->end() && m_id==(*ih).CSCid() && m_barcode==(*ih).truthBarcode() ) {
221 // merge energy deposits and move on
222 m_edeposit += (*ih).energyDeposit();
223 ++ih;
224 }
225 m_t_simHits->Fill();
226 }
227 }
228
229 for (int ipileup=0;ipileup<2;ipileup++) {
230 m_pileup = ipileup;
231
233 if (pixHits.isValid() && pixHits->size()) {
234 SiHitCollection::const_iterator ih=pixHits->begin();
235 while (ih!=pixHits->end()) {
236 m_type = 5;
237 m_id = (*ih).identify();
238 HepMcParticleLink HMPL = (*ih).particleLink();
239 if (HMPL.isValid()) m_mother = (HMPL.cptr())->pdg_id();
240 else m_mother=0;
241 m_time = (*ih).meanTime();
242 m_drift = 0.;
243 m_edeposit = (*ih).energyLoss();
244 m_barcode = (*ih).truthBarcode();
246
247 ++ih;
248 while (ih!=pixHits->end() && ((unsigned int)m_id)==(*ih).identify() && m_barcode==(*ih).truthBarcode() ) {
249 // merge energy deposits and move on
250 m_edeposit += (*ih).energyLoss();
251 ++ih;
252 }
253 m_t_simHits->Fill();
254 }
255 }
256
258 if (sctHits.isValid() && sctHits->size()) {
259 SiHitCollection::const_iterator ih=sctHits->begin();
260 while (ih!=sctHits->end()) {
261 m_type = 6;
262 m_id = (*ih).identify();
263 HepMcParticleLink HMPL = (*ih).particleLink();
264 if (HMPL.isValid()) m_mother = (HMPL.cptr())->pdg_id();
265 else m_mother=0;
266 m_time = (*ih).meanTime();
267 m_drift = 0.;
268 m_edeposit = (*ih).energyLoss();
269 m_barcode = (*ih).truthBarcode();
271
272 ++ih;
273 while (ih!=sctHits->end() && ((unsigned int)m_id)==(*ih).identify() && m_barcode==(*ih).truthBarcode() ) {
274 // merge energy deposits and move on
275 m_edeposit += (*ih).energyLoss();
276 ++ih;
277 }
278 m_t_simHits->Fill();
279 }
280 }
281
283 if (trtHits.isValid() && trtHits->size()) {
285 while ( ih!=trtHits->end()) {
286 m_type = 7;
287 m_id = (*ih).GetHitID();
288 m_mother = (*ih).GetParticleEncoding();
289 HepMcParticleLink HMPL = (*ih).particleLink();
290 m_time = (*ih).GetGlobalTime();
291 m_drift = 0.;
292 m_edeposit = (*ih).GetEnergyDeposit();
293 m_barcode = (*ih).truthBarcode();
295
296 ++ih;
297 while (ih!=trtHits->end() && m_id==(*ih).GetHitID() && m_barcode==(*ih).truthBarcode() ) {
298 // merge energy deposits and move on
299 m_edeposit += (*ih).GetEnergyDeposit();
300 ++ih;
301 }
302 m_t_simHits->Fill();
303 }
304 }
305 }
306 return StatusCode::SUCCESS;
307}
308
309//** Add information from HepMcParticleLink to TTree - common for all SimHit types */
311 if (HMPL.isValid()) {
312 auto t_mom=(HMPL.cptr())->momentum();
313 m_momentum = std::sqrt( t_mom.x()*t_mom.x()+t_mom.y()*t_mom.y()+t_mom.z()*t_mom.z());
314 m_eta= (HMPL.cptr())->momentum().eta();
315 m_theta= (HMPL.cptr())->momentum().theta();
316 m_phi = (HMPL.cptr())->momentum().phi();
317 }
318 else {
319 m_momentum=-1.;
320 m_theta=-1.;
321 m_eta=-10.;
322 m_phi = -10.;
323 }
324 return;
325}
#define ATH_CHECK
Evaluate an expression and check for errors.
Handle class for reading from StoreGate.
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
CONT::const_iterator const_iterator
virtual ~SimHitTreeCreator()
Destructor.
StatusCode fillSimHitsTree()
Fill the simhits tree - validation mode only.
SG::ReadHandleKey< SiHitCollection > m_bcmHits
virtual StatusCode execute() override final
Athena algorithm's interface method execute()
SG::ReadHandleKey< SiHitCollection > m_blmHits
SG::ReadHandleKey< TGCSimHitCollection > m_tgcHits
virtual StatusCode initialize() override final
Athena algorithm's interface method initialize()
SG::ReadHandleKey< RPCSimHitCollection > m_rpcHits
SG::ReadHandleKey< SiHitCollection > m_sctHits
std::string m_validationStream
validation THist stream name
SG::ReadHandleKey< MDTSimHitCollection > m_mdtHits
SG::ReadHandleKey< TRTUncompressedHitCollection > m_trtPileupHits
SimHitTreeCreator(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters.
SG::ReadHandleKey< SiHitCollection > m_pixPileupHits
SG::ReadHandleKey< SiHitCollection > m_pixHits
SG::ReadHandleKey< TRTUncompressedHitCollection > m_trtHits
ServiceHandle< ITHistSvc > m_thistSvc
the histogram service
SG::ReadHandleKey< SiHitCollection > m_sctPileupHits
SG::ReadHandleKey< CSCSimHitCollection > m_cscHits
StatusCode createSimHitsTree()
Create the simhits tree - validation mode only.
TTree * m_t_simHits
Validation output TTree (+variables)
void addHepMcParticleLinkInfoToTree(HepMcParticleLink &HMPL)
virtual bool isValid() override final
Can the handle be successfully dereferenced?