CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 G4Envelope * fRegion
 
G4ThreeVector fTailPos
 
const TrackingActionfTrackingAction
 

Detailed Description

Definition at line 16 of file LowEnergyFastSimModel.h.

Constructor & Destructor Documentation

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

Definition at line 16 of file LowEnergyFastSimModel.cc.

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

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

Member Function Documentation

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

Definition at line 43 of file LowEnergyFastSimModel.cc.

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

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

Definition at line 25 of file LowEnergyFastSimModel.cc.

References funct::abs().

25  {
26  return (11 == std::abs(particle.GetPDGEncoding()));
27 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
G4bool LowEnergyFastSimModel::ModelTrigger ( const G4FastTrack &  fastTrack)
override

Definition at line 29 of file LowEnergyFastSimModel.cc.

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

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

Member Data Documentation

G4bool LowEnergyFastSimModel::fCheck
private

Definition at line 28 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

G4double LowEnergyFastSimModel::fEmax
private

Definition at line 25 of file LowEnergyFastSimModel.h.

Referenced by LowEnergyFastSimModel(), and ModelTrigger().

GFlashHitMaker LowEnergyFastSimModel::fHitMaker
private

Definition at line 30 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

LowEnergyFastSimParam LowEnergyFastSimModel::fParam
private

Definition at line 31 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

const G4Envelope* LowEnergyFastSimModel::fRegion
private

Definition at line 26 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().

G4ThreeVector LowEnergyFastSimModel::fTailPos
private

Definition at line 29 of file LowEnergyFastSimModel.h.

Referenced by DoIt().

const TrackingAction* LowEnergyFastSimModel::fTrackingAction
private

Definition at line 27 of file LowEnergyFastSimModel.h.

Referenced by ModelTrigger().