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