ATLAS Offline Software
Loading...
Searching...
No Matches
HGTDSensorGmxSD.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// HGTD Sensitive Detector.
6// The Hits are processed here. For every hit, the position and information
7// on the sensor in which the interaction happened are obtained.
8
9// Class header
10#include "HGTDSensorGmxSD.h"
11
12// Athena headers
13#include "MCTruth/TrackHelper.h"
14
15// Geant4 headers
16#include "G4ChargedGeantino.hh"
17#include "G4Geantino.hh"
18#include "G4SDManager.hh"
19#include "G4Step.hh"
20#include "G4VisExtent.hh" // Added to get extent of sensors
21#include "G4VSolid.hh" // Added to get extent of sensors
22
23// GeoModel headers
24#include <GeoModelKernel/GeoFullPhysVol.h>
25#include <GeoModelRead/ReadGeoModel.h>
26
28
29// CLHEP headers
30#include "CLHEP/Geometry/Transform3D.h"
31#include "CLHEP/Units/SystemOfUnits.h"
32#include <stdint.h>
33#include <string.h>
34
35HGTDSensorGmxSD::HGTDSensorGmxSD(const std::string& name, const std::string& hitCollectionName, GeoModelIO::ReadGeoModel * sqlreader)
36 : G4VSensitiveDetector( name ),
37 m_HitColl( hitCollectionName ),
38 m_sqlreader( sqlreader )
39{
40
41}
42
43// Initialize from G4
44void HGTDSensorGmxSD::Initialize(G4HCofThisEvent *)
45{
46 if (!m_HitColl.isValid()) m_HitColl = std::make_unique<SiHitCollection>();
47}
48
49G4bool HGTDSensorGmxSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
50{
51 if (verboseLevel>5) G4cout << "Process Hit" << G4endl;
52
53 G4double edep = aStep->GetTotalEnergyDeposit();
54 edep *= CLHEP::MeV;
55
56 if (edep==0.) {
57 if (aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
58 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
59 return false;
60 }
61
62 //
63 // Get the Touchable History:
64 //
65 const G4TouchableHistory* myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
66 if (not myTouch){
67 G4cout<<"HGTDSensorGmxSD::ProcessHits: dynamic cast failed"<<G4endl;
68 return false;
69 }
70 if(verboseLevel>5){
71 for (int i=0;i<myTouch->GetHistoryDepth();i++){
72 std::string detname = myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
73 int copyno = myTouch->GetVolume(i)->GetCopyNo();
74 G4cout << "Volume " << detname << " Copy Nr. " << copyno << G4endl;
75 }
76 }
77
78 //
79 // Get the hit coordinates. Start and End Point
80 //
81 G4ThreeVector startCoord = aStep->GetPreStepPoint()->GetPosition();
82 G4ThreeVector endCoord = aStep->GetPostStepPoint()->GetPosition();
83
84 // Create the SiHits
85
86 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
87
88 G4ThreeVector localPosition1 = transformation.TransformPoint(startCoord);
89 G4ThreeVector localPosition2 = transformation.TransformPoint(endCoord);
90
91 HepGeom::Point3D<double> lP1,lP2;
92 //TODO need to check
93 lP1[SiHit::xEta] = localPosition1[1]*CLHEP::mm; //long edge of the module
94 lP1[SiHit::xPhi] = localPosition1[0]*CLHEP::mm; //short edge of the module
95 lP1[SiHit::xDep] = localPosition1[2]*CLHEP::mm; //depth (z)
96
97 lP2[SiHit::xEta] = localPosition2[1]*CLHEP::mm;
98 lP2[SiHit::xPhi] = localPosition2[0]*CLHEP::mm;
99 lP2[SiHit::xDep] = localPosition2[2]*CLHEP::mm;
100
101 // get the HepMcParticleLink from the TrackHelper
102 TrackHelper trHelp(aStep->GetTrack());
103
104 if(m_sqlreader){
105 // if sqlite inputs, Identifier indices come from PhysVol Name
106 std::string physVolName = myTouch->GetVolume()->GetName();
107
108 int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringHGTD(2,physVolName);
109
110 m_HitColl->Emplace(lP1,
111 lP2,
112 edep,
113 aStep->GetPreStepPoint()->GetGlobalTime(),
114 trHelp.GenerateParticleLink(),
115 hitIdOfWafer);
116
117 return true;
118 }
119
120 // if not from SQLite, we assume that the Identifier has already been written in as the copy number
121 // (it should hsave done if GeoModel building ran within Athena)
122 //
123 // Get the indexes of which detector the hit is in
124 //
125 const int id = myTouch->GetVolume()->GetCopyNo();
126
127 m_HitColl->Emplace(lP1,
128 lP2,
129 edep,
130 aStep->GetPreStepPoint()->GetGlobalTime(),
131 trHelp.GenerateParticleLink(),
132 id);
133
134 return true;
135}
GeoModelIO::ReadGeoModel * m_sqlreader
SG::WriteHandle< SiHitCollection > m_HitColl
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
HGTDSensorGmxSD(const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
void Initialize(G4HCofThisEvent *) override final
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
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52