ATLAS Offline Software
Loading...
Searching...
No Matches
PixelSensorGmxSD.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5//
6// Pixel Sensitive Detector - specialisation for GeoModelXml
7// The Hits are processed here. For every hit I get the position and
8// an information on the sensor in which the interaction happened
9//
10
11// Class header
12#include "PixelSensorGmxSD.h"
13
14// Athena headers
16#include "MCTruth/TrackHelper.h"
17
18// Geant4 headers
19#include "G4ChargedGeantino.hh"
20#include <G4EventManager.hh>
21#include "G4Geantino.hh"
22#include "G4SDManager.hh"
23#include "G4Step.hh"
24#include "G4ThreeVector.hh"
25
26// CLHEP headers
27#include "CLHEP/Geometry/Transform3D.h"
28#include "CLHEP/Units/SystemOfUnits.h"
29
30#include <GeoModelKernel/GeoFullPhysVol.h>
31#include <GeoModelRead/ReadGeoModel.h>
32
34
35
36PixelSensorGmxSD::PixelSensorGmxSD(const std::string& name, const std::string& hitCollectionName,GeoModelIO::ReadGeoModel * sqlreader)
37 : G4VSensitiveDetector( name )
38 , m_HitCollName( hitCollectionName )
39{
40 m_sqlreader = sqlreader;
41}
42
43// Initialize from G4 - cache the hit collection for the current event
44void PixelSensorGmxSD::Initialize(G4HCofThisEvent *)
45{
46 // ISF calls G4SDManager::PrepareNewEvent() before the Geant4 event loop starts...
47 if(auto* eventManger = G4EventManager::GetEventManager()){
48 if(auto* eventInfo = static_cast<AtlasG4EventUserInfo*>(eventManger->GetUserInformation())){
49 m_HitColl = eventInfo->GetHitCollectionMap()->Find<SiHitCollection>(m_HitCollName);
50 m_g4UserEventInfo = eventInfo;
51 }
52 }
53
54}
55
56
57G4bool PixelSensorGmxSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
58{
59 if (verboseLevel>5) G4cout << "Process Hit" << G4endl;
60
61 G4double edep = aStep->GetTotalEnergyDeposit();
62 edep *= CLHEP::MeV;
63 if(edep==0.) {
64 if(aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
65 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
66 return false;
67 }
68
69 //use the global time. i.e. the time from the beginning of the event
70 //
71 // Get the Touchable History:
72 //
73 const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
74 if (not myTouch) {
75 G4cout << "PixelSensorGmxSD::ProcessHits bad dynamic_cast" << G4endl;
76 return false;
77 }
78 if(verboseLevel>5){
79 for (int i=0;i<myTouch->GetHistoryDepth();i++){
80 std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
81 int copyno=myTouch->GetVolume(i)->GetCopyNo();
82 G4cout << "Volume " <<detname <<" Copy Nr. " << copyno << G4endl;
83 }
84 }
85 //
86 // Get the hit coordinates. Start and End Point
87 //
88 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
89 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
90
91 // Calculate the local step begin and end position.
92 // From a G4 FAQ:
93 // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
94 //
95 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
96 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
97 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
98
99 HepGeom::Point3D<double> lP1,lP2;
100 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
101 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
102 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
103
104 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
105 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
106 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
107
108 TrackHelper trHelp(aStep->GetTrack());
109 auto mcParticleLink = trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr);
110
111 if(m_sqlreader){
112 //if sqlite inputs, Identifier indices come from PhysVol Name
113 std::string physVolName = myTouch->GetVolume(0)->GetName();
114
115 int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringITk(0,physVolName);
116
117
118 m_HitColl->Emplace(lP1,
119 lP2,
120 edep,
121 aStep->GetPreStepPoint()->GetGlobalTime(),//use the global time. i.e. the time from the beginning of the event
122 std::move(mcParticleLink),
123 hitIdOfWafer);
124 return true;
125
126 }
127 // if not from SQLite, we assume that the Identifier has already been written in as the copy number
128 // (it should have done if GeoModel building ran within Athena)
129
130 int id = myTouch->GetVolume()->GetCopyNo();
131
132 m_HitColl->Emplace(lP1,
133 lP2,
134 edep,
135 aStep->GetPreStepPoint()->GetGlobalTime(),
136 std::move(mcParticleLink),
137 id);
138 return true;
139}
AtlasHitsVector< SiHit > SiHitCollection
This class is attached to G4Event objects as UserInformation.
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
PixelSensorGmxSD(const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
std::string m_HitCollName
virtual void Initialize(G4HCofThisEvent *) override final
AtlasG4EventUserInfo * m_g4UserEventInfo
SiHitCollection * m_HitColl
GeoModelIO::ReadGeoModel * m_sqlreader
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
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52