ATLAS Offline Software
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
18 #include "MCTruth/TrackHelper.h"
19 
20 // Geant4 includes
21 #include "G4Step.hh"
22 #include "G4ThreeVector.hh"
23 #include "G4SDManager.hh"
24 #include "G4Geantino.hh"
25 #include "G4ChargedGeantino.hh"
26 
27 #include <GeoModelKernel/GeoFullPhysVol.h>
28 #include <GeoModelRead/ReadGeoModel.h>
29 
31 
32 SctSensorGmxSD::SctSensorGmxSD(const std::string& name, const std::string& hitCollectionName,GeoModelIO::ReadGeoModel * sqlreader)
33  : SctSensorSD( name, hitCollectionName) {
34  m_sqlreader = sqlreader;
35  }
36 
38 
39 G4bool SctSensorGmxSD::ProcessHits(G4Step *aStep, G4TouchableHistory * /* not used */) {
40  double edep = aStep->GetTotalEnergyDeposit();
41  if (edep == 0.) {
42  if (aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
43  aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
44  return false;
45  }
46  edep *= CLHEP::MeV;
47 
48  const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
49  if (not myTouch) {
50  G4cout << "SctSensorGmxSD::ProcessHits bad dynamic_cast" << G4endl;
51  return false;
52  }
53  //
54  // Get the hit start and end point local coordinates
55  //
56  G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
57  G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
58  const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
59  G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
60  G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
61  HepGeom::Point3D<double> lP1, lP2;
62 
63 
64  lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
65  lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
66  lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
67 
68  lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
69  lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
70  lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
71 
72 
73  // get the HepMcParticleLink from the TrackHelper
74  TrackHelper trHelp(aStep->GetTrack());
75 
76  if(m_sqlreader){
77  //if sqlite inputs, Identifier indices come from PhysVol Name
78  std::string physVolName = myTouch->GetVolume(0)->GetName();
79 
80  int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringITk(1,physVolName);
81 
82  m_HitColl->Emplace(lP1,
83  lP2,
84  edep,
85  aStep->GetPreStepPoint()->GetGlobalTime(),//use the global time. i.e. the time from the beginning of the event
86  trHelp.GenerateParticleLink(),
87  hitIdOfWafer);
88  return true;
89  }
90 
91  // if not from SQLite, we assume that the Identifier has already been written in as the copy number
92  // (it should hsave done if GeoModel building ran within Athena)
93  //
94  // Get the indexes of which detector the hit is in
95  //
96  const int id = myTouch->GetVolume(0)->GetCopyNo();
97 
98  m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(), trHelp.GenerateParticleLink(), id);
99 
100  return true;
101 }
SiHit::xDep
@ xDep
Definition: SiHit.h:162
TrackHelper.h
SiHit::xPhi
@ xPhi
Definition: SiHit.h:162
python.SystemOfUnits.MeV
int MeV
Definition: SystemOfUnits.py:154
SctSensorSD::m_HitColl
SG::WriteHandle< SiHitCollection > m_HitColl
Definition: SctSensorSD.h:52
SctSensorGmxSD::~SctSensorGmxSD
~SctSensorGmxSD()
Definition: SctSensorGmxSD.cxx:37
SctSensorGmxSD.h
SiHitIdHelper.h
TrackHelper
Definition: TrackHelper.h:14
SctSensorGmxSD::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
Definition: SctSensorGmxSD.cxx:39
SctSensorGmxSD::SctSensorGmxSD
SctSensorGmxSD(const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
Definition: SctSensorGmxSD.cxx:32
SctSensorGmxSD::m_sqlreader
GeoModelIO::ReadGeoModel * m_sqlreader
Definition: SctSensorGmxSD.h:37
SiHitIdHelper::buildHitIdFromStringITk
int buildHitIdFromStringITk(int part, const std::string &) const
Definition: SiHitIdHelper.cxx:131
SctSensorSD
Definition: SctSensorSD.h:25
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
SiHit::xEta
@ xEta
Definition: SiHit.h:162
python.SystemOfUnits.mm
int mm
Definition: SystemOfUnits.py:83
SiHitIdHelper::GetHelper
static const SiHitIdHelper * GetHelper()
Definition: SiHitIdHelper.cxx:19
TrackHelper::GenerateParticleLink
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition: TrackHelper.h:35