Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 56 of file GFlashShowerModel.cc.

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

◆ GFlashShowerModel() [2/2]

GFlashShowerModel::GFlashShowerModel ( G4String modelName)

Definition at line 68 of file GFlashShowerModel.cc.

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

◆ ~GFlashShowerModel()

GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 79 of file GFlashShowerModel.cc.

80{
81 delete Messenger;
82}

Member Function Documentation

◆ DoIt()

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

Implements G4VFastSimulationModel.

Definition at line 181 of file GFlashShowerModel.cc.

182{
183 // parametrise electrons
184 if(fastTrack.GetPrimaryTrack()->GetDefinition()
186 fastTrack.GetPrimaryTrack()->GetDefinition()
188 ElectronDoIt(fastTrack,fastStep);
189}
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 90 of file GFlashShowerModel.hh.

91 { return FlagParamType; }

Referenced by GFlashShowerModelMessenger::GetCurrentValue().

◆ GetFlagParticleContainment()

G4int GFlashShowerModel::GetFlagParticleContainment ( )
inline

Definition at line 92 of file GFlashShowerModel.hh.

93 { return FlagParticleContainment; }

◆ GetStepInX0()

G4double GFlashShowerModel::GetStepInX0 ( )
inline

Definition at line 94 of file GFlashShowerModel.hh.

95 { return StepInX0; }

◆ IsApplicable()

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

Implements G4VFastSimulationModel.

Definition at line 85 of file GFlashShowerModel.cc.

86{
87 return
88 &particleType == G4Electron::ElectronDefinition() ||
89 &particleType == G4Positron::PositronDefinition();
90}

◆ ModelTrigger()

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

Implements G4VFastSimulationModel.

Definition at line 96 of file GFlashShowerModel.cc.

98{
99 G4bool select = false;
100 if(FlagParamType != 0)
101 {
102 G4double ParticleEnergy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
104 *(fastTrack.GetPrimaryTrack()->GetDefinition());
105 if(ParticleEnergy > PBound->GetMinEneToParametrise(ParticleType) &&
106 ParticleEnergy < PBound->GetMaxEneToParametrise(ParticleType) )
107 {
108 // check conditions depending on particle flavour
109 // performance to be optimized @@@@@@@
111 select = CheckParticleDefAndContainment(fastTrack);
112 if (select) EnergyStop= PBound->GetEneToKill(ParticleType);
113 }
114 }
115
116 return select;
117}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4double GetKineticEnergy() const
G4double GetMinEneToParametrise(G4ParticleDefinition &particleType)
G4double GetEneToKill(G4ParticleDefinition &particleType)
virtual void GenerateLongitudinalProfile(G4double Energy)=0

◆ SetFlagParamType()

void GFlashShowerModel::SetFlagParamType ( G4int I)
inline

Definition at line 75 of file GFlashShowerModel.hh.

76 { FlagParamType = I; }

Referenced by GFlashShowerModelMessenger::SetNewValue().

◆ SetFlagParticleContainment()

void GFlashShowerModel::SetFlagParticleContainment ( G4int I)
inline

Definition at line 77 of file GFlashShowerModel.hh.

78 { FlagParticleContainment = I; }

Referenced by GFlashShowerModelMessenger::SetNewValue().

◆ SetHitMaker()

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker & Maker)
inline

Definition at line 83 of file GFlashShowerModel.hh.

84 { HMaker=&Maker; }

◆ SetParameterisation()

void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation & DP)
inline

Definition at line 81 of file GFlashShowerModel.hh.

82 { Parameterisation=&DP;}

◆ SetParticleBounds()

void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds & SpecificBound)
inline

Definition at line 85 of file GFlashShowerModel.hh.

86 { PBound =&SpecificBound; }

◆ SetStepInX0()

void GFlashShowerModel::SetStepInX0 ( G4double Lenght)
inline

Definition at line 79 of file GFlashShowerModel.hh.

80 { StepInX0=Lenght; }

Referenced by GFlashShowerModelMessenger::SetNewValue().

Member Data Documentation

◆ Parameterisation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

Definition at line 101 of file GFlashShowerModel.hh.

Referenced by ModelTrigger(), and SetParameterisation().

◆ PBound


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