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
 
const G4Material * fMaterial
 
LowEnergyFastSimParam fParam
 
const G4ParticleDefinition * fPositron
 
const G4Envelope * fRegion
 
G4ThreeVector fTailPos
 
const TrackingActionfTrackingAction
 

Detailed Description

Definition at line 18 of file LowEnergyFastSimModel.h.

Constructor & Destructor Documentation

◆ LowEnergyFastSimModel()

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

Definition at line 20 of file LowEnergyFastSimModel.cc.

References fEmax, fMaterial, fPositron, edm::ParameterSet::getParameter(), dqmiodumpmetadata::n, and TableParser::table.

21  : G4VFastSimulationModel(name, region),
22  fRegion(region),
23  fTrackingAction(nullptr),
24  fCheck(false),
25  fTailPos(0., 0., 0.) {
26  fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
27  fPositron = G4Positron::Positron();
28  fMaterial = nullptr;
29  auto table = G4Material::GetMaterialTable();
30  for (auto& mat : *table) {
31  G4String nam = mat->GetName();
32  size_t n = nam.size();
33  if (n > 4) {
34  G4String sn = nam.substr(n - 5, 5);
35  if (sn == "PbWO4") {
36  fMaterial = mat;
37  break;
38  }
39  }
40  }
41  G4String nm = (nullptr == fMaterial) ? "not found!" : fMaterial->GetName();
42  edm::LogVerbatim("LowEnergyFastSimModel") << "LowEGFlash material: <" << nm << ">";
43 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const G4ParticleDefinition * fPositron
const G4Envelope * fRegion
const G4Material * fMaterial
const TrackingAction * fTrackingAction

Member Function Documentation

◆ DoIt()

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

Definition at line 68 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, scaleFactor, funct::sin(), HLT_2022v14_cff::track, twomass, and z.

68  {
69  fastStep.KillPrimaryTrack();
70  fastStep.SetPrimaryTrackPathLength(0.0);
71  auto track = fastTrack.GetPrimaryTrack();
72  G4double energy = track->GetKineticEnergy() * scaleFactor;
73 
74  const G4ThreeVector& pos = track->GetPosition();
75 
76  G4double inPointEnergy = fParam.GetInPointEnergyFraction(energy) * energy;
77 
78  // take into account positron annihilation (not included in in-point)
79  if (fPositron == track->GetDefinition())
80  energy += twomass;
81 
82  const G4ThreeVector& momDir = track->GetMomentumDirection();
83 
84  // in point energy deposition
85  GFlashEnergySpot spot;
86  spot.SetEnergy(inPointEnergy);
87  spot.SetPosition(pos);
88  fHitMaker.make(&spot, &fastTrack);
89 
90  // tail energy deposition
91  const G4double etail = energy - inPointEnergy;
92  const G4int nspots = etail;
93  const G4double tailEnergy = etail / (nspots + 1);
94  /*
95  edm::LogVerbatim("LowEnergyFastSimModel") << track->GetDefinition()->GetParticleName()
96  << " Ekin(MeV)=" << energy << " material: <"
97  << track->GetMaterial()->GetName()
98  << "> Elocal=" << inPointEnergy
99  << " Etail=" << tailEnergy
100  << " Nspots=" << nspots+1;
101  */
102  for (G4int i = 0; i <= nspots; ++i) {
103  const G4double r = fParam.GetRadius(energy);
104  const G4double z = fParam.GetZ();
105 
106  const G4double phi = CLHEP::twopi * G4UniformRand();
107  fTailPos.set(r * std::cos(phi), r * std::sin(phi), z);
108  fTailPos.rotateUz(momDir);
109  fTailPos += pos;
110 
111  spot.SetEnergy(tailEnergy);
112  spot.SetPosition(fTailPos);
113  fHitMaker.make(&spot, &fastTrack);
114  }
115 }
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
constexpr G4double scaleFactor
G4double GetRadius(G4double energy) const

◆ IsApplicable()

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

Definition at line 45 of file LowEnergyFastSimModel.cc.

References funct::abs().

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

◆ ModelTrigger()

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

Definition at line 49 of file LowEnergyFastSimModel.cc.

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

49  {
50  const G4Track* track = fastTrack.GetPrimaryTrack();
51  if (fCheck) {
52  if (nullptr == fTrackingAction) {
53  fTrackingAction = static_cast<const TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
54  }
55  int pdgMother = std::abs(fTrackingAction->geant4Track()->GetDefinition()->GetPDGEncoding());
56  if (pdgMother == 11 || pdgMother == 22)
57  return false;
58  }
59  G4double energy = track->GetKineticEnergy();
60  /*
61  edm::LogVerbatim("LowEnergyFastSimModel") << track->GetDefinition()->GetParticleName()
62  << " Ekin(MeV)=" << energy << " material: <"
63  << track->GetMaterial()->GetName() << ">";
64  */
65  return (energy < fEmax && fMaterial == track->GetMaterial());
66 }
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 32 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

◆ fEmax

G4double LowEnergyFastSimModel::fEmax
private

Definition at line 27 of file LowEnergyFastSimModel.h.

Referenced by LowEnergyFastSimModel().

◆ fHitMaker

GFlashHitMaker LowEnergyFastSimModel::fHitMaker
private

Definition at line 34 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fMaterial

const G4Material* LowEnergyFastSimModel::fMaterial
private

Definition at line 31 of file LowEnergyFastSimModel.h.

Referenced by LowEnergyFastSimModel().

◆ fParam

LowEnergyFastSimParam LowEnergyFastSimModel::fParam
private

Definition at line 35 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fPositron

const G4ParticleDefinition* LowEnergyFastSimModel::fPositron
private

Definition at line 30 of file LowEnergyFastSimModel.h.

Referenced by DoIt(), and LowEnergyFastSimModel().

◆ fRegion

const G4Envelope* LowEnergyFastSimModel::fRegion
private

Definition at line 28 of file LowEnergyFastSimModel.h.

◆ fTailPos

G4ThreeVector LowEnergyFastSimModel::fTailPos
private

Definition at line 33 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

◆ fTrackingAction

const TrackingAction* LowEnergyFastSimModel::fTrackingAction
private

Definition at line 29 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().