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

#include <HGTDSensorGmxSD.h>

Inheritance diagram for HGTDSensorGmxSD:
Collaboration diagram for HGTDSensorGmxSD:

Public Member Functions

 HGTDSensorGmxSD (const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
virtual ~HGTDSensorGmxSD ()
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
void Initialize (G4HCofThisEvent *) override final

Private Member Functions

SiHitCollectiongetHitCollection () const

Private Attributes

std::string m_hitCollectionName
SiHitCollectionm_hitColl {}
GeoModelIO::ReadGeoModel * m_sqlreader

Detailed Description

Definition at line 32 of file HGTDSensorGmxSD.h.

Constructor & Destructor Documentation

◆ HGTDSensorGmxSD()

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

Definition at line 37 of file HGTDSensorGmxSD.cxx.

38 : G4VSensitiveDetector( name ),
39 m_hitCollectionName( hitCollectionName ),
40 m_sqlreader( sqlreader )
41{
42
43}
GeoModelIO::ReadGeoModel * m_sqlreader
std::string m_hitCollectionName

◆ ~HGTDSensorGmxSD()

virtual HGTDSensorGmxSD::~HGTDSensorGmxSD ( )
inlinevirtual

Definition at line 41 of file HGTDSensorGmxSD.h.

41{}

Member Function Documentation

◆ getHitCollection()

SiHitCollection * HGTDSensorGmxSD::getHitCollection ( ) const
private

Definition at line 146 of file HGTDSensorGmxSD.cxx.

147{
148 auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo();
149 if (!eventInfo) {
150 return nullptr;
151 }
152
153 auto hitCollections = eventInfo->GetHitCollectionMap();
154 return hitCollections ? hitCollections->Find<SiHitCollection>(m_hitCollectionName) : nullptr;
155}
AtlasHitsVector< SiHit > SiHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()

◆ Initialize()

void HGTDSensorGmxSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 46 of file HGTDSensorGmxSD.cxx.

47{
49}
SiHitCollection * m_hitColl
SiHitCollection * getHitCollection() const

◆ ProcessHits()

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

Definition at line 51 of file HGTDSensorGmxSD.cxx.

52{
53 if (!m_hitColl) {
55 if (!m_hitColl) {
56 return false;
57 }
58 }
59
60 if (verboseLevel>5) G4cout << "Process Hit" << G4endl;
61
62 G4double edep = aStep->GetTotalEnergyDeposit();
63 edep *= CLHEP::MeV;
64
65 if (edep==0.) {
66 if (aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
67 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
68 return false;
69 }
70
71 //
72 // Get the Touchable History:
73 //
74 const G4TouchableHistory* myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
75 if (not myTouch){
76 G4cout<<"HGTDSensorGmxSD::ProcessHits: dynamic cast failed"<<G4endl;
77 return false;
78 }
79 if(verboseLevel>5){
80 for (int i=0;i<myTouch->GetHistoryDepth();i++){
81 std::string detname = myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
82 int copyno = myTouch->GetVolume(i)->GetCopyNo();
83 G4cout << "Volume " << detname << " Copy Nr. " << copyno << G4endl;
84 }
85 }
86
87 //
88 // Get the hit coordinates. Start and End Point
89 //
90 G4ThreeVector startCoord = aStep->GetPreStepPoint()->GetPosition();
91 G4ThreeVector endCoord = aStep->GetPostStepPoint()->GetPosition();
92
93 // Create the SiHits
94
95 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
96
97 G4ThreeVector localPosition1 = transformation.TransformPoint(startCoord);
98 G4ThreeVector localPosition2 = transformation.TransformPoint(endCoord);
99
100 HepGeom::Point3D<double> lP1,lP2;
101 //TODO need to check
102 lP1[SiHit::xEta] = localPosition1[1]*CLHEP::mm; //long edge of the module
103 lP1[SiHit::xPhi] = localPosition1[0]*CLHEP::mm; //short edge of the module
104 lP1[SiHit::xDep] = localPosition1[2]*CLHEP::mm; //depth (z)
105
106 lP2[SiHit::xEta] = localPosition2[1]*CLHEP::mm;
107 lP2[SiHit::xPhi] = localPosition2[0]*CLHEP::mm;
108 lP2[SiHit::xDep] = localPosition2[2]*CLHEP::mm;
109
110 // get the HepMcParticleLink from the TrackHelper
111 TrackHelper trHelp(aStep->GetTrack());
112
113 if(m_sqlreader){
114 // if sqlite inputs, Identifier indices come from PhysVol Name
115 std::string physVolName = myTouch->GetVolume()->GetName();
116
117 int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringHGTD(2,physVolName);
118
119 m_hitColl->Emplace(lP1,
120 lP2,
121 edep,
122 aStep->GetPreStepPoint()->GetGlobalTime(),
123 trHelp.GenerateParticleLink(),
124 hitIdOfWafer);
125
126 return true;
127 }
128
129 // if not from SQLite, we assume that the Identifier has already been written in as the copy number
130 // (it should hsave done if GeoModel building ran within Athena)
131 //
132 // Get the indexes of which detector the hit is in
133 //
134 const int id = myTouch->GetVolume()->GetCopyNo();
135
136 m_hitColl->Emplace(lP1,
137 lP2,
138 edep,
139 aStep->GetPreStepPoint()->GetGlobalTime(),
140 trHelp.GenerateParticleLink(),
141 id);
142
143 return true;
144}
int buildHitIdFromStringHGTD(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_hitColl

SiHitCollection* HGTDSensorGmxSD::m_hitColl {}
private

Definition at line 55 of file HGTDSensorGmxSD.h.

55{};

◆ m_hitCollectionName

std::string HGTDSensorGmxSD::m_hitCollectionName
private

Definition at line 53 of file HGTDSensorGmxSD.h.

◆ m_sqlreader

GeoModelIO::ReadGeoModel* HGTDSensorGmxSD::m_sqlreader
private

Definition at line 56 of file HGTDSensorGmxSD.h.


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