Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GFlashShowerModel Class Reference

#include <GFlashShowerModel.hh>

+ Inheritance diagram for GFlashShowerModel:

Public Member Functions

 GFlashShowerModel (G4String, G4Envelope *)
 
 GFlashShowerModel (G4String)
 
 ~GFlashShowerModel ()
 
G4bool ModelTrigger (const G4FastTrack &)
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void DoIt (const G4FastTrack &, G4FastStep &)
 
void SetFlagParamType (G4int I)
 
void SetFlagParticleContainment (G4int I)
 
void SetStepInX0 (G4double Lenght)
 
void SetParameterisation (GVFlashShowerParameterisation &DP)
 
void SetHitMaker (GFlashHitMaker &Maker)
 
void SetParticleBounds (GFlashParticleBounds &SpecificBound)
 
G4int GetFlagParamType ()
 
G4int GetFlagParticleContainment ()
 
G4double GetStepInX0 ()
 
- Public Member Functions inherited from G4VFastSimulationModel
 G4VFastSimulationModel (const G4String &aName)
 
 G4VFastSimulationModel (const G4String &aName, G4Envelope *, G4bool IsUnique=FALSE)
 
virtual ~G4VFastSimulationModel ()=default
 
virtual G4bool AtRestModelTrigger (const G4FastTrack &)
 
virtual void AtRestDoIt (const G4FastTrack &, G4FastStep &)
 
virtual void Flush ()
 
const G4String GetName () const
 
G4bool operator== (const G4VFastSimulationModel &) const
 

Public Attributes

GFlashParticleBoundsPBound
 
GVFlashShowerParameterisationParameterisation
 

Detailed Description

Definition at line 59 of file GFlashShowerModel.hh.

Constructor & Destructor Documentation

◆ GFlashShowerModel() [1/2]

GFlashShowerModel::GFlashShowerModel ( G4String modelName,
G4Envelope * envelope )

Definition at line 55 of file GFlashShowerModel.cc.

56 : G4VFastSimulationModel(modelName, envelope), PBound(0), Parameterisation(0), HMaker(0)
57{
58 FlagParamType = 0;
59 FlagParticleContainment = 1;
60 StepInX0 = 0.1;
61 EnergyStop = 0.0;
62 Messenger = new GFlashShowerModelMessenger(this);
63}
G4VFastSimulationModel(const G4String &aName)
GFlashParticleBounds * PBound
GVFlashShowerParameterisation * Parameterisation

◆ GFlashShowerModel() [2/2]

GFlashShowerModel::GFlashShowerModel ( G4String modelName)

Definition at line 65 of file GFlashShowerModel.cc.

66 : G4VFastSimulationModel(modelName), PBound(0), Parameterisation(0), HMaker(0)
67{
68 FlagParamType = 1;
69 FlagParticleContainment = 1;
70 StepInX0 = 0.1;
71 EnergyStop = 0.0;
72 Messenger = new GFlashShowerModelMessenger(this);
73}

◆ ~GFlashShowerModel()

GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 75 of file GFlashShowerModel.cc.

76{
77 delete Messenger;
78}

Member Function Documentation

◆ DoIt()

void GFlashShowerModel::DoIt ( const G4FastTrack & fastTrack,
G4FastStep & fastStep )
virtual

Implements G4VFastSimulationModel.

Definition at line 163 of file GFlashShowerModel.cc.

164{
165 // parametrise electrons
168 ElectronDoIt(fastTrack, fastStep);
169}
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
const G4Track * GetPrimaryTrack() const
static G4Positron * PositronDefinition()
Definition G4Positron.cc:85
G4ParticleDefinition * GetDefinition() const

◆ GetFlagParamType()

G4int GFlashShowerModel::GetFlagParamType ( )
inline

Definition at line 83 of file GFlashShowerModel.hh.

83{ return FlagParamType; }

◆ GetFlagParticleContainment()

G4int GFlashShowerModel::GetFlagParticleContainment ( )
inline

Definition at line 84 of file GFlashShowerModel.hh.

84{ return FlagParticleContainment; }

◆ GetStepInX0()

G4double GFlashShowerModel::GetStepInX0 ( )
inline

Definition at line 85 of file GFlashShowerModel.hh.

85{ return StepInX0; }

◆ IsApplicable()

G4bool GFlashShowerModel::IsApplicable ( const G4ParticleDefinition & particleType)
virtual

Implements G4VFastSimulationModel.

Definition at line 80 of file GFlashShowerModel.cc.

81{
82 return &particleType == G4Electron::ElectronDefinition()
83 || &particleType == G4Positron::PositronDefinition();
84}

◆ ModelTrigger()

G4bool GFlashShowerModel::ModelTrigger ( const G4FastTrack & fastTrack)
virtual

Implements G4VFastSimulationModel.

Definition at line 90 of file GFlashShowerModel.cc.

92{
93 G4bool select = false;
94 if (FlagParamType != 0) {
95 G4double ParticleEnergy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
96 G4ParticleDefinition& ParticleType = *(fastTrack.GetPrimaryTrack()->GetDefinition());
97 if (ParticleEnergy > PBound->GetMinEneToParametrise(ParticleType)
98 && ParticleEnergy < PBound->GetMaxEneToParametrise(ParticleType))
99 {
100 // check conditions depending on particle flavour
101 // performance to be optimized @@@@@@@
102 Parameterisation->GenerateLongitudinalProfile(ParticleEnergy);
103 select = CheckParticleDefAndContainment(fastTrack);
104 if (select) EnergyStop = PBound->GetEneToKill(ParticleType);
105 }
106 }
107
108 return select;
109}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4double GetKineticEnergy() const

◆ SetFlagParamType()

void GFlashShowerModel::SetFlagParamType ( G4int I)
inline

Definition at line 74 of file GFlashShowerModel.hh.

74{ FlagParamType = I; }

◆ SetFlagParticleContainment()

void GFlashShowerModel::SetFlagParticleContainment ( G4int I)
inline

Definition at line 75 of file GFlashShowerModel.hh.

75{ FlagParticleContainment = I; }

◆ SetHitMaker()

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker & Maker)
inline

Definition at line 78 of file GFlashShowerModel.hh.

78{ HMaker = &Maker; }

◆ SetParameterisation()

void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation & DP)
inline

Definition at line 77 of file GFlashShowerModel.hh.

77{ Parameterisation = &DP; }

◆ SetParticleBounds()

void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds & SpecificBound)
inline

Definition at line 79 of file GFlashShowerModel.hh.

79{ PBound = &SpecificBound; }

◆ SetStepInX0()

void GFlashShowerModel::SetStepInX0 ( G4double Lenght)
inline

Definition at line 76 of file GFlashShowerModel.hh.

76{ StepInX0 = Lenght; }

Member Data Documentation

◆ Parameterisation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

◆ PBound

GFlashParticleBounds* GFlashShowerModel::PBound

The documentation for this class was generated from the following files: