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
14#include "MCTruth/TrackHelper.h"
15
16
17// Geant4 headers
18#include "G4Step.hh"
19#include "G4Track.hh"
20#include "G4VProcess.hh"
21#include "G4OpticalPhoton.hh"
22
23// CLHEP headers
24#include "CLHEP/Units/PhysicalConstants.h"
25
26//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
27
28LUCID_SensitiveDetector::LUCID_SensitiveDetector(const std::string& name, const std::string& hitCollectionName)
29 : G4VSensitiveDetector( name )
30 , m_HitColl( hitCollectionName )
31{
32 m_hit = new LUCID_HitHelper();
33}
34
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
36// Initialize from G4 - necessary to new the write handle for now
38{
39 if (!m_HitColl.isValid()) m_HitColl = std::make_unique<LUCID_SimHitCollection>(m_HitColl.name());
40}
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44bool LUCID_SensitiveDetector::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
45
46 if (verboseLevel>5)
47 {
48 G4cout << "LUCID_SensitiveDetector::ProcessHits - Begin" << G4endl;
49 }
50 G4Track* aTrack = aStep->GetTrack();
51
52 if (aTrack->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition()) return false;
53
54 if (verboseLevel>5)
55 {
56 G4cout << "LUCID_SensitiveDetector::ProcessHits(): There is an OpticalPhoton " << G4endl;
57 }
58
59 aTrack->SetTrackStatus(fKillTrackAndSecondaries);
60
61 if (aTrack->GetCreatorProcess()->GetProcessName() != "Cerenkov") return false;
62
63 if (verboseLevel>5)
64 {
65 G4cout << "LUCID_SensitiveDetector::ProcessHits(): It is from a Cerenkov process " << G4endl;
66 }
67
68 TrackHelper trHelp(aTrack);
69 double energy = aTrack->GetKineticEnergy()/CLHEP::eV;
70 double lambda = m_hit->GetWaveLength(energy);
71
72 m_HitColl->Emplace(m_hit->GetTubNumber(aStep),
73 aTrack->GetDefinition()->GetPDGEncoding(),
74 trHelp.GenerateParticleLink(),
75 LUCID_HitHelper::GetVolNumber (aTrack->GetLogicalVolumeAtVertex()->GetName()),
76 m_hit->GetPreStepPoint (aStep).x(),
77 m_hit->GetPreStepPoint (aStep).y(),
78 m_hit->GetPreStepPoint (aStep).z(),
79 m_hit->GetPostStepPoint(aStep).x(),
80 m_hit->GetPostStepPoint(aStep).y(),
81 m_hit->GetPostStepPoint(aStep).z(),
82 m_hit->GetPreStepTime (aStep),
83 m_hit->GetPostStepTime (aStep),
84 lambda,
85 energy);
86 return true;
87}
static int GetVolNumber(const G4String &)
LUCID_SensitiveDetector(const std::string &name, const std::string &hitCollectionName)
SG::WriteHandle< LUCID_SimHitCollection > m_HitColl
bool ProcessHits(G4Step *, G4TouchableHistory *) override final
void Initialize(G4HCofThisEvent *) override final
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52