Geant4 9.6.0
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 ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)=0
 
virtual G4bool ModelTrigger (const G4FastTrack &)=0
 
virtual void DoIt (const G4FastTrack &, G4FastStep &)=0
 
virtual G4bool AtRestModelTrigger (const G4FastTrack &)
 
virtual void AtRestDoIt (const G4FastTrack &, G4FastStep &)
 
const G4String GetName () const
 
G4bool operator== (const G4VFastSimulationModel &) const
 

Public Attributes

GFlashParticleBoundsPBound
 
GVFlashShowerParameterisationParameterisation
 

Detailed Description

Definition at line 60 of file GFlashShowerModel.hh.

Constructor & Destructor Documentation

◆ GFlashShowerModel() [1/2]

GFlashShowerModel::GFlashShowerModel ( G4String  modelName,
G4Envelope envelope 
)

Definition at line 57 of file GFlashShowerModel.cc.

59 : G4VFastSimulationModel(modelName, envelope),
60 PBound(0), Parameterisation(0), HMaker(0)
61{
62 FlagParamType = 0;
63 FlagParticleContainment = 1;
64 StepInX0 = 0.1;
65 Messenger = new GFlashShowerModelMessenger(this);
66}
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 Messenger = new GFlashShowerModelMessenger(this);
76}

◆ ~GFlashShowerModel()

GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 78 of file GFlashShowerModel.cc.

79{
80 delete Messenger;
81}

Member Function Documentation

◆ DoIt()

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

Implements G4VFastSimulationModel.

Definition at line 180 of file GFlashShowerModel.cc.

181{
182 // parametrise electrons
183 if(fastTrack.GetPrimaryTrack()->GetDefinition()
185 fastTrack.GetPrimaryTrack()->GetDefinition()
187 ElectronDoIt(fastTrack,fastStep);
188}
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
G4ParticleDefinition * GetDefinition() const

◆ GetFlagParamType()

G4int GFlashShowerModel::GetFlagParamType ( )
inline

Definition at line 91 of file GFlashShowerModel.hh.

92 { return FlagParamType; }

Referenced by GFlashShowerModelMessenger::GetCurrentValue().

◆ GetFlagParticleContainment()

G4int GFlashShowerModel::GetFlagParticleContainment ( )
inline

Definition at line 93 of file GFlashShowerModel.hh.

94 { return FlagParticleContainment; }

◆ GetStepInX0()

G4double GFlashShowerModel::GetStepInX0 ( )
inline

Definition at line 95 of file GFlashShowerModel.hh.

96 { return StepInX0; }

◆ IsApplicable()

G4bool GFlashShowerModel::IsApplicable ( const G4ParticleDefinition particleType)
virtual

Implements G4VFastSimulationModel.

Definition at line 84 of file GFlashShowerModel.cc.

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

◆ ModelTrigger()

G4bool GFlashShowerModel::ModelTrigger ( const G4FastTrack fastTrack)
virtual

Implements G4VFastSimulationModel.

Definition at line 95 of file GFlashShowerModel.cc.

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

◆ SetFlagParamType()

void GFlashShowerModel::SetFlagParamType ( G4int  I)
inline

Definition at line 76 of file GFlashShowerModel.hh.

77 { FlagParamType = I; }

Referenced by GFlashShowerModelMessenger::SetNewValue().

◆ SetFlagParticleContainment()

void GFlashShowerModel::SetFlagParticleContainment ( G4int  I)
inline

Definition at line 78 of file GFlashShowerModel.hh.

79 { FlagParticleContainment = I; }

Referenced by GFlashShowerModelMessenger::SetNewValue().

◆ SetHitMaker()

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker Maker)
inline

Definition at line 84 of file GFlashShowerModel.hh.

85 { HMaker=&Maker; }

◆ SetParameterisation()

void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation DP)
inline

Definition at line 82 of file GFlashShowerModel.hh.

83 { Parameterisation=&DP;}

◆ SetParticleBounds()

void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds SpecificBound)
inline

Definition at line 86 of file GFlashShowerModel.hh.

87 { PBound =&SpecificBound; }

◆ SetStepInX0()

void GFlashShowerModel::SetStepInX0 ( G4double  Lenght)
inline

Definition at line 80 of file GFlashShowerModel.hh.

81 { StepInX0=Lenght; }

Referenced by GFlashShowerModelMessenger::SetNewValue().

Member Data Documentation

◆ Parameterisation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

Definition at line 102 of file GFlashShowerModel.hh.

Referenced by ModelTrigger(), and SetParameterisation().

◆ PBound


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