ATLAS Offline Software
Loading...
Searching...
No Matches
SctSensorGmxSD Class Reference

#include <SctSensorGmxSD.h>

Inheritance diagram for SctSensorGmxSD:
Collaboration diagram for SctSensorGmxSD:

Public Member Functions

 SctSensorGmxSD (const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
 ~SctSensorGmxSD ()
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
void Initialize (G4HCofThisEvent *) override final
template<class... Args>
void AddHit (Args &&... args)
 Templated method to stuff a single hit into the sensitive detector class.

Protected Attributes

std::string m_HitCollName
SiHitCollectionm_HitColl {nullptr}
AtlasG4EventUserInfom_g4UserEventInfo {nullptr}

Private Member Functions

 FRIEND_TEST (SctSensorSDtest, ProcessHits)
 FRIEND_TEST (SctSensorSDtest, indexMethod)
 FRIEND_TEST (SctSensorSDtest, AddHit)
void indexMethod (const G4TouchableHistory *myTouch, double coord1z, int &brlEcap, int &layerDisk, int &etaMod, int &phiMod, int &side)

Private Attributes

GeoModelIO::ReadGeoModel * m_sqlreader {nullptr}

Detailed Description

Definition at line 24 of file SctSensorGmxSD.h.

Constructor & Destructor Documentation

◆ SctSensorGmxSD()

SctSensorGmxSD::SctSensorGmxSD ( const std::string & name,
const std::string & hitCollectionName,
GeoModelIO::ReadGeoModel * sqlreader = nullptr )

Definition at line 33 of file SctSensorGmxSD.cxx.

34 : SctSensorSD( name, hitCollectionName) {
35 m_sqlreader = sqlreader;
36 }
GeoModelIO::ReadGeoModel * m_sqlreader
SctSensorSD(const std::string &name, const std::string &hitCollectionName)

◆ ~SctSensorGmxSD()

SctSensorGmxSD::~SctSensorGmxSD ( )

Definition at line 38 of file SctSensorGmxSD.cxx.

38{ }

Member Function Documentation

◆ AddHit()

template<class... Args>
void SctSensorSD::AddHit ( Args &&... args)
inlineinherited

Templated method to stuff a single hit into the sensitive detector class.

This could get rather tricky, but the idea is to allow fast simulations to use the very same SD classes as the standard simulation.

Definition at line 45 of file SctSensorSD.h.

45{ m_HitColl->Emplace( args... ); }
SiHitCollection * m_HitColl
Definition SctSensorSD.h:52

◆ FRIEND_TEST() [1/3]

SctSensorSD::FRIEND_TEST ( SctSensorSDtest ,
AddHit  )
privateinherited

◆ FRIEND_TEST() [2/3]

SctSensorSD::FRIEND_TEST ( SctSensorSDtest ,
indexMethod  )
privateinherited

◆ FRIEND_TEST() [3/3]

SctSensorSD::FRIEND_TEST ( SctSensorSDtest ,
ProcessHits  )
privateinherited

◆ indexMethod()

void SctSensorSD::indexMethod ( const G4TouchableHistory * myTouch,
double coord1z,
int & brlEcap,
int & layerDisk,
int & etaMod,
int & phiMod,
int & side )
privateinherited

Definition at line 119 of file SctSensorSD.cxx.

120 {
121
122
123 //
124 // In the case of the TB the positioning integers won't be initialized
125 // and the identifying integer will be zero. There is no need to do
126 // anything else
127
128 const int sensorCopyNo = myTouch->GetVolume()->GetCopyNo();
129
130 if (verboseLevel>5){
131 for (int i=0;i<myTouch->GetHistoryDepth();i++){
132 std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
133 int copyno=myTouch->GetVolume(i)->GetCopyNo();
134 G4cout << "Volume "<<detname << " Copy Nr. " << copyno << G4endl;
135 }
136 }
137
138 // Barrel geometry from GeoModel. Flag it with a 1000 in the IdTag
139 //
140 if(sensorCopyNo/1000) {
141 // Barrel
142 if(sensorCopyNo == 1000) {
143 side = myTouch->GetVolume(1)->GetCopyNo();
144 etaMod = myTouch->GetVolume(2)->GetCopyNo();
145 phiMod = myTouch->GetVolume(3)->GetCopyNo();
146 layerDisk = myTouch->GetVolume(5)->GetCopyNo();
147 } else if(sensorCopyNo == 1100) {
148 // SLHC stave layout
149 int etaModTmp = myTouch->GetVolume(1)->GetCopyNo();
150 int sign = (etaModTmp>=0)? +1 : -1;
151 side = std::abs(etaModTmp)%2;
152 etaMod = sign * std::abs(etaModTmp)/2;
153 phiMod = myTouch->GetVolume(2)->GetCopyNo();
154 layerDisk = myTouch->GetVolume(4)->GetCopyNo();
155 } else if(sensorCopyNo/100 == 12) {
156 // Segmented sensor (SLHC only)
157 int iSegment = sensorCopyNo%100;
158 int sensorEnvCopyNo = myTouch->GetVolume(1)->GetCopyNo();
159 if (sensorEnvCopyNo == 1000) {
160 side = myTouch->GetVolume(2)->GetCopyNo();
161 etaMod = myTouch->GetVolume(3)->GetCopyNo() + iSegment;
162 phiMod = myTouch->GetVolume(4)->GetCopyNo();
163 layerDisk = myTouch->GetVolume(6)->GetCopyNo();
164 } else if (sensorEnvCopyNo == 1100) {
165 // Stave layout
166 int etaModTmp = myTouch->GetVolume(2)->GetCopyNo();
167 int sign = (etaModTmp>=0)? +1 : -1;
168 side = std::abs(etaModTmp)%2;
169 etaMod = sign * std::abs(etaModTmp)/2 + iSegment;
170 phiMod = myTouch->GetVolume(3)->GetCopyNo();
171 layerDisk = myTouch->GetVolume(5)->GetCopyNo();
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", "UnrecognizedSCTGeometry1", FatalException, description);
176 abort();
177 }
178 } else {
179 G4ExceptionDescription description;
180 description << "ProcessHits: Unrecognized geometry in SCT sensitive detector. Please contact the maintainer of the SCT Detector Description.";
181 G4Exception("SctSensorSD", "UnrecognizedSCTGeometry2", FatalException, description);
182 abort();
183 }
184 if (verboseLevel>5){
185 G4cout << "In the SCT Barrel" << G4endl;
186 G4cout << "----- side : " << side << G4endl;
187 G4cout << "----- eta_module : " << etaMod << G4endl;
188 G4cout << "----- phi_module : " << phiMod << G4endl;
189 G4cout << "----- layer : " << layerDisk << G4endl;
190 G4cout << "----- barrel_ec : " << brlEcap<< G4endl;
191 }
192 } else {
193 // Endcap geometry from GeoModel. Flag it with a 500 or 600 in the IdTag
194 // The level in the hierarchy is different for different geometries
195 // 500 For pre DC3
196 // 600 For DC3 and later (including SLHC)
197 int sensorCopyNoTest = sensorCopyNo/100;
198 if(sensorCopyNoTest == 5 || sensorCopyNoTest == 6) {
199 int signZ = (coord1z > 0) ? 1 : -1;
200 brlEcap = 2 * signZ;
201 side = myTouch->GetVolume(0)->GetCopyNo() % 2;
202 if (sensorCopyNoTest == 5) { // DC2 and Rome
203 etaMod = myTouch->GetVolume(3)->GetCopyNo();
204 layerDisk = myTouch->GetVolume(4)->GetCopyNo();
205 } else { // DC3 and later and SLHC
206 etaMod = myTouch->GetVolume(2)->GetCopyNo();
207 layerDisk = myTouch->GetVolume(3)->GetCopyNo();
208 }
209 int phiNumber = myTouch->GetVolume(1)->GetCopyNo();
210 // Number of modules in ring and module which becomes 0 is encoded in the copy number.
211 // This is in order to recreate the correct id in the negative endcap.
212 int maxModules = (phiNumber & 0x00ff0000) >> 16;
213 int moduleZero = (phiNumber & 0xff000000) >> 24; // When phiMod = moduleZero -> 0 after inverting
214 phiMod = phiNumber & 0x0000ffff;
215 // If -ve endcap then invert numbering
216 // maxModules is non-zero, but check for maxModules != 0 safeguards
217 // against divide by zero (in case we create a geometry without the encoding).
218 if (maxModules && signZ < 0) phiMod = (maxModules + moduleZero - phiMod) % maxModules;
219
220 // Some printout
221 if (verboseLevel>5){
222 G4cout << "In the SCT EndCap" << G4endl;
223 G4cout << "----- side : " << side << G4endl;
224 G4cout << "----- phi_module : " << phiMod<< G4endl;
225 G4cout << "----- eta_module : " << etaMod << G4endl;
226 G4cout << "----- layerDisk : " << layerDisk << G4endl;
227 G4cout << "----- barrel_ec : " << brlEcap << G4endl;
228 }
229 } else {
230 // Do not expect other numbers. Need to fix SCT_GeoModel or SCT_SLHC_GeoModel if this occurs.
231 G4ExceptionDescription description;
232 description << "ProcessHits: Unrecognized geometry in SCT sensitive detector. Please contact the maintainer of the SCT Detector Description.";
233 G4Exception("SctSensorSD", "UnrecognizedSCTGeometry3", FatalException, description);
234 abort();
235 }
236 }
237 return;
238}
int sign(int a)
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91

◆ Initialize()

void SctSensorSD::Initialize ( G4HCofThisEvent * )
finaloverrideinherited

Definition at line 39 of file SctSensorSD.cxx.

40{
41 // ISF calls G4SDManager::PrepareNewEvent() before the Geant4 event loop starts...
42 if(auto* eventManger = G4EventManager::GetEventManager()){
43 if(auto* eventInfo = static_cast<AtlasG4EventUserInfo*>(eventManger->GetUserInformation())){
44 m_HitColl = eventInfo->GetHitCollectionMap()->Find<SiHitCollection>(m_HitCollName);
45 m_g4UserEventInfo = eventInfo;
46 }
47 }
48}
AtlasHitsVector< SiHit > SiHitCollection
std::string m_HitCollName
Definition SctSensorSD.h:51
AtlasG4EventUserInfo * m_g4UserEventInfo
Definition SctSensorSD.h:53

◆ ProcessHits()

G4bool SctSensorGmxSD::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
finaloverride

Definition at line 40 of file SctSensorGmxSD.cxx.

40 {
41 double edep = aStep->GetTotalEnergyDeposit();
42 if (edep == 0.) {
43 if (aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
44 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
45 return false;
46 }
47 edep *= CLHEP::MeV;
48
49 const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
50 if (not myTouch) {
51 G4cout << "SctSensorGmxSD::ProcessHits bad dynamic_cast" << G4endl;
52 return false;
53 }
54 //
55 // Get the hit start and end point local coordinates
56 //
57 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
58 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
59 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
60 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
61 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
62 HepGeom::Point3D<double> lP1, lP2;
63
64
65 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
66 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
67 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
68
69 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
70 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
71 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
72
73
74 // get the HepMcParticleLink from the TrackHelper
75 TrackHelper trHelp(aStep->GetTrack());
76 auto particleLink = trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr);
77
78 if(m_sqlreader){
79 //if sqlite inputs, Identifier indices come from PhysVol Name
80 std::string physVolName = myTouch->GetVolume(0)->GetName();
81
82 int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringITk(1,physVolName);
83
84 m_HitColl->Emplace(lP1,
85 lP2,
86 edep,
87 aStep->GetPreStepPoint()->GetGlobalTime(),//use the global time. i.e. the time from the beginning of the event
88 std::move(particleLink),
89 hitIdOfWafer);
90 return true;
91 }
92
93 // if not from SQLite, we assume that the Identifier has already been written in as the copy number
94 // (it should hsave done if GeoModel building ran within Athena)
95 //
96 // Get the indexes of which detector the hit is in
97 //
98 const int id = myTouch->GetVolume(0)->GetCopyNo();
99
100 m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(), std::move(particleLink), id);
101
102 return true;
103}
int buildHitIdFromStringITk(int part, const std::string &) const
static const SiHitIdHelper * GetHelper()
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162

Member Data Documentation

◆ m_g4UserEventInfo

AtlasG4EventUserInfo* SctSensorSD::m_g4UserEventInfo {nullptr}
protectedinherited

Definition at line 53 of file SctSensorSD.h.

53{nullptr};

◆ m_HitColl

SiHitCollection* SctSensorSD::m_HitColl {nullptr}
protectedinherited

Definition at line 52 of file SctSensorSD.h.

52{nullptr};

◆ m_HitCollName

std::string SctSensorSD::m_HitCollName
protectedinherited

Definition at line 51 of file SctSensorSD.h.

◆ m_sqlreader

GeoModelIO::ReadGeoModel* SctSensorGmxSD::m_sqlreader {nullptr}
private

Definition at line 37 of file SctSensorGmxSD.h.

37{nullptr};

The documentation for this class was generated from the following files: