ATLAS Offline Software
Loading...
Searching...
No Matches
BCMSensorSD Class Reference

#include <BCMSensorSD.h>

Inheritance diagram for BCMSensorSD:
Collaboration diagram for BCMSensorSD:

Public Member Functions

 BCMSensorSD (const std::string &name, const std::string &hitCollectionName)
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
void Initialize (G4HCofThisEvent *) override final
template<class... Args>
void AddHit (Args &&... args)
 Templated method to stuff a single hit into the sensitive detector class.

Private Member Functions

 FRIEND_TEST (BCMSensorSDtest, ProcessHits)
 FRIEND_TEST (BCMSensorSDtest, AddHit)

Private Attributes

std::string m_HitCollName
 Name of the hit collection.
SiHitCollectionm_HitColl {nullptr}
 Pointer to the hit collection.
AtlasG4EventUserInfom_g4UserEventInfo {nullptr}

Detailed Description

Definition at line 24 of file BCMSensorSD.h.

Constructor & Destructor Documentation

◆ BCMSensorSD()

BCMSensorSD::BCMSensorSD ( const std::string & name,
const std::string & hitCollectionName )

Definition at line 31 of file BCMSensorSD.cxx.

32 : G4VSensitiveDetector( name )
33 , m_HitCollName( hitCollectionName )
34{
35}
std::string m_HitCollName
Name of the hit collection.
Definition BCMSensorSD.h:43

Member Function Documentation

◆ AddHit()

template<class... Args>
void BCMSensorSD::AddHit ( Args &&... args)
inline

Templated method to stuff a single hit into the sensitive detector class.

This could get rather tricky, but the idea is to allow fast simulations to use the very same SD classes as the standard simulation.

Definition at line 41 of file BCMSensorSD.h.

41{ m_HitColl->Emplace( args... ); }
SiHitCollection * m_HitColl
Pointer to the hit collection.
Definition BCMSensorSD.h:44

◆ FRIEND_TEST() [1/2]

BCMSensorSD::FRIEND_TEST ( BCMSensorSDtest ,
AddHit  )
private

◆ FRIEND_TEST() [2/2]

BCMSensorSD::FRIEND_TEST ( BCMSensorSDtest ,
ProcessHits  )
private

◆ Initialize()

void BCMSensorSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 39 of file BCMSensorSD.cxx.

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}
AtlasHitsVector< SiHit > SiHitCollection
AtlasG4EventUserInfo * m_g4UserEventInfo
Definition BCMSensorSD.h:45

◆ ProcessHits()

G4bool BCMSensorSD::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
finaloverride

Definition at line 50 of file BCMSensorSD.cxx.

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}
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162

Member Data Documentation

◆ m_g4UserEventInfo

AtlasG4EventUserInfo* BCMSensorSD::m_g4UserEventInfo {nullptr}
private

Definition at line 45 of file BCMSensorSD.h.

45{nullptr};

◆ m_HitColl

SiHitCollection* BCMSensorSD::m_HitColl {nullptr}
private

Pointer to the hit collection.

Definition at line 44 of file BCMSensorSD.h.

44{nullptr};

◆ m_HitCollName

std::string BCMSensorSD::m_HitCollName
private

Name of the hit collection.

Definition at line 43 of file BCMSensorSD.h.


The documentation for this class was generated from the following files: