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

G4bool fCheck
 
G4double fEmax
 
GFlashHitMaker fHitMaker
 
LowEnergyFastSimParam fParam
 
const G4ParticleDefinition * fPositron
 
const G4Envelope * fRegion
 
G4ThreeVector fTailPos
 
const TrackingActionfTrackingAction
 

Detailed Description

Definition at line 17 of file LowEnergyFastSimModel.h.

Constructor & Destructor Documentation

◆ LowEnergyFastSimModel()

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

Definition at line 18 of file LowEnergyFastSimModel.cc.

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

19  : G4VFastSimulationModel(name, region),
20  fRegion(region),
21  fTrackingAction(nullptr),
22  fCheck(false),
23  fTailPos(0., 0., 0.) {
24  fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
25  fPositron = G4Positron::Positron();
26 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const G4ParticleDefinition * fPositron
const G4Envelope * fRegion
const TrackingAction * fTrackingAction

Member Function Documentation

◆ DoIt()

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

Definition at line 46 of file LowEnergyFastSimModel.cc.

References funct::cos(), HCALHighEnergyHPDFilter_cfi::energy, fHitMaker, fParam, fPositron, fTailPos, LowEnergyFastSimParam::GetInPointEnergyFraction(), LowEnergyFastSimParam::GetRadius(), LowEnergyFastSimParam::GetZ(), mps_fire::i, phi, alignCSCRings::r, funct::sin(), HLT_2022v12_cff::track, twomass, and z.

46  {
47  fastStep.KillPrimaryTrack();
48  fastStep.SetPrimaryTrackPathLength(0.0);
49  auto track = fastTrack.GetPrimaryTrack();
50  G4double energy = track->GetKineticEnergy();
51 
52  const G4ThreeVector& pos = track->GetPosition();
53 
54  G4double inPointEnergy = fParam.GetInPointEnergyFraction(energy) * energy;
55 
56  // take into account positron annihilation (not included in in-point)
57  if (fPositron == track->GetDefinition())
58  energy += twomass;
59 
60  const G4ThreeVector& momDir = track->GetMomentumDirection();
61 
62  // in point energy deposition
63  GFlashEnergySpot spot;
64  spot.SetEnergy(inPointEnergy);
65  spot.SetPosition(pos);
66  fHitMaker.make(&spot, &fastTrack);
67 
68  // tail energy deposition
69  G4double etail = energy - inPointEnergy;
70  const G4int nspots = (G4int)(etail) + 1;
71  const G4double tailEnergy = etail / (G4double)nspots;
72  for (G4int i = 0; i < nspots; ++i) {
73  const G4double r = fParam.GetRadius(energy);
74  const G4double z = fParam.GetZ();
75 
76  const G4double phi = CLHEP::twopi * G4UniformRand();
77  fTailPos.set(r * std::cos(phi), r * std::sin(phi), z);
78  fTailPos.rotateUz(momDir);
79  fTailPos += pos;
80 
81  spot.SetEnergy(tailEnergy);
82  spot.SetPosition(fTailPos);
83  fHitMaker.make(&spot, &fastTrack);
84  }
85 }
const G4ParticleDefinition * fPositron
constexpr G4double twomass
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
LowEnergyFastSimParam fParam
G4double GetInPointEnergyFraction(G4double energy) const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
G4double GetRadius(G4double energy) const

◆ IsApplicable()

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

Definition at line 28 of file LowEnergyFastSimModel.cc.

References funct::abs().

28  {
29  return (11 == std::abs(particle.GetPDGEncoding()));
30 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ ModelTrigger()

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

Definition at line 32 of file LowEnergyFastSimModel.cc.

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

32  {
33  const G4Track* track = fastTrack.GetPrimaryTrack();
34  if (fCheck) {
35  if (nullptr == fTrackingAction) {
36  fTrackingAction = static_cast<const TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
37  }
38  int pdgMother = std::abs(fTrackingAction->geant4Track()->GetDefinition()->GetPDGEncoding());
39  if (pdgMother == 11 || pdgMother == 22)
40  return false;
41  }
42  G4double energy = track->GetKineticEnergy();
43  return (energy < fEmax && fRegion == fastTrack.GetEnvelope());
44 }
const G4Envelope * fRegion
const G4Track * geant4Track() const
const TrackingAction * fTrackingAction
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

Member Data Documentation

◆ fCheck

G4bool LowEnergyFastSimModel::fCheck
private

Definition at line 30 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

◆ fEmax

G4double LowEnergyFastSimModel::fEmax
private

Definition at line 26 of file LowEnergyFastSimModel.h.

Referenced by LowEnergyFastSimModel(), and ModelTrigger().

◆ fHitMaker

GFlashHitMaker LowEnergyFastSimModel::fHitMaker
private

Definition at line 32 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fParam

LowEnergyFastSimParam LowEnergyFastSimModel::fParam
private

Definition at line 33 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fPositron

const G4ParticleDefinition* LowEnergyFastSimModel::fPositron
private

Definition at line 29 of file LowEnergyFastSimModel.h.

Referenced by DoIt(), and LowEnergyFastSimModel().

◆ fRegion

const G4Envelope* LowEnergyFastSimModel::fRegion
private

Definition at line 27 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

◆ fTailPos

G4ThreeVector LowEnergyFastSimModel::fTailPos
private

Definition at line 31 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fTrackingAction

const TrackingAction* LowEnergyFastSimModel::fTrackingAction
private

Definition at line 28 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().