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 21 of file LowEnergyFastSimModel.cc.

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

22  : G4VFastSimulationModel(name, region),
23  fRegion(region),
24  fTrackingAction(nullptr),
25  fCheck(true),
26  fTailPos(0., 0., 0.) {
27  fEmax = parSet.getParameter<double>("LowEnergyGflashEcalEmax") * CLHEP::GeV;
28  fPositron = G4Positron::Positron();
29  fMaterial = nullptr;
30  auto table = G4Material::GetMaterialTable();
31  for (auto& mat : *table) {
32  G4String nam = mat->GetName();
33  size_t n = nam.size();
34  if (n > 4) {
35  G4String sn = nam.substr(n - 5, 5);
36  if (sn == "PbWO4") {
37  fMaterial = mat;
38  break;
39  }
40  }
41  }
42  G4String nm = (nullptr == fMaterial) ? "not found!" : fMaterial->GetName();
43  edm::LogVerbatim("LowEnergyFastSimModel") << "LowEGFlash material: <" << nm << ">";
44 }
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 76 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_2022v15_cff::track, twomass, and z.

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

References funct::abs().

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

◆ ModelTrigger()

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

Definition at line 50 of file LowEnergyFastSimModel.cc.

References funct::abs(), HCALHighEnergyHPDFilter_cfi::energy, fCheck, fEmax, fMaterial, fTrackingAction, TrackingAction::geant4Track(), TrackInformation::isPrimary(), and HLT_2022v15_cff::track.

50  {
51  const G4Track* track = fastTrack.GetPrimaryTrack();
52  G4double energy = track->GetKineticEnergy();
53  if (fMaterial != track->GetMaterial() || energy >= fEmax)
54  return false;
55 
56  /*
57  edm::LogVerbatim("LowEnergyFastSimModel") << track->GetDefinition()->GetParticleName()
58  << " Ekin(MeV)=" << energy << " material: <"
59  << track->GetMaterial()->GetName() << ">";
60  */
61  if (fCheck) {
62  if (nullptr == fTrackingAction) {
63  fTrackingAction = static_cast<const TrackingAction*>(G4EventManager::GetEventManager()->GetUserTrackingAction());
64  }
65  const G4Track* mother = fTrackingAction->geant4Track();
66  const TrackInformation* ptr = static_cast<TrackInformation*>(mother->GetUserInformation());
67  if (ptr->isPrimary()) {
68  int pdgMother = mother->GetDefinition()->GetPDGEncoding();
69  if (std::abs(pdgMother) == 11 || pdgMother == 22)
70  return false;
71  }
72  }
73  return true;
74 }
bool isPrimary() const
const G4Material * fMaterial
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(), and ModelTrigger().

◆ 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(), and ModelTrigger().

◆ 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().