Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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}
43
44// --------------------------------------------------------------------
46{
47}
48
49// --------------------------------------------------------------------
53{
55 isMomentumChanged = r.isMomentumChanged;
56 theMaterialChange = r.theMaterialChange;
57 theMaterialCutsCoupleChange = r.theMaterialCutsCoupleChange;
58 theSensitiveDetectorChange = r.theSensitiveDetectorChange;
59}
60
61// --------------------------------------------------------------------
64{
65 if(this != &r)
66 {
72 theMaterialChange = r.theMaterialChange;
73 theMaterialCutsCoupleChange = r.theMaterialCutsCoupleChange;
74 theSensitiveDetectorChange = r.theSensitiveDetectorChange;
84 }
85 return *this;
86}
87
88// --------------------------------------------------------------------
90{
91 // Update the G4Step specific attributes
92 return UpdateStepInfo(pStep);
93}
94
95// --------------------------------------------------------------------
97{
98 // Smooth curved tajectory representation: let the Step know about
99 // the auxiliary trajectory points (jacek 30/10/2002)
100 pStep->SetPointerToVectorOfAuxiliaryPoints(fpVectorOfAuxiliaryPointsPointer);
101
102 // copy of G4ParticleChange::UpdateStepForAlongStep
103 // i.e. no effect for touchable
104
105 // A physics process always calculates the final state of the
106 // particle relative to the initial state at the beginning
107 // of the Step, i.e., based on information of G4Track (or
108 // equivalently the PreStepPoint).
109 // So, the differences (delta) between these two states have to be
110 // calculated and be accumulated in PostStepPoint.
111
112 // Take note that the return type of GetMomentumChange is a
113 // pointer to G4ThreeVector. Also it is a normalized
114 // momentum vector.
115
116 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
117 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
118 G4Track* aTrack = pStep->GetTrack();
119 G4double mass = aTrack->GetDynamicParticle()->GetMass();
120
121 // update kinetic energy
122 // now assume that no energy change in transportation
123 // However it is not true in electric fields
124 // Case for changing energy will be implemented in future
125
126 // update momentum direction and energy
127 if(isMomentumChanged)
128 {
129 G4double energy;
130 energy = pPostStepPoint->GetKineticEnergy() +
131 (theEnergyChange - pPreStepPoint->GetKineticEnergy());
132
133 // calculate new momentum
134 G4ThreeVector pMomentum =
135 pPostStepPoint->GetMomentum() +
137 pPreStepPoint->GetMomentum());
138 G4double tMomentum = pMomentum.mag();
139 G4ThreeVector direction(1.0, 0.0, 0.0);
140 if(tMomentum > 0.)
141 {
142 G4double inv_Momentum = 1.0 / tMomentum;
143 direction = pMomentum * inv_Momentum;
144 }
145 pPostStepPoint->SetMomentumDirection(direction);
146 pPostStepPoint->SetKineticEnergy(energy);
147 }
149 pPostStepPoint->SetVelocity(theVelocityChange);
150
151 // stop case should not occur
152 // pPostStepPoint->SetMomentumDirection(G4ThreeVector(1., 0., 0.));
153
154 // update polarization
155 pPostStepPoint->AddPolarization(thePolarizationChange -
156 pPreStepPoint->GetPolarization());
157
158 // update position and time
159 pPostStepPoint->AddPosition(thePositionChange - pPreStepPoint->GetPosition());
160 pPostStepPoint->AddGlobalTime(theTimeChange - pPreStepPoint->GetLocalTime());
161 pPostStepPoint->AddLocalTime(theTimeChange - pPreStepPoint->GetLocalTime());
162 pPostStepPoint->AddProperTime(theProperTimeChange -
163 pPreStepPoint->GetProperTime());
164
165#ifdef G4VERBOSE
166 if(debugFlag) { CheckIt(*aTrack); }
167#endif
168
169 // Update the G4Step specific attributes
170 // pStep->SetStepLength( theTrueStepLength );
171 // pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
173
174 return pStep;
175}
176
177// --------------------------------------------------------------------
179{
180 // A physics process always calculates the final state of the particle
181
182 // Change volume only if some kinetic energy remains
183 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
184 if(pPostStepPoint->GetKineticEnergy() > 0.0)
185 {
186 // update next touchable
187 // (touchable can be changed only at PostStepDoIt)
188 pPostStepPoint->SetTouchableHandle(theTouchableHandle);
189
190 pPostStepPoint->SetMaterial(theMaterialChange);
191 pPostStepPoint->SetMaterialCutsCouple(theMaterialCutsCoupleChange);
192 pPostStepPoint->SetSensitiveDetector(theSensitiveDetectorChange);
193 }
194 if(this->GetFirstStepInVolume())
195 {
196 pStep->SetFirstStepFlag();
197 }
198 else
199 {
200 pStep->ClearFirstStepFlag();
201 }
202 if(this->GetLastStepInVolume())
203 {
204 pStep->SetLastStepFlag();
205 }
206 else
207 {
208 pStep->ClearLastStepFlag();
209 }
210 // It used to call base class's method
211 // - but this would copy uninitialised data members
212 // return G4ParticleChange::UpdateStepForPostStep(pStep);
213
214 // Copying what the base class does would instead
215 // - also not useful
216 // return G4VParticleChange::UpdateStepInfo(pStep);
217
218 return pStep;
219}
220
221// --------------------------------------------------------------------
223{
224 // use base-class DumpInfo
226
227 G4int oldprc = G4cout.precision(3);
228 G4cout << " Touchable (pointer) : " << std::setw(20)
230 G4cout.precision(oldprc);
231}
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)
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)
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