ATLAS Offline Software
NeutronFastSim.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
7 
8 #include "G4Event.hh"
9 #include "G4Neutron.hh"
10 #include "G4SDManager.hh"
11 #include "G4VSensitiveDetector.hh"
12 #include "G4EventManager.hh"
14 
16 
17 NeutronFastSim::NeutronFastSim(const std::string& name, const std::string& fsSDname, const double etaCut, const double timeCut)
18  : G4VFastSimulationModel(name)
19  , m_Energy(5)
20  , m_fsSD(0)
21  , m_init(false)
22  , m_fsSDname(fsSDname)
23  , m_etaCut(etaCut)
24  , m_timeCut(timeCut)
25 {
26 }
27 
28 G4bool NeutronFastSim::IsApplicable(const G4ParticleDefinition&)
29 {
30  if (!m_init){
31  m_init = true;
32 
33  G4SDManager *sdm = G4SDManager::GetSDMpointer();
34  G4VSensitiveDetector * vsd = sdm->FindSensitiveDetector( m_fsSDname );
35  if (!vsd) {
36  G4cout << "NeutronFastSim::IsApplicable WARNING Could not get TrackFastSimSD sensitive detector. If you are not writing track records this is expected." << G4endl;
37  } else {
38  m_fsSD = dynamic_cast<TrackFastSimSD*>(vsd);
39  if (!m_fsSD) {
40  G4cout << "NeutronFastSim::IsApplicable WARNING Could not cast the SD. If you are not writing track records this is expected." << G4endl;
41  }
42  } // found the SD
43  } // End of lazy init
44  return true;
45 }
46 
47 G4bool NeutronFastSim::ModelTrigger(const G4FastTrack& fastTrack)
48 {
49  // Trigger if the neutron energy is below our threshold or if the time is over 150 ns
50  if (fastTrack.GetPrimaryTrack()->GetDefinition() == G4Neutron::NeutronDefinition() ){
51  return (m_Energy<0?true:fastTrack.GetPrimaryTrack()->GetKineticEnergy()<m_Energy) || fastTrack.GetPrimaryTrack()->GetGlobalTime()>m_timeCut;
52  }
53 
54  // Not a neutron... Pick it up if the primary had eta>6.0
55  AtlasG4EventUserInfo *atlasG4EvtUserInfo=static_cast<AtlasG4EventUserInfo*>(G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetUserInformation());
56  HepMC::ConstGenParticlePtr primaryGenParticle = atlasG4EvtUserInfo->GetCurrentPrimaryGenParticle();
57  if (std::abs(primaryGenParticle->momentum().eta())>m_etaCut && !HepMC::is_simulation_particle(primaryGenParticle)){
58  return true;
59  } else {
60  return false;
61  }
62 }
63 
64 void NeutronFastSim::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep)
65 {
66  if (m_fsSD) m_fsSD->WriteTrack( fastTrack.GetPrimaryTrack() , false , false );
67  fastStep.KillPrimaryTrack();
68 }
AtlasG4EventUserInfo
This class is attached to G4Event objects as UserInformation. It holds a pointer to the HepMC::GenEve...
Definition: AtlasG4EventUserInfo.h:21
TrackFastSimSD
Definition: TrackFastSimSD.h:24
TrackFastSimSD.h
NeutronFastSim::ModelTrigger
virtual G4bool ModelTrigger(const G4FastTrack &) override final
Definition: NeutronFastSim.cxx:47
NeutronFastSim::m_fsSDname
std::string m_fsSDname
Definition: NeutronFastSim.h:33
NeutronFastSim::m_Energy
G4double m_Energy
Definition: NeutronFastSim.h:30
NeutronFastSim::DoIt
void DoIt(const G4FastTrack &, G4FastStep &) override final
Definition: NeutronFastSim.cxx:64
NeutronFastSim.h
HepMC::is_simulation_particle
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
Definition: MagicNumbers.h:342
NeutronFastSim::m_timeCut
double m_timeCut
Definition: NeutronFastSim.h:35
HepMC::ConstGenParticlePtr
const GenParticle * ConstGenParticlePtr
Definition: GenParticle.h:38
NeutronFastSim::NeutronFastSim
NeutronFastSim(const std::string &name, const std::string &fsSDname, const double etaCut, const double timeCut)
Definition: NeutronFastSim.cxx:17
NeutronFastSim::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &) override final
Definition: NeutronFastSim.cxx:28
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
MagicNumbers.h
AtlasG4EventUserInfo::GetCurrentPrimaryGenParticle
HepMC::ConstGenParticlePtr GetCurrentPrimaryGenParticle() const
return a pointer to the HepMC::GenParticle used to create the current G4PrimaryParticle.
Definition: AtlasG4EventUserInfo.h:44
NeutronFastSim::m_etaCut
double m_etaCut
Definition: NeutronFastSim.h:34
AtlasG4EventUserInfo.h
TrackFastSimSD::WriteTrack
void WriteTrack(const G4Track *, const bool, const bool)
Definition: TrackFastSimSD.cxx:97
NeutronFastSim::m_fsSD
TrackFastSimSD * m_fsSD
Definition: NeutronFastSim.h:31
NeutronFastSim::m_init
bool m_init
Definition: NeutronFastSim.h:32