GFlashShowerModel Class Reference

#include <GFlashShowerModel.hh>

Inheritance diagram for GFlashShowerModel:

G4VFastSimulationModel

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 ()

Data Fields

GFlashParticleBoundsPBound
GVFlashShowerParameterisationParameterisation

Detailed Description

Definition at line 60 of file GFlashShowerModel.hh.


Constructor & Destructor Documentation

GFlashShowerModel::GFlashShowerModel ( G4String  ,
G4Envelope  
)

Definition at line 57 of file GFlashShowerModel.cc.

00059   : G4VFastSimulationModel(modelName, envelope),
00060     PBound(0), Parameterisation(0), HMaker(0)
00061 {
00062   FlagParamType           = 0;
00063   FlagParticleContainment = 1;  
00064   StepInX0 = 0.1;
00065   Messenger       = new GFlashShowerModelMessenger(this);
00066 }

GFlashShowerModel::GFlashShowerModel ( G4String   ) 

Definition at line 68 of file GFlashShowerModel.cc.

00069   : G4VFastSimulationModel(modelName),
00070     PBound(0), Parameterisation(0), HMaker(0)
00071 {
00072   FlagParamType           =1;
00073   FlagParticleContainment = 1;  
00074   StepInX0 = 0.1; 
00075   Messenger       = new GFlashShowerModelMessenger(this); 
00076 }

GFlashShowerModel::~GFlashShowerModel (  ) 

Definition at line 78 of file GFlashShowerModel.cc.

00079 {
00080   delete Messenger;
00081 }


Member Function Documentation

void GFlashShowerModel::DoIt ( const G4FastTrack ,
G4FastStep  
) [virtual]

Implements G4VFastSimulationModel.

Definition at line 180 of file GFlashShowerModel.cc.

References G4Electron::ElectronDefinition(), G4Track::GetDefinition(), G4FastTrack::GetPrimaryTrack(), and G4Positron::PositronDefinition().

00181 {
00182   // parametrise electrons
00183   if(fastTrack.GetPrimaryTrack()->GetDefinition()
00184      == G4Electron::ElectronDefinition() || 
00185      fastTrack.GetPrimaryTrack()->GetDefinition()
00186      == G4Positron::PositronDefinition() ) 
00187   ElectronDoIt(fastTrack,fastStep);
00188 }

G4int GFlashShowerModel::GetFlagParamType (  )  [inline]

Definition at line 91 of file GFlashShowerModel.hh.

Referenced by GFlashShowerModelMessenger::GetCurrentValue().

00092       { return FlagParamType; }

G4int GFlashShowerModel::GetFlagParticleContainment (  )  [inline]

Definition at line 93 of file GFlashShowerModel.hh.

00094       { return FlagParticleContainment; }  

G4double GFlashShowerModel::GetStepInX0 (  )  [inline]

Definition at line 95 of file GFlashShowerModel.hh.

00096       { return StepInX0; }

G4bool GFlashShowerModel::IsApplicable ( const G4ParticleDefinition  )  [virtual]

Implements G4VFastSimulationModel.

Definition at line 84 of file GFlashShowerModel.cc.

References G4Electron::ElectronDefinition(), and G4Positron::PositronDefinition().

00085 { 
00086   return 
00087   &particleType == G4Electron::ElectronDefinition() ||
00088   &particleType == G4Positron::PositronDefinition(); 
00089 }

G4bool GFlashShowerModel::ModelTrigger ( const G4FastTrack  )  [virtual]

Implements G4VFastSimulationModel.

Definition at line 95 of file GFlashShowerModel.cc.

References GVFlashShowerParameterisation::GenerateLongitudinalProfile(), G4Track::GetDefinition(), GFlashParticleBounds::GetEneToKill(), G4Track::GetKineticEnergy(), GFlashParticleBounds::GetMaxEneToParametrise(), GFlashParticleBounds::GetMinEneToParametrise(), G4FastTrack::GetPrimaryTrack(), Parameterisation, and PBound.

00097 {
00098   G4bool select = false;
00099   if(FlagParamType != 0)                  
00100   {
00101     G4double  ParticleEnergy = fastTrack.GetPrimaryTrack()->GetKineticEnergy(); 
00102     G4ParticleDefinition &ParticleType =
00103       *(fastTrack.GetPrimaryTrack()->GetDefinition()); 
00104     if(ParticleEnergy > PBound->GetMinEneToParametrise(ParticleType) &&
00105        ParticleEnergy < PBound->GetMaxEneToParametrise(ParticleType) )
00106     {
00107       // check conditions depending on particle flavour
00108       // performance to be optimized @@@@@@@
00109       Parameterisation->GenerateLongitudinalProfile(ParticleEnergy);
00110       select     = CheckParticleDefAndContainment(fastTrack);  
00111       if (select) EnergyStop= PBound->GetEneToKill(ParticleType);
00112     }
00113   }
00114 
00115   return select; 
00116 }

void GFlashShowerModel::SetFlagParamType ( G4int  I  )  [inline]

Definition at line 76 of file GFlashShowerModel.hh.

Referenced by GFlashShowerModelMessenger::SetNewValue().

00077       { FlagParamType = I; }

void GFlashShowerModel::SetFlagParticleContainment ( G4int  I  )  [inline]

Definition at line 78 of file GFlashShowerModel.hh.

Referenced by GFlashShowerModelMessenger::SetNewValue().

00079       { FlagParticleContainment = I; }

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker Maker  )  [inline]

Definition at line 84 of file GFlashShowerModel.hh.

00085       { HMaker=&Maker; }

void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation DP  )  [inline]

Definition at line 82 of file GFlashShowerModel.hh.

References Parameterisation.

00083       { Parameterisation=&DP;}

void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds SpecificBound  )  [inline]

Definition at line 86 of file GFlashShowerModel.hh.

References PBound.

00087       { PBound =&SpecificBound; }

void GFlashShowerModel::SetStepInX0 ( G4double  Lenght  )  [inline]

Definition at line 80 of file GFlashShowerModel.hh.

Referenced by GFlashShowerModelMessenger::SetNewValue().

00081       { StepInX0=Lenght; } 


Field Documentation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

Definition at line 102 of file GFlashShowerModel.hh.

Referenced by ModelTrigger(), and SetParameterisation().

GFlashParticleBounds* GFlashShowerModel::PBound

Definition at line 101 of file GFlashShowerModel.hh.

Referenced by GFlashShowerModelMessenger::GetCurrentValue(), ModelTrigger(), GFlashShowerModelMessenger::SetNewValue(), and SetParticleBounds().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:54:02 2013 for Geant4 by  doxygen 1.4.7