CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
LowEnergyFastSimModel Class Reference

#include <LowEnergyFastSimModel.h>

Inheritance diagram for LowEnergyFastSimModel:

Public Member Functions

void DoIt (const G4FastTrack &fastTrack, G4FastStep &fastStep) override
 
G4bool IsApplicable (const G4ParticleDefinition &particle) override
 
 LowEnergyFastSimModel (const G4String &name, G4Region *region, const edm::ParameterSet &parSet)
 
G4bool ModelTrigger (const G4FastTrack &fastTrack) override
 

Private Attributes

G4double fEmax
 
GFlashHitMaker fHitMaker
 
const G4Envelope * fRegion
 
const TrackingActionfTrackingAction
 
LowEnergyFastSimParam param
 

Detailed Description

Definition at line 15 of file LowEnergyFastSimModel.h.

Constructor & Destructor Documentation

◆ LowEnergyFastSimModel()

LowEnergyFastSimModel::LowEnergyFastSimModel ( const G4String &  name,
G4Region *  region,
const edm::ParameterSet parSet 
)

Definition at line 16 of file LowEnergyFastSimModel.cc.

17  : G4VFastSimulationModel(name, region), fRegion(region), fTrackingAction(nullptr) {
18  fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
19 }

References fEmax, edm::ParameterSet::getParameter(), and GeV.

Member Function Documentation

◆ DoIt()

void LowEnergyFastSimModel::DoIt ( const G4FastTrack &  fastTrack,
G4FastStep &  fastStep 
)
override

Definition at line 37 of file LowEnergyFastSimModel.cc.

37  {
38  fastStep.KillPrimaryTrack();
39  fastStep.SetPrimaryTrackPathLength(0.0);
40  double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
41 
42  const G4ThreeVector& pos = fastTrack.GetPrimaryTrack()->GetPosition();
43 
44  G4double inPointEnergy = param.GetInPointEnergyFraction(energy) * energy;
45 
46  // take into account positron annihilation (not included in in-point)
47  if (-11 == fastTrack.GetPrimaryTrack()->GetDefinition()->GetPDGEncoding())
48  energy += twomass;
49 
50  const G4ThreeVector& momDir = fastTrack.GetPrimaryTrack()->GetMomentumDirection();
51  const G4ThreeVector& ortho = momDir.orthogonal();
52  const G4ThreeVector& cross = momDir.cross(ortho);
53 
54  // in point energy deposition
55  GFlashEnergySpot spot;
56  spot.SetEnergy(inPointEnergy);
57  spot.SetPosition(pos);
58  fHitMaker.make(&spot, &fastTrack);
59 
60  // tail energy deposition
61  G4double etail = energy - inPointEnergy;
62  const G4int nspots = int(etail) + 1;
63  const G4double tailEnergy = etail / (G4double)nspots;
64  for (G4int i = 0; i < nspots; ++i) {
65  const G4double radius = param.GetRadius(energy);
66  const G4double z = param.GetZ();
67 
68  const G4double phi = CLHEP::twopi * G4UniformRand();
69  const G4ThreeVector tailPos = pos + z * momDir + radius * std::cos(phi) * ortho + radius * std::sin(phi) * cross;
70 
71  spot.SetEnergy(tailEnergy);
72  spot.SetPosition(tailPos);
73  fHitMaker.make(&spot, &fastTrack);
74  }
75 }

References funct::cos(), Basic3DVector< T >::cross(), cross(), HCALHighEnergyHPDFilter_cfi::energy, fHitMaker, LowEnergyFastSimParam::GetInPointEnergyFraction(), LowEnergyFastSimParam::GetRadius(), LowEnergyFastSimParam::GetZ(), mps_fire::i, createfilelist::int, param, phi, CosmicsPD_Skims::radius, funct::sin(), twomass, and z.

◆ IsApplicable()

G4bool LowEnergyFastSimModel::IsApplicable ( const G4ParticleDefinition &  particle)
override

Definition at line 21 of file LowEnergyFastSimModel.cc.

21  {
22  return (11 == std::abs(particle.GetPDGEncoding()));
23 }

References funct::abs().

◆ ModelTrigger()

G4bool LowEnergyFastSimModel::ModelTrigger ( const G4FastTrack &  fastTrack)
override

Definition at line 25 of file LowEnergyFastSimModel.cc.

25  {
26  const G4Track* track = fastTrack.GetPrimaryTrack();
27  if (nullptr == fTrackingAction) {
28  fTrackingAction = static_cast<const TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
29  }
30  int pdgMother = std::abs(fTrackingAction->geant4Track()->GetDefinition()->GetPDGEncoding());
31  if (pdgMother == 11 || pdgMother == 22)
32  return false;
33  G4double energy = track->GetKineticEnergy();
34  return energy < fEmax && fRegion == fastTrack.GetEnvelope();
35 }

References funct::abs(), HCALHighEnergyHPDFilter_cfi::energy, fEmax, fRegion, fTrackingAction, TrackingAction::geant4Track(), and HLT_FULL_cff::track.

Member Data Documentation

◆ fEmax

G4double LowEnergyFastSimModel::fEmax
private

Definition at line 24 of file LowEnergyFastSimModel.h.

Referenced by LowEnergyFastSimModel(), and ModelTrigger().

◆ fHitMaker

GFlashHitMaker LowEnergyFastSimModel::fHitMaker
private

Definition at line 27 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fRegion

const G4Envelope* LowEnergyFastSimModel::fRegion
private

Definition at line 25 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

◆ fTrackingAction

const TrackingAction* LowEnergyFastSimModel::fTrackingAction
private

Definition at line 26 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

◆ param

LowEnergyFastSimParam LowEnergyFastSimModel::param
private

Definition at line 28 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

mps_fire.i
i
Definition: mps_fire.py:428
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
twomass
constexpr double twomass
Definition: LowEnergyFastSimModel.cc:14
pos
Definition: PixelAliasList.h:18
TrackingAction::geant4Track
const G4Track * geant4Track() const
Definition: TrackingAction.h:25
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
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::fEmax
G4double fEmax
Definition: LowEnergyFastSimModel.h:24
GeV
const double GeV
Definition: MathUtil.h:16
LowEnergyFastSimModel::param
LowEnergyFastSimParam param
Definition: LowEnergyFastSimModel.h:28
HLT_FULL_cff.region
region
Definition: HLT_FULL_cff.py:88271
createfilelist.int
int
Definition: createfilelist.py:10
LowEnergyFastSimParam::GetZ
G4double GetZ() const
Definition: LowEnergyFastSimParam.h:24
LowEnergyFastSimModel::fRegion
const G4Envelope * fRegion
Definition: LowEnergyFastSimModel.h:25
DDAxes::phi
LowEnergyFastSimModel::fHitMaker
GFlashHitMaker fHitMaker
Definition: LowEnergyFastSimModel.h:27
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LowEnergyFastSimModel::fTrackingAction
const TrackingAction * fTrackingAction
Definition: LowEnergyFastSimModel.h:26