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