Geant4-11
Public Member Functions | Data Fields | Private Member Functions | Private Attributes
GFlashShowerModel Class Reference

#include <GFlashShowerModel.hh>

Inheritance diagram for GFlashShowerModel:
G4VFastSimulationModel

Public Member Functions

virtual void AtRestDoIt (const G4FastTrack &, G4FastStep &)
 
virtual G4bool AtRestModelTrigger (const G4FastTrack &)
 
void DoIt (const G4FastTrack &, G4FastStep &)
 
G4int GetFlagParamType ()
 
G4int GetFlagParticleContainment ()
 
const G4String GetName () const
 
G4double GetStepInX0 ()
 
 GFlashShowerModel (G4String)
 
 GFlashShowerModel (G4String, G4Envelope *)
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
G4bool ModelTrigger (const G4FastTrack &)
 
G4bool operator== (const G4VFastSimulationModel &) const
 
void SetFlagParamType (G4int I)
 
void SetFlagParticleContainment (G4int I)
 
void SetHitMaker (GFlashHitMaker &Maker)
 
void SetParameterisation (GVFlashShowerParameterisation &DP)
 
void SetParticleBounds (GFlashParticleBounds &SpecificBound)
 
void SetStepInX0 (G4double Lenght)
 
 ~GFlashShowerModel ()
 

Data Fields

GVFlashShowerParameterisationParameterisation
 
GFlashParticleBoundsPBound
 

Private Member Functions

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

Private Attributes

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

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

References FlagParamType, FlagParticleContainment, Messenger, and StepInX0.

◆ GFlashShowerModel() [2/2]

GFlashShowerModel::GFlashShowerModel ( G4String  modelName)

Definition at line 67 of file GFlashShowerModel.cc.

68 : G4VFastSimulationModel(modelName),
70{
73 StepInX0 = 0.1;
75}

References FlagParamType, FlagParticleContainment, Messenger, and StepInX0.

◆ ~GFlashShowerModel()

GFlashShowerModel::~GFlashShowerModel ( )

Definition at line 77 of file GFlashShowerModel.cc.

78{
79 delete Messenger;
80}

References Messenger.

Member Function Documentation

◆ AtRestDoIt()

virtual void G4VFastSimulationModel::AtRestDoIt ( const G4FastTrack ,
G4FastStep  
)
inlinevirtualinherited

Definition at line 126 of file G4VFastSimulationModel.hh.

126{}

Referenced by G4FastSimulationManager::InvokeAtRestDoIt().

◆ AtRestModelTrigger()

virtual G4bool G4VFastSimulationModel::AtRestModelTrigger ( const G4FastTrack )
inlinevirtualinherited

Definition at line 115 of file G4VFastSimulationModel.hh.

115{return false;}

◆ CheckContainment()

G4bool GFlashShowerModel::CheckContainment ( const G4FastTrack fastTrack)
private

Definition at line 137 of file GFlashShowerModel.cc.

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

References CLHEP::Hep3Vector::cross(), GVFlashShowerParameterisation::GetAveR90(), GVFlashShowerParameterisation::GetAveT90(), G4FastTrack::GetEnvelopeSolid(), G4FastTrack::GetPrimaryTrackLocalDirection(), G4FastTrack::GetPrimaryTrackLocalPosition(), G4VSolid::Inside(), kOutside, CLHEP::Hep3Vector::orthogonal(), Parameterisation, and Z.

Referenced by CheckParticleDefAndContainment().

◆ CheckParticleDefAndContainment()

G4bool GFlashShowerModel::CheckParticleDefAndContainment ( const G4FastTrack fastTrack)
private

Definition at line 119 of file GFlashShowerModel.cc.

120{
121 G4bool filter=false;
123 fastTrack.GetPrimaryTrack()->GetDefinition();
124
127 {
128 filter=true;
130 {
131 filter=CheckContainment(fastTrack);
132 }
133 }
134 return filter;
135}
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
G4bool CheckContainment(const G4FastTrack &fastTrack)

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

Referenced by ModelTrigger().

◆ 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}
void ElectronDoIt(const G4FastTrack &, G4FastStep &)

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

◆ 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 190 of file GFlashShowerModel.cc.

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

References CLHEP::Hep3Vector::cross(), G4VSolid::DistanceToOut(), EnergyStop, GVFlashShowerParameterisation::GenerateLongitudinalProfile(), GVFlashShowerParameterisation::GeneratePhi(), GVFlashShowerParameterisation::GenerateRadius(), G4FastTrack::GetEnvelopeSolid(), G4Track::GetKineticEnergy(), G4Track::GetMomentumDirection(), GVFlashShowerParameterisation::GetNspot(), G4Track::GetPosition(), G4FastTrack::GetPrimaryTrack(), G4FastTrack::GetPrimaryTrackLocalDirection(), G4FastTrack::GetPrimaryTrackLocalPosition(), GVFlashShowerParameterisation::GetX0(), HMaker, GVFlashShowerParameterisation::IntegrateEneLongitudinal(), GVFlashShowerParameterisation::IntegrateNspLongitudinal(), G4FastStep::KillPrimaryTrack(), GFlashHitMaker::make(), G4INCL::Math::max(), G4INCL::Math::min(), CLHEP::Hep3Vector::orthogonal(), Parameterisation, pos, GFlashEnergySpot::SetEnergy(), GFlashEnergySpot::SetPosition(), G4FastStep::SetPrimaryTrackPathLength(), G4FastStep::SetTotalEnergyDeposited(), G4InuclParticleNames::sp, and StepInX0.

