ATLAS Offline Software
TRTProcessingOfBarrelHits.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // Class header
7 
8 // Package headers
9 #include "TRTSensitiveDetector.h"
10 
11 // Athena headers
12 #include "MCTruth/TrackHelper.h"
13 #include "TRT_G4Utilities/TRTParameters.hh"
14 
15 // Geant4 headers
16 #include "G4AffineTransform.hh"
17 #include "G4NavigationHistory.hh"
18 #include "G4Step.hh"
19 #include "G4StepPoint.hh"
20 #include "G4ThreeVector.hh"
21 #include "G4TouchableHandle.hh"
22 #include "G4TouchableHistory.hh"
23 #include "G4Track.hh"
24 
25 // STL headers
26 #include <cmath>
27 
28 // Called by TRTSensitiveDetector::InitializeHitProcessing
29 
31 (TRTSensitiveDetector* pSensitiveDet):
32  m_printMessages(0), //FIXME obsolete?
33  m_barrelIdentifier(0),
34  m_testLocalCoordinatesOfHits(0),
35  m_numberOfStrawLayersA(0),
36  m_numberOfStrawLayersB(0),
37  m_numberOfStrawLayersC(0),
38  m_numberOfStrawsA(0),
39  m_numberOfStrawsB(0),
40  m_numberOfStrawsC(0),
41  m_integralDistributionOfStrawsA(nullptr),
42  m_integralDistributionOfStrawsB(nullptr),
43  m_integralDistributionOfStrawsC(nullptr),
44  m_strawIDToLayerIDA(nullptr),
45  m_strawIDToLayerIDB(nullptr),
46  m_strawIDToLayerIDC(nullptr),
47  m_pParameters(nullptr),
48  m_pSensitiveDetector(pSensitiveDet),
49  m_verboseLevel(pSensitiveDet->verboseLevel)
50 {
51  m_pParameters = TRTParameters::GetPointer();
52 
53  m_printMessages = m_pParameters->GetInteger("PrintMessages"); //FIXME obsolete?
54 
55  if (m_verboseLevel>5) { G4cout << "##### Constructor TRTProcessingOfBarrelHits" << G4endl; }
56 
57  this->Initialize();
58 
59  m_pParameters = nullptr;
60 
61  if (m_verboseLevel>5) { G4cout << "##### Constructor TRTProcessingOfBarrelHits done" << G4endl; }
62 }
63 
64 
65 // Called by TRTSensitiveDetector::DeleteObjects
66 
68 {
69 
70  if (m_verboseLevel>5) { G4cout << "##### Destructor TRTProcessingOfBarrelHits" << G4endl; }
71 
72  this->DeleteArrays();
73 
74  if (m_verboseLevel>5) { G4cout << "##### Destructor TRTProcessingOfBarrelHits done" << G4endl; }
75 
76 }
77 
78 
79 // Called by TRTProcessingOfBarrelHits
80 
82 {
83  if (m_verboseLevel>5) { G4cout << "######### Method TRTProcessingOfBarrelHits::Initialize" << G4endl; }
84 
85 
86  m_barrelIdentifier = m_pParameters->GetInteger("BarrelIdentifier");
88  m_pParameters->GetInteger("TestLocalCoordinatesOfHits");
89 
90  if (m_barrelIdentifier == 1)
91  {
92  m_strawIDToLayerIDA = nullptr;
93  m_strawIDToLayerIDB = nullptr;
94  m_strawIDToLayerIDC = nullptr;
95  }
96 
97 
98  if (m_barrelIdentifier == 2)
99  {
100  m_numberOfStrawsA = m_pParameters->GetInteger("NumberOfStrawsA");
101  m_numberOfStrawsB = m_pParameters->GetInteger("NumberOfStrawsB");
102  m_numberOfStrawsC = m_pParameters->GetInteger("NumberOfStrawsC");
103 
107  }
108 
109  m_numberOfStrawLayersA = m_pParameters->GetInteger("NumberOfStrawLayersA");
110  m_numberOfStrawLayersB = m_pParameters->GetInteger("NumberOfStrawLayersB");
111  m_numberOfStrawLayersC = m_pParameters->GetInteger("NumberOfStrawLayersC");
112 
116 
117  if (m_verboseLevel>5) { G4cout << "######### Method TRTProcessingOfBarrelHits::Initialize done" << G4endl; }
118 }
119 
120 
121 // Called by Geant4
122 
124 {
125  G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
126 
127  G4TouchableHandle pTouchableHandle = pPreStepPoint->GetTouchableHandle();
128  int depth = 2;
129 
130  int gasID = pTouchableHandle->GetReplicaNumber();
131  int strawID = pTouchableHandle->GetReplicaNumber(depth++);
132  int ringID = pTouchableHandle->GetReplicaNumber(depth++);
133  int moduleID = pTouchableHandle->GetReplicaNumber(depth);
134 
135  int hitID = gasID % 2;
136 
137  G4Track* pTrack = pStep->GetTrack();
138  // get the HepMC barcode using the track helper
139  TrackHelper trHelp(pTrack);
140 
141  G4ThreeVector globalPreStepPoint = pPreStepPoint->GetPosition();
142 
143  G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
144  G4ThreeVector globalPostStepPoint = pPostStepPoint->GetPosition();
145 
146 
147  const G4TouchableHistory* pTouchableHistory =
148  static_cast<const G4TouchableHistory*>(pPreStepPoint->GetTouchable());
149 
150  const G4AffineTransform& topTransform = pTouchableHistory->GetHistory()->
151  GetTopTransform();
152 
153  G4ThreeVector localPreStepPoint = topTransform.TransformPoint(globalPreStepPoint);
154  G4ThreeVector localPostStepPoint = topTransform.TransformPoint(globalPostStepPoint);
155 
156  double preStepZ = localPreStepPoint.z();
157  double postStepZ = localPostStepPoint.z();
158 
159 
160  double preStepX = localPreStepPoint.x();
161  double preStepY = localPreStepPoint.y();
162 
163  double postStepX = localPostStepPoint.x();
164  double postStepY = localPostStepPoint.y();
165 
167  {
168  double preStepR = std::sqrt(preStepX * preStepX + preStepY * preStepY);
169  double postStepR = std::sqrt(postStepX * postStepX + postStepY * postStepY);
170 
171  if (preStepR > 2.0000001 || postStepR > 2.0000001)
172  {
173  G4int particleEncoding = m_pSensitiveDetector->m_particleEncoding;
174  G4double kineticEnergy = m_pSensitiveDetector->m_kineticEnergy;
176 
177  std::cout << "!!!!! Barrel. Error in local coordinates of hits!" << std::endl;
178  std::cout << " barrelID=" << hitID << " ringID=" << ringID
179  << " moduleID=" << moduleID << " strawID=" << strawID
180  << " trackID=" << trHelp.GetUniqueID() << std::endl;
181  std::cout << " particleEncoding=" << particleEncoding;
182 
183  if (kineticEnergy < 0.0001)
184  std::cout << " kineticEnergy=" << kineticEnergy * 1000000. << " eV";
185  else if (kineticEnergy < 0.1)
186  std::cout << " kineticEnergy=" << kineticEnergy * 1000. << " keV";
187  else if (kineticEnergy >= 100.)
188  std::cout << " kineticEnergy=" << kineticEnergy / 1000. << " GeV";
189  else if (kineticEnergy >= 100000.)
190  std::cout << " kineticEnergy=" << kineticEnergy / 1000000. << " TeV";
191  else
192  std::cout << " kineticEnergy=" << kineticEnergy << " MeV";
193 
194  if (energyDeposit < 0.1)
195  std::cout << " energyDeposit=" << energyDeposit * 1000. << " eV"
196  << std::endl;
197  else if (energyDeposit >= 100.)
198  std::cout << " energyDeposit=" << energyDeposit / 1000. << " MeV"
199  << std::endl;
200  else
201  std::cout << " energyDeposit=" << energyDeposit << " keV" << std::endl;
202 
203  std::cout << " preStepX=" << preStepX << " mm preStepY=" << preStepY
204  << " mm preStepR=" << preStepR << " mm" << std::endl;
205  std::cout << " postStepX=" << postStepX << " mm postStepY=" << postStepY
206  << " mm postStepR=" << postStepR << " mm" << std::endl << std::endl;
207  }
208  }
209 
210  double preStepGlobalTime = pPreStepPoint->GetGlobalTime();
211  double postStepGlobalTime = pPostStepPoint->GetGlobalTime();
212  double globalTime = (preStepGlobalTime + postStepGlobalTime) / 2.;
213 
214  int layerID;
215 
216  if (m_barrelIdentifier == 1)
217  layerID = GetLayerID1(ringID, strawID);
218  else
219  layerID = GetLayerID2(ringID, strawID);
220 
221  hitID <<= 20;
222  hitID += (ringID << 15);
223  hitID += (moduleID << 10);
224  hitID += (layerID << 5);
225  hitID += strawID;
226 
227  m_pSensitiveDetector->m_hitID = hitID;
229  m_pSensitiveDetector->m_preStepX = preStepX;
230  m_pSensitiveDetector->m_preStepY = preStepY;
231  m_pSensitiveDetector->m_preStepZ = preStepZ;
232  m_pSensitiveDetector->m_postStepX = postStepX;
233  m_pSensitiveDetector->m_postStepY = postStepY;
234  m_pSensitiveDetector->m_postStepZ = postStepZ;
235  m_pSensitiveDetector->m_globalTime = globalTime;
236 
237  return true;
238 }
239 
240 
241 // Called by ProcessHit
242 
244  int& strawID) const
245 {
246  ++strawID;
247  int layerID = 15;
248 
249  if (ringID == 0)
250  {
251  if (strawID <= m_integralDistributionOfStrawsA[layerID])
252  {
253  layerID -= 8;
254  if (strawID > m_integralDistributionOfStrawsA[layerID])
255  layerID += 4;
256  else
257  layerID -= 4;
258  }
259  if (strawID > m_integralDistributionOfStrawsA[layerID])
260  layerID += 2;
261  else
262  layerID -= 2;
263  if (strawID > m_integralDistributionOfStrawsA[layerID])
264  ++layerID;
265  else
266  --layerID;
267  if (strawID > m_integralDistributionOfStrawsA[layerID])
268  ++layerID;
269 
270  if (layerID > 0)
271  strawID -= m_integralDistributionOfStrawsA[layerID - 1];
272  }
273  else if (ringID == 1)
274  {
275  if (strawID > m_integralDistributionOfStrawsB[layerID])
276  layerID += 8;
277  else
278  layerID -= 8;
279  if (strawID > m_integralDistributionOfStrawsB[layerID])
280  layerID += 4;
281  else
282  layerID -= 4;
283  if (strawID > m_integralDistributionOfStrawsB[layerID])
284  layerID += 2;
285  else
286  layerID -= 2;
287  if (strawID > m_integralDistributionOfStrawsB[layerID])
288  ++layerID;
289  else
290  --layerID;
291  if (strawID > m_integralDistributionOfStrawsB[layerID])
292  ++layerID;
293 
294  if (layerID > 0)
295  strawID -= m_integralDistributionOfStrawsB[layerID - 1];
296  }
297  else
298  {
299  if (strawID > m_integralDistributionOfStrawsC[layerID])
300  layerID += 8;
301  else
302  layerID -= 8;
303  if (strawID > m_integralDistributionOfStrawsC[layerID])
304  layerID += 4;
305  else
306  layerID -= 4;
307  if (strawID > m_integralDistributionOfStrawsC[layerID])
308  layerID += 2;
309  else
310  layerID -= 2;
311  if (strawID > m_integralDistributionOfStrawsC[layerID])
312  ++layerID;
313  else
314  --layerID;
315  if (strawID > m_integralDistributionOfStrawsC[layerID])
316  ++layerID;
317 
318  if (layerID > 0)
319  strawID -= m_integralDistributionOfStrawsC[layerID - 1];
320  }
321 
322  --strawID;
323  return layerID;
324 }
325 
326 
327 // Called by ProcessHit
328 
330  int& strawID) const
331 {
332  int layerID;
333  if (ringID == 0)
334  {
335  layerID = m_strawIDToLayerIDA[strawID];
336  if (layerID > 0)
337  strawID -= m_integralDistributionOfStrawsA[layerID - 1];
338  }
339  else if (ringID == 1)
340  {
341  layerID = m_strawIDToLayerIDB[strawID];
342  if (layerID > 0)
343  strawID -= m_integralDistributionOfStrawsB[layerID - 1];
344  }
345  else
346  {
347  layerID = m_strawIDToLayerIDC[strawID];
348  if (layerID > 0)
349  strawID -= m_integralDistributionOfStrawsC[layerID - 1];
350  }
351 
352  return layerID;
353 }
354 
355 
356 // Called by ~TRTProcessingOfBarrelHits
357 
359 {
360 
361  if (m_verboseLevel>5) { G4cout << "######### Method TRTProcessingOfBarrelHits::DeleteArrays" << G4endl; }
362 
366 
367  if (m_barrelIdentifier == 2)
368  {
369  delete [] m_strawIDToLayerIDA;
370  delete [] m_strawIDToLayerIDB;
371  delete [] m_strawIDToLayerIDC;
372  }
373 
374  if (m_verboseLevel>5) { G4cout << "######### Method TRTProcessingOfBarrelHits::DeleteArrays done" << G4endl; }
375 }
TRTProcessingOfBarrelHits::m_printMessages
int m_printMessages
Definition: TRTProcessingOfBarrelHits.h:33
TRTProcessingOfBarrelHits::m_integralDistributionOfStrawsB
int * m_integralDistributionOfStrawsB
Definition: TRTProcessingOfBarrelHits.h:47
egammaParameters::depth
@ depth
pointing depth of the shower as calculated in egammaqgcld
Definition: egammaParamDefs.h:276
TRTProcessingOfBarrelHits::m_numberOfStrawsC
int m_numberOfStrawsC
Definition: TRTProcessingOfBarrelHits.h:44
TRTProcessingOfBarrelHits::m_numberOfStrawLayersA
int m_numberOfStrawLayersA
Definition: TRTProcessingOfBarrelHits.h:38
TRTProcessingOfBarrelHits::m_integralDistributionOfStrawsC
int * m_integralDistributionOfStrawsC
Definition: TRTProcessingOfBarrelHits.h:48
TRTProcessingOfBarrelHits::m_strawIDToLayerIDB
int * m_strawIDToLayerIDB
Definition: TRTProcessingOfBarrelHits.h:51
TRTProcessingOfBarrelHits::TRTProcessingOfBarrelHits
TRTProcessingOfBarrelHits(TRTSensitiveDetector *)
Definition: TRTProcessingOfBarrelHits.cxx:31
TrackHelper.h
TRTSensitiveDetector::m_energyDeposit
double m_energyDeposit
Definition: TRTSensitiveDetector.h:73
TRTProcessingOfBarrelHits::m_numberOfStrawsA
int m_numberOfStrawsA
Definition: TRTProcessingOfBarrelHits.h:42
TRTProcessingOfBarrelHits::Initialize
void Initialize()
Definition: TRTProcessingOfBarrelHits.cxx:81
TRTSensitiveDetector
Definition: TRTSensitiveDetector.h:25
Trk::energyDeposit
@ energyDeposit
Definition: MeasurementType.h:32
TRTProcessingOfBarrelHits::m_pSensitiveDetector
TRTSensitiveDetector * m_pSensitiveDetector
Definition: TRTProcessingOfBarrelHits.h:55
TRTSensitiveDetector::m_partLink
HepMcParticleLink m_partLink
Definition: TRTSensitiveDetector.h:70
TrackHelper
Definition: TrackHelper.h:14
TRTProcessingOfBarrelHits::m_numberOfStrawLayersC
int m_numberOfStrawLayersC
Definition: TRTProcessingOfBarrelHits.h:40
TRTSensitiveDetector::m_postStepX
double m_postStepX
Definition: TRTSensitiveDetector.h:78
TRTProcessingOfBarrelHits::GetLayerID1
int GetLayerID1(const int &, int &) const
Definition: TRTProcessingOfBarrelHits.cxx:243
TRTProcessingOfBarrelHits::m_pParameters
const TRTParameters * m_pParameters
Definition: TRTProcessingOfBarrelHits.h:54
TRTSensitiveDetector::m_globalTime
double m_globalTime
Definition: TRTSensitiveDetector.h:81
TRTSensitiveDetector::m_kineticEnergy
double m_kineticEnergy
Definition: TRTSensitiveDetector.h:72
TRTProcessingOfBarrelHits::m_strawIDToLayerIDA
int * m_strawIDToLayerIDA
Definition: TRTProcessingOfBarrelHits.h:50
TRTProcessingOfBarrelHits::m_numberOfStrawLayersB
int m_numberOfStrawLayersB
Definition: TRTProcessingOfBarrelHits.h:39
TRTProcessingOfBarrelHits::m_numberOfStrawsB
int m_numberOfStrawsB
Definition: TRTProcessingOfBarrelHits.h:43
TRTSensitiveDetector.h
TRTSensitiveDetector::m_postStepY
double m_postStepY
Definition: TRTSensitiveDetector.h:79
TRTSensitiveDetector::m_preStepX
double m_preStepX
Definition: TRTSensitiveDetector.h:75
TRTProcessingOfBarrelHits.h
TRTProcessingOfBarrelHits::GetLayerID2
int GetLayerID2(const int &, int &) const
Definition: TRTProcessingOfBarrelHits.cxx:329
TRTSensitiveDetector::m_preStepY
double m_preStepY
Definition: TRTSensitiveDetector.h:76
TRTProcessingOfBarrelHits::DeleteArrays
void DeleteArrays()
Definition: TRTProcessingOfBarrelHits.cxx:358
TRTSensitiveDetector::m_postStepZ
double m_postStepZ
Definition: TRTSensitiveDetector.h:80
TRTProcessingOfBarrelHits::ProcessHit
bool ProcessHit(G4Step *)
Definition: TRTProcessingOfBarrelHits.cxx:123
TRTSensitiveDetector::m_preStepZ
double m_preStepZ
Definition: TRTSensitiveDetector.h:77
TRTProcessingOfBarrelHits::m_testLocalCoordinatesOfHits
int m_testLocalCoordinatesOfHits
Definition: TRTProcessingOfBarrelHits.h:36
TRTProcessingOfBarrelHits::m_integralDistributionOfStrawsA
int * m_integralDistributionOfStrawsA
Definition: TRTProcessingOfBarrelHits.h:46
TRTProcessingOfBarrelHits::m_verboseLevel
int m_verboseLevel
Definition: TRTProcessingOfBarrelHits.h:57
TRTSensitiveDetector::m_particleEncoding
int m_particleEncoding
Definition: TRTSensitiveDetector.h:71
TrackHelper::GenerateParticleLink
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition: TrackHelper.h:35
TrackHelper::GetUniqueID
int GetUniqueID() const
Definition: TrackHelper.cxx:41
TRTSensitiveDetector::m_hitID
int m_hitID
Properties of current TRTUncompressedHit, set by TRTProcessingOfBarrelHits and TRTProcessingOfEndCapH...
Definition: TRTSensitiveDetector.h:69
TRTProcessingOfBarrelHits::m_barrelIdentifier
int m_barrelIdentifier
Definition: TRTProcessingOfBarrelHits.h:35
TRTProcessingOfBarrelHits::m_strawIDToLayerIDC
int * m_strawIDToLayerIDC
Definition: TRTProcessingOfBarrelHits.h:52
TRTProcessingOfBarrelHits::~TRTProcessingOfBarrelHits
~TRTProcessingOfBarrelHits()
Definition: TRTProcessingOfBarrelHits.cxx:67