Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
28//
29//
30// --------------------------------------------------------------
31// GEANT 4 class implementation file
32//
33//
34//
35// ------------------------------------------------------------
36// Implemented for the new scheme 10 May. 1998 H.Kurahige
37// Correct tratment of fpNextTouchable 12 May. 1998 H.Kurashige
38// Change to the next volume only if energy>0 19 Jan. 2004 V.Ivanchenko
39// --------------------------------------------------------------
40
42#include "G4TouchableHandle.hh"
43#include "G4Track.hh"
44#include "G4Step.hh"
45#include "G4TrackFastVector.hh"
46#include "G4DynamicParticle.hh"
47
49 : G4ParticleChange(), isMomentumChanged(false), theMaterialChange(0),
50 theMaterialCutsCoupleChange(0), theSensitiveDetectorChange(0),
51 fpVectorOfAuxiliaryPointsPointer(0)
52{
53 if (verboseLevel>2) {
54 G4cout << "G4ParticleChangeForTransport::G4ParticleChangeForTransport() "
55 << G4endl;
56 }
57}
58
60{
61 if (verboseLevel>2) {
62 G4cout << "G4ParticleChangeForTransport::~G4ParticleChangeForTransport() "
63 << G4endl;
64 }
65}
66
70 fpVectorOfAuxiliaryPointsPointer(0)
71{
72 if (verboseLevel>0) {
73 G4cout << "G4ParticleChangeForTransport:: copy constructor is called "
74 << G4endl;
75 }
77 isMomentumChanged = r.isMomentumChanged;
78 theMaterialChange = r.theMaterialChange;
79 theMaterialCutsCoupleChange = r.theMaterialCutsCoupleChange;
80 theSensitiveDetectorChange = r.theSensitiveDetectorChange;
81}
82
83// assignemnt operator
86{
87 if (verboseLevel>1) {
88 G4cout << "G4ParticleChangeForTransport:: assignment operator is called "
89 << G4endl;
90 }
91 if (this != &r)
92 {
98 theMaterialChange = r.theMaterialChange;
99 theMaterialCutsCoupleChange = r.theMaterialCutsCoupleChange;
100 theSensitiveDetectorChange = r.theSensitiveDetectorChange;
110 }
111 return *this;
112}
113
114//----------------------------------------------------------------
115// methods for updating G4Step
116//
117
119{
120 // Nothing happens for AtRestDoIt
121 if (verboseLevel>0) {
122 G4cout << "G4ParticleChangeForTransport::UpdateStepForAtRest() is called"
123 << G4endl;
124 G4cout << " Nothing happens for this method " << G4endl;
125 }
126 // Update the G4Step specific attributes
127 return UpdateStepInfo(pStep);
128}
129
130
132{
133 // Smooth curved tajectory representation: let the Step know about
134 // the auxiliary trajectory points (jacek 30/10/2002)
135 pStep->SetPointerToVectorOfAuxiliaryPoints(fpVectorOfAuxiliaryPointsPointer);
136
137 // copy of G4ParticleChange::UpdateStepForAlongStep
138 // i.e. no effect for touchable
139
140 // A physics process always calculates the final state of the
141 // particle relative to the initial state at the beginning
142 // of the Step, i.e., based on information of G4Track (or
143 // equivalently the PreStepPoint).
144 // So, the differences (delta) between these two states have to be
145 // calculated and be accumulated in PostStepPoint.
146
147 // Take note that the return type of GetMomentumChange is a
148 // pointer to G4ThreeVector. Also it is a normalized
149 // momentum vector.
150
151 G4StepPoint* pPreStepPoint = pStep->GetPreStepPoint();
152 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
153 G4Track* aTrack = pStep->GetTrack();
154 G4double mass = aTrack->GetDynamicParticle()->GetMass();
155
156 // uodate kinetic energy
157 // now assume that no energy change in transportation
158 // However it is not true in electric fields
159 // Case for changing energy will be implemented in future
160
161
162 // update momentum direction and energy
163 if (isMomentumChanged) {
164 G4double energy;
165 energy= pPostStepPoint->GetKineticEnergy()
166 + (theEnergyChange - pPreStepPoint->GetKineticEnergy());
167
168 // calculate new momentum
169 G4ThreeVector pMomentum = pPostStepPoint->GetMomentum()
171 - pPreStepPoint->GetMomentum());
172 G4double tMomentum = pMomentum.mag();
173 G4ThreeVector direction(1.0,0.0,0.0);
174 if( tMomentum > 0. ){
175 G4double inv_Momentum= 1.0 / tMomentum;
176 direction= pMomentum * inv_Momentum;
177 }
178 pPostStepPoint->SetMomentumDirection(direction);
179 pPostStepPoint->SetKineticEnergy( energy );
180 }
181 if (isVelocityChanged) pPostStepPoint->SetVelocity(theVelocityChange);
182
183 // stop case should not occur
184 //pPostStepPoint->SetMomentumDirection(G4ThreeVector(1., 0., 0.));
185
186
187 // update polarization
189 - pPreStepPoint->GetPolarization());
190
191 // update position and time
192 pPostStepPoint->AddPosition( thePositionChange
193 - pPreStepPoint->GetPosition() );
194 pPostStepPoint->AddGlobalTime( theTimeChange
195 - pPreStepPoint->GetLocalTime());
196 pPostStepPoint->AddLocalTime( theTimeChange
197 - pPreStepPoint->GetLocalTime());
198 pPostStepPoint->AddProperTime( theProperTimeChange
199 - pPreStepPoint->GetProperTime());
200
201#ifdef G4VERBOSE
202 if (debugFlag) CheckIt(*aTrack);
203#endif
204
205 // Update the G4Step specific attributes
206 //pStep->SetStepLength( theTrueStepLength );
207 // pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
209 return pStep;
210 // return UpdateStepInfo(pStep);
211}
212
214{
215 // A physics process always calculates the final state of the particle
216
217 // Change volume only if some kinetic energy remains
218 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
219 if(pPostStepPoint->GetKineticEnergy() > 0.0) {
220
221 // update next touchable
222 // (touchable can be changed only at PostStepDoIt)
223 pPostStepPoint->SetTouchableHandle( theTouchableHandle );
224
225 pPostStepPoint->SetMaterial( theMaterialChange );
226 pPostStepPoint->SetMaterialCutsCouple( theMaterialCutsCoupleChange );
227 pPostStepPoint->SetSensitiveDetector( theSensitiveDetectorChange );
228 }
229 if( this->GetLastStepInVolume() ){
230 pStep->SetLastStepFlag();
231 }else{
232 pStep->ClearLastStepFlag();
233 }
234 // It used to call base class's method
235 // - but this would copy uninitialised data members
236 // return G4ParticleChange::UpdateStepForPostStep(pStep);
237
238 // Copying what the base class does would instead
239 // - also not useful
240 // return G4VParticleChange::UpdateStepInfo(pStep);
241
242 return pStep;
243}
244
245
246//----------------------------------------------------------------
247// methods for printing messages
248//
249
251{
252// use base-class DumpInfo
254
255 G4int oldprc = G4cout.precision(3);
256 G4cout << " Touchable (pointer) : "
257 << std::setw(20) << theTouchableHandle() << G4endl;
258 G4cout.precision(oldprc);
259}
260
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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:78
G4Track * GetTrack() const
void SetLastStepFlag()
G4StepPoint * GetPreStepPoint() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *theNewVectorPointer)
Definition: G4Step.hh:237
void SetControlFlag(G4SteppingControl StepControlFlag)
void ClearLastStepFlag()
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
G4TrackStatus theStatusChange
G4TrackFastVector * theListOfSecondaries
G4SteppingControl theSteppingControlFlag
G4bool GetLastStepInVolume() const