Referenced by DoIt().

◆ GetFlagParamType()

G4int GFlashShowerModel::GetFlagParamType ( )
inline

Definition at line 90 of file GFlashShowerModel.hh.

91 { return FlagParamType; }

References FlagParamType.

Referenced by GFlashShowerModelMessenger::GetCurrentValue().

◆ GetFlagParticleContainment()

G4int GFlashShowerModel::GetFlagParticleContainment ( )
inline

Definition at line 92 of file GFlashShowerModel.hh.

93 { return FlagParticleContainment; }

References FlagParticleContainment.

◆ GetName()

const G4String G4VFastSimulationModel::GetName ( ) const
inlineinherited

Definition at line 146 of file G4VFastSimulationModel.hh.

147{
148 return theModelName;
149}

References G4VFastSimulationModel::theModelName.

◆ GetStepInX0()

G4double GFlashShowerModel::GetStepInX0 ( )
inline

Definition at line 94 of file GFlashShowerModel.hh.

95 { return StepInX0; }

References 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}

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

◆ 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);
111 }
112 }
113
114 return select;
115}
G4double GetMinEneToParametrise(G4ParticleDefinition &particleType)
G4double GetMaxEneToParametrise(G4ParticleDefinition &particleType)
G4double GetEneToKill(G4ParticleDefinition &particleType)
G4bool CheckParticleDefAndContainment(const G4FastTrack &fastTrack)

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

◆ operator==()

G4bool G4VFastSimulationModel::operator== ( const G4VFastSimulationModel fsm) const
inlineinherited

Definition at line 152 of file G4VFastSimulationModel.hh.

153{
154 return (this==&fsm) ? true : false;
155}

◆ SetFlagParamType()

void GFlashShowerModel::SetFlagParamType ( G4int  I)
inline

Definition at line 75 of file GFlashShowerModel.hh.

76 { FlagParamType = I; }

References FlagParamType.

Referenced by GFlashShowerModelMessenger::SetNewValue().

◆ SetFlagParticleContainment()

void GFlashShowerModel::SetFlagParticleContainment ( G4int  I)
inline

Definition at line 77 of file GFlashShowerModel.hh.

References FlagParticleContainment.

Referenced by GFlashShowerModelMessenger::SetNewValue().

◆ SetHitMaker()

void GFlashShowerModel::SetHitMaker ( GFlashHitMaker Maker)
inline

Definition at line 83 of file GFlashShowerModel.hh.

84 { HMaker=&Maker; }

References HMaker.

◆ SetParameterisation()

void GFlashShowerModel::SetParameterisation ( GVFlashShowerParameterisation DP)
inline

Definition at line 81 of file GFlashShowerModel.hh.

82 { Parameterisation=&DP;}

References Parameterisation.

◆ SetParticleBounds()

void GFlashShowerModel::SetParticleBounds ( GFlashParticleBounds SpecificBound)
inline

Definition at line 85 of file GFlashShowerModel.hh.

86 { PBound =&SpecificBound; }

References PBound.

◆ SetStepInX0()

void GFlashShowerModel::SetStepInX0 ( G4double  Lenght)
inline

Definition at line 79 of file GFlashShowerModel.hh.

80 { StepInX0=Lenght; }

References StepInX0.

Referenced by GFlashShowerModelMessenger::SetNewValue().

Field Documentation

◆ EnergyStop

G4double GFlashShowerModel::EnergyStop
private

Definition at line 120 of file GFlashShowerModel.hh.

Referenced by ElectronDoIt(), and ModelTrigger().

◆ FlagParamType

G4int GFlashShowerModel::FlagParamType
private

◆ FlagParticleContainment

G4int GFlashShowerModel::FlagParticleContainment
private

0=no GFlash 1=only em showers parametrized

Definition at line 118 of file GFlashShowerModel.hh.

Referenced by CheckParticleDefAndContainment(), GetFlagParticleContainment(), GFlashShowerModel(), and SetFlagParticleContainment().

◆ HMaker

GFlashHitMaker* GFlashShowerModel::HMaker
private

Definition at line 113 of file GFlashShowerModel.hh.

Referenced by ElectronDoIt(), and SetHitMaker().

◆ Messenger

GFlashShowerModelMessenger* GFlashShowerModel::Messenger
private

Definition at line 114 of file GFlashShowerModel.hh.

Referenced by GFlashShowerModel(), and ~GFlashShowerModel().

◆ Parameterisation

GVFlashShowerParameterisation* GFlashShowerModel::Parameterisation

◆ PBound

GFlashParticleBounds* GFlashShowerModel::PBound

◆ StepInX0

G4double GFlashShowerModel::StepInX0
private

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

Definition at line 119 of file GFlashShowerModel.hh.

Referenced by ElectronDoIt(), GetStepInX0(), GFlashShowerModel(), and SetStepInX0().

◆ theModelName

G4String G4VFastSimulationModel::theModelName
privateinherited

Definition at line 143 of file G4VFastSimulationModel.hh.

Referenced by G4VFastSimulationModel::GetName().


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