CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
LowEnergyFastSimModel.cc
Go to the documentation of this file.
3 
6 
7 #include "G4VFastSimulationModel.hh"
8 #include "G4EventManager.hh"
9 #include "G4Electron.hh"
10 #include "GFlashHitMaker.hh"
11 #include "G4Region.hh"
12 #include "G4PhysicalConstants.hh"
13 
14 constexpr G4double twomass = 2 * CLHEP::electron_mass_c2;
15 
17  : G4VFastSimulationModel(name, region),
18  fRegion(region),
19  fTrackingAction(nullptr),
20  fCheck(false),
21  fTailPos(0., 0., 0.) {
22  fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
23 }
24 
25 G4bool LowEnergyFastSimModel::IsApplicable(const G4ParticleDefinition& particle) {
26  return (11 == std::abs(particle.GetPDGEncoding()));
27 }
28 
29 G4bool LowEnergyFastSimModel::ModelTrigger(const G4FastTrack& fastTrack) {
30  const G4Track* track = fastTrack.GetPrimaryTrack();
31  if (fCheck) {
32  if (nullptr == fTrackingAction) {
33  fTrackingAction = static_cast<const TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
34  }
35  int pdgMother = std::abs(fTrackingAction->geant4Track()->GetDefinition()->GetPDGEncoding());
36  if (pdgMother == 11 || pdgMother == 22)
37  return false;
38  }
39  G4double energy = track->GetKineticEnergy();
40  return (energy < fEmax && fRegion == fastTrack.GetEnvelope());
41 }
42 
43 void LowEnergyFastSimModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
44  fastStep.KillPrimaryTrack();
45  fastStep.SetPrimaryTrackPathLength(0.0);
46  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
47 
48  const G4ThreeVector& pos = fastTrack.GetPrimaryTrack()->GetPosition();
49 
50  G4double inPointEnergy = fParam.GetInPointEnergyFraction(energy) * energy;
51 
52  // take into account positron annihilation (not included in in-point)
53  if (-11 == fastTrack.GetPrimaryTrack()->GetDefinition()->GetPDGEncoding())
54  energy += twomass;
55 
56  const G4ThreeVector& momDir = fastTrack.GetPrimaryTrack()->GetMomentumDirection();
57 
58  // in point energy deposition
59  GFlashEnergySpot spot;
60  spot.SetEnergy(inPointEnergy);
61  spot.SetPosition(pos);
62  fHitMaker.make(&spot, &fastTrack);
63 
64  // tail energy deposition
65  G4double etail = energy - inPointEnergy;
66  const G4int nspots = (G4int)(etail) + 1;
67  const G4double tailEnergy = etail / (G4double)nspots;
68  for (G4int i = 0; i < nspots; ++i) {
69  const G4double r = fParam.GetRadius(energy);
70  const G4double z = fParam.GetZ();
71 
72  const G4double phi = CLHEP::twopi * G4UniformRand();
73  fTailPos.set(r * std::cos(phi), r * std::sin(phi), z);
74  fTailPos.rotateUz(momDir);
75  fTailPos += pos;
76 
77  spot.SetEnergy(tailEnergy);
78  spot.SetPosition(fTailPos);
79  fHitMaker.make(&spot, &fastTrack);
80  }
81 }
const double GeV
Definition: MathUtil.h:16
G4double GetRadius(G4double energy) const
G4bool IsApplicable(const G4ParticleDefinition &particle) override
constexpr G4double twomass
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const G4Envelope * fRegion
LowEnergyFastSimParam fParam
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
G4double GetInPointEnergyFraction(G4double energy) const
void DoIt(const G4FastTrack &fastTrack, G4FastStep &fastStep) override
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
G4bool ModelTrigger(const G4FastTrack &fastTrack) override
const G4Track * geant4Track() const