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

#include <LUCID_SensitiveDetector.h>

Inheritance diagram for LUCID_SensitiveDetector:
Collaboration diagram for LUCID_SensitiveDetector:

Public Member Functions

 LUCID_SensitiveDetector (const std::string &name, const std::string &hitCollectionName)
 ~LUCID_SensitiveDetector ()
void Initialize (G4HCofThisEvent *) override final
bool ProcessHits (G4Step *, G4TouchableHistory *) 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 (LUCID_SensitiveDetectortest, ProcessHits)
 FRIEND_TEST (LUCID_SensitiveDetectortest, AddHit)
LUCID_SimHitCollectiongetHitCollection () const
 LUCID_SensitiveDetector (const LUCID_SensitiveDetector &)
LUCID_SensitiveDetectoroperator= (const LUCID_SensitiveDetector &)

Private Attributes

std::string m_hitCollectionName
LUCID_SimHitCollectionm_HitColl {}
LUCID_HitHelperm_hit

Detailed Description

Definition at line 22 of file LUCID_SensitiveDetector.h.

Constructor & Destructor Documentation

◆ LUCID_SensitiveDetector() [1/2]

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

Definition at line 30 of file LUCID_SensitiveDetector.cxx.

31 : G4VSensitiveDetector( name )
32 , m_hitCollectionName( hitCollectionName )
33{
34 m_hit = new LUCID_HitHelper();
35}

◆ ~LUCID_SensitiveDetector()

LUCID_SensitiveDetector::~LUCID_SensitiveDetector ( )
inline

Definition at line 31 of file LUCID_SensitiveDetector.h.

31{ /* I don't own myHitColl if all has gone well */ }

◆ LUCID_SensitiveDetector() [2/2]

LUCID_SensitiveDetector::LUCID_SensitiveDetector ( const LUCID_SensitiveDetector & )
private

Member Function Documentation

◆ AddHit()

template<class... Args>
void LUCID_SensitiveDetector::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 42 of file LUCID_SensitiveDetector.h.

43 {
44 if (m_HitColl) {
45 m_HitColl->Emplace(std::forward<Args>(args)...);
46 }
47 }
LUCID_SimHitCollection * m_HitColl

◆ FRIEND_TEST() [1/2]

LUCID_SensitiveDetector::FRIEND_TEST ( LUCID_SensitiveDetectortest ,
AddHit  )
private

◆ FRIEND_TEST() [2/2]

LUCID_SensitiveDetector::FRIEND_TEST ( LUCID_SensitiveDetectortest ,
ProcessHits  )
private

◆ getHitCollection()

LUCID_SimHitCollection * LUCID_SensitiveDetector::getHitCollection ( ) const
private

Definition at line 97 of file LUCID_SensitiveDetector.cxx.

98{
100 if (!eventInfo) {
101 return nullptr;
102 }
103 auto hitCollections = eventInfo->GetHitCollectionMap();
104 return hitCollections ? hitCollections->Find<LUCID_SimHitCollection>(m_hitCollectionName) : nullptr;
105}
AtlasHitsVector< LUCID_SimHit > LUCID_SimHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()

◆ Initialize()

void LUCID_SensitiveDetector::Initialize ( G4HCofThisEvent * )
finaloverride

Definition at line 39 of file LUCID_SensitiveDetector.cxx.

40{
42}
LUCID_SimHitCollection * getHitCollection() const

◆ operator=()

LUCID_SensitiveDetector & LUCID_SensitiveDetector::operator= ( const LUCID_SensitiveDetector & )
private

◆ ProcessHits()

bool LUCID_SensitiveDetector::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )
finaloverride

Definition at line 46 of file LUCID_SensitiveDetector.cxx.

46 {
47 if (!m_HitColl) {
49 if (!m_HitColl) {
50 return false;
51 }
52 }
53
54 if (verboseLevel>5)
55 {
56 G4cout << "LUCID_SensitiveDetector::ProcessHits - Begin" << G4endl;
57 }
58 G4Track* aTrack = aStep->GetTrack();
59
60 if (aTrack->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition()) return false;
61
62 if (verboseLevel>5)
63 {
64 G4cout << "LUCID_SensitiveDetector::ProcessHits(): There is an OpticalPhoton " << G4endl;
65 }
66
67 aTrack->SetTrackStatus(fKillTrackAndSecondaries);
68
69 if (aTrack->GetCreatorProcess()->GetProcessName() != "Cerenkov") return false;
70
71 if (verboseLevel>5)
72 {
73 G4cout << "LUCID_SensitiveDetector::ProcessHits(): It is from a Cerenkov process " << G4endl;
74 }
75
76 TrackHelper trHelp(aTrack);
77 double energy = aTrack->GetKineticEnergy()/CLHEP::eV;
78 double lambda = m_hit->GetWaveLength(energy);
79
80 m_HitColl->Emplace(m_hit->GetTubNumber(aStep),
81 aTrack->GetDefinition()->GetPDGEncoding(),
82 trHelp.GenerateParticleLink(),
83 LUCID_HitHelper::GetVolNumber (aTrack->GetLogicalVolumeAtVertex()->GetName()),
84 m_hit->GetPreStepPoint (aStep).x(),
85 m_hit->GetPreStepPoint (aStep).y(),
86 m_hit->GetPreStepPoint (aStep).z(),
87 m_hit->GetPostStepPoint(aStep).x(),
88 m_hit->GetPostStepPoint(aStep).y(),
89 m_hit->GetPostStepPoint(aStep).z(),
90 m_hit->GetPreStepTime (aStep),
91 m_hit->GetPostStepTime (aStep),
92 lambda,
93 energy);
94 return true;
95}
static int GetVolNumber(const G4String &)

Member Data Documentation

◆ m_hit

LUCID_HitHelper* LUCID_SensitiveDetector::m_hit
private

Definition at line 58 of file LUCID_SensitiveDetector.h.

◆ m_HitColl

LUCID_SimHitCollection* LUCID_SensitiveDetector::m_HitColl {}
private

Definition at line 57 of file LUCID_SensitiveDetector.h.

57{};

◆ m_hitCollectionName

std::string LUCID_SensitiveDetector::m_hitCollectionName
private

Definition at line 55 of file LUCID_SensitiveDetector.h.


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