Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Track.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 //
32 // G4Track.cc
33 //
34 //---------------------------------------------------------------
35 // Add copy constructor Hisaya Feb. 07 01
36 // Fix GetVelocity Hisaya Feb. 17 01
37 // Modification for G4TouchableHandle 22 Oct. 2001 R.Chytracek//
38 // Fix GetVelocity (bug report #741) Horton-Smith Apr 14 2005
39 // Remove massless check in GetVelocity 02 Apr. 09 H.Kurashige
40 // Use G4VelocityTable 17 AUg. 2011 H.Kurashige
41 
42 #include "G4Track.hh"
43 #include "G4PhysicalConstants.hh"
44 #include "G4ParticleTable.hh"
45 #include "G4VelocityTable.hh"
46 
47 #include <iostream>
48 #include <iomanip>
49 
51 
52 G4VelocityTable* G4Track::velTable=0;
53 
55 G4Track::G4Track(G4DynamicParticle* apValueDynamicParticle,
56  G4double aValueTime,
57  const G4ThreeVector& aValuePosition)
59  : fCurrentStepNumber(0), fPosition(aValuePosition),
60  fGlobalTime(aValueTime), fLocalTime(0.),
61  fTrackLength(0.),
62  fParentID(0), fTrackID(0),
63  fVelocity(c_light),
64  fpDynamicParticle(apValueDynamicParticle),
65  fTrackStatus(fAlive),
66  fBelowThreshold(false), fGoodForTracking(false),
67  fStepLength(0.0), fWeight(1.0),
68  fpStep(0),
69  fVtxKineticEnergy(0.0),
70  fpLVAtVertex(0), fpCreatorProcess(0),
71  fpUserInformation(0),
72  prev_mat(0), groupvel(0),
73  prev_velocity(0.0), prev_momentum(0.0),
74  is_OpticalPhoton(false),
75  useGivenVelocity(false)
76 {
77  static G4bool isFirstTime = true;
78  static G4ParticleDefinition* fOpticalPhoton =0;
79  if ( isFirstTime ) {
80  isFirstTime = false;
81  // set fOpticalPhoton
82  fOpticalPhoton = G4ParticleTable::GetParticleTable()->FindParticle("opticalphoton");
83  }
84  // check if the particle type is Optical Photon
85  is_OpticalPhoton = (fpDynamicParticle->GetDefinition() == fOpticalPhoton);
86 
87  if (velTable ==0 ) velTable = G4VelocityTable::GetVelocityTable();
88 
89  fVelocity = CalculateVelocity();
90 
91 }
92 
96  : fCurrentStepNumber(0),
97  fGlobalTime(0), fLocalTime(0.),
98  fTrackLength(0.),
99  fParentID(0), fTrackID(0),
100  fVelocity(c_light),
101  fpDynamicParticle(0),
102  fTrackStatus(fAlive),
103  fBelowThreshold(false), fGoodForTracking(false),
104  fStepLength(0.0), fWeight(1.0),
105  fpStep(0),
106  fVtxKineticEnergy(0.0),
107  fpLVAtVertex(0), fpCreatorProcess(0),
108  fpUserInformation(0),
109  prev_mat(0), groupvel(0),
110  prev_velocity(0.0), prev_momentum(0.0),
111  is_OpticalPhoton(false),
112  useGivenVelocity(false)
113 {
114 }
118  : fCurrentStepNumber(0),
119  fGlobalTime(0), fLocalTime(0.),
120  fTrackLength(0.),
121  fParentID(0), fTrackID(0),
122  fVelocity(c_light),
123  fpDynamicParticle(0),
124  fTrackStatus(fAlive),
125  fBelowThreshold(false), fGoodForTracking(false),
126  fStepLength(0.0), fWeight(1.0),
127  fpStep(0),
128  fVtxKineticEnergy(0.0),
129  fpLVAtVertex(0), fpCreatorProcess(0),
130  fpUserInformation(0),
131  prev_mat(0), groupvel(0),
132  prev_velocity(0.0), prev_momentum(0.0),
133  is_OpticalPhoton(false),
134  useGivenVelocity(false)
135 {
136  *this = right;
137 }
138 
142 {
143  delete fpDynamicParticle;
144  delete fpUserInformation;
145 }
146 
148 G4Track & G4Track::operator=(const G4Track &right)
150 {
151  if (this != &right) {
152  fPosition = right.fPosition;
153  fGlobalTime = right.fGlobalTime;
154  fLocalTime = right.fLocalTime;
155  fTrackLength = right.fTrackLength;
156  fWeight = right.fWeight;
157  fStepLength = right.fStepLength;
158 
159  // Track ID (and Parent ID) is not copied and set to zero for new track
160  fTrackID = 0;
161  fParentID =0;
162 
163  // CurrentStepNumber is set to be 0
164  fCurrentStepNumber = 0;
165 
166  // velocity information
167  fVelocity = right.fVelocity;
168 
169  // dynamic particle information
170  fpDynamicParticle = new G4DynamicParticle(*(right.fpDynamicParticle));
171 
172  // track status and flags for tracking
173  fTrackStatus = right.fTrackStatus;
174  fBelowThreshold = right.fBelowThreshold;
175  fGoodForTracking = right.fGoodForTracking;
176 
177  // Step information (Step Length, Step Number, pointer to the Step,)
178  // are not copied
179  fpStep=0;
180 
181  // vertex information
182  fVtxPosition = right.fVtxPosition;
183  fpLVAtVertex = right.fpLVAtVertex;
184  fVtxKineticEnergy = right.fVtxKineticEnergy;
185  fVtxMomentumDirection = right.fVtxMomentumDirection;
186 
187  // CreatorProcess and UserInformation are not copied
188  fpCreatorProcess = 0;
189  fpUserInformation = 0;
190 
191  prev_mat = right.prev_mat;
192  groupvel = right.groupvel;
193  prev_velocity = right.prev_velocity;
194  prev_momentum = right.prev_momentum;
195 
196  is_OpticalPhoton = right.is_OpticalPhoton;
197  useGivenVelocity = right.useGivenVelocity;
198  }
199  return *this;
200 }
201 
203 void G4Track::CopyTrackInfo(const G4Track& right)
205 {
206  *this = right;
207 }
208 
212 {
213  if (useGivenVelocity) return fVelocity;
214 
215  G4double velocity = c_light ;
216 
217  G4double mass = fpDynamicParticle->GetMass();
218 
219  // special case for photons
220  if ( is_OpticalPhoton ) return CalculateVelocityForOpticalPhoton();
221 
222  // particles other than optical photon
223  if (mass<DBL_MIN) {
224  // Zero Mass
225  velocity = c_light;
226  } else {
227  G4double T = (fpDynamicParticle->GetKineticEnergy())/mass;
228  if (T > GetMaxTOfVelocityTable()) {
229  velocity = c_light;
230  } else if (T<DBL_MIN) {
231  velocity =0.;
232  } else if (T<GetMinTOfVelocityTable()) {
233  velocity = c_light*std::sqrt(T*(T+2.))/(T+1.0);
234  } else {
235  velocity = velTable->Value(T);
236  }
237 
238  }
239  return velocity ;
240 }
241 
245 {
246 
247  G4double velocity = c_light ;
248 
249 
250  G4Material* mat=0;
251  G4bool update_groupvel = false;
252  if ( fpStep !=0 ){
253  mat= this->GetMaterial(); // Fix for repeated volumes
254  }else{
255  if (fpTouchable!=0){
256  mat=fpTouchable->GetVolume()->GetLogicalVolume()->GetMaterial();
257  }
258  }
259  // check if previous step is in the same volume
260  // and get new GROUPVELOCITY table if necessary
261  if ((mat != 0) && ((mat != prev_mat)||(groupvel==0))) {
262  groupvel = 0;
263  if(mat->GetMaterialPropertiesTable() != 0)
264  groupvel = mat->GetMaterialPropertiesTable()->GetProperty("GROUPVEL");
265  update_groupvel = true;
266  }
267  prev_mat = mat;
268 
269  if (groupvel != 0 ) {
270  // light velocity = c/(rindex+d(rindex)/d(log(E_phot)))
271  // values stored in GROUPVEL material properties vector
272  velocity = prev_velocity;
273 
274  // check if momentum is same as in the previous step
275  // and calculate group velocity if necessary
276  G4double current_momentum = fpDynamicParticle->GetTotalMomentum();
277  if( update_groupvel || (current_momentum != prev_momentum) ) {
278  velocity =
279  groupvel->Value(current_momentum);
280  prev_velocity = velocity;
281  prev_momentum = current_momentum;
282  }
283  }
284 
285  return velocity ;
286 }
287 
291 {
294 }
295 
300 
305 
310