ATLAS Offline Software
Loading...
Searching...
No Matches
PixelSensorSD.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
5//
6// Pixel Sensitive Detector.
7// The Hits are processed here. For every hit I get the position and
8// an information on the sensor in which the interaction happened
9//
10
11// Class header
12#include "PixelSensorSD.h"
13#include <GaudiKernel/EventContext.h>
14
15// Athena headers
18#include "MCTruth/TrackHelper.h"
19
20// Geant4 headers
21#include "G4ChargedGeantino.hh"
22#include <G4EventManager.hh>
23#include "G4Geantino.hh"
24#include "G4SDManager.hh"
25#include "G4Step.hh"
26#include "G4ThreeVector.hh"
27
28// CLHEP headers
29#include "CLHEP/Geometry/Transform3D.h"
30#include "CLHEP/Units/SystemOfUnits.h"
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34PixelSensorSD::PixelSensorSD(const std::string& name, const std::string& hitCollectionName)
35 : G4VSensitiveDetector( name )
36 , m_HitCollName( hitCollectionName )
37{
38}
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41// Initialize from G4 - cache the hit collection for the current event
42void PixelSensorSD::Initialize(G4HCofThisEvent *)
43{
44 // ISF calls G4SDManager::PrepareNewEvent() before the Geant4 event loop starts...
45 if(auto* eventManger = G4EventManager::GetEventManager()){
46 if(auto* eventInfo = static_cast<AtlasG4EventUserInfo*>(eventManger->GetUserInformation())){
47 m_HitColl = eventInfo->GetHitCollectionMap()->Find<SiHitCollection>(m_HitCollName);
48 m_g4UserEventInfo = eventInfo;
49 }
50 }
51}
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54
55G4bool PixelSensorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
56{
57 if (verboseLevel>5) G4cout << "Process Hit" << G4endl;
58
59 G4double edep = aStep->GetTotalEnergyDeposit();
60 edep *= CLHEP::MeV;
61 if(edep==0.) {
62 if(aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
63 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
64 return false;
65 }
66
67 //use the global time. i.e. the time from the beginning of the event
68 //
69 // Get the Touchable History:
70 //
71 const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
72 if (not myTouch) {
73 G4cout << "PixelSensorSD::ProcessHits bad dynamic_cast" << G4endl;
74 return false;
75 }
76 if(verboseLevel>5){
77 for (int i=0;i<myTouch->GetHistoryDepth();i++){
78 std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
79 int copyno=myTouch->GetVolume(i)->GetCopyNo();
80 G4cout << "Volume " <<detname <<" Copy Nr. " << copyno << G4endl;
81 }
82 }
83 //
84 // Get the hit coordinates. Start and End Point
85 //
86 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
87 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
88
89 // Calculate the local step begin and end position.
90 // From a G4 FAQ:
91 // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
92 //
93 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
94 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
95 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
96
97 HepGeom::Point3D<double> lP1,lP2;
98 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
99 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
100 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
101
102 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
103 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
104 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
105
106 //
107 // Get it into a vector in local coords and with the right units:
108
109
110 // Now Navigate the history to know in what detector the step is:
111 // and finally set the ID of det element in which the hit is.
112 //
113 //G4int History;
114 //
115 // Is it TB, barrel or endcap?
116 //
117 int BrlEcap = 0;
118 int LayerDisk = 0;
119 int etaMod = 0;
120 int phiMod = 0;
121 int side = 0;
122
123 // G4VPhysicalVolume* TheVolume = myTouch->GetVolume();
124 //
125 // In the case of the TB the positioning integers won't be initialized
126 // and the identifying integer will be zero. There is no need to do
127 // anything else
128
129 int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
130
131 // Barrel
132 if(BEcopyNo == 100) {
133 // Move up of one level in the tree to get the module number:
134 //
135 etaMod = myTouch->GetVolume(1)->GetCopyNo();
136 //
137 // Another level for the ladder
138 //
139 phiMod = myTouch->GetVolume(2)->GetCopyNo();
140 //
141 // Another one for the layer
142 //
143 LayerDisk = myTouch->GetVolume(3)->GetCopyNo();
144 //
145 // Some printout
146 //
147 if (verboseLevel>5){
148 G4cout << "In the Pixel Barrel" << G4endl;
149 G4cout << "----- Phi Module # " << phiMod << G4endl;
150 G4cout << "----- Eta Ladder # " << etaMod << G4endl;
151 G4cout << "----- Layer # " << LayerDisk << G4endl;
152 }
153
154 // Standard ATLAS EndCap
155 } else if (BEcopyNo == 200) {
156
157 // Use the copy no. to get the id:
158 phiMod = myTouch->GetVolume(1)->GetCopyNo();
159 if (phiMod == 48) phiMod = 0;
160
161 // etaMod always 0.
162
163 // Move up of one level in the tree to get the disk number:
164 LayerDisk = myTouch->GetVolume(2)->GetCopyNo();
165
166 // Move up of one level in the tree to get the EC number:
167 BrlEcap = myTouch->GetVolume(3)->GetCopyNo();
168 //The following is no longer necessary
169 //if(BrlEcap == -2) BrlEcap = 1;
170
171 // Workaround for bug in cosmic setup
172 if (BrlEcap == 0) {
173 BrlEcap = (coord1.z() > 0) ? 2 : -2;
174 }
175
176 // Some printout
177 if (verboseLevel>5){
178 G4cout << "In the Pixel EndCap" << G4endl;
179 G4cout << "----- PhiModule # " << phiMod << G4endl;
180 G4cout << "----- Disk # " << LayerDisk << G4endl;
181 G4cout << "----- Endcap # " << BrlEcap << G4endl;
182 }
183
184 // SLHC EndCap
185 } else if (BEcopyNo == 300) {
186
187 // Use the copy no. to get the id:
188 phiMod = myTouch->GetVolume(1)->GetCopyNo();
189
190 // Note: copyNo is 2*ieta for front face, 2*ieta+1 for back face
191 etaMod = myTouch->GetVolume(2)->GetCopyNo()/2;
192
193 // Move up of one level in the tree to get the disk number:
194 LayerDisk = myTouch->GetVolume(3)->GetCopyNo();
195
196 // Move up of one level in the tree to get the EC number:
197 BrlEcap = myTouch->GetVolume(4)->GetCopyNo();
198 //if(BrlEcap == -2) BrlEcap = 1;
199
200 // Some printout
201 if (verboseLevel>5){
202 G4cout << "In the SLHC Pixel EndCap" << G4endl;
203 G4cout << "----- PhiModule # " << phiMod << G4endl;
204 G4cout << "----- Ring/Eta # " << etaMod << G4endl;
205 G4cout << "----- Disk # " << LayerDisk << G4endl;
206 G4cout << "----- Endcap # " << BrlEcap << G4endl;
207 }
208
209 } else if(BEcopyNo == 400) { //DBM
210
211 // Note: there is one volume (3-layers unit) contain the diamond+chip volume
212 // so one increment is needed after 'layerDisk'
213
214 // only one eta module for the DBM
215 etaMod = 0;
216
217 // Move up one level to get the layer number, 3 layers per telescope
218 LayerDisk = myTouch->GetVolume(1)->GetCopyNo();
219
220 // Move up two level to get the phi module, azimuthal position of telescope
221 phiMod = myTouch->GetVolume(3)->GetCopyNo();
222
223 // Move up one level to get the DBM side, A side or C side
224 BrlEcap = myTouch->GetVolume(4)->GetCopyNo();
225
226 // Some printout
227 if (verboseLevel>5){
228 G4cout << "In the DBM" << G4endl;
229 G4cout << "----- PhiModule # " << phiMod << G4endl;
230 G4cout << "----- Ring/Eta # " << etaMod << G4endl;
231 G4cout << "----- Disk # " << LayerDisk << G4endl;
232 G4cout << "----- Endcap # " << BrlEcap << G4endl;
233 }
234
235 } else if(BEcopyNo == 500) { //Inclined
236
237 std::string volName = myTouch->GetVolume()->GetName();
238
239 // split the volname vs '_'
240 std::replace(volName.begin(),volName.end(),'_',' ');
241
242 std::vector<std::string> v;
243 std::istringstream s(volName);
244 std::string tmp;
245 while (s >> tmp) {
246 v.push_back(tmp);
247 }
248
249 BrlEcap = 0; // no endcap defined in alpine
250 LayerDisk = atoi(v[1].c_str());
251 phiMod = atoi(v[2].c_str());
252 etaMod = atoi(v[3].c_str());
253
254 if (verboseLevel>5){
255 G4cout << "Volume name " << volName <<G4endl;
256 double xpos = coord1.x();
257 double ypos = coord1.y();
258 double zpos = coord1.z();
259 double r = sqrt(xpos*xpos+ypos*ypos);
260 G4cout << "In the Alpine " << G4endl;
261 G4cout << "----- PhiModule # " << phiMod << G4endl;
262 G4cout << "----- Ring/Eta # " << etaMod << G4endl;
263 G4cout << "----- Disk # " << LayerDisk << G4endl;
264 G4cout << "----- Endcap # " << BrlEcap << G4endl;
265 G4cout << "----- Pos # " << r<<" "<<zpos << G4endl;
266 G4cout << "----- volume " << myTouch->GetVolume()->GetName()<< G4endl;
267 G4cout << " " << myTouch->GetVolume(1)->GetName()<<" "<<
268 myTouch->GetVolume(2)->GetName()<< " " << myTouch->GetVolume(3)->GetName()<<G4endl;
269 }
270
271 } else if(BEcopyNo == 600) { //ITk ECring
272
273 std::string volName = myTouch->GetVolume()->GetName();
274
275 // split the volname vs '_'
276 std::replace(volName.begin(),volName.end(),'_',' ');
277
278 std::vector<std::string> v;
279 std::istringstream s(volName);
280 std::string tmp;
281 while (s >> tmp) {
282 v.push_back(tmp);
283 }
284
285 BrlEcap = atoi(v[1].c_str());
286 LayerDisk = atoi(v[2].c_str());
287 phiMod = atoi(v[3].c_str());
288 etaMod = atoi(v[4].c_str());
289
290 if (verboseLevel>5){
291 double xpos = coord1.x();
292 double ypos = coord1.y();
293 double zpos = coord1.z();
294 double r = sqrt(xpos*xpos+ypos*ypos);
295 G4cout << "Volume name " << volName <<G4endl;
296 G4cout << "In the ITk EC ring " << G4endl;
297 G4cout << "----- PhiModule # " << phiMod << G4endl;
298 G4cout << "----- Ring/Eta # " << etaMod << G4endl;
299 G4cout << "----- Disk # " << LayerDisk << G4endl;
300 G4cout << "----- Endcap # " << BrlEcap << G4endl;
301 G4cout << "----- Pos # " << r<<" "<<zpos << G4endl;
302 G4cout << "----- volume " << myTouch->GetVolume()->GetName()<< G4endl;
303 G4cout << " " << myTouch->GetVolume(1)->GetName()<<" "<<
304 myTouch->GetVolume(2)->GetName()<< " " << myTouch->GetVolume(3)->GetName()<<G4endl;
305 }
306 }
307 else {
308 // Do not expect other numbers. Need to fix PixelGeoModel if this occurs.
309 G4ExceptionDescription description;
310 description << "ProcessHits: Unrecognized geometry in Pixel sensitive detector. Please contact the maintainer of the Pixel Detector Description.";
311 description << "If you are processing a GeoModelXML geometry, please ensure a PixelSensorGmxSD instance is used instead of this (PixelSensorSD).";
312 G4Exception("PixelSensorSD", "UnrecognizedPixelGeometry", FatalException, description);
313 abort();
314 }
315
316 // get the HepMcParticleLink from the TrackHelper
317
318 TrackHelper trHelp(aStep->GetTrack());
319 // Temporary solution to get EventContext from a Geant4 thread
320 m_HitColl->Emplace(lP1,
321 lP2,
322 edep,
323 aStep->GetPreStepPoint()->GetGlobalTime(),
324 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
325 0,BrlEcap,LayerDisk,etaMod,phiMod,side);
326 return true;
327}
AtlasHitsVector< SiHit > SiHitCollection
This class is attached to G4Event objects as UserInformation.
void Initialize(G4HCofThisEvent *) override final
PixelSensorSD(const std::string &name, const std::string &hitCollectionName)
AtlasG4EventUserInfo * m_g4UserEventInfo
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
SiHitCollection * m_HitColl
std::string m_HitCollName
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
int r
Definition globals.cxx:22
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91