ATLAS Offline Software
Loading...
Searching...
No Matches
SctSensorGmxSD.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// SCT Sensitive Detectors where the sensitive volume identifier is already the required identifier
7// (previous versions built up the identifier from the parent volume numeric identifiers)
8// GeoModelXml introduced this way of doing things.
9//
10// The Geant 4 Hits are processed here. For every hit: get the local start and end points,
11// energy deposit, hit time, barcode, and identifier; put these in a list to be persistified.
12//
13
14// class headers
15#include "SctSensorGmxSD.h"
16
17// athena includes
19#include "MCTruth/TrackHelper.h"
20
21// Geant4 includes
22#include "G4Step.hh"
23#include "G4ThreeVector.hh"
24#include "G4SDManager.hh"
25#include "G4Geantino.hh"
26#include "G4ChargedGeantino.hh"
27
28#include <GeoModelKernel/GeoFullPhysVol.h>
29#include <GeoModelRead/ReadGeoModel.h>
30
32
33SctSensorGmxSD::SctSensorGmxSD(const std::string& name, const std::string& hitCollectionName,GeoModelIO::ReadGeoModel * sqlreader)
34 : SctSensorSD( name, hitCollectionName) {
35 m_sqlreader = sqlreader;
36 }
37
39
40G4bool SctSensorGmxSD::ProcessHits(G4Step *aStep, G4TouchableHistory * /* not used */) {
41 double edep = aStep->GetTotalEnergyDeposit();
42 if (edep == 0.) {
43 if (aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
44 aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
45 return false;
46 }
47 edep *= CLHEP::MeV;
48
49 const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
50 if (not myTouch) {
51 G4cout << "SctSensorGmxSD::ProcessHits bad dynamic_cast" << G4endl;
52 return false;
53 }
54 //
55 // Get the hit start and end point local coordinates
56 //
57 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
58 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
59 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
60 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
61 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
62 HepGeom::Point3D<double> lP1, lP2;
63
64
65 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
66 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
67 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
68
69 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
70 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
71 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
72
73
74 // get the HepMcParticleLink from the TrackHelper
75 TrackHelper trHelp(aStep->GetTrack());
76 auto particleLink = trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr);
77
78 if(m_sqlreader){
79 //if sqlite inputs, Identifier indices come from PhysVol Name
80 std::string physVolName = myTouch->GetVolume(0)->GetName();
81
82 int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringITk(1,physVolName);
83
84 m_HitColl->Emplace(lP1,
85 lP2,
86 edep,
87 aStep->GetPreStepPoint()->GetGlobalTime(),//use the global time. i.e. the time from the beginning of the event
88 std::move(particleLink),
89 hitIdOfWafer);
90 return true;
91 }
92
93 // if not from SQLite, we assume that the Identifier has already been written in as the copy number
94 // (it should hsave done if GeoModel building ran within Athena)
95 //
96 // Get the indexes of which detector the hit is in
97 //
98 const int id = myTouch->GetVolume(0)->GetCopyNo();
99
100 m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(), std::move(particleLink), id);
101
102 return true;
103}
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
GeoModelIO::ReadGeoModel * m_sqlreader
SctSensorGmxSD(const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
SctSensorSD(const std::string &name, const std::string &hitCollectionName)
SiHitCollection * m_HitColl
Definition SctSensorSD.h:52
AtlasG4EventUserInfo * m_g4UserEventInfo
Definition SctSensorSD.h:53
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