ATLAS Offline Software
Loading...
Searching...
No Matches
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
31BCMSensorSD::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
39void 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 m_g4UserEventInfo = eventInfo;
46 }
47 }
48}
49
50G4bool BCMSensorSD::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
51{
52 G4double edep = aStep->GetTotalEnergyDeposit();
53 edep *= CLHEP::MeV;
54
55 //if there is no energy deposition skip everything
56 if(edep==0.)
57 {
58 if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() && aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition()) return false;
59 }
60
61 //Get the Touchable History:
62 const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
63 if (not myTouch) {
64 G4cout << "BCMSensorSD::ProcessHits bad dynamic_cast" << G4endl;
65 return false;
66 }
67
68 int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
69
70 // Get the hit coordinates. Start and End Point
71 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
72 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
73
74 // Calculate the local step begin and end position.
75 // From a G4 FAQ:
76 // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
77 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
78 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
79 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
80
81 HepGeom::Point3D<double> lP1,lP2;
82 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
83 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
84 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
85
86 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
87 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
88 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
89
90 //BCM hit stuff
91 if(BEcopyNo == 11950 || BEcopyNo == 11951)
92 {
93 TrackHelper trHelp(aStep->GetTrack());
94 //primary or not
95 int primaren = 0;
96 if(trHelp.IsPrimary())
97 primaren = 1;
98 else if(trHelp.IsRegeneratedPrimary())
99 primaren = 2;
100 else if(trHelp.IsSecondary())
101 primaren = 3;
102 else if(trHelp.IsRegisteredSecondary())
103 primaren = 4;
104
105 int produced_in_diamond = 0;
106 if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::bcmDiamondLog")
107 produced_in_diamond = 1;
108 else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::bcmModLog")
109 produced_in_diamond = 2;
110 else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::bcmWallLog")
111 produced_in_diamond = 3;
112 // Temporary solution to get EventContext from a Geant4 thread
113 m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(),
114 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
115 0, 0, myTouch->GetVolume(1)->GetCopyNo()-951, BEcopyNo - 11950, primaren, produced_in_diamond);
116 }
117 return true;
118}
AtlasHitsVector< SiHit > SiHitCollection
This class is attached to G4Event objects as UserInformation.
AtlasG4EventUserInfo * m_g4UserEventInfo
Definition BCMSensorSD.h:45
std::string m_HitCollName
Name of the hit collection.
Definition BCMSensorSD.h:43
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
SiHitCollection * m_HitColl
Pointer to the hit collection.
Definition BCMSensorSD.h:44
BCMSensorSD(const std::string &name, const std::string &hitCollectionName)
void Initialize(G4HCofThisEvent *) override final
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162
bool IsPrimary() const
bool IsSecondary() const
bool IsRegisteredSecondary() const
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
bool IsRegeneratedPrimary() const