Geant4 11.1.1
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 &)
 
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 Messenger = new GFlashShowerModelMessenger(this);
65}
GFlashParticleBounds * PBound
GVFlashShowerParameterisation * Parameterisation

◆ GFlashShowerModel() [2/2]

GFlashShowerModel::GFlashShowerModel ( G4String  modelName)

Definition at line 67 of file GFlashShowerModel.cc.

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

◆ ~GFlashShowerModel()

GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 77 of file GFlashShowerModel.cc.

78{
79 delete Messenger;
80}

Member Function Documentation

◆ DoIt()

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

Implements G4VFastSimulationModel.

Definition at line 179 of file GFlashShowerModel.cc.

180{
181 // parametrise electrons
182 if(fastTrack.GetPrimaryTrack()->GetDefinition()
184 fastTrack.GetPrimaryTrack()->GetDefinition()
186 ElectronDoIt(fastTrack,fastStep);
187}
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:88
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 83 of file GFlashShowerModel.cc.

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

◆ ModelTrigger()

G4bool GFlashShowerModel::ModelTrigger ( const G4FastTrack fastTrack)
virtual

Implements G4VFastSimulationModel.

Definition at line 94 of file GFlashShowerModel.cc.

96{
97 G4bool select = false;
98 if(FlagParamType != 0)
99 {
100 G4double ParticleEnergy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
102 *(fastTrack.GetPrimaryTrack()->GetDefinition());
103 if(ParticleEnergy > PBound->GetMinEneToParametrise(ParticleType) &&
104 ParticleEnergy < PBound->GetMaxEneToParametrise(ParticleType) )
105 {
106 // check conditions depending on particle flavour
107 // performance to be optimized @@@@@@@
109 select = CheckParticleDefAndContainment(fastTrack);
110 if (select) EnergyStop= PBound->GetEneToKill(ParticleType);
111 }
112 }
113
114 return select;
115}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
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 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: