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

#include <HGTDSensorSD.h>

Inheritance diagram for HGTDSensorSD:
Collaboration diagram for HGTDSensorSD:

Public Member Functions

 HGTDSensorSD (const std::string &name, const std::string &hitCollectionName)
virtual ~HGTDSensorSD ()
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
void Initialize (G4HCofThisEvent *) override final

Private Member Functions

SiHitCollectiongetHitCollection () const

Private Attributes

std::string m_hitCollectionName
SiHitCollectionm_hitColl {}

Detailed Description

Definition at line 28 of file HGTDSensorSD.h.

Constructor & Destructor Documentation

◆ HGTDSensorSD()

HGTDSensorSD::HGTDSensorSD ( const std::string & name,
const std::string & hitCollectionName )

Definition at line 31 of file HGTDSensorSD.cxx.

32 : G4VSensitiveDetector( name ),
33 m_hitCollectionName( hitCollectionName )
34{
35
36}
std::string m_hitCollectionName

◆ ~HGTDSensorSD()

virtual HGTDSensorSD::~HGTDSensorSD ( )
inlinevirtual

Definition at line 37 of file HGTDSensorSD.h.

37{}

Member Function Documentation

◆ getHitCollection()

SiHitCollection * HGTDSensorSD::getHitCollection ( ) const
private

Definition at line 193 of file HGTDSensorSD.cxx.

194{
195 auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo();
196 if (!eventInfo) {
197 return nullptr;
198 }
199
200 auto hitCollections = eventInfo->GetHitCollectionMap();
201 return hitCollections ? hitCollections->Find<SiHitCollection>(m_hitCollectionName) : nullptr;
202}
AtlasHitsVector< SiHit > SiHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()

◆ Initialize()

void HGTDSensorSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 39 of file HGTDSensorSD.cxx.

40{
42}
SiHitCollection * m_hitColl
SiHitCollection * getHitCollection() const

◆ ProcessHits()

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

Definition at line 44 of file HGTDSensorSD.cxx.

