source: trunk/source/track/include/G4ParticleChangeForLoss.hh@ 1340

Last change on this file since 1340 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 10.4 KB
Line 
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,v 1.22 2009/06/17 17:25:57 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
29//
30//
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34//
35// ------------------------------------------------------------
36// Implemented for the new scheme 23 Mar. 1998 H.Kurahige
37//
38// Modified:
39// 16.01.04 V.Ivanchenko update for model variant of energy loss
40// 15.04.05 V.Ivanchenko inline update methods
41// 30.01.06 V.Ivanchenko add ProposedMomentumDirection for AlongStep
42// and ProposeWeight for PostStep
43// 07.06.06 V.Ivanchenko RemoveProposedMomentumDirection from AlongStep
44// 28.08.06 V.Ivanchenko Added access to current track and polarizaion
45// 17.06.09 V.Ivanchenko Added SetLowEnergyLimit method
46//
47// ------------------------------------------------------------
48//
49// Class Description
50// This class is a concrete class for ParticleChange for EnergyLoss
51//
52#ifndef G4ParticleChangeForLoss_h
53#define G4ParticleChangeForLoss_h 1
54
55#include "globals.hh"
56#include "G4ios.hh"
57#include "G4VParticleChange.hh"
58
59class G4DynamicParticle;
60
61class G4ParticleChangeForLoss: public G4VParticleChange
62{
63public:
64 // default constructor
65 G4ParticleChangeForLoss();
66
67 // destructor
68 virtual ~G4ParticleChangeForLoss();
69
70 // with description
71 // ----------------------------------------------------
72 // --- the following methods are for updating G4Step -----
73
74 G4Step* UpdateStepForAlongStep(G4Step* Step);
75 G4Step* UpdateStepForPostStep(G4Step* Step);
76 // A physics process gives the final state of the particle
77 // based on information of G4Track
78
79 void InitializeForAlongStep(const G4Track&);
80 void InitializeForPostStep(const G4Track&);
81 //Initialize all propoerties by using G4Track information
82
83 void AddSecondary(G4DynamicParticle* aParticle);
84 // Add next secondary
85
86 G4double GetProposedCharge() const;
87 void SetProposedCharge(G4double theCharge);
88 // Get/Set theCharge
89
90 G4double GetCharge() const;
91 void ProposeCharge(G4double finalCharge);
92 // Get/Propose the final dynamical Charge in G4DynamicParticle
93
94 G4double GetProposedKineticEnergy() const;
95 void SetProposedKineticEnergy(G4double proposedKinEnergy);
96 // Get/Set the final kinetic energy of the current particle.
97
98 const G4ThreeVector& GetProposedMomentumDirection() const;
99 void SetProposedMomentumDirection(const G4ThreeVector& dir);
100 const G4ThreeVector& GetMomentumDirection() const;
101 void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz);
102 void ProposeMomentumDirection(const G4ThreeVector& Pfinal);
103 // Get/Propose the MomentumDirection vector: it is the final momentum direction.
104
105 const G4ThreeVector& GetProposedPolarization() const;
106 void ProposePolarization(const G4ThreeVector& dir);
107 void ProposePolarization(G4double Px, G4double Py, G4double Pz);
108
109 const G4Track* GetCurrentTrack() const;
110
111 void SetLowEnergyLimit(G4double elimit);
112
113 virtual void DumpInfo() const;
114
115 // for Debug
116 virtual G4bool CheckIt(const G4Track&);
117
118protected:
119 // hide copy constructor and assignment operaor as protected
120 G4ParticleChangeForLoss(const G4ParticleChangeForLoss &right);
121 G4ParticleChangeForLoss & operator=(const G4ParticleChangeForLoss &right);
122
123private:
124
125 const G4Track* currentTrack;
126 // The pointer to G4Track
127
128 G4double proposedKinEnergy;
129 // The final kinetic energy of the current particle.
130
131 G4double lowEnergyLimit;
132 // The limit kinetic energy below which particle is stopped
133
134 G4double currentCharge;
135 // The final charge of the current particle.
136
137 G4ThreeVector proposedMomentumDirection;
138 // The final momentum direction of the current particle.
139
140 G4ThreeVector proposedPolarization;
141 // The final polarization of the current particle.
142};
143
144// ------------------------------------------------------------
145
146inline G4double G4ParticleChangeForLoss::GetProposedKineticEnergy() const
147{
148 return proposedKinEnergy;
149}
150
151inline void G4ParticleChangeForLoss::SetProposedKineticEnergy(G4double energy)
152{
153 proposedKinEnergy = energy;
154}
155
156inline G4double G4ParticleChangeForLoss::GetProposedCharge() const
157{
158 return currentCharge;
159}
160
161inline G4double G4ParticleChangeForLoss::GetCharge() const
162{
163 return currentCharge;
164}
165
166inline void G4ParticleChangeForLoss::SetProposedCharge(G4double theCharge)
167{
168 currentCharge = theCharge;
169}
170
171inline void G4ParticleChangeForLoss::ProposeCharge(G4double theCharge)
172{
173 currentCharge = theCharge;
174}
175
176inline
177 const G4ThreeVector& G4ParticleChangeForLoss::GetProposedMomentumDirection() const
178{
179 return proposedMomentumDirection;
180}
181
182inline
183 const G4ThreeVector& G4ParticleChangeForLoss::GetMomentumDirection() const
184{
185 return proposedMomentumDirection;
186}
187
188inline
189 void G4ParticleChangeForLoss::ProposeMomentumDirection(const G4ThreeVector& dir)
190{
191 proposedMomentumDirection = dir;
192}
193
194inline
195 void G4ParticleChangeForLoss::SetProposedMomentumDirection(const G4ThreeVector& dir)
196{
197 proposedMomentumDirection = dir;
198}
199
200inline
201 void G4ParticleChangeForLoss::ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
202{
203 proposedMomentumDirection.setX(Px);
204 proposedMomentumDirection.setY(Py);
205 proposedMomentumDirection.setZ(Pz);
206}
207
208inline const G4Track* G4ParticleChangeForLoss::GetCurrentTrack() const
209{
210 return currentTrack;
211}
212
213inline
214 const G4ThreeVector& G4ParticleChangeForLoss::GetProposedPolarization() const
215{
216 return proposedPolarization;
217}
218
219inline
220 void G4ParticleChangeForLoss::ProposePolarization(const G4ThreeVector& dir)
221{
222 proposedPolarization = dir;
223}
224
225inline
226 void G4ParticleChangeForLoss::ProposePolarization(G4double Px, G4double Py, G4double Pz)
227{
228 proposedPolarization.setX(Px);
229 proposedPolarization.setY(Py);
230 proposedPolarization.setZ(Pz);
231}
232
233inline void G4ParticleChangeForLoss::InitializeForAlongStep(const G4Track& track)
234{
235 theStatusChange = track.GetTrackStatus();
236 theLocalEnergyDeposit = 0.0;
237 theNonIonizingEnergyDeposit = 0.0;
238 InitializeSecondaries(track);
239 theParentWeight = track.GetWeight();
240 proposedKinEnergy = track.GetKineticEnergy();
241 currentCharge = track.GetDynamicParticle()->GetCharge();
242}
243
244inline void G4ParticleChangeForLoss::InitializeForPostStep(const G4Track& track)
245{
246 theStatusChange = track.GetTrackStatus();
247 theLocalEnergyDeposit = 0.0;
248 theNonIonizingEnergyDeposit = 0.0;
249 InitializeSecondaries(track);
250 theParentWeight = track.GetWeight();
251 proposedKinEnergy = track.GetKineticEnergy();
252 currentCharge = track.GetDynamicParticle()->GetCharge();
253 proposedMomentumDirection = track.GetMomentumDirection();
254 proposedPolarization = track.GetPolarization();
255 currentTrack = &track;
256}
257
258//----------------------------------------------------------------
259// methods for updating G4Step
260//
261
262inline G4Step* G4ParticleChangeForLoss::UpdateStepForAlongStep(G4Step* pStep)
263{
264 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
265
266 // accumulate change of the kinetic energy
267 G4double kinEnergy = pPostStepPoint->GetKineticEnergy() +
268 (proposedKinEnergy - pStep->GetPreStepPoint()->GetKineticEnergy());
269
270 // update kinetic energy and charge
271 if (kinEnergy < lowEnergyLimit) {
272 theLocalEnergyDeposit += kinEnergy;
273 kinEnergy = 0.0;
274 } else {
275 pPostStepPoint->SetCharge( currentCharge );
276 }
277 pPostStepPoint->SetKineticEnergy( kinEnergy );
278
279 // update weight
280 // this feature is commented out, it should be overwritten in case
281 // if energy loss processes will use biasing
282 // G4double newWeight = theParentWeight*(pPostStepPoint->GetWeight())
283 // /(pPreStepPoint->GetWeight());
284 // pPostStepPoint->SetWeight( newWeight );
285 pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
286 pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
287 return pStep;
288}
289
290inline G4Step* G4ParticleChangeForLoss::UpdateStepForPostStep(G4Step* pStep)
291{
292 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
293 pPostStepPoint->SetCharge( currentCharge );
294 pPostStepPoint->SetMomentumDirection( proposedMomentumDirection );
295 pPostStepPoint->SetKineticEnergy( proposedKinEnergy );
296 pPostStepPoint->SetPolarization( proposedPolarization );
297 // update weight if process cannot do that
298 if (!fSetParentWeightByProcess)
299 pPostStepPoint->SetWeight( theParentWeight );
300
301 pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
302 pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
303 return pStep;
304}
305
306inline void G4ParticleChangeForLoss::AddSecondary(G4DynamicParticle* aParticle)
307{
308 // create track
309 G4Track* aTrack = new G4Track(aParticle, currentTrack->GetGlobalTime(),
310 currentTrack->GetPosition());
311
312 // Touchable handle is copied to keep the pointer
313 aTrack->SetTouchableHandle(currentTrack->GetTouchableHandle());
314
315 // add a secondary
316 G4VParticleChange::AddSecondary(aTrack);
317}
318
319inline void G4ParticleChangeForLoss::SetLowEnergyLimit(G4double elimit)
320{
321 lowEnergyLimit = elimit;
322}
323
324#endif
325
Note: See TracBrowser for help on using the repository browser.