ATLAS Offline Software
StoppedParticleFastSim.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3 */
4 
7 
8 #include "G4FastTrack.hh"
9 #include "G4Track.hh"
10 #include "G4DynamicParticle.hh"
11 #include "G4ParticleDefinition.hh"
12 #include "G4Material.hh"
13 #include "G4ElementVector.hh"
14 #include "G4SDManager.hh"
15 #include "G4VSensitiveDetector.hh"
16 
17 #include "CLHEP/Units/PhysicalConstants.h"
18 
19 #include <cmath>
20 
21 StoppedParticleFastSim::StoppedParticleFastSim(const std::string& name, const std::string& fsSDname)
22  : G4VFastSimulationModel(name)
23  , m_fsSD(nullptr)
24  , m_init(false)
25  , m_fsSDname(fsSDname)
26 {
27 }
28 
29 G4bool StoppedParticleFastSim::IsApplicable(const G4ParticleDefinition&)
30 {
31  return true;
32 }
33 
34 G4bool StoppedParticleFastSim::ModelTrigger(const G4FastTrack& fastTrack)
35 {
36  // Trigger if the energy is below our threshold or if the time is over 150 ns
37  int id = fastTrack.GetPrimaryTrack()->GetDynamicParticle()->GetDefinition()->GetPDGEncoding();
38  if (id<1000000 || id>1100000) return true;
39  else if (isSUSYParticle(id)){
40  G4Material * mat = fastTrack.GetPrimaryTrack()->GetMaterial();
41  double minA=1500000.;
42  for (unsigned int i=0;i<mat->GetNumberOfElements();++i){
43  if (mat->GetElement(i) &&
44  minA>mat->GetElement(i)->GetN()){
45  minA=mat->GetElement(i)->GetN();
46  }
47  }
48  if (fastTrack.GetPrimaryTrack()->GetVelocity()<0.15*std::pow(minA,-2./3.)*CLHEP::c_light) return true;
49  return false;
50  }
51  return true;
52 }
53 
54 void StoppedParticleFastSim::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep)
55 {
56  if (!m_init){
57  m_init = true;
58 
59  G4SDManager *sdm = G4SDManager::GetSDMpointer();
60  G4VSensitiveDetector * vsd = sdm->FindSensitiveDetector( m_fsSDname );
61  if (!vsd) {
62  G4cout << "StoppedParticleFastSim::DoIt WARNING Could not get TrackFastSimSD sensitive detector. If you are not writing track records this is expected." << G4endl;
63  m_fsSD = dynamic_cast<TrackFastSimSD*>(vsd);
64  if (!m_fsSD) {
65  G4cout << "StoppedParticleFastSim::DoIt WARNING Could not cast the SD. If you are not writing track records this is expected." << G4endl;
66  }
67  } // found the SD
68  } // End of lazy init
69 
70  if (isSUSYParticle(fastTrack.GetPrimaryTrack()->GetDynamicParticle()->GetDefinition()->GetPDGEncoding()) &&
71  m_fsSD) {
72  m_fsSD->WriteTrack( fastTrack.GetPrimaryTrack() , false , true );
73  }
74  fastStep.KillPrimaryTrack();
75 }
76 
78 {
79  if (id==1000021 || id==1000005 || id==1000006 || id==1000512 || id==1000522 || id==1000991 || id==1000993 ||
80  id==1000612 || id==1000622 || id==1000632 || id==1000642 || id==1000652 || id==1005211 ||
81  id==1006113 || id==1006211 || id==1006213 || id==1006223 || id==1006311 ||
82  id==1006313 || id==1006321 || id==1006323 || id==1006333 ||
83  id==1009111 || id==1009113 || id==1009211 || id==1009213 || id==1009311 ||
84  id==1009313 || id==1009321 || id==1009323 || id==1009223 || id==1009333 ||
85  id==1092112 || id==1091114 || id==1092114 || id==1092212 || id==1092214 || id==1092224 ||
86  id==1093114 || id==1093122 || id==1093214 || id==1093224 || id==1093314 || id==1093324 || id==1093334)
87  return true;
88  return false;
89 }
StoppedParticleFastSim::m_fsSD
TrackFastSimSD * m_fsSD
Definition: StoppedParticleFastSim.h:31
TrackFastSimSD
Definition: TrackFastSimSD.h:24
TrackFastSimSD.h
mat
GeoMaterial * mat
Definition: LArDetectorConstructionTBEC.cxx:53
conifer::pow
constexpr int pow(int x)
Definition: conifer.h:20
StoppedParticleFastSim.h
StoppedParticleFastSim::m_fsSDname
std::string m_fsSDname
Definition: StoppedParticleFastSim.h:33
StoppedParticleFastSim::StoppedParticleFastSim
StoppedParticleFastSim(const std::string &name, const std::string &fsSDname)
Definition: StoppedParticleFastSim.cxx:21
StoppedParticleFastSim::m_init
bool m_init
Definition: StoppedParticleFastSim.h:32
lumiFormat.i
int i
Definition: lumiFormat.py:92
StoppedParticleFastSim::isSUSYParticle
bool isSUSYParticle(const int) const
Definition: StoppedParticleFastSim.cxx:77
StoppedParticleFastSim::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &) override final
Definition: StoppedParticleFastSim.cxx:29
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
python.PhysicalConstants.c_light
float c_light
Definition: PhysicalConstants.py:63
StoppedParticleFastSim::ModelTrigger
virtual G4bool ModelTrigger(const G4FastTrack &) override final
Definition: StoppedParticleFastSim.cxx:34
TrackFastSimSD::WriteTrack
void WriteTrack(const G4Track *, const bool, const bool)
Definition: TrackFastSimSD.cxx:97
StoppedParticleFastSim::DoIt
void DoIt(const G4FastTrack &, G4FastStep &) override final
Definition: StoppedParticleFastSim.cxx:54