45{
46 if (!m_hitColl) {
48 if (!m_hitColl) {
49 return false;
50 }
51 }
52
53 if (verboseLevel>5) G4cout << "Process Hit" << G4endl;
54
55 G4double edep = aStep->GetTotalEnergyDeposit();
56 edep *= CLHEP::MeV;
57
58 if (edep==0.) {
59 if (aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
60 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
61 return false;
62 }
63
64 //
65 // Get the Touchable History:
66 //
67 const G4TouchableHistory* myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
68 if (not myTouch){
69 G4cout << "HGTDSensorSD::ProcessHits: Dynamic cast failed"<<G4endl;
70 return false;
71 }
72 if(verboseLevel>5){
73 for (int i=0;i<myTouch->GetHistoryDepth();i++){
74 std::string detname = myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
75 int copyno = myTouch->GetVolume(i)->GetCopyNo();
76 G4cout << "Volume " << detname << " Copy Nr. " << copyno << G4endl;
77 }
78 }
79
80 //
81 // Get the hit coordinates. Start and End Point
82 //
83 G4ThreeVector startCoord = aStep->GetPreStepPoint()->GetPosition();
84 G4ThreeVector endCoord = aStep->GetPostStepPoint()->GetPosition();
85 // get the average position for writing
86 G4ThreeVector globalPosition = (startCoord + endCoord)/2;
87
88 // now deal with the Identifier
89 int posNegEndcap = ( globalPosition.z() > 0. ? 1 : -1 );
90
91 // make sure that it's a Sensor:
92 std::string detname_layer=myTouch->GetVolume(0)->GetLogicalVolume()->GetName();
93 if ( detname_layer.find("Sensor") == std::string::npos ) {
94 // Do not expect other names. Need to fix HGTDGeoModel if this occurs.
95 G4ExceptionDescription description;
96 description << "ProcessHits: No HGTD sensitive detector with substring Sensor found. Check HGTD Detector Description.";
97 G4Exception("HGTDSensorSD", "UnrecognizedHGTDGeometry", FatalException, description);
98 abort();
99 }
100
101 if(verboseLevel>5){
102 // The following code has been used to check the transformations between the detector elements
103 // TRANSFORM from local to global / step by step
104 for ( int i = myTouch->GetHistory()->GetDepth(); i >= 0 ; i-- ) {
105 std::string detname = myTouch->GetHistory()->GetVolume(i)->GetLogicalVolume()->GetName();
106 G4VisExtent extent = myTouch->GetHistory()->GetVolume(i)->GetLogicalVolume()->GetSolid()->GetExtent();
107 int copyno=myTouch->GetHistory()->GetVolume(i)->GetCopyNo();
108 const G4AffineTransform transformation1 = myTouch->GetHistory()->GetTransform( i ); // Transformation from global to current
109 const G4AffineTransform transformationInverse = transformation1.Inverse(); // Transformation from current to global
110 G4ThreeVector pos_center_local(0.0, 0.0, 0.0);
111 G4ThreeVector pos_center_global = transformationInverse.TransformPoint( pos_center_local );
112 G4ThreeVector pos_current = transformation1.TransformPoint( globalPosition );
113 G4cout << "DEBUG HGTDG4SD : "
114 << "VOLUME: " << i << " detname: " << detname
115 << ", center of element : " << pos_center_global.x()*CLHEP::mm << ", y: " << pos_center_global.y()*CLHEP::mm << ", z: " << pos_center_global.z()*CLHEP::mm
116 << ", extent: x: " << extent.GetXmax() - extent.GetXmin() << ", y: " << extent.GetYmax() - extent.GetYmin() << ", z: " << extent.GetZmax() - extent.GetZmin()
117 << ", copyno: " << copyno << G4endl;
118 G4cout << "DEBUG HGTDG4SD : LOCAL: x: " << pos_current.x()*CLHEP::mm << ", y: " << pos_current.y()*CLHEP::mm << ", z: " << pos_current.z()*CLHEP::mm << G4endl;
119
120 if ( i >= 1) {
121 const G4AffineTransform transformation2 = myTouch->GetHistory()->GetTransform( i-1 ); // Transformation from global to 1 up
122 G4AffineTransform transformation_up; // Transformation from current to 1 up
123 transformation_up.Product( transformationInverse, transformation2 );
124 G4ThreeVector pos_up = transformation_up.TransformPoint( pos_current );
125 G4RotationMatrix rotmat = transformation_up.NetRotation(); // https://www-zeuthen.desy.de/geant4/clhep-2.0.4.3/classCLHEP_1_1HepRotation.html
126 G4ThreeVector translation = transformation_up.NetTranslation(); // https://www-zeuthen.desy.de/ILC/geant4/clhep-2.0.4.3/classCLHEP_1_1Hep3Vector.html
127 G4cout << "DEBUG HGTDG4SD : Rotation:"
128 << "| xx:" << rotmat.xx() << ", xy: " << rotmat.xy() << ", xz: " << rotmat.xz()
129 << "| yx:" << rotmat.yx() << ", yy: " << rotmat.yy() << ", yz: " << rotmat.yz()
130 << "| zx:" << rotmat.zx() << ", zy: " << rotmat.zy() << ", zz: " << rotmat.zz() << " | " << G4endl;
131 G4cout << "DEBUG HGTDG4SD : Translation: x: " << translation.x() << ", y:" << translation.y() << ", z:" << translation.z() << G4endl;
132 G4cout << "DEBUG HGTDG4SD : TRANSFORMED: x:" << pos_up.x()*CLHEP::mm << ", y:" << pos_up.y()*CLHEP::mm << ", z:" << pos_up.z()*CLHEP::mm << G4endl;
133 }
134 } // end loop stackdepth
135 }
136
137 // Create the SiHits
138
139 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
140
141 G4ThreeVector localPosition1 = transformation.TransformPoint(startCoord);
142 G4ThreeVector localPosition2 = transformation.TransformPoint(endCoord);
143
144 if(verboseLevel>5){
145 G4cout << " PreStepPoint " << G4endl;
146 G4cout << " x (global/local) " << startCoord.x()*CLHEP::mm << " " << localPosition1[0]*CLHEP::mm << G4endl;
147 G4cout << " y (global/local) " << startCoord.y()*CLHEP::mm << " " << localPosition1[1]*CLHEP::mm << G4endl;
148 G4cout << " z (global/local) " << startCoord.z()*CLHEP::mm << " " << localPosition1[2]*CLHEP::mm << G4endl;
149
150 G4cout << " PostStepPoint: " << G4endl;
151 G4cout << " x (global/local) " << endCoord.x()*CLHEP::mm << " " << localPosition2[0]*CLHEP::mm << G4endl;
152 G4cout << " y (global/local) " << endCoord.y()*CLHEP::mm << " " << localPosition2[1]*CLHEP::mm << G4endl;
153 G4cout << " z (global/local) " << endCoord.z()*CLHEP::mm << " " << localPosition2[2]*CLHEP::mm << G4endl;
154 }
155
156 HepGeom::Point3D<double> lP1,lP2;
157 lP1[SiHit::xEta] = localPosition1[1]*CLHEP::mm; //long edge of the module
158 lP1[SiHit::xPhi] = localPosition1[0]*CLHEP::mm; //short edge of the module
159 lP1[SiHit::xDep] = localPosition1[2]*CLHEP::mm; //depth (z)
160
161 lP2[SiHit::xEta] = localPosition2[1]*CLHEP::mm;
162 lP2[SiHit::xPhi] = localPosition2[0]*CLHEP::mm;
163 lP2[SiHit::xDep] = localPosition2[2]*CLHEP::mm;
164
165 std::string module_indices = myTouch->GetVolume(1)->GetLogicalVolume()->GetName();
166 std::size_t found = module_indices.find_last_of("_");
167
168 // get indices from the volume name
169 // nomenclature is expected to be e.g. "HGTDModule0_layer_0_1_2"
170 // for layer=0, phi=1, eta=2 (defined from HGTD_DetectorFactory)
171 int eta = atoi((module_indices.substr(found+1)).c_str());
172 module_indices.erase(found);
173 found = module_indices.find_last_of("_");
174 int phi = atoi((module_indices.substr(found+1)).c_str());
175 module_indices.erase(found);
176 found = module_indices.find_last_of("_");
177 int layer = atoi((module_indices.substr(found+1)).c_str());
178
179 int endcap_side = 2*posNegEndcap;
180
181 TrackHelper trHelp(aStep->GetTrack());
182 m_hitColl->Emplace(lP1,
183 lP2,
184 edep,
185 aStep->GetPreStepPoint()->GetGlobalTime(),
186 trHelp.GenerateParticleLink(),
187 // this is hgtd, endcap_barrel, layer_disk, eta_module, phi_module, side
188 2,endcap_side,layer,eta,phi,0);
189
190 return true;
191}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:93
int atoi(std::string_view str)
Helper functions to unpack numbers decoded in string into integers and doubles The strings are requir...
@ layer
Definition HitInfo.h:79

Member Data Documentation

◆ m_hitColl

SiHitCollection* HGTDSensorSD::m_hitColl {}
private

Definition at line 51 of file HGTDSensorSD.h.

51{};

◆ m_hitCollectionName

std::string HGTDSensorSD::m_hitCollectionName
private

Definition at line 49 of file HGTDSensorSD.h.


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