Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
GFlashShowerModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28// ------------------------------------------------------------
29// GEANT 4 class implementation
30//
31// ---------------- GFlashShowerModel ----------------
32//
33// Authors: E.Barberio & Joanna Weng - 9.11.2004
34// ------------------------------------------------------------
35
36#include "G4Electron.hh"
37#include "G4Positron.hh"
38#include "G4NeutrinoE.hh"
39#include "G4NeutrinoMu.hh"
40#include "G4NeutrinoTau.hh"
41#include "G4AntiNeutrinoE.hh"
42#include "G4AntiNeutrinoMu.hh"
43#include "G4AntiNeutrinoTau.hh"
44#include "G4PionZero.hh"
45#include "G4VProcess.hh"
46#include "G4ios.hh"
47#include "G4LogicalVolume.hh"
48#include "geomdefs.hh"
49
50#include "GFlashShowerModel.hh"
53#include "GFlashEnergySpot.hh"
54
55
57 G4Envelope* envelope)
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}
66
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}
76
78{
79 delete Messenger;
80}
81
84{
85 return
86 &particleType == G4Electron::ElectronDefinition() ||
87 &particleType == G4Positron::PositronDefinition();
88}
89
90/**********************************************************************/
91/* Checks whether conditions of fast parameterisation are fullfilled */
92/**********************************************************************/
93
95
96{
97 G4bool select = false;
98 if(FlagParamType != 0)
99 {
100 G4double ParticleEnergy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
101 G4ParticleDefinition &ParticleType =
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}
116
117
118G4bool
119GFlashShowerModel::CheckParticleDefAndContainment(const G4FastTrack& fastTrack)
120{
121 G4bool filter=false;
122 G4ParticleDefinition * ParticleType =
123 fastTrack.GetPrimaryTrack()->GetDefinition();
124
125 if( ParticleType == G4Electron::ElectronDefinition() ||
126 ParticleType == G4Positron::PositronDefinition() )
127 {
128 filter=true;
129 if(FlagParticleContainment == 1)
130 {
131 filter=CheckContainment(fastTrack);
132 }
133 }
134 return filter;
135}
136
137G4bool GFlashShowerModel::CheckContainment(const G4FastTrack& fastTrack)
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}
176
177
178void
180{
181 // parametrise electrons
182 if(fastTrack.GetPrimaryTrack()->GetDefinition()
184 fastTrack.GetPrimaryTrack()->GetDefinition()
186 ElectronDoIt(fastTrack,fastStep);
187}
188
189void
190GFlashShowerModel::ElectronDoIt(const G4FastTrack& fastTrack,
191 G4FastStep& fastStep)
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 //--------------------------------
214 ///Generate longitudinal profile
215 //--------------------------------
217 // performance iteration @@@@@@@
218
219 ///Initialisation of long. loop variables
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 //--------------------------
243 /// Begin Longitudinal Loop
244 //-------------------------
245
246 do
247 {
248 //determine step size=min(1Xo,next boundary)
249 G4double stepLength = StepInX0*Parameterisation->GetX0();
250 if(Bound < stepLength)
251 {
252 Dz = Bound;
253 Bound = 0.00;
254 }
255 else
256 {
257 Dz = stepLength;
258 Bound = Bound-Dz;
259 }
260 ZEndStep=ZEndStep+Dz;
261
262 // Determine Energy Release in Step
263 if(EnergyNow > EnergyStop)
264 {
265 LastEneIntegral = EneIntegral;
266 EneIntegral = Parameterisation->IntegrateEneLongitudinal(ZEndStep);
267 DEne = std::min( EnergyNow,
268 (EneIntegral-LastEneIntegral)*Energy);
269 LastNspIntegral = NspIntegral;
270 NspIntegral = Parameterisation->IntegrateNspLongitudinal(ZEndStep);
271 DNsp = std::max(1., std::floor( (NspIntegral-LastNspIntegral)
273 }
274 // end of the shower
275 else
276 {
277 DEne = EnergyNow;
278 DNsp = std::max(1., std::floor( (1.- NspIntegral)
280 }
281 EnergyNow = EnergyNow - DEne;
282
283 // Apply sampling fluctuation - only in sampling calorimeters
284 //
287 if (sp)
288 {
289 G4double DEneSampling = sp->ApplySampling(DEne,Energy);
290 DEne = DEneSampling;
291 }
292
293 //move particle in the middle of the step
294 StepLenght = StepLenght + Dz/2.00;
295 NewPositionShower = NewPositionShower +
296 StepLenght*DirectionShower;
297 StepLenght = Dz/2.00;
298
299 //generate spots & hits:
300 for (G4int i = 0; i < DNsp; ++i)
301 {
302 GFlashEnergySpot Spot;
303
304 //Spot energy: the same for all spots
305 Spot.SetEnergy( DEne / DNsp );
306 G4double PhiSpot = Parameterisation->GeneratePhi(); // phi of spot
307 G4double RSpot = Parameterisation // radius of spot
308 ->GenerateRadius(i,Energy,ZEndStep-Dz/2.);
309
310 // check reference-> may be need to introduce rot matrix @@@
311 // Position: equally spaced in z
312
313 G4ThreeVector SpotPosition = NewPositionShower +
314 Dz/DNsp*DirectionShower*(i+1/2.-DNsp/2.) +
315 RSpot*std::cos(PhiSpot)*OrthoShower +
316 RSpot*std::sin(PhiSpot)*CrossShower;
317 Spot.SetPosition(SpotPosition);
318
319 //Generate Hits of this spot
320 HMaker->make(&Spot, &fastTrack);
321 }
322 }
323 while(EnergyNow > 0.0 && Bound> 0.0);
324
325 //---------------
326 /// End Loop
327 //---------------
328}
329
330/*
331
332void
333GFlashShowerModel::GammaDoIt(const G4FastTrack& fastTrack,
334 G4FastStep& fastStep)
335{
336
337 if( fastTrack.GetPrimaryTrack()->GetKineticEnergy() > EnergyStop )
338 return;
339
340 //deposita in uno spot unico l'energia
341 //con andamento exp decrescente.
342
343 // Kill the particle to be parametrised
344 fastStep.KillPrimaryTrack();
345 fastStep.SetPrimaryTrackPathLength(0.0);
346 fastStep.SetTotalEnergyDeposited(fastTrack.GetPrimaryTrack()
347 ->GetKineticEnergy());
348 // other settings????
349 feSpotList.clear();
350
351 //-----------------------------
352 // Get track parameters
353 //-----------------------------
354
355 // E,vect{p} and t,vec(x)
356 G4double Energy =
357 fastTrack.GetPrimaryTrack()->GetKineticEnergy();
358 // axis of the shower, in global reference frame:
359 G4ThreeVector DirectionShower =
360 fastTrack.GetPrimaryTrack()->GetMomentumDirection();
361 // starting point of the shower:
362 G4ThreeVector PositionShower =
363 fastTrack.GetPrimaryTrack()->GetPosition();
364
365 //G4double DEneSampling = Parameterisation->ApplySampling(Energy,Energy);
366 //if(DEneSampling <= 0.00) DEneSampling=Energy;
367
368 if(Energy > 0.0)
369 {
370 G4double dist = Parameterisation->GenerateExponential(Energy);
371
372 GFlashEnergySpot Spot;
373 Spot.SetEnergy( Energy );
374 G4ThreeVector SpotPosition = PositionShower + dist*DirectionShower;
375 Spot.SetPosition(SpotPosition);
376
377 // Record the Spot:
378 feSpotList.push_back(Spot);
379
380 //Generate Hits of this spot
381 HMaker->make(Spot);
382 }
383}
384
385*/
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
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
void SetTotalEnergyDeposited(G4double anEnergyPart)
void KillPrimaryTrack()
Definition: G4FastStep.cc:87
void SetPrimaryTrackPathLength(G4double)
G4ThreeVector GetPrimaryTrackLocalPosition() const
Definition: G4FastTrack.hh:211
const G4Track * GetPrimaryTrack() const
Definition: G4FastTrack.hh:206
G4ThreeVector GetPrimaryTrackLocalDirection() const
Definition: G4FastTrack.hh:221
G4VSolid * GetEnvelopeSolid() const
Definition: G4FastTrack.hh:201
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:88
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetDefinition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
virtual EInside Inside(const G4ThreeVector &p) const =0
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)
G4double GetMinEneToParametrise(G4ParticleDefinition &particleType)
G4double GetMaxEneToParametrise(G4ParticleDefinition &particleType)
G4double GetEneToKill(G4ParticleDefinition &particleType)
void DoIt(const G4FastTrack &, G4FastStep &)
G4bool ModelTrigger(const G4FastTrack &)
GFlashShowerModel(G4String, G4Envelope *)
G4bool IsApplicable(const G4ParticleDefinition &)
GFlashParticleBounds * PBound
GVFlashShowerParameterisation * Parameterisation
virtual G4double IntegrateEneLongitudinal(G4double LongitudinalStep)=0
virtual G4double GetAveR90()=0
virtual void GenerateLongitudinalProfile(G4double Energy)=0
virtual G4double GenerateRadius(G4int ispot, G4double Energy, G4double LongitudinalPosition)=0
virtual G4double GetAveT90()=0
virtual G4double GetX0()=0
virtual G4double IntegrateNspLongitudinal(G4double LongitudinalStep)=0
virtual G4double GetNspot()=0
@ kOutside
Definition: geomdefs.hh:68