Geant4_10
G4ParticleChangeForLoss.hh
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 // $Id: G4ParticleChangeForLoss.hh 68795 2013-04-05 13:24:46Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 //
34 // ------------------------------------------------------------
35 // Implemented for the new scheme 23 Mar. 1998 H.Kurahige
36 //
37 // Modified:
38 // 16.01.04 V.Ivanchenko update for model variant of energy loss
39 // 15.04.05 V.Ivanchenko inline update methods
40 // 30.01.06 V.Ivanchenko add ProposedMomentumDirection for AlongStep
41 // and ProposeWeight for PostStep
42 // 07.06.06 V.Ivanchenko RemoveProposedMomentumDirection from AlongStep
43 // 28.08.06 V.Ivanchenko Added access to current track and polarizaion
44 // 17.06.09 V.Ivanchenko Added SetLowEnergyLimit method
45 //
46 // ------------------------------------------------------------
47 //
48 // Class Description
49 // This class is a concrete class for ParticleChange for EnergyLoss
50 //
51 #ifndef G4ParticleChangeForLoss_h
52 #define G4ParticleChangeForLoss_h 1
53 
54 #include "globals.hh"
55 #include "G4ios.hh"
56 #include "G4VParticleChange.hh"
57 
58 class G4DynamicParticle;
59 
61 {
62 public:
63  // default constructor
65 
66  // destructor
67  virtual ~G4ParticleChangeForLoss();
68 
69  // with description
70  // ----------------------------------------------------
71  // --- the following methods are for updating G4Step -----
72 
75  // A physics process gives the final state of the particle
76  // based on information of G4Track
77 
78  void InitializeForAlongStep(const G4Track&);
79  void InitializeForPostStep(const G4Track&);
80  //Initialize all propoerties by using G4Track information
81 
82  // void AddSecondary(G4DynamicParticle* aParticle);
83  // Add next secondary
84 
85  inline G4double GetProposedCharge() const;
86  inline void SetProposedCharge(G4double theCharge);
87  // Get/Set theCharge
88 
89  inline G4double GetCharge() const;
90  inline void ProposeCharge(G4double finalCharge);
91  // Get/Propose the final dynamical Charge in G4DynamicParticle
92 
93  inline G4double GetProposedKineticEnergy() const;
94  inline void SetProposedKineticEnergy(G4double proposedKinEnergy);
95  // Get/Set the final kinetic energy of the current particle.
96 
97  inline const G4ThreeVector& GetProposedMomentumDirection() const;
99  inline const G4ThreeVector& GetMomentumDirection() const;
100  inline void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz);
101  inline void ProposeMomentumDirection(const G4ThreeVector& Pfinal);
102  // Get/Propose the MomentumDirection vector: it is the final momentum direction.
103 
104  inline const G4ThreeVector& GetProposedPolarization() const;
105  inline void ProposePolarization(const G4ThreeVector& dir);
106  inline void ProposePolarization(G4double Px, G4double Py, G4double Pz);
107 
108  inline const G4Track* GetCurrentTrack() const;
109 
110  inline void SetLowEnergyLimit(G4double elimit);
111 
112  virtual void DumpInfo() const;
113 
114  // for Debug
115  virtual G4bool CheckIt(const G4Track&);
116 
117 protected:
118  // hide copy constructor and assignment operaor as protected
121 
122 private:
123 
124  const G4Track* currentTrack;
125  // The pointer to G4Track
126 
127  G4double proposedKinEnergy;
128  // The final kinetic energy of the current particle.
129 
130  G4double lowEnergyLimit;
131  // The limit kinetic energy below which particle is stopped
132 
133  G4double currentCharge;
134  // The final charge of the current particle.
135 
136  G4ThreeVector proposedMomentumDirection;
137  // The final momentum direction of the current particle.
138 
139  G4ThreeVector proposedPolarization;
140  // The final polarization of the current particle.
141 };
142 
143 // ------------------------------------------------------------
144 
146 {
147  return proposedKinEnergy;
148 }
149 
151 {
152  proposedKinEnergy = energy;
153 }
154 
156 {
157  return currentCharge;
158 }
159 
161 {
162  return currentCharge;
163 }
164 
166 {
167  currentCharge = theCharge;
168 }
169 
171 {
172  currentCharge = theCharge;
173 }
174 
175 inline
177 {
178  return proposedMomentumDirection;
179 }
180 
181 inline
183 {
184  return proposedMomentumDirection;
185 }
186 
187 inline
189 {
190  proposedMomentumDirection = dir;
191 }
192 
193 inline
195 {
196  proposedMomentumDirection = dir;
197 }
198 
199 inline
201 {
202  proposedMomentumDirection.setX(Px);
203  proposedMomentumDirection.setY(Py);
204  proposedMomentumDirection.setZ(Pz);
205 }
206 
208 {
209  return currentTrack;
210 }
211 
212 inline
214 {
215  return proposedPolarization;
216 }
217 
218 inline
220 {
221  proposedPolarization = dir;
222 }
223 
224 inline
226 {
227  proposedPolarization.setX(Px);
228  proposedPolarization.setY(Py);
229  proposedPolarization.setZ(Pz);
230 }
231 
233 {
235  theLocalEnergyDeposit = 0.0;
237  InitializeSecondaries(track);
238  theParentWeight = track.GetWeight();
239  // isParentWeightProposed = false;
240  proposedKinEnergy = track.GetKineticEnergy();
241  currentCharge = track.GetDynamicParticle()->GetCharge();
242 }
243 
245 {
247  theLocalEnergyDeposit = 0.0;
249  InitializeSecondaries(track);
250  theParentWeight = track.GetWeight();
251  // isParentWeightProposed = false;
252  proposedKinEnergy = track.GetKineticEnergy();
253  currentCharge = track.GetDynamicParticle()->GetCharge();
254  proposedMomentumDirection = track.GetMomentumDirection();
255  proposedPolarization = track.GetPolarization();
256  currentTrack = &track;
257 }
258 
259 /*
260 inline void G4ParticleChangeForLoss::AddSecondary(G4DynamicParticle* aParticle)
261 {
262  // create track
263  G4Track* aTrack = new G4Track(aParticle, currentTrack->GetGlobalTime(),
264  currentTrack->GetPosition());
265 
266  // Touchable handle is copied to keep the pointer
267  aTrack->SetTouchableHandle(currentTrack->GetTouchableHandle());
268 
269  // add a secondary
270  G4VParticleChange::AddSecondary(aTrack);
271 }
272 */
274 {
275  lowEnergyLimit = elimit;
276 }
277 
278 #endif
279 
void ProposePolarization(const G4ThreeVector &dir)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4ThreeVector & GetPolarization() const
virtual void DumpInfo() const
virtual G4bool CheckIt(const G4Track &)
void InitializeForPostStep(const G4Track &)
G4Step * UpdateStepForAlongStep(G4Step *Step)
const G4DynamicParticle * GetDynamicParticle() const
void SetLowEnergyLimit(G4double elimit)
const G4ThreeVector & GetProposedPolarization() const
G4TrackStatus GetTrackStatus() const
const G4ThreeVector & GetMomentumDirection() const
void ProposeCharge(G4double finalCharge)
const G4ThreeVector & GetProposedMomentumDirection() const
void InitializeForAlongStep(const G4Track &)
G4double GetProposedKineticEnergy() const
void setY(double)
void setZ(double)
void setX(double)
G4Step * UpdateStepForPostStep(G4Step *Step)
G4ParticleChangeForLoss & operator=(const G4ParticleChangeForLoss &right)
G4double GetKineticEnergy() const
double energy
Definition: plottest35.C:25
bool G4bool
Definition: G4Types.hh:79
G4double GetCharge() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
Definition: G4Step.hh:76
void SetProposedMomentumDirection(const G4ThreeVector &dir)
void InitializeSecondaries(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
void SetProposedCharge(G4double theCharge)
Definition: Step.hh:41
G4double GetWeight() const
const G4Track * GetCurrentTrack() const
G4TrackStatus theStatusChange
double G4double
Definition: G4Types.hh:76
G4double theNonIonizingEnergyDeposit
TDirectory * dir
Definition: macro.C:5