Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4SpecialCuts.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Geant4 Virtual Monte Carlo package
3// Copyright (C) 2007 - 2014 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
14
15#include "TG4SpecialCuts.h"
16#include "TG4Limits.h"
17
18//
19// Class TG4SpecialCutsForChargedHadron implementation
20//
21
22//_____________________________________________________________________________
24 const G4String& processName)
25 : TG4VSpecialCuts(processName)
26{
28}
29
30//_____________________________________________________________________________
35
36//_____________________________________________________________________________
38 const TG4Limits& limits, const G4Track& track) const
39{
41 // ---
42
43 return limits.GetMinEkineForChargedHadron(track);
44}
45
46//
47// Class TG4SpecialCutsForElectron implementation
48//
49
50//_____________________________________________________________________________
52 const G4String& processName)
53 : TG4VSpecialCuts(processName)
54{
56}
57
58//_____________________________________________________________________________
63
64//_____________________________________________________________________________
66 const TG4Limits& limits, const G4Track& track) const
67{
69
70 return limits.GetMinEkineForElectron(track);
71}
72
73//
74// Class TG4SpecialCutsForGamma implementation
75//
76
77//_____________________________________________________________________________
79 : TG4VSpecialCuts(processName)
80{
82}
83
84//_____________________________________________________________________________
89
90//_____________________________________________________________________________
92 const TG4Limits& limits, const G4Track& track) const
93{
95
96 return limits.GetMinEkineForGamma(track);
97}
98
99//
100// Class TG4SpecialCutsForMuon implementation
101//
102
103//_____________________________________________________________________________
105 : TG4VSpecialCuts(processName)
106{
108}
109
110//_____________________________________________________________________________
115
116//_____________________________________________________________________________
118 const TG4Limits& limits, const G4Track& track) const
119{
121
122 return limits.GetMinEkineForMuon(track);
123}
124
125//
126// Class TG4SpecialCutsForNeutralHadron implementation
127//
128
129//_____________________________________________________________________________
131 const G4String& processName)
132 : TG4VSpecialCuts(processName)
133{
135}
136
137//_____________________________________________________________________________
142
143//_____________________________________________________________________________
145 const TG4Limits& limits, const G4Track& track) const
146{
148
149 return limits.GetMinEkineForNeutralHadron(track);
150}
151
152//
153// Class TG4SpecialCutsForNeutralHadron implementation
154//
155
156//_____________________________________________________________________________
158 : TG4VSpecialCuts(processName)
159{
161}
162
163//_____________________________________________________________________________
168
169//_____________________________________________________________________________
171 const TG4Limits& limits, const G4Track& track) const
172{
174
175 return limits.GetMinEkineForNeutralHadron(track);
176}
177
178//_____________________________________________________________________________
180 const G4Track& track, const G4Step& /*step*/)
181{
184
185 aParticleChange.Initialize(track);
186 aParticleChange.ProposeEnergy(0.);
187 aParticleChange.ProposeLocalEnergyDeposit(track.GetKineticEnergy());
188 aParticleChange.ProposeTrackStatus(fStopAndKill);
189
190 return &aParticleChange;
191}
Definition of the TG4Limits class.
Definition of the TG4SpecialCutsFor* classes.
Extended G4UserLimits class.
Definition TG4Limits.h:38
G4double GetMinEkineForNeutralHadron(const G4Track &track) const
G4double GetMinEkineForMuon(const G4Track &track) const
G4double GetMinEkineForElectron(const G4Track &track) const
G4double GetMinEkineForGamma(const G4Track &track) const
G4double GetMinEkineForChargedHadron(const G4Track &track) const
TG4SpecialCutsForChargedHadron(const G4String &processName="specialCutForChargedHadron")
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const
Return the kinetic energy limit.
TG4SpecialCutsForElectron(const G4String &processName="specialCutForElectron")
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const
Return the kinetic energy limit.
TG4SpecialCutsForGamma(const G4String &processName="specialCutForGamma")
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const
Return the kinetic energy limit.
TG4SpecialCutsForMuon(const G4String &processName="specialCutForMuon")
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const
Return the kinetic energy limit.
TG4SpecialCutsForNeutralHadron(const G4String &processName="specialCutForNeutralHadron")
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const
Return the kinetic energy limit.
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
virtual G4double GetMinEkine(const TG4Limits &limits, const G4Track &track) const
Return the kinetic energy limit.
TG4SpecialCutsForNeutron(const G4String &processName="specialCutForNeutron")
Abstract base class for a special process that activates kinetic energy cuts.