ATLAS Offline Software
Loading...
Searching...
No Matches
TRTProcessingOfEndCapHits.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
10
11// Athena headers
13#include "MCTruth/TrackHelper.h"
14#include "TRT_G4Utilities/TRTParameters.hh"
15
16// Geant4 headers
17#include "G4AffineTransform.hh"
18#include "G4NavigationHistory.hh"
19#include "G4Step.hh"
20#include "G4StepPoint.hh"
21#include "G4ThreeVector.hh"
22#include "G4TouchableHandle.hh"
23#include "G4TouchableHistory.hh"
24#include "G4Track.hh"
25
26// STL headers
27#include <cmath>
28
29// Called by TRTSensitiveDetector::InitializeHitProcessing
30
32(TRTSensitiveDetector* pSensitiveDet):
33 m_printMessages(0), //FIXME obsolete?
34 m_sectorsABC(0),
45 m_pParameters(nullptr),
46 m_pSensitiveDetector(pSensitiveDet),
47 m_verboseLevel(pSensitiveDet->verboseLevel)
48{
49 m_pParameters = TRTParameters::GetPointer();
50
51 m_printMessages = m_pParameters->GetInteger("PrintMessages"); //FIXME Obsolete?
52
53 if (m_verboseLevel>5) { G4cout << "##### Constructor TRTProcessingOfEndCapHits" << G4endl; }
54
55 this->Initialize();
56
57 m_pParameters = nullptr;
58
59 if (m_verboseLevel>5) { G4cout << "##### Constructor TRTProcessingOfEndCapHits done" << G4endl; }
60}
61
62
63// Called by TRTSensitiveDetector::DeleteObjects
64
68
69
70// Called by TRTProcessingOfEndCapHits
71
73{
74 if (m_verboseLevel>5) { G4cout << "######### Method TRTProcessingOfEndCapHits::Initialize" << G4endl; }
75
76 m_sectorsABC = m_pParameters->GetInteger("SectorsABC");
77
79 m_pParameters->GetInteger("TestLocalCoordinatesOfHits");
80
81 if (m_verboseLevel>5) { G4cout << "######### Method TRTProcessingOfEndCapHits::Initialize done" << G4endl; }
82
83}
84
85
86// Called by Geant4
87
89{
90 G4Track* pTrack = pStep->GetTrack();
91 // get the HepMC barcode using the track helper
92 TrackHelper trHelp(pTrack);
93
94 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
95 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
96
97 G4TouchableHandle pTouchableHandle = pPreStepPoint->GetTouchableHandle();
98 int depth = 1;
99
100 int strawID = pTouchableHandle->GetReplicaNumber(depth++);
101
102 int sectorID = 0;
103
104 if (m_sectorsABC)
105 sectorID = pTouchableHandle->GetReplicaNumber(depth++);
106
107 int planeID = pTouchableHandle->GetReplicaNumber(depth++);
108
109 int wheelID = pTouchableHandle->GetReplicaNumber(depth);
110
111 G4ThreeVector globalPreStepPoint = pPreStepPoint->GetPosition();
112
113 int endCapID;
114
115 if (globalPreStepPoint.z() > 0.)
116 endCapID = 0;
117 else
118 endCapID = 1;
119
120 if (m_sectorsABC)
121 {
122 if (wheelID < m_numberOfWheelsAB)
123 {
124 strawID += m_numberOfStrawsInSectorsAB * sectorID -
126 if (strawID < 0)
128 }
129 else
130 {
132 if (strawID < 0)
133 strawID += m_numberOfStrawsInPlaneC;
134 }
135 }
136
137 G4ThreeVector globalPostStepPoint = pPostStepPoint->GetPosition();
138
139 const G4TouchableHistory* pTouchableHistory =
140 static_cast<const G4TouchableHistory*>(pPreStepPoint->GetTouchable());
141
142 const G4AffineTransform& topTransform = pTouchableHistory->GetHistory()->
143 GetTopTransform();
144
145 G4ThreeVector localPreStepPoint = topTransform.TransformPoint(globalPreStepPoint);
146 G4ThreeVector localPostStepPoint = topTransform.TransformPoint(globalPostStepPoint);
147
148 double preStepX = localPreStepPoint.x();
149 double preStepY = localPreStepPoint.y();
150 double preStepZ = localPreStepPoint.z();
151 double postStepX = localPostStepPoint.x();
152 double postStepY = localPostStepPoint.y();
153 double postStepZ = localPostStepPoint.z();
154
155
157 {
158 double preStepR = std::sqrt(preStepX * preStepX + preStepY * preStepY);
159 double postStepR = std::sqrt(postStepX * postStepX + postStepY * postStepY);
160
161 if (preStepR > 2.0000001 || postStepR > 2.0000001)
162 {
163 G4int particleEncoding = m_pSensitiveDetector->m_particleEncoding;
164 G4double kineticEnergy = m_pSensitiveDetector->m_kineticEnergy;
165 G4double energyDeposit = m_pSensitiveDetector->m_energyDeposit;
166
167 std::cout << "!!!!! End-caps. Error in local coordinates of hits!" << std::endl;
168 std::cout << " endCapID=" << endCapID << " wheelID=" << wheelID
169 << " planeID=" << planeID << " strawID=" << strawID
170 << " trackID=" << trHelp.GetUniqueID() << std::endl;
171 std::cout << " particleEncoding=" << particleEncoding;
172
173 if (kineticEnergy < 0.0001)
174 std::cout << " kineticEnergy=" << kineticEnergy * 1000000. << " eV";
175 else if (kineticEnergy < 0.1)
176 std::cout << " kineticEnergy=" << kineticEnergy * 1000. << " keV";
177 else if (kineticEnergy >= 100.)
178 std::cout << " kineticEnergy=" << kineticEnergy / 1000. << " GeV";
179 else if (kineticEnergy >= 100000.)
180 std::cout << " kineticEnergy=" << kineticEnergy / 1000000. << " TeV";
181 else
182 std::cout << " kineticEnergy=" << kineticEnergy << " MeV";
183
184 if (energyDeposit < 0.1)
185 std::cout << " energyDeposit=" << energyDeposit * 1000. << " eV"
186 << std::endl;
187 else if (energyDeposit >= 100.)
188 std::cout << " energyDeposit=" << energyDeposit / 1000. << " MeV"
189 << std::endl;
190 else
191 std::cout << " energyDeposit=" << energyDeposit << " keV" << std::endl;
192
193 std::cout << " preStepX=" << preStepX << " mm preStepY=" << preStepY
194 << " mm preStepR=" << preStepR << " mm" << std::endl;
195 std::cout << " postStepX=" << postStepX << " mm postStepY=" << postStepY
196 << " mm postStepR=" << postStepR << " mm" << std::endl << std::endl;
197 }
198 }
199
200 double preStepGlobalTime = pPreStepPoint->GetGlobalTime();
201 double postStepGlobalTime = pPostStepPoint->GetGlobalTime();
202 double globalTime = (preStepGlobalTime + postStepGlobalTime) / 2.;
203
204
205 int numberOfStrawsInPlane;
206 int numberOfStrawsInIDSector;
207
208 if (wheelID < m_numberOfWheelsAB) {
209 numberOfStrawsInPlane = m_numberOfStrawsInPlanesAB;
210 numberOfStrawsInIDSector = m_numberOfStrawsInIDSectorsAB;
211 } else {
212 numberOfStrawsInPlane = m_numberOfStrawsInPlaneC;
213 numberOfStrawsInIDSector = m_numberOfStrawsInIDSectorC;
214 }
215
216 // Mapping for negative endcap is
217 // nSectors = 32
218 // nStraws = num straws in sector
219 // straw = straw number within sector
220 // sector -> (nSectors + nSectors/2 - sector - 1) % nSectors
221 // straw -> nStraws - 1 - straw
222 //
223 // Since we calculate sector number below from strawID (where
224 // strawID is straw number in plane), its equivalent to map strawID
225 // first and then calculate sector number and straw number within sector.
226
227 // Mapping straw number in plane for negative endcap
228 if (endCapID > 0) {
229 strawID = (3*numberOfStrawsInPlane/2 - strawID - 1) % numberOfStrawsInPlane;
230 }
231
232 // Calculate sector number and straw number within sector.
233 sectorID = strawID / numberOfStrawsInIDSector;
234 strawID %= numberOfStrawsInIDSector;
235
236
237 int hitID;
238 if (endCapID == 0)
239 hitID = 0x00200000;
240 else
241 hitID = 0x00300000;
242
243 hitID += (wheelID << 15);
244 hitID += (sectorID << 10);
245 hitID += (planeID << 5);
246 hitID += strawID;
247
248 m_pSensitiveDetector->m_hitID = hitID;
249 auto* eventInfo = m_pSensitiveDetector->m_g4UserEventInfo;
250 m_pSensitiveDetector->m_partLink = trHelp.GenerateParticleLink(eventInfo ? eventInfo->GetEventStore() : nullptr);
251 m_pSensitiveDetector->m_preStepX = preStepX;
252 m_pSensitiveDetector->m_preStepY = preStepY;
253 m_pSensitiveDetector->m_preStepZ = preStepZ;
254 m_pSensitiveDetector->m_postStepX = postStepX;
255 m_pSensitiveDetector->m_postStepY = postStepY;
256 m_pSensitiveDetector->m_postStepZ = postStepZ;
257 m_pSensitiveDetector->m_globalTime = globalTime;
258
259 return true;
260}
TRTProcessingOfEndCapHits(TRTSensitiveDetector *)
TRTSensitiveDetector * m_pSensitiveDetector
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
int GetUniqueID() const
std::string depth
tag string for intendation
Definition fastadd.cxx:46