CMS 3D CMS Logo

LowEnergyFastSimModel.cc
Go to the documentation of this file.
3 
5 
6 #include "G4VFastSimulationModel.hh"
7 #include "G4Electron.hh"
8 #include "GFlashHitMaker.hh"
9 #include "G4Region.hh"
10 
12  : G4VFastSimulationModel(name, region),
13  fEmax(parSet.getParameter<double>("LowEnergyGflashEcalEmax")),
14  fRegion(region) {}
15 
16 G4bool LowEnergyFastSimModel::IsApplicable(const G4ParticleDefinition& particle) {
17  return &particle == G4Electron::Definition();
18 }
19 
20 G4bool LowEnergyFastSimModel::ModelTrigger(const G4FastTrack& fastTrack) {
21  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
22  return energy < fEmax && fRegion == fastTrack.GetEnvelope();
23 }
24 
25 void LowEnergyFastSimModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
26  fastStep.KillPrimaryTrack();
27  fastStep.SetPrimaryTrackPathLength(0.0);
28  fastStep.SetTotalEnergyDeposited(fastTrack.GetPrimaryTrack()->GetKineticEnergy());
29 
30  const G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
31  const G4ThreeVector& pos = fastTrack.GetPrimaryTrack()->GetPosition();
32 
33  G4double inPointEnergy = param.GetInPointEnergyFraction(energy) * energy;
34 
35  const G4ThreeVector& momDir = fastTrack.GetPrimaryTrack()->GetMomentumDirection();
36  const G4ThreeVector& ortho = momDir.orthogonal();
37  const G4ThreeVector& cross = momDir.cross(ortho);
38 
39  // in point energy deposition
40  GFlashEnergySpot spot;
41  spot.SetEnergy(inPointEnergy);
42  spot.SetPosition(pos);
43  fHitMaker.make(&spot, &fastTrack);
44 
45  // tail energy deposition
46  G4double etail = energy - inPointEnergy;
47  const G4int nspots = int(etail) + 1;
48  for (G4int i = 0; i < nspots; ++i) {
49  const G4double radius = param.GetRadius(energy);
50  const G4double z = param.GetZ();
51 
52  const G4double phi = CLHEP::twopi * G4UniformRand();
53  const G4ThreeVector tailPos = pos + z * momDir + radius * std::cos(phi) * ortho + radius * std::sin(phi) * cross;
54 
55  const G4double tailEnergy = etail / nspots;
56 
57  spot.SetEnergy(tailEnergy);
58  spot.SetPosition(tailPos);
59  fHitMaker.make(&spot, &fastTrack);
60  }
61 }
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
pos
Definition: PixelAliasList.h:18
LowEnergyFastSimModel.h
LowEnergyFastSimParam::GetInPointEnergyFraction
G4double GetInPointEnergyFraction(G4double energy) const
Definition: LowEnergyFastSimParam.h:9
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Basic3DVector::cross
Basic3DVector cross(const Basic3DVector &lh) const
Vector product, or "cross" product, with a vector of same type.
Definition: extBasic3DVector.h:203
DDAxes::z
LowEnergyFastSimParam::GetRadius
G4double GetRadius(G4double energy) const
Definition: LowEnergyFastSimParam.h:15
LowEnergyFastSimModel::fEmax
const G4double fEmax
Definition: LowEnergyFastSimModel.h:22
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
cross
Basic3DVector cross(const Basic3DVector &v) const
Vector product, or "cross" product, with a vector of same type.
Definition: Basic3DVectorLD.h:225
LowEnergyFastSimModel::DoIt
void DoIt(const G4FastTrack &fastTrack, G4FastStep &fastStep) override
Definition: LowEnergyFastSimModel.cc:25
edm::ParameterSet
Definition: ParameterSet.h:36
LowEnergyFastSimModel::param
LowEnergyFastSimParam param
Definition: LowEnergyFastSimModel.h:25
createfilelist.int
int
Definition: createfilelist.py:10
LowEnergyFastSimParam::GetZ
G4double GetZ() const
Definition: LowEnergyFastSimParam.h:24
DDAxes::phi
LowEnergyFastSimModel::fHitMaker
GFlashHitMaker fHitMaker
Definition: LowEnergyFastSimModel.h:24
LowEnergyFastSimModel::fRegion
const G4Envelope *const fRegion
Definition: LowEnergyFastSimModel.h:23
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
HLT_2018_cff.region
region
Definition: HLT_2018_cff.py:81479
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
LowEnergyFastSimModel::ModelTrigger
G4bool ModelTrigger(const G4FastTrack &fastTrack) override
Definition: LowEnergyFastSimModel.cc:20
LowEnergyFastSimModel::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &particle) override
Definition: LowEnergyFastSimModel.cc:16
ParameterSet.h
LowEnergyFastSimModel::LowEnergyFastSimModel
LowEnergyFastSimModel(const G4String &name, G4Region *region, const edm::ParameterSet &parSet)
Definition: LowEnergyFastSimModel.cc:11