Geant4  10.02.p03
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
 

Private Member Functions

void ElectronDoIt (const G4FastTrack &, G4FastStep &)
 
G4bool CheckParticleDefAndContainment (const G4FastTrack &fastTrack)
 
G4bool CheckContainment (const G4FastTrack &fastTrack)
 

Private Attributes

GFlashHitMakerHMaker
 
GFlashShowerModelMessengerMessenger
 
G4int FlagParamType
 
G4int FlagParticleContainment
 0=no GFlash 1=only em showers parametrized More...
 
G4double StepInX0
 0=no check ///1=only fully contained... More...
 
G4double EnergyStop
 

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;
64  StepInX0 = 0.1;
66 }
GFlashParticleBounds * PBound
G4int FlagParticleContainment
0=no GFlash 1=only em showers parametrized
GFlashShowerModelMessenger * Messenger
GFlashHitMaker * HMaker
G4double StepInX0
0=no check ///1=only fully contained...
G4VFastSimulationModel(const G4String &aName)
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;
74  StepInX0 = 0.1;
76 }
GFlashParticleBounds * PBound
G4int FlagParticleContainment
0=no GFlash 1=only em showers parametrized
GFlashShowerModelMessenger * Messenger
GFlashHitMaker * HMaker
G4double StepInX0
0=no check ///1=only fully contained...
G4VFastSimulationModel(const G4String &aName)
GVFlashShowerParameterisation * Parameterisation

◆ ~GFlashShowerModel()

GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 78 of file GFlashShowerModel.cc.

79 {
80  delete Messenger;
81 }
GFlashShowerModelMessenger * Messenger

Member Function Documentation

◆ CheckContainment()

G4bool GFlashShowerModel::CheckContainment ( const G4FastTrack fastTrack)
private

Definition at line 138 of file GFlashShowerModel.cc.

139 {
140  G4bool filter=false;
141  // track informations
142  G4ThreeVector DirectionShower=fastTrack.GetPrimaryTrackLocalDirection();
143  G4ThreeVector InitialPositionShower=fastTrack.GetPrimaryTrackLocalPosition();
144 
145  G4ThreeVector OrthoShower, CrossShower;
146  // Returns orthogonal vector
147  OrthoShower = DirectionShower.orthogonal();
148  // Shower in direction perpendicular to OrthoShower and DirectionShower
149  CrossShower = DirectionShower.cross(OrthoShower);
150 
153  G4int CosPhi[4] = {1,0,-1,0};
154  G4int SinPhi[4] = {0,1,0,-1};
155 
156  G4ThreeVector Position;
157  G4int NlateralInside=0;
158  // pointer to solid we're in
159  G4VSolid *SolidCalo = fastTrack.GetEnvelopeSolid();
160  for(int i=0; i<4 ;i++)
161  {
162  // polar coordinates
163  Position = InitialPositionShower +
164  Z*DirectionShower +
165  R*CosPhi[i]*OrthoShower +
166  R*SinPhi[i]*CrossShower ;
167 
168  if(SolidCalo->Inside(Position) != kOutside)
169  NlateralInside++;
170  }
171 
172  // choose to parameterise or flag when all inetc...
173  if(NlateralInside==4) filter=true;
174  // std::cout << " points = " <<NlateralInside << std::endl;
175  return filter;
176 }
G4VSolid * GetEnvelopeSolid() const
Definition: G4FastTrack.hh:203
int G4int
Definition: G4Types.hh:78
Hep3Vector cross(const Hep3Vector &) const
virtual EInside Inside(const G4ThreeVector &p) const =0
G4ThreeVector GetPrimaryTrackLocalDirection() const
Definition: G4FastTrack.hh:223
Float_t Z
bool G4bool
Definition: G4Types.hh:79
Hep3Vector orthogonal() const
virtual G4double GetAveR90()=0
virtual G4double GetAveT90()=0
G4ThreeVector GetPrimaryTrackLocalPosition() const
Definition: G4FastTrack.hh:213
double G4double
Definition: G4Types.hh:76
GVFlashShowerParameterisation * Parameterisation
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CheckParticleDefAndContainment()

