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
14 #include "MCTruth/TrackHelper.h"
15 
16 // Geant4 headers
17 #include "G4Step.hh"
18 #include "G4ThreeVector.hh"
19 #include "G4Geantino.hh"
20 #include "G4ChargedGeantino.hh"
21 
22 // CLHEP headers
23 #include "CLHEP/Geometry/Transform3D.h"
24 #include "CLHEP/Units/SystemOfUnits.h"
25 
26 #include <memory> // For make unique
27 
28 BLMSensorSD::BLMSensorSD(const std::string& name, const std::string& hitCollectionName)
29  : G4VSensitiveDetector( name )
30  , m_HitColl( hitCollectionName )
31 {
32 }
33 
34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35 // Initialize from G4 - necessary to new the write handle for now
36 void BLMSensorSD::Initialize(G4HCofThisEvent *)
37 {
38  if (!m_HitColl.isValid()) m_HitColl = std::make_unique<SiHitCollection>();
39 }
40 
41 G4bool BLMSensorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
42 {
43  G4double edep = aStep->GetTotalEnergyDeposit();
44  edep *= CLHEP::MeV;
45 
46  //if there is no energy deposition skip everything
47  if(edep==0.)
48  {
49  if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() && aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition()) return false;
50  }
51 
52  //Get the Touchable History:
53  const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
54  if (not myTouch) {
55  G4cout << "BLMSensorSD::ProcessHits bad dynamic_cast" << G4endl;
56  return false;
57  }
58  int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
59 
60  // Get the hit coordinates. Start and End Point
61  G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
62  G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
63 
64  // Calculate the local step begin and end position.
65  // From a G4 FAQ:
66  // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
67  const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
68  G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
69  G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
70 
71  HepGeom::Point3D<double> lP1,lP2;
72  lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
73  lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
74  lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
75 
76  lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
77  lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
78  lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
79 
80  //BLM hit stuff
81  if(BEcopyNo == 2009)
82  {
83  TrackHelper trHelp(aStep->GetTrack());
84  //primary or not
85  int primaren = 0;
86  if(trHelp.IsPrimary())
87  primaren = 1;
88  else if(trHelp.IsRegeneratedPrimary())
89  primaren = 2;
90  else if(trHelp.IsSecondary())
91  primaren = 3;
92  else if(trHelp.IsRegisteredSecondary())
93  primaren = 4;
94 
95  int produced_in_diamond = 0;
96  if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::blmDiamondLog")
97  produced_in_diamond = 1;
98  else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::blmModLog")
99  produced_in_diamond = 2;
100  else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::blmWallLog")
101  produced_in_diamond = 3;
102 
103  m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(), trHelp.GenerateParticleLink(),
104  0, 0, myTouch->GetVolume(1)->GetCopyNo()-222, 0, primaren, produced_in_diamond);
105  }
106  return true;
107 }
SiHit::xDep
@ xDep
Definition: SiHit.h:162
TrackHelper::IsRegeneratedPrimary
bool IsRegeneratedPrimary() const
Definition: TrackHelper.cxx:20
TrackHelper.h
SiHit::xPhi
@ xPhi
Definition: SiHit.h:162
python.SystemOfUnits.MeV
int MeV
Definition: SystemOfUnits.py:154
BLMSensorSD::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
Definition: BLMSensorSD.cxx:41
TrackHelper
Definition: TrackHelper.h:14
BLMSensorSD::BLMSensorSD
BLMSensorSD(const std::string &name, const std::string &hitCollectionName)
Definition: BLMSensorSD.cxx:28
TrackHelper::IsSecondary
bool IsSecondary() const
Definition: TrackHelper.cxx:30
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
BLMSensorSD.h
SiHit::xEta
@ xEta
Definition: SiHit.h:162
python.SystemOfUnits.mm
int mm
Definition: SystemOfUnits.py:83
TrackHelper::IsRegisteredSecondary
bool IsRegisteredSecondary() const
Definition: TrackHelper.cxx:25
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::m_HitColl
SG::WriteHandle< SiHitCollection > m_HitColl
Definition: BLMSensorSD.h:49
BLMSensorSD::Initialize
void Initialize(G4HCofThisEvent *) override final
Definition: BLMSensorSD.cxx:36