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