CMS 3D CMS Logo

LowEnergyFastSimModel.cc
Go to the documentation of this file.
3 
7 
8 #include "G4VFastSimulationModel.hh"
9 #include "G4EventManager.hh"
10 #include "G4Electron.hh"
11 #include "GFlashHitMaker.hh"
12 #include "G4Region.hh"
13 #include "G4Material.hh"
14 #include "G4Positron.hh"
15 #include "G4ParticleDefinition.hh"
16 #include "G4PhysicalConstants.hh"
17 
18 constexpr G4double twomass = 2 * CLHEP::electron_mass_c2;
19 constexpr G4double scaleFactor = 1.025;
20 
22  : G4VFastSimulationModel(name, region),
23  fRegion(region),
24  fTrackingAction(nullptr),
25  fCheck(true),
26  fTailPos(0., 0., 0.) {
27  fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
28  fPositron = G4Positron::Positron();
29  fMaterial = nullptr;
30  auto table = G4Material::GetMaterialTable();
31  for (auto& mat : *table) {
32  G4String nam = mat->GetName();
33  size_t n = nam.size();
34  if (n > 4) {
35  G4String sn = nam.substr(n - 5, 5);
36  if (sn == "PbWO4") {
37  fMaterial = mat;
38  break;
39  }
40  }
41  }
42  G4String nm = (nullptr == fMaterial) ? "not found!" : fMaterial->GetName();
43  edm::LogVerbatim("LowEnergyFastSimModel") << "LowEGFlash material: <" << nm << ">";
44 }
45 
46 G4bool LowEnergyFastSimModel::IsApplicable(const G4ParticleDefinition& particle) {
47  return (11 == std::abs(particle.GetPDGEncoding()));
48 }
49 
50 G4bool LowEnergyFastSimModel::ModelTrigger(const G4FastTrack& fastTrack) {
51  const G4Track* track = fastTrack.GetPrimaryTrack();
52  G4double energy = track->GetKineticEnergy();
53  if (fMaterial != track->GetMaterial() || energy >= fEmax)
54  return false;
55 
56  /*
57  edm::LogVerbatim("LowEnergyFastSimModel") << track->GetDefinition()->GetParticleName()
58  << " Ekin(MeV)=" << energy << " material: <"
59  << track->GetMaterial()->GetName() << ">";
60  */
61  if (fCheck) {
62  if (nullptr == fTrackingAction) {
63  fTrackingAction = static_cast<const TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
64  }
65  const TrackInformation* ptrti = static_cast<TrackInformation*>(track->GetUserInformation());
66  int pdg = ptrti->genParticlePID();
67  if (std::abs(pdg) == 11 || pdg == 22)
68  return false;
69  }
70  return true;
71 }
72 
73 void LowEnergyFastSimModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
74  fastStep.KillPrimaryTrack();
75  fastStep.SetPrimaryTrackPathLength(0.0);
76  auto track = fastTrack.GetPrimaryTrack();
77  G4double energy = track->GetKineticEnergy() * scaleFactor;
78 
79  const G4ThreeVector& pos = track->GetPosition();
80 
81  G4double inPointEnergy = fParam.GetInPointEnergyFraction(energy) * energy;
82 
83  // take into account positron annihilation (not included in in-point)
84  if (fPositron == track->GetDefinition())
85  energy += twomass;
86 
87  const G4ThreeVector& momDir = track->GetMomentumDirection();
88 
89  // in point energy deposition
90  GFlashEnergySpot spot;
91  spot.SetEnergy(inPointEnergy);
92  spot.SetPosition(pos);
93  fHitMaker.make(&spot, &fastTrack);
94 
95  // tail energy deposition
96  const G4double etail = energy - inPointEnergy;
97  const G4int nspots = etail;
98  const G4double tailEnergy = etail / (nspots + 1);
99  /*
100  edm::LogVerbatim("LowEnergyFastSimModel") << track->GetDefinition()->GetParticleName()
101  << " Ekin(MeV)=" << energy << " material: <"
102  << track->GetMaterial()->GetName()
103  << "> Elocal=" << inPointEnergy
104  << " Etail=" << tailEnergy
105  << " Nspots=" << nspots+1;
106  */
107  for (G4int i = 0; i <= nspots; ++i) {
108  const G4double r = fParam.GetRadius(energy);
109  const G4double z = fParam.GetZ();
110 
111  const G4double phi = CLHEP::twopi * G4UniformRand();
112  fTailPos.set(r * std::cos(phi), r * std::sin(phi), z);
113  fTailPos.rotateUz(momDir);
114  fTailPos += pos;
115 
116  spot.SetEnergy(tailEnergy);
117  spot.SetPosition(fTailPos);
118  fHitMaker.make(&spot, &fastTrack);
119  }
120 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const G4ParticleDefinition * fPositron
G4bool IsApplicable(const G4ParticleDefinition &particle) override
constexpr G4double twomass
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
LowEnergyFastSimParam fParam
const G4Material * fMaterial
G4double GetInPointEnergyFraction(G4double energy) const
const TrackingAction * fTrackingAction
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
LowEnergyFastSimModel(const G4String &name, G4Region *region, const edm::ParameterSet &parSet)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr G4double scaleFactor
int genParticlePID() const
void DoIt(const G4FastTrack &fastTrack, G4FastStep &fastStep) override
G4bool ModelTrigger(const G4FastTrack &fastTrack) override
G4double GetRadius(G4double energy) const