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

#include <BLMSensorSD.h>

Inheritance diagram for BLMSensorSD:
Collaboration diagram for BLMSensorSD:

Public Member Functions

 BLMSensorSD (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 (BLMSensorSDtest, ProcessHits)
 FRIEND_TEST (BLMSensorSDtest, AddHit)

Private Attributes

std::string m_HitCollName
SiHitCollectionm_HitColl {nullptr}
AtlasG4EventUserInfom_g4UserEventInfo {nullptr}

Detailed Description

Definition at line 24 of file BLMSensorSD.h.

Constructor & Destructor Documentation

◆ BLMSensorSD()

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

Definition at line 29 of file BLMSensorSD.cxx.

30 : G4VSensitiveDetector( name )
31 , m_HitCollName( hitCollectionName )
32{
33}
std::string m_HitCollName
Definition BLMSensorSD.h:44

Member Function Documentation

◆ AddHit()

template<class... Args>
void BLMSensorSD::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 BLMSensorSD.h.

41{ m_HitColl->Emplace( args... ); }
SiHitCollection * m_HitColl
Definition BLMSensorSD.h:45

◆ FRIEND_TEST() [1/2]

BLMSensorSD::FRIEND_TEST ( BLMSensorSDtest ,
AddHit  )
private

◆ FRIEND_TEST() [2/2]

BLMSensorSD::FRIEND_TEST ( BLMSensorSDtest ,
ProcessHits  )
private

◆ Initialize()

void BLMSensorSD::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 37 of file BLMSensorSD.cxx.

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 m_g4UserEventInfo = eventInfo;
44 }
45 }
46}
AtlasHitsVector< SiHit > SiHitCollection
AtlasG4EventUserInfo * m_g4UserEventInfo
Definition BLMSensorSD.h:46

◆ ProcessHits()

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

Definition at line 48 of file BLMSensorSD.cxx.

49{
50 G4double edep = aStep->GetTotalEnergyDeposit();
51 edep *= CLHEP::MeV;
52
53 //if there is no energy deposition skip everything
54 if(edep==0.)
55 {
56 if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() && aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition()) return false;
57 }
58
59 //Get the Touchable History:
60 const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
61 if (not myTouch) {
62 G4cout << "BLMSensorSD::ProcessHits bad dynamic_cast" << G4endl;
63 return false;
64 }
65 int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
66
67 // Get the hit coordinates. Start and End Point
68 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
69 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
70
71 // Calculate the local step begin and end position.
72 // From a G4 FAQ:
73 // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
74 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
75 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
76 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
77
78 HepGeom::Point3D<double> lP1,lP2;
79 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
80 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
81 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
82
83 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
84 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
85 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
86
87 //BLM hit stuff
88 if(BEcopyNo == 2009)
89 {
90 TrackHelper trHelp(aStep->GetTrack());
91 //primary or not
92 int primaren = 0;
93 if(trHelp.IsPrimary())
94 primaren = 1;
95 else if(trHelp.IsRegeneratedPrimary())
96 primaren = 2;
97 else if(trHelp.IsSecondary())
98 primaren = 3;
99 else if(trHelp.IsRegisteredSecondary())
100 primaren = 4;
101
102 int produced_in_diamond = 0;
103 if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::blmDiamondLog")
104 produced_in_diamond = 1;
105 else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::blmModLog")
106 produced_in_diamond = 2;
107 else if(aStep->GetTrack()->GetLogicalVolumeAtVertex()->GetName() == "Pixel::blmWallLog")
108 produced_in_diamond = 3;
109
110 m_HitColl->Emplace(lP1, lP2, edep, aStep->GetPreStepPoint()->GetGlobalTime(),
111 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
112 0, 0, myTouch->GetVolume(1)->GetCopyNo()-222, 0, primaren, produced_in_diamond);
113 }
114 return true;
115}
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162

Member Data Documentation

◆ m_g4UserEventInfo

AtlasG4EventUserInfo* BLMSensorSD::m_g4UserEventInfo {nullptr}
private

Definition at line 46 of file BLMSensorSD.h.

46{nullptr};

◆ m_HitColl

SiHitCollection* BLMSensorSD::m_HitColl {nullptr}
private

Definition at line 45 of file BLMSensorSD.h.

45{nullptr};

◆ m_HitCollName

std::string BLMSensorSD::m_HitCollName
private

Definition at line 44 of file BLMSensorSD.h.


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