G4bool GFlashShowerModel::CheckParticleDefAndContainment ( const G4FastTrack fastTrack)
private

Definition at line 120 of file GFlashShowerModel.cc.

121 {
122  G4bool filter=false;
124  fastTrack.GetPrimaryTrack()->GetDefinition();
125 
126  if( ParticleType == G4Electron::ElectronDefinition() ||
127  ParticleType == G4Positron::PositronDefinition() )
128  {
129  filter=true;
130  if(FlagParticleContainment == 1)
131  {
132  filter=CheckContainment(fastTrack);
133  }
134  }
135  return filter;
136 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4int FlagParticleContainment
0=no GFlash 1=only em showers parametrized
bool G4bool
Definition: G4Types.hh:79
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
G4bool CheckContainment(const G4FastTrack &fastTrack)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
void ElectronDoIt(const G4FastTrack &, G4FastStep &)
Here is the call graph for this function:

◆ ElectronDoIt()

void GFlashShowerModel::ElectronDoIt ( const G4FastTrack fastTrack,
G4FastStep fastStep 
)
private

Generate longitudinal profile

Initialisation of long. loop variables

Begin Longitudinal Loop

End Loop

Definition at line 191 of file GFlashShowerModel.cc.

193 {
194  // std::cout<<"--- ElectronDoit --- "<<std::endl;
195 
196  fastStep.KillPrimaryTrack();
197  fastStep.SetPrimaryTrackPathLength(0.0);
198  fastStep.SetTotalEnergyDeposited(fastTrack.GetPrimaryTrack()->
199  GetKineticEnergy());
200 
201  //-----------------------------
202  // Get track parameters
203  //-----------------------------
204  //E,vect{p} and t,vec(x)
205  G4double Energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
206 
207  // axis of the shower, in global reference frame:
208  G4ThreeVector DirectionShower =
209  fastTrack.GetPrimaryTrack()->GetMomentumDirection();
210  G4ThreeVector OrthoShower, CrossShower;
211  OrthoShower = DirectionShower.orthogonal();
212  CrossShower = DirectionShower.cross(OrthoShower);
213 
214  //--------------------------------
216  //--------------------------------
218  // performance iteration @@@@@@@
219 
221  G4VSolid *SolidCalo = fastTrack.GetEnvelopeSolid();
224  G4double Bound = SolidCalo->DistanceToOut(pos,dir);
225 
226  G4double Dz = 0.00;
227  G4double ZEndStep = 0.00;
228 
229  G4double EnergyNow = Energy;
230  G4double EneIntegral = 0.00;
231  G4double LastEneIntegral = 0.00;
232  G4double DEne = 0.00;
233 
234  G4double NspIntegral = 0.00;
235  G4double LastNspIntegral = 0.00;
236  G4double DNsp = 0.00;
237 
238  // starting point of the shower:
239  G4ThreeVector PositionShower = fastTrack.GetPrimaryTrack()->GetPosition();
240  G4ThreeVector NewPositionShower = PositionShower;
241  G4double StepLenght = 0.00;
242 
243  G4int NSpotDeposited =0;
244 
245  //--------------------------
247  //-------------------------
248 
249  do
250  {
251  //determine step size=min(1Xo,next boundary)
252  G4double stepLength = StepInX0*Parameterisation->GetX0();
253  if(Bound < stepLength)
254  {
255  Dz = Bound;
256  Bound = 0.00;
257  }
258  else
259  {
260  Dz = stepLength;
261  Bound = Bound-Dz;
262  }
263  ZEndStep=ZEndStep+Dz;
264 
265  // Determine Energy Release in Step
266  if(EnergyNow > EnergyStop)
267  {
268  LastEneIntegral = EneIntegral;
269  EneIntegral = Parameterisation->IntegrateEneLongitudinal(ZEndStep);
270  DEne = std::min( EnergyNow,
271  (EneIntegral-LastEneIntegral)*Energy);
272  LastNspIntegral = NspIntegral;
273  NspIntegral = Parameterisation->IntegrateNspLongitudinal(ZEndStep);
274  DNsp = std::max(1., std::floor( (NspIntegral-LastNspIntegral)
275  *Parameterisation->GetNspot() ));
276  }
277  // end of the shower
278  else
279  {
280  DEne = EnergyNow;
281  DNsp = std::max(1., std::floor( (1.- NspIntegral)
282  *Parameterisation->GetNspot() ));
283  }
284  EnergyNow = EnergyNow - DEne;
285 
286  // Apply sampling fluctuation - only in sampling calorimeters
287  //
290  if (sp)
291  {
292  G4double DEneSampling = sp->ApplySampling(DEne,Energy);
293  DEne = DEneSampling;
294  }
295 
296  //move particle in the middle of the step
297  StepLenght = StepLenght + Dz/2.00;
298  NewPositionShower = NewPositionShower +
299  StepLenght*DirectionShower;
300  StepLenght = Dz/2.00;
301 
302  //generate spots & hits:
303  for (int i = 0; i < DNsp; i++)
304  {
305  NSpotDeposited++;
306  GFlashEnergySpot Spot;
307 
308  //Spot energy: the same for all spots
309  Spot.SetEnergy( DEne / DNsp );
310  G4double PhiSpot = Parameterisation->GeneratePhi(); // phi of spot
311  G4double RSpot = Parameterisation // radius of spot
312  ->GenerateRadius(i,Energy,ZEndStep-Dz/2.);
313 
314  // check reference-> may be need to introduce rot matrix @@@
315  // Position: equally spaced in z
316 
317  G4ThreeVector SpotPosition = NewPositionShower +
318  Dz/DNsp*DirectionShower*(i+1/2.-DNsp/2.) +
319  RSpot*std::cos(PhiSpot)*OrthoShower +
320  RSpot*std::sin(PhiSpot)*CrossShower;
321  Spot.SetPosition(SpotPosition);
322 
323  //Generate Hits of this spot
324  HMaker->make(&Spot, &fastTrack);
325  }
326  }
327  while(EnergyNow > 0.0 && Bound> 0.0);
328 
329  //---------------
331  //---------------
332 }
virtual G4double GetNspot()=0
G4VSolid * GetEnvelopeSolid() const
Definition: G4FastTrack.hh:203
virtual G4double GetX0()=0
void SetPosition(const G4ThreeVector &point)
void make(GFlashEnergySpot *aSpot, const G4FastTrack *aT)
TDirectory * dir
int G4int
Definition: G4Types.hh:78
Hep3Vector cross(const Hep3Vector &) const
G4ThreeVector GetPrimaryTrackLocalDirection() const
Definition: G4FastTrack.hh:223
virtual void GenerateLongitudinalProfile(G4double Energy)=0
void SetEnergy(const G4double &E)
Hep3Vector orthogonal() const
G4double ApplySampling(const G4double DEne, const G4double Energy)
void SetPrimaryTrackPathLength(G4double)
GFlashHitMaker * HMaker
G4ThreeVector GetPrimaryTrackLocalPosition() const
Definition: G4FastTrack.hh:213
void KillPrimaryTrack()
Definition: G4FastStep.cc:87
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
virtual G4double IntegrateEneLongitudinal(G4double LongitudinalStep)=0
G4double StepInX0
0=no check ///1=only fully contained...
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
virtual G4double GenerateRadius(G4int ispot, G4double Energy, G4double LongitudinalPosition)=0
void SetTotalEnergyDeposited(G4double anEnergyPart)
static const G4double pos
GVFlashShowerParameterisation * Parameterisation
virtual G4double IntegrateNspLongitudinal(G4double LongitudinalStep)=0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetFlagParamType()

G4int GFlashShowerModel::GetFlagParamType ( )
inline

Definition at line 91 of file GFlashShowerModel.hh.

92  { return FlagParamType; }
Here is the caller graph for this function:

◆ GetFlagParticleContainment()

G4int GFlashShowerModel::GetFlagParticleContainment ( )
inline

Definition at line 93 of file GFlashShowerModel.hh.

94  { return FlagParticleContainment; }
G4int FlagParticleContainment
0=no GFlash 1=only em showers parametrized

◆ GetStepInX0()

G4double GFlashShowerModel::GetStepInX0 ( )
inline

Definition at line 95 of file GFlashShowerModel.hh.

96  { return StepInX0; }
G4double StepInX0
0=no check ///1=only fully contained...

◆ 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 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
Here is the call graph for this function:

◆ 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 }
GFlashParticleBounds * PBound
G4bool CheckParticleDefAndContainment(const G4FastTrack &fastTrack)
G4double GetMaxEneToParametrise(G4ParticleDefinition &particleType)
G4double GetMinEneToParametrise(G4ParticleDefinition &particleType)
virtual void GenerateLongitudinalProfile(G4double Energy)=0
bool G4bool
Definition: G4Types.hh:79
G4double GetEneToKill(G4ParticleDefinition &particleType)
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:208
double G4double
Definition: G4Types.hh:76
GVFlashShowerParameterisation * Parameterisation
Here is the call graph for this function:

