ATLAS Offline Software
SctSensorSD.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 //
6 // SCT 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 headers
12 #include "SctSensorSD.h"
13 
14 // athena includes
15 #include "MCTruth/TrackHelper.h"
16 
17 // Geant4 includes
18 #include "G4Step.hh"
19 #include "G4ThreeVector.hh"
20 #include "G4SDManager.hh"
21 #include "G4Geantino.hh"
22 #include "G4ChargedGeantino.hh"
23 
24 // CLHEP transform
25 #include "CLHEP/Geometry/Transform3D.h"
26 
27 #include <memory> // For make unique
28 
29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30 
31 SctSensorSD::SctSensorSD( const std::string& name, const std::string& hitCollectionName )
32  : G4VSensitiveDetector( name )
33  , m_HitColl( hitCollectionName )
34 {
35 }
36 
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
38 
39 void SctSensorSD::Initialize(G4HCofThisEvent *)
40 {
41  if (!m_HitColl.isValid()) m_HitColl = std::make_unique<SiHitCollection>();
42 }
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45 
46 G4bool SctSensorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
47 {
48  double edep = aStep->GetTotalEnergyDeposit();
49  if(edep==0.) {
50  if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() &&
51  aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition())
52  return false;
53  }
54  edep *= CLHEP::MeV;
55  //
56  // Get the Touchable History:
57  //
58  const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
59  if (not myTouch) {
60  G4cout << "SctSensorSD::ProcessHits bad dynamic_cast" << G4endl;
61  return false;
62  }
63  //
64  // Get the hit coordinates. Start and End Point
65  //
66  G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
67  G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
68  //
69  // Calculate the local step begin and end position.
70  // From a G4 FAQ:
71  // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
72  //
73  const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
74  G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
75  G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
76  //
77  // Get it into a vector in local coords and with the right units:
78  //
79  HepGeom::Point3D<double> lP1,lP2;
80 
81  lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
82  lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
83  lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
84 
85  lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
86  lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
87  lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
88 
89  // Now Navigate the history to know in what detector the step is:
90  // and finally set the ID of det element in which the hit is.
91  //
92  //G4int History;
93  //
94  // Is it TB, barrel or endcap?
95  //
96  int brlEcap = 0;
97  int layerDisk = 0;
98  int etaMod = 0;
99  int phiMod = 0;
100  int side = 0;
101  this->indexMethod(myTouch, coord1.z(), brlEcap, layerDisk, etaMod, phiMod, side);
102  // get the HepMcParticleLink from the TrackHelper
103  TrackHelper trHelp(aStep->GetTrack());
104  m_HitColl->Emplace(lP1,
105  lP2,
106  edep,
107  aStep->GetPreStepPoint()->GetGlobalTime(),//use the global time. i.e. the time from the beginning of the event
108  trHelp.GenerateParticleLink(),
109  1,brlEcap,layerDisk,etaMod,phiMod,side);
110  return true;
111 }
112 
113 void SctSensorSD::indexMethod(const G4TouchableHistory *myTouch, double coord1z,
114  int &brlEcap, int &layerDisk, int &etaMod, int &phiMod, int &side) {
115 
116 
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  const int sensorCopyNo = myTouch->GetVolume()->GetCopyNo();
123 
124  if (verboseLevel>5){
125  for (int i=0;i<myTouch->GetHistoryDepth();i++){
126  std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
127  int copyno=myTouch->GetVolume(i)->GetCopyNo();
128  G4cout << "Volume "<<detname << " Copy Nr. " << copyno << G4endl;
129  }
130  }
131 
132  // Barrel geometry from GeoModel. Flag it with a 1000 in the IdTag
133  //
134  if(sensorCopyNo/1000) {
135  // Barrel
136  if(sensorCopyNo == 1000) {
137  side = myTouch->GetVolume(1)->GetCopyNo();
138  etaMod = myTouch->GetVolume(2)->GetCopyNo();
139  phiMod = myTouch->GetVolume(3)->GetCopyNo();
140  layerDisk = myTouch->GetVolume(5)->GetCopyNo();
141  } else if(sensorCopyNo == 1100) {
142  // SLHC stave layout
143  int etaModTmp = myTouch->GetVolume(1)->GetCopyNo();
144  int sign = (etaModTmp>=0)? +1 : -1;
145  side = std::abs(etaModTmp)%2;
146  etaMod = sign * std::abs(etaModTmp)/2;
147  phiMod = myTouch->GetVolume(2)->GetCopyNo();
148  layerDisk = myTouch->GetVolume(4)->GetCopyNo();
149  } else if(sensorCopyNo/100 == 12) {
150  // Segmented sensor (SLHC only)
151  int iSegment = sensorCopyNo%100;
152  int sensorEnvCopyNo = myTouch->GetVolume(1)->GetCopyNo();
153  if (sensorEnvCopyNo == 1000) {
154  side = myTouch->GetVolume(2)->GetCopyNo();
155  etaMod = myTouch->GetVolume(3)->GetCopyNo() + iSegment;
156  phiMod = myTouch->GetVolume(4)->GetCopyNo();
157  layerDisk = myTouch->GetVolume(6)->GetCopyNo();
158  } else if (sensorEnvCopyNo == 1100) {
159  // Stave layout
160  int etaModTmp = myTouch->GetVolume(2)->GetCopyNo();
161  int sign = (etaModTmp>=0)? +1 : -1;
162  side = std::abs(etaModTmp)%2;
163  etaMod = sign * std::abs(etaModTmp)/2 + iSegment;
164  phiMod = myTouch->GetVolume(3)->GetCopyNo();
165  layerDisk = myTouch->GetVolume(5)->GetCopyNo();
166  } else {
167  G4ExceptionDescription description;
168  description << "ProcessHits: Unrecognized geometry in SCT sensitive detector. Please contact the maintainer of the SCT Detector Description.";
169  G4Exception("SctSensorSD", "UnrecognizedSCTGeometry1", FatalException, description);
170  abort();
171  }
172  } else {
173  G4ExceptionDescription description;
174  description << "ProcessHits: Unrecognized geometry in SCT sensitive detector. Please contact the maintainer of the SCT Detector Description.";
175  G4Exception("SctSensorSD", "UnrecognizedSCTGeometry2", FatalException, description);
176  abort();
177  }
178  if (verboseLevel>5){
179  G4cout << "In the SCT Barrel" << G4endl;
180  G4cout << "----- side : " << side << G4endl;
181  G4cout << "----- eta_module : " << etaMod << G4endl;
182  G4cout << "----- phi_module : " << phiMod << G4endl;
183  G4cout << "----- layer : " << layerDisk << G4endl;
184  G4cout << "----- barrel_ec : " << brlEcap<< G4endl;
185  }
186  } else {
187  // Endcap geometry from GeoModel. Flag it with a 500 or 600 in the IdTag
188  // The level in the hierarchy is different for different geometries
189  // 500 For pre DC3
190  // 600 For DC3 and later (including SLHC)
191  int sensorCopyNoTest = sensorCopyNo/100;
192  if(sensorCopyNoTest == 5 || sensorCopyNoTest == 6) {
193  int signZ = (coord1z > 0) ? 1 : -1;
194  brlEcap = 2 * signZ;
195  side = myTouch->GetVolume(0)->GetCopyNo() % 2;
196  if (sensorCopyNoTest == 5) { // DC2 and Rome
197  etaMod = myTouch->GetVolume(3)->GetCopyNo();
198  layerDisk = myTouch->GetVolume(4)->GetCopyNo();
199  } else { // DC3 and later and SLHC
200  etaMod = myTouch->GetVolume(2)->GetCopyNo();
201  layerDisk = myTouch->GetVolume(3)->GetCopyNo();
202  }
203  int phiNumber = myTouch->GetVolume(1)->GetCopyNo();
204  // Number of modules in ring and module which becomes 0 is encoded in the copy number.
205  // This is in order to recreate the correct id in the negative endcap.
206  int maxModules = (phiNumber & 0x00ff0000) >> 16;
207  int moduleZero = (phiNumber & 0xff000000) >> 24; // When phiMod = moduleZero -> 0 after inverting
208  phiMod = phiNumber & 0x0000ffff;
209  // If -ve endcap then invert numbering
210  // maxModules is non-zero, but check for maxModules != 0 safeguards
211  // against divide by zero (in case we create a geometry without the encoding).
212  if (maxModules && signZ < 0) phiMod = (maxModules + moduleZero - phiMod) % maxModules;
213 
214  // Some printout
215  if (verboseLevel>5){
216  G4cout << "In the SCT EndCap" << G4endl;
217  G4cout << "----- side : " << side << G4endl;
218  G4cout << "----- phi_module : " << phiMod<< G4endl;
219  G4cout << "----- eta_module : " << etaMod << G4endl;
220  G4cout << "----- layerDisk : " << layerDisk << G4endl;
221  G4cout << "----- barrel_ec : " << brlEcap << G4endl;
222  }
223  } else {
224  // Do not expect other numbers. Need to fix SCT_GeoModel or SCT_SLHC_GeoModel if this occurs.
225  G4ExceptionDescription description;
226  description << "ProcessHits: Unrecognized geometry in SCT sensitive detector. Please contact the maintainer of the SCT Detector Description.";
227  G4Exception("SctSensorSD", "UnrecognizedSCTGeometry3", FatalException, description);
228  abort();
229  }
230  }
231  return;
232 }
SctSensorSD::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: SctSensorSD.cxx:46
SiHit::xDep
@ xDep
Definition: SiHit.h:162
SctSensorSD::SctSensorSD
SctSensorSD(const std::string &name, const std::string &hitCollectionName)
Definition: SctSensorSD.cxx:31
TrackHelper.h
SiHit::xPhi
@ xPhi
Definition: SiHit.h:162
python.SystemOfUnits.MeV
int MeV
Definition: SystemOfUnits.py:154
SctSensorSD::m_HitColl
SG::WriteHandle< SiHitCollection > m_HitColl
Definition: SctSensorSD.h:52
SctSensorSD.h
TRT::Hit::side
@ side
Definition: HitInfo.h:83
TrackHelper
Definition: TrackHelper.h:14
lumiFormat.i
int i
Definition: lumiFormat.py:85
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:107
SctSensorSD::Initialize
void Initialize(G4HCofThisEvent *) override final
Definition: SctSensorSD.cxx:39
SctSensorSD::indexMethod
void indexMethod(const G4TouchableHistory *myTouch, double coord1z, int &brlEcap, int &layerDisk, int &etaMod, int &phiMod, int &side)
Definition: SctSensorSD.cxx:113
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
SiHit::xEta
@ xEta
Definition: SiHit.h:162
python.SystemOfUnits.mm
int mm
Definition: SystemOfUnits.py:83
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