Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GFlashShowerModel Class Reference

#include <GFlashShowerModel.hh>

Inheritance diagram for GFlashShowerModel:
Collaboration 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 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::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
G4VFastSimulationModel(const G4String &aName)
GVFlashShowerParameterisation * Parameterisation
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 }
GFlashParticleBounds * PBound
G4VFastSimulationModel(const G4String &aName)
GVFlashShowerParameterisation * Parameterisation
GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 78 of file GFlashShowerModel.cc.

79 {
80  delete Messenger;
81 }

Member Function Documentation

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 }
G4ParticleDefinition * GetDefinition() const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89

Here is the call graph for this function:

G4int GFlashShowerModel::GetFlagParamType ( )
inline

Definition at line 91 of file GFlashShowerModel.hh.

92  { return FlagParamType; }

Here is the caller graph for this function:

G4int GFlashShowerModel::GetFlagParticleContainment ( )
inline

Definition at line 93 of file GFlashShowerModel.hh.

94  { return FlagParticleContainment; }
G4double GFlashShowerModel::GetStepInX0 ( )
inline

Definition at line 95 of file GFlashShowerModel.hh.

96  { return StepInX0; }
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 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89

Here is the call graph for this function:

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 }
G4ParticleDefinition * GetDefinition() const
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
GFlashParticleBounds * PBound
G4double GetMaxEneToParametrise(G4ParticleDefinition &particleType)
G4double GetKineticEnergy() const
G4double GetMinEneToParametrise(G4ParticleDefinition &particleType)
virtual void GenerateLongitudinalProfile(G4double Energy)=0
bool G4bool
Definition: G4Types.hh:79
G4double GetEneToKill(G4ParticleDefinition &particleType)
double G4double
Definition: G4Types.hh:76
GVFlashShowerParameterisation * Parameterisation

Here is the call graph for this function:

void GFlashShowerModel::SetFlagParamType ( G4int  I)
inline

Definition at line 76 of file GFlashShowerModel.hh.

77  { FlagParamType = I; }

Here is the caller graph for this function:

void GFlashShowerModel::SetFlagParticleContainment ( G4int  I)
inline

Definition at line 78 of file GFlashShowerModel.hh.

79  { FlagParticleContainment = I; }

Here is the caller graph for this function:

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker Maker)
inline

Definition at line 84 of file GFlashShowerModel.hh.

85  { HMaker=&Maker; }
void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation DP)
inline

Definition at line 82 of file GFlashShowerModel.hh.

83  { Parameterisation=&DP;}
GVFlashShowerParameterisation * Parameterisation
void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds SpecificBound)
inline

Definition at line 86 of file GFlashShowerModel.hh.

87  { PBound =&SpecificBound; }
GFlashParticleBounds * PBound
void GFlashShowerModel::SetStepInX0 ( G4double  Lenght)
inline

Definition at line 80 of file GFlashShowerModel.hh.

81  { StepInX0=Lenght; }

Here is the caller graph for this function:

Member Data Documentation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

Definition at line 102 of file GFlashShowerModel.hh.

GFlashParticleBounds* GFlashShowerModel::PBound

Definition at line 101 of file GFlashShowerModel.hh.


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