Geant4-11
G4ParticleChangeForTransport.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// G4ParticleChangeForTransport class implementation
27//
28// Author: Hisaya Kurashige, 10 May 1998
29// --------------------------------------------------------------------
30
32#include "G4TouchableHandle.hh"
33#include "G4Track.hh"
34#include "G4Step.hh"
35#include "G4TrackFastVector.hh"
36#include "G4DynamicParticle.hh"
37
38// --------------------------------------------------------------------
41{
42 // Disable flag that is enabled in G4VParticleChange if G4VERBOSE.
43 debugFlag = false;
44}
45
46// --------------------------------------------------------------------
48{
49}
50
51// --------------------------------------------------------------------
55{
61}
62
63// --------------------------------------------------------------------
66{
67 if(this != &r)
68 {
86 }
87 return *this;
88}
89
90// --------------------------------------------------------------------
92{
93 // Update the G4Step specific attributes
94 return UpdateStepInfo(pStep);
95}
96
97// --------------------------------------------------------------------
99{
100 // Smooth curved tajectory representation: let the Step know about
101 // the auxiliary trajectory points (jacek 30/10/2002)
103
104 // copy of G4ParticleChange::UpdateStepForAlongStep
105 // i.e. no effect for touchable
106
107 // A physics process always calculates the final state of the
108 // particle relative to the initial state at the beginning
109 // of the Step, i.e., based on information of G4Track (or
110 // equivalently the PreStepPoint).
111 // So, the differences (delta) between these two states have to be
112 // calculated and be accumulated in PostStepPoint.
113
114 // Take note that the return type of GetMomentumChange is a
115 // pointer to G4ThreeVector. Also it is a normalized
116 // momentum vector.
117
118 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
119 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
120 G4Track* aTrack = pStep->GetTrack();
121 G4double mass = aTrack->GetDynamicParticle()->GetMass();
122
123 // update kinetic energy
124 // now assume that no energy change in transportation
125 // However it is not true in electric fields
126 // Case for changing energy will be implemented in future
127
128 // update momentum direction and energy
130 {
132 energy = pPostStepPoint->GetKineticEnergy() +
133 (theEnergyChange - pPreStepPoint->GetKineticEnergy());
134
135 // calculate new momentum
136 G4ThreeVector pMomentum =
137 pPostStepPoint->GetMomentum() +
139 pPreStepPoint->GetMomentum());
140 G4double tMomentum = pMomentum.mag();
141 G4ThreeVector direction(1.0, 0.0, 0.0);
142 if(tMomentum > 0.)
143 {
144 G4double inv_Momentum = 1.0 / tMomentum;
145 direction = pMomentum * inv_Momentum;
146 }
147 pPostStepPoint->SetMomentumDirection(direction);
148 pPostStepPoint->SetKineticEnergy(energy);
149 }
151 pPostStepPoint->SetVelocity(theVelocityChange);
152
153 // stop case should not occur
154 // pPostStepPoint->SetMomentumDirection(G4ThreeVector(1., 0., 0.));
155
156 // update polarization
157 pPostStepPoint->AddPolarization(thePolarizationChange -
158 pPreStepPoint->GetPolarization());
159
160 // update position and time
161 pPostStepPoint->AddPosition(thePositionChange - pPreStepPoint->GetPosition());
162 pPostStepPoint->AddGlobalTime(theTimeChange - pPreStepPoint->GetLocalTime());
163 pPostStepPoint->AddLocalTime(theTimeChange - pPreStepPoint->GetLocalTime());
164 pPostStepPoint->AddProperTime(theProperTimeChange -
165 pPreStepPoint->GetProperTime());
166
167#ifdef G4VERBOSE
168 if(debugFlag) { CheckIt(*aTrack); }
169#endif
170
171 // Update the G4Step specific attributes
172 // pStep->SetStepLength( theTrueStepLength );
173 // pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
175
176 return pStep;
177}
178
179// --------------------------------------------------------------------
181{
182 // A physics process always calculates the final state of the particle
183
184 // Change volume only if some kinetic energy remains
185 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
186 if(pPostStepPoint->GetKineticEnergy() > 0.0)
187 {
188 // update next touchable
189 // (touchable can be changed only at PostStepDoIt)
190 pPostStepPoint->SetTouchableHandle(theTouchableHandle);
191
192 pPostStepPoint->SetMaterial(theMaterialChange);
195 }
196 if(this->GetFirstStepInVolume())
197 {
198 pStep->SetFirstStepFlag();
199 }
200 else
201 {
202 pStep->ClearFirstStepFlag();
203 }
204 if(this->GetLastStepInVolume())
205 {
206 pStep->SetLastStepFlag();
207 }
208 else
209 {
210 pStep->ClearLastStepFlag();
211 }
212 // It used to call base class's method
213 // - but this would copy uninitialised data members
214 // return G4ParticleChange::UpdateStepForPostStep(pStep);
215
216 // Copying what the base class does would instead
217 // - also not useful
218 // return G4VParticleChange::UpdateStepInfo(pStep);
219
220 return pStep;
221}
222
223// --------------------------------------------------------------------
225{
226 // use base-class DumpInfo
228
229 G4int oldprc = G4cout.precision(3);
230 G4cout << " Touchable (pointer) : " << std::setw(20)
232 G4cout.precision(oldprc);
233}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag() const
G4double GetMass() const
G4ParticleChangeForTransport & operator=(const G4ParticleChangeForTransport &right)
const G4MaterialCutsCouple * theMaterialCutsCoupleChange
G4VSensitiveDetector * theSensitiveDetectorChange
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)
std::vector< G4ThreeVector > * fpVectorOfAuxiliaryPointsPointer
G4ThreeVector CalcMomentum(G4double energy, G4ThreeVector direction, G4double mass) const
virtual void DumpInfo() const
G4ThreeVector thePositionChange
G4ThreeVector theMomentumDirectionChange
G4Step * UpdateStepInfo(G4Step *Step)
G4double theProperTimeChange
virtual G4bool CheckIt(const G4Track &)
G4ThreeVector thePolarizationChange
void AddPolarization(const G4ThreeVector &aValue)
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetKineticEnergy(const G4double aValue)
void SetMaterial(G4Material *)
G4double GetProperTime() const
void AddProperTime(const G4double aValue)
void SetVelocity(G4double v)
void AddPosition(const G4ThreeVector &aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPosition() const
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
void AddGlobalTime(const G4double aValue)
G4double GetLocalTime() const
const G4ThreeVector & GetPolarization() const
G4double GetKineticEnergy() const
void AddLocalTime(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
Definition: G4Step.hh:62
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *vec)
G4Track * GetTrack() const
void SetLastStepFlag()
G4StepPoint * GetPreStepPoint() const
void SetControlFlag(G4SteppingControl StepControlFlag)
void ClearLastStepFlag()
void SetFirstStepFlag()
void ClearFirstStepFlag()
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
G4bool GetFirstStepInVolume() const
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
G4bool GetLastStepInVolume() const
G4double energy(const ThreeVector &p, const G4double m)