Geant4 11.2.2
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 EnergyStop = 0.0;
65 Messenger = new GFlashShowerModelMessenger(this);
66}
67
69 : G4VFastSimulationModel(modelName),
70 PBound(0), Parameterisation(0), HMaker(0)
71{
72 FlagParamType =1;
73 FlagParticleContainment = 1;
74 StepInX0 = 0.1;
75 EnergyStop = 0.0;
76 Messenger = new GFlashShowerModelMessenger(this);
77}
78
80{
81 delete Messenger;
82}
83
86{
87 return
88 &particleType == G4Electron::ElectronDefinition() ||
89 &particleType == G4Positron::PositronDefinition();
90}
91
92/**********************************************************************/
93/* Checks whether conditions of fast parameterisation are fullfilled */
94/**********************************************************************/
95
97
98{
99 G4bool select = false;
100 if(FlagParamType != 0)
101 {
102 G4double ParticleEnergy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
103 G4ParticleDefinition &ParticleType =
104 *(fastTrack.GetPrimaryTrack()->GetDefinition());
105 if(ParticleEnergy > PBound->GetMinEneToParametrise(ParticleType) &&
106 ParticleEnergy < PBound->GetMaxEneToParametrise(ParticleType) )
107 {
108 // check conditions depending on particle flavour
109 // performance to be optimized @@@@@@@
111 select = CheckParticleDefAndContainment(fastTrack);
112 if (select) EnergyStop= PBound->GetEneToKill(ParticleType);
113 }
114 }
115
116 return select;
117}
118
119
120G4bool
121GFlashShowerModel::CheckParticleDefAndContainment(const G4FastTrack& fastTrack)
122{
123 G4bool filter=false;
124 G4ParticleDefinition * ParticleType =
125 fastTrack.GetPrimaryTrack()->GetDefinition();
126
127 if( ParticleType == G4Electron::ElectronDefinition() ||
128 ParticleType == G4Positron::PositronDefinition() )
129 {
130 filter=true;
131 if(FlagParticleContainment == 1)
132 {
133 filter=CheckContainment(fastTrack);
134 }
135 }
136 return filter;
137}
138
139G4bool GFlashShowerModel::CheckContainment(const G4FastTrack& fastTrack)
140{
141 G4bool filter=false;
142 // track informations
143 G4ThreeVector DirectionShower=fastTrack.GetPrimaryTrackLocalDirection();
144 G4ThreeVector InitialPositionShower=fastTrack.GetPrimaryTrackLocalPosition();
145
146 G4ThreeVector OrthoShower, CrossShower;
147 // Returns orthogonal vector
148 OrthoShower = DirectionShower.orthogonal();
149 // Shower in direction perpendicular to OrthoShower and DirectionShower
150 CrossShower = DirectionShower.cross(OrthoShower);
151
154 G4int CosPhi[4] = {1,0,-1,0};
155 G4int SinPhi[4] = {0,1,0,-1};
156
158 G4int NlateralInside=0;
159 // pointer to solid we're in
160 G4VSolid *SolidCalo = fastTrack.GetEnvelopeSolid();
161 for(int i=0; i<4 ;i++)
162 {
163 // polar coordinates
164 Position = InitialPositionShower +
165 Z*DirectionShower +
166 R*CosPhi[i]*OrthoShower +
167 R*SinPhi[i]*CrossShower ;
168
169 if(SolidCalo->Inside(Position) != kOutside)
170 NlateralInside++;
171 }
172
173 // choose to parameterise or flag when all inetc...
174 if(NlateralInside==4) filter=true;
175 // std::cout << " points = " <<NlateralInside << std::endl;
176 return filter;
177}
178
179
180void
182{
183 // parametrise electrons
184 if(fastTrack.GetPrimaryTrack()->GetDefinition()
186 fastTrack.GetPrimaryTrack()->GetDefinition()
188 ElectronDoIt(fastTrack,fastStep);
189}
190
191void
192GFlashShowerModel::ElectronDoIt(const G4FastTrack& fastTrack,
193 G4FastStep& fastStep)
194{
195 // std::cout<<"--- ElectronDoit --- "<<std::endl;
196
197 fastStep.KillPrimaryTrack();
198 fastStep.ProposePrimaryTrackPathLength(0.0);
199 fastStep.ProposeTotalEnergyDeposited(fastTrack.GetPrimaryTrack()->
200 GetKineticEnergy());
201
202 //-----------------------------
203 // Get track parameters
204 //-----------------------------
205 //E,vect{p} and t,vec(x)
206 G4double Energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy();
207
208 // axis of the shower, in global reference frame:
209 G4ThreeVector DirectionShower =
211 G4ThreeVector OrthoShower, CrossShower;
212 OrthoShower = DirectionShower.orthogonal();
213 CrossShower = DirectionShower.cross(OrthoShower);
214
215 //--------------------------------
216 ///Generate longitudinal profile
217 //--------------------------------
219 // performance iteration @@@@@@@
220
221 ///Initialisation of long. loop variables
222 G4VSolid *SolidCalo = fastTrack.GetEnvelopeSolid();
225 G4double Bound = SolidCalo->DistanceToOut(pos,dir);
226
227 G4double Dz = 0.00;
228 G4double ZEndStep = 0.00;
229
230 G4double EnergyNow = Energy;
231 G4double EneIntegral = 0.00;
232 G4double LastEneIntegral = 0.00;
233 G4double DEne = 0.00;
234
235 G4double NspIntegral = 0.00;
236 G4double LastNspIntegral = 0.00;
237 G4double DNsp = 0.00;
238
239 // starting point of the shower:
240 G4ThreeVector PositionShower = fastTrack.GetPrimaryTrack()->GetPosition();
241 G4ThreeVector NewPositionShower = PositionShower;
242 G4double StepLenght = 0.00;
243
244 //--------------------------
245 /// Begin Longitudinal Loop
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 (G4int i = 0; i < DNsp; ++i)
303 {
304 GFlashEnergySpot Spot;
305
306 //Spot energy: the same for all spots
307 Spot.SetEnergy( DEne / DNsp );
308 G4double PhiSpot = Parameterisation->GeneratePhi(); // phi of spot
309 G4double RSpot = Parameterisation // radius of spot
310 ->GenerateRadius(i,Energy,ZEndStep-Dz/2.);
311
312 // check reference-> may be need to introduce rot matrix @@@
313 // Position: equally spaced in z
314
315 G4ThreeVector SpotPosition = NewPositionShower +
316 Dz/DNsp*DirectionShower*(i+1/2.-DNsp/2.) +
317 RSpot*std::cos(PhiSpot)*OrthoShower +
318 RSpot*std::sin(PhiSpot)*CrossShower;
319 Spot.SetPosition(SpotPosition);
320
321 //Generate Hits of this spot
322 HMaker->make(&Spot, &fastTrack);
323 }
324 }
325 while(EnergyNow > 0.0 && Bound> 0.0);
326
327 //---------------
328 /// End Loop
329 //---------------
330}
331
332/*
333
334void
335GFlashShowerModel::GammaDoIt(const G4FastTrack& fastTrack,
336 G4FastStep& fastStep)
337{
338
339 if( fastTrack.GetPrimaryTrack()->GetKineticEnergy() > EnergyStop )
340 return;
341
342 //deposita in uno spot unico l'energia
343 //con andamento exp decrescente.
344
345 // Kill the particle to be parametrised
346 fastStep.KillPrimaryTrack();
347 fastStep.SetPrimaryTrackPathLength(0.0);
348 fastStep.SetTotalEnergyDeposited(fastTrack.GetPrimaryTrack()
349 ->GetKineticEnergy());
350 // other settings????
351 feSpotList.clear();
352
353 //-----------------------------
354 // Get track parameters
355 //-----------------------------
356
357 // E,vect{p} and t,vec(x)
358 G4double Energy =
359 fastTrack.GetPrimaryTrack()->GetKineticEnergy();
360 // axis of the shower, in global reference frame:
361 G4ThreeVector DirectionShower =
362 fastTrack.GetPrimaryTrack()->GetMomentumDirection();
363 // starting point of the shower:
364 G4ThreeVector PositionShower =
365 fastTrack.GetPrimaryTrack()->GetPosition();
366
367 //G4double DEneSampling = Parameterisation->ApplySampling(Energy,Energy);
368 //if(DEneSampling <= 0.00) DEneSampling=Energy;
369
370 if(Energy > 0.0)
371 {
372 G4double dist = Parameterisation->GenerateExponential(Energy);
373
374 GFlashEnergySpot Spot;
375 Spot.SetEnergy( Energy );
376 G4ThreeVector SpotPosition = PositionShower + dist*DirectionShower;
377 Spot.SetPosition(SpotPosition);
378
379 // Record the Spot:
380 feSpotList.push_back(Spot);
381
382 //Generate Hits of this spot
383 HMaker->make(Spot);
384 }
385}
386
387*/
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
void KillPrimaryTrack()
Definition G4FastStep.cc:87
void ProposeTotalEnergyDeposited(G4double anEnergyPart)
void ProposePrimaryTrackPathLength(G4double)
G4ThreeVector GetPrimaryTrackLocalPosition() const
const G4Track * GetPrimaryTrack() const
G4ThreeVector GetPrimaryTrackLocalDirection() const
G4VSolid * GetEnvelopeSolid() const
static G4Positron * PositronDefinition()
Definition G4Positron.cc:85
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 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