ATLAS Offline Software
BCMSensorSD.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 //###############################################
6 // BCM Sensitive Detector class
7 // Bostjan Macek 15.may.2007
8 //###############################################
9 
10 // Class header
11 #include "BCMSensorSD.h"
12 
13 // Package headers
14 #include "BCMExtra.h"
15 
16 // Athena headers
18 #include "MCTruth/TrackHelper.h"
19 
20 // Geant4 headers
21 #include <G4EventManager.hh>
22 #include "G4Step.hh"
23 #include "G4ThreeVector.hh"
24 #include "G4Geantino.hh"
25 #include "G4ChargedGeantino.hh"
26 
27 // CLHEP headers
28 #include "CLHEP/Geometry/Transform3D.h"
29 #include "CLHEP/Units/SystemOfUnits.h"
30 
31 BCMSensorSD::BCMSensorSD(const std::string& name, const std::string& hitCollectionName)
32  : G4VSensitiveDetector( name )
33  , m_HitCollName( hitCollectionName )
34 {
35 }
36 
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
38 // Initialize from G4 - cache the hit collection for the current event
39 void BCMSensorSD::Initialize(G4HCofThisEvent *)
40 {
41  // ISF calls G4SDManager::PrepareNewEvent() before the Geant4 event loop starts...
42  if(auto* eventManger = G4EventManager::GetEventManager()){
43  if(auto* eventInfo = static_cast<AtlasG4EventUserInfo*>(eventManger->GetUserInformation())){
44  m_HitColl = eventInfo->GetHitCollectionMap()->Find<SiHitCollection>(m_HitCollName);
45  }
46  }
47 }
48 
49 G4bool BCMSensorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
50 {
51  G4double edep = aStep->GetTotalEnergyDeposit();
52  edep *= CLHEP::MeV;
53 
54  //if there is no energy deposition skip everything
55  if(edep==0.)
56  {
57  if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() && aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition()) return false;
58  }
59 
60  //Get the Touchable History:
61  const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
62  if (not myTouch) {
63  G4cout << "BCMSensorSD::ProcessHits bad dynamic_cast" << G4endl;
64  return false;
65  }
66 
67  int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
68 
69  // Get the hit coordinates. Start and End Point
70  G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
71  G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
72 
73  // Calculate the local step begin and end position.
74  // From a G4 FAQ:
75  // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
76  const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
77  G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
78  G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
79 
80  HepGeom::Point3D<double> lP1,lP2;
81  lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
82  lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
83  lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
84 
85  lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
86  lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
87  lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
88 
89  //BCM hit stuff
90  if(BEcopyNo == 11950 || BEcopyNo == 11951)
91  {
92  TrackHelper trHelp(aStep->GetTrack());
93  //primary or not
94  int primaren = 0;
95  if(trHelp.IsPrimary())
96  primaren = 1;
97  else if(trHelp.IsRegeneratedPrimary())
98  primaren = 2;
99  else if(trHelp.IsSecondary())
100  primaren = 3;
101  else if(trHelp.IsRegisteredSecondary())
102  primaren = 4;
103 
104  int produced_in_diamond = 0;
105  if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::bcmDiamondLog")
106  produced_in_diamond = 1;
107  else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::bcmModLog")
108  produced_in_diamond = 2;
109  else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::bcmWallLog")
110  produced_in_diamond = 3;
111 
112  m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(), trHelp.GenerateParticleLink(),
113  0, 0, myTouch->GetVolume(1)->GetCopyNo()-951, BEcopyNo - 11950, primaren, produced_in_diamond);
114  }
115  return true;
116 }
AtlasG4EventUserInfo
This class is attached to G4Event objects as UserInformation. It holds a pointer to the HepMC::GenEve...
Definition: AtlasG4EventUserInfo.h:23
SiHit::xDep
@ xDep
Definition: SiHit.h:162
python.SystemOfUnits.mm
float mm
Definition: SystemOfUnits.py:98
TrackHelper::IsRegeneratedPrimary
bool IsRegeneratedPrimary() const
Definition: TrackHelper.cxx:20
TrackHelper.h
SiHit::xPhi
@ xPhi
Definition: SiHit.h:162
BCMSensorSD::Initialize
void Initialize(G4HCofThisEvent *) override final
Definition: BCMSensorSD.cxx:39
BCMSensorSD::BCMSensorSD
BCMSensorSD(const std::string &name, const std::string &hitCollectionName)
Definition: BCMSensorSD.cxx:31
AtlasHitsVector< SiHit >
BCMSensorSD::m_HitColl
SiHitCollection * m_HitColl
Pointer to the hit collection.
Definition: BCMSensorSD.h:43
BCMSensorSD.h
python.SystemOfUnits.MeV
float MeV
Definition: SystemOfUnits.py:172
AtlasHitsVector::Emplace
void Emplace(Args &&... args)
Definition: AtlasHitsVector.h:80
TrackHelper
Definition: TrackHelper.h:14
BCMExtra.h
TrackHelper::IsSecondary
bool IsSecondary() const
Definition: TrackHelper.cxx:30
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
SiHit::xEta
@ xEta
Definition: SiHit.h:162
TrackHelper::IsRegisteredSecondary
bool IsRegisteredSecondary() const
Definition: TrackHelper.cxx:25
AtlasG4EventUserInfo.h
TrackHelper::IsPrimary
bool IsPrimary() const
Definition: TrackHelper.cxx:15
TrackHelper::GenerateParticleLink
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition: TrackHelper.h:35
BCMSensorSD::m_HitCollName
std::string m_HitCollName
Name of the hit collection.
Definition: BCMSensorSD.h:42
BCMSensorSD::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
Definition: BCMSensorSD.cxx:49