source: trunk/source/track/include/G4Track.icc@ 1204

Last change on this file since 1204 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 10.1 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: G4Track.icc,v 1.16 2008/10/24 08:22:20 kurasige Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//-----------------------------------------------------------------
31// Definitions of inline functions
32//-----------------------------------------------------------------
33// change GetMaterial 16 Feb. 2000 H.Kurashige
34
35#if defined G4TRACK_ALLOC_EXPORT
36 extern G4DLLEXPORT G4Allocator<G4Track> aTrackAllocator;
37#else
38 extern G4DLLIMPORT G4Allocator<G4Track> aTrackAllocator;
39#endif
40
41// Operators
42
43 inline void* G4Track::operator new(size_t)
44 { void *aTrack;
45 aTrack = (void *) aTrackAllocator.MallocSingle();
46 return aTrack;
47 }
48 // Override "new" for "G4Allocator".
49
50 inline void G4Track::operator delete(void *aTrack)
51 { aTrackAllocator.FreeSingle((G4Track *) aTrack);}
52 // Override "delete" for "G4Allocator".
53
54 inline G4bool G4Track::operator==( const G4Track& s)
55 { return (this==&s); }
56 // Define "==" operator because "G4TrackVector" uses
57 // "RWPtrOrderdVector" which requires this.
58
59// Get/Set functions
60// dynamic particle
61 inline const G4DynamicParticle* G4Track::GetDynamicParticle() const
62 { return fpDynamicParticle; }
63
64// particle definition
65 inline G4ParticleDefinition* G4Track::GetDefinition() const
66 { return fpDynamicParticle->GetDefinition(); }
67
68// parent track ID
69 inline G4int G4Track::GetParentID() const
70 { return fParentID; }
71
72 inline void G4Track::SetParentID(const G4int aValue)
73 { fParentID = aValue; }
74
75// current track ID
76 inline G4int G4Track::GetTrackID() const
77 { return fTrackID; }
78
79 inline void G4Track::SetTrackID(const G4int aValue)
80 { fTrackID = aValue; }
81
82// position
83 inline const G4ThreeVector& G4Track::GetPosition() const
84 { return fPosition; }
85
86 inline void G4Track::SetPosition(const G4ThreeVector& aValue)
87 { fPosition = aValue; }
88
89// global time
90 inline G4double G4Track::GetGlobalTime() const
91 { return fGlobalTime; }
92
93 inline void G4Track::SetGlobalTime(const G4double aValue)
94 { fGlobalTime = aValue; }
95 // Time since the event in which the track belongs is created.
96
97// local time
98 inline G4double G4Track::GetLocalTime() const
99 { return fLocalTime; }
100
101 inline void G4Track::SetLocalTime(const G4double aValue)
102 { fLocalTime = aValue; }
103 // Time since the current track is created.
104
105// proper time
106 inline G4double G4Track::GetProperTime() const
107 { return fpDynamicParticle->GetProperTime(); }
108
109 inline void G4Track::SetProperTime(const G4double aValue)
110 { fpDynamicParticle->SetProperTime(aValue); }
111 // Proper time of the current track
112
113// volume
114 inline G4VPhysicalVolume* G4Track::GetVolume() const
115 { if ( fpTouchable ==0 ) return 0;
116 return fpTouchable->GetVolume(); }
117
118 inline G4VPhysicalVolume* G4Track::GetNextVolume() const
119 { if ( fpNextTouchable ==0 ) return 0;
120 return fpNextTouchable->GetVolume(); }
121
122// material
123 inline
124 const G4MaterialCutsCouple* G4Track::GetMaterialCutsCouple() const
125 { return fpStep->GetPreStepPoint()->GetMaterialCutsCouple(); }
126
127 inline
128 const G4MaterialCutsCouple* G4Track::GetNextMaterialCutsCouple() const
129 { return fpStep->GetPostStepPoint()->GetMaterialCutsCouple(); }
130
131// material
132 inline G4Material* G4Track::GetMaterial() const
133 { return fpStep->GetPreStepPoint()->GetMaterial(); }
134
135 inline G4Material* G4Track::GetNextMaterial() const
136 { return fpStep->GetPostStepPoint()->GetMaterial(); }
137
138
139// touchable
140 inline const G4VTouchable* G4Track::GetTouchable() const
141 {
142 return fpTouchable();
143 }
144
145 inline const G4TouchableHandle& G4Track::GetTouchableHandle() const
146 {
147 return fpTouchable;
148 }
149
150 inline void G4Track::SetTouchableHandle( const G4TouchableHandle& apValue)
151 {
152 fpTouchable = apValue;
153 }
154
155 inline const G4VTouchable* G4Track::GetNextTouchable() const
156 {
157 return fpNextTouchable();
158 }
159
160 inline const G4TouchableHandle& G4Track::GetNextTouchableHandle() const
161 {
162 return fpNextTouchable;
163 }
164
165 inline void G4Track::SetNextTouchableHandle( const G4TouchableHandle& apValue)
166 {
167 fpNextTouchable = apValue;
168 }
169
170
171// kinetic energy
172 inline G4double G4Track::GetKineticEnergy() const
173 { return fpDynamicParticle->GetKineticEnergy(); }
174
175 inline void G4Track::SetKineticEnergy(const G4double aValue)
176 { fpDynamicParticle->SetKineticEnergy(aValue); }
177
178// total energy
179 inline G4double G4Track::GetTotalEnergy() const
180 { return fpDynamicParticle->GetTotalEnergy(); }
181
182// momentum
183 inline G4ThreeVector G4Track::GetMomentum() const
184 { return fpDynamicParticle->GetMomentum(); }
185
186// momentum (direction)
187 inline const G4ThreeVector& G4Track::GetMomentumDirection() const
188 { return fpDynamicParticle->GetMomentumDirection(); }
189
190 inline void G4Track::SetMomentumDirection(const G4ThreeVector& aValue)
191 { fpDynamicParticle->SetMomentumDirection(aValue) ;}
192
193// polarization
194 inline const G4ThreeVector& G4Track::GetPolarization() const
195 { return fpDynamicParticle->GetPolarization(); }
196
197 inline void G4Track::SetPolarization(const G4ThreeVector& aValue)
198 { fpDynamicParticle->SetPolarization(aValue.x(),
199 aValue.y(),
200 aValue.z()); }
201
202// track status
203 inline G4TrackStatus G4Track::GetTrackStatus() const
204 { return fTrackStatus; }
205
206 inline void G4Track::SetTrackStatus(const G4TrackStatus aTrackStatus)
207 { fTrackStatus = aTrackStatus; }
208
209// track length
210 inline G4double G4Track::GetTrackLength() const
211 { return fTrackLength; }
212
213 inline void G4Track::AddTrackLength(const G4double aValue)
214 { fTrackLength += aValue; }
215 // Accumulated track length
216
217// step number
218 inline G4int G4Track::GetCurrentStepNumber() const
219 { return fCurrentStepNumber; }
220
221 inline void G4Track::IncrementCurrentStepNumber()
222 { fCurrentStepNumber++; }
223
224// step length
225 inline G4double G4Track::GetStepLength() const
226 { return fStepLength; }
227
228 inline void G4Track::SetStepLength(G4double value)
229 { fStepLength = value; }
230
231// vertex (where this track was created) information
232 inline const G4ThreeVector& G4Track::GetVertexPosition() const
233 { return fVtxPosition; }
234
235 inline void G4Track::SetVertexPosition(const G4ThreeVector& aValue)
236 { fVtxPosition = aValue; }
237
238 inline const G4ThreeVector& G4Track::GetVertexMomentumDirection() const
239 { return fVtxMomentumDirection; }
240 inline void G4Track::SetVertexMomentumDirection(const G4ThreeVector& aValue)
241 { fVtxMomentumDirection = aValue ;}
242
243 inline G4double G4Track::GetVertexKineticEnergy() const
244 { return fVtxKineticEnergy; }
245
246 inline void G4Track::SetVertexKineticEnergy(const G4double aValue)
247 { fVtxKineticEnergy = aValue; }
248
249 inline const G4LogicalVolume* G4Track::GetLogicalVolumeAtVertex() const
250 { return fpLVAtVertex; }
251
252 inline void G4Track::SetLogicalVolumeAtVertex(const G4LogicalVolume* aValue)
253 { fpLVAtVertex = aValue; }
254
255 inline const G4VProcess* G4Track::GetCreatorProcess() const
256 { return fpCreatorProcess; }
257 // If the pointer is 0, this means the track is created
258 // by the event generator, i.e. the primary track.If it is not
259 // 0, it points to the process which created this track.
260
261 inline void G4Track::SetCreatorProcess(const G4VProcess* aValue)
262 { fpCreatorProcess = aValue; }
263
264// flag for "Below Threshold"
265 inline G4bool G4Track::IsBelowThreshold() const
266 { return fBelowThreshold; }
267
268 inline void G4Track::SetBelowThresholdFlag(G4bool value)
269 { fBelowThreshold = value; }
270
271// flag for " Good for Tracking"
272 inline G4bool G4Track::IsGoodForTracking() const
273 { return fGoodForTracking; }
274
275 inline void G4Track::SetGoodForTrackingFlag(G4bool value)
276 { fGoodForTracking = value; }
277
278// track weight
279 inline void G4Track::SetWeight(G4double aValue)
280 { fWeight = aValue; }
281
282 inline G4double G4Track::GetWeight() const
283 { return fWeight; }
284
285// user information
286 inline G4VUserTrackInformation* G4Track::GetUserInformation() const
287 { return fpUserInformation; }
288 inline void G4Track::SetUserInformation(G4VUserTrackInformation* aValue)
289 { fpUserInformation = aValue; }
290
291
292
293//-------------------------------------------------------------
294// To implement bi-directional association between G4Step and
295// and G4Track, a combined usage of 'forward declaration' and
296// 'include' is necessary.
297//-------------------------------------------------------------
298#include "G4Step.hh"
299
300 inline const G4Step* G4Track::GetStep() const
301 { return fpStep; }
302
303 inline void G4Track::SetStep(const G4Step* aValue)
304 { fpStep = (aValue); }
305
306
307
308
309
Note: See TracBrowser for help on using the repository browser.