ATLAS Offline Software
Loading...
Searching...
No Matches
LUCID_SensitiveDetector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6// Class header
8
9// Package headers
10#include "LUCID_HitHelper.h"
11
12// Athena headers
16#include "MCTruth/TrackHelper.h"
17
18
19// Geant4 headers
20#include "G4Step.hh"
21#include "G4Track.hh"
22#include "G4VProcess.hh"
23#include "G4OpticalPhoton.hh"
24
25// CLHEP headers
26#include "CLHEP/Units/PhysicalConstants.h"
27
28//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
29
30LUCID_SensitiveDetector::LUCID_SensitiveDetector(const std::string& name, const std::string& hitCollectionName)
31 : G4VSensitiveDetector( name )
32 , m_hitCollectionName( hitCollectionName )
33{
34 m_hit = new LUCID_HitHelper();
35}
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
38// Initialize from G4.
40{
42}
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
45
46bool LUCID_SensitiveDetector::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
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}
96
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()
static int GetVolNumber(const G4String &)
LUCID_SensitiveDetector(const std::string &name, const std::string &hitCollectionName)
LUCID_SimHitCollection * m_HitColl
bool ProcessHits(G4Step *, G4TouchableHistory *) override final
LUCID_SimHitCollection * getHitCollection() const
void Initialize(G4HCofThisEvent *) override final
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52