◆ SetFlagParamType()

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:

◆ SetFlagParticleContainment()

void GFlashShowerModel::SetFlagParticleContainment ( G4int  I)
inline

Definition at line 78 of file GFlashShowerModel.hh.

G4int FlagParticleContainment
0=no GFlash 1=only em showers parametrized
Here is the caller graph for this function:

◆ SetHitMaker()

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker Maker)
inline

Definition at line 84 of file GFlashShowerModel.hh.

85  { HMaker=&Maker; }
GFlashHitMaker * HMaker
Here is the caller graph for this function:

◆ SetParameterisation()

void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation DP)
inline

Definition at line 82 of file GFlashShowerModel.hh.

83  { Parameterisation=&DP;}
GVFlashShowerParameterisation * Parameterisation
Here is the caller graph for this function:

◆ SetParticleBounds()

void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds SpecificBound)
inline

Definition at line 86 of file GFlashShowerModel.hh.

87  { PBound =&SpecificBound; }
GFlashParticleBounds * PBound
Here is the caller graph for this function:

◆ SetStepInX0()

void GFlashShowerModel::SetStepInX0 ( G4double  Lenght)
inline

Definition at line 80 of file GFlashShowerModel.hh.

81  { StepInX0=Lenght; }
G4double StepInX0
0=no check ///1=only fully contained...
Here is the caller graph for this function:

Member Data Documentation

◆ EnergyStop

G4double GFlashShowerModel::EnergyStop
private

Definition at line 121 of file GFlashShowerModel.hh.

◆ FlagParamType

G4int GFlashShowerModel::FlagParamType
private

Definition at line 118 of file GFlashShowerModel.hh.

◆ FlagParticleContainment

G4int GFlashShowerModel::FlagParticleContainment
private

0=no GFlash 1=only em showers parametrized

Definition at line 119 of file GFlashShowerModel.hh.

◆ HMaker

GFlashHitMaker* GFlashShowerModel::HMaker
private

Definition at line 114 of file GFlashShowerModel.hh.

◆ Messenger

GFlashShowerModelMessenger* GFlashShowerModel::Messenger
private

Definition at line 115 of file GFlashShowerModel.hh.

◆ Parameterisation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

Definition at line 102 of file GFlashShowerModel.hh.

◆ PBound

GFlashParticleBounds* GFlashShowerModel::PBound

Definition at line 101 of file GFlashShowerModel.hh.

◆ StepInX0

G4double GFlashShowerModel::StepInX0
private

0=no check ///1=only fully contained...

Definition at line 120 of file GFlashShowerModel.hh.


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