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 Attributes

SG::WriteHandle< SiHitCollectionm_HitColl

Detailed Description

Definition at line 29 of file HGTDSensorSD.h.

Constructor & Destructor Documentation

◆ HGTDSensorSD()

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

Definition at line 29 of file HGTDSensorSD.cxx.

30 : G4VSensitiveDetector( name ),
31 m_HitColl( hitCollectionName )
32{
33
34}
SG::WriteHandle< SiHitCollection > m_HitColl

◆ ~HGTDSensorSD()

virtual HGTDSensorSD::~HGTDSensorSD ( )
inlinevirtual

Definition at line 38 of file HGTDSensorSD.h.

38{}

Member Function Documentation

◆ Initialize()

void HGTDSensorSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 37 of file HGTDSensorSD.cxx.

38{
39 if (!m_HitColl.isValid()) m_HitColl = std::make_unique<SiHitCollection>();
40}

◆ ProcessHits()

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

Definition at line 42 of file HGTDSensorSD.cxx.

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

SG::WriteHandle<SiHitCollection> HGTDSensorSD::m_HitColl
private

Definition at line 48 of file HGTDSensorSD.h.


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