Geant4  10.01.p03
G4TrackingInformation.hh
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 // $Id: G4TrackingInformation.hh 90769 2015-06-09 10:33:41Z gcosmo $
27 //
28 // Author: Mathieu Karamitros, kara@cenbg.in2p3.fr
29 
30 // The code is developed in the framework of the ESA AO7146
31 //
32 // We would be very happy hearing from you, send us your feedback! :)
33 //
34 // In order for Geant4-DNA to be maintained and still open-source,
35 // article citations are crucial.
36 // If you use Geant4-DNA chemistry and you publish papers about your software,
37 // in addition to the general paper on Geant4-DNA:
38 //
39 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
40 //
41 // we would be very happy if you could please also cite the following
42 // reference papers on chemistry:
43 //
44 // J. Comput. Phys. 274 (2014) 841-882
45 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508
46 
47 #ifndef G4TRACKINGINFORMATION_HH
48 #define G4TRACKINGINFORMATION_HH
49 
50 #include "globals.hh"
51 #include <vector>
52 #include <map>
53 #include "G4StepStatus.hh"
54 #include "G4ThreeVector.hh"
55 #include "G4TouchableHandle.hh"
56 #include "G4TrackState.hh"
57 #include "G4memory.hh"
58 
60 
61 typedef std::vector<G4int> G4SelectedAtRestDoItVector;
62 typedef std::vector<G4int> G4SelectedAlongStepDoItVector;
63 typedef std::vector<G4int> G4SelectedPostStepDoItVector;
64 typedef std::vector<G4int> G4SelectedPostStepAtTimeDoItVector;
65 
66 class G4Trajectory_Lock;
67 class G4Track;
68 struct G4ProcessState_Lock;
70 class G4SaveNavigatorState_Lock;
71 struct G4ITNavigatorState_Lock;
72 
74 {
75  friend class G4TrackingInformation;
76 protected:
78  {
79  ;
80  }
81 };
82 
90 {
91 public:
94 
95  //________________________________________________
100  inline bool IsLeadingStep()
101  {
102  return fStepLeader;
103  }
104  inline void SetLeadingStep(bool value)
105  {
106  fStepLeader = value;
107  }
108 
109  //________________________________________________
114  G4shared_ptr<G4ProcessState_Lock> GetProcessState(size_t index);
115 
116  inline void RecordProcessState(G4shared_ptr<G4ProcessState_Lock>,
117  size_t index);
118 
119  //___________________________________________________
120 
123 
124  /*
125  std::map<int,G4VTrackStateHandle> fTrackStates;
126  std::map<void*,G4VTrackStateHandle> fMultipleTrackStates;
127 
128  void SetTrackState(void* adress, G4VTrackStateHandle state)
129  {
130  fMultipleTrackStates[adress] = state;
131  }
132  G4VTrackStateHandle GetTrackState(void* adress)
133  {
134  return fMultipleTrackStates[adress];
135  }
136 
137  void SetTrackState(G4VTrackStateHandle state)
138  {
139  fTrackStates[state->GetID()] = state;
140  }
141  template<typename T> G4VTrackStateHandle GetTrackState()
142  {
143  return fTrackStates[G4TrackStateID<T>::GetID()] ;
144  }
145  */
146 
148 
150  {
151  return fTrackStateManager;
152  }
153  /*
154  G4TrackStateManager& GetTrackStateManager() const
155  {
156  return fTrackStateManager;
157  }
158  */
160  {
161  return fpTrajectory_Lock;
162  }
163 
164  inline void SetTrajectory_Lock(G4Trajectory_Lock* trajLock)
165  {
166  fpTrajectory_Lock = trajLock;
167  }
168 
170  inline const G4ThreeVector& GetPreStepPosition() const;
171  inline G4double GetPreStepLocalTime() const;
172  inline G4double GetPreStepGlobalTime() const;
173 
174  inline void SetNavigatorState(G4ITNavigatorState_Lock *);
175  inline G4ITNavigatorState_Lock* GetNavigatorState() const;
176 
177  //-------------
178 protected:
179  //-------------
180  friend class G4ITStepProcessor;
181  //_______________________________________________________
183  //_______________________________________________________
185 
186  //_______________________________________________________
190 
191  //_______________________________________________________
192  G4ITNavigatorState_Lock* fNavigatorState;
193 // G4SaveNavigatorState_Lock* fNavigatorState;
194 
195 //_______________________________________________________
200 // std::vector<G4ProcessState_Lock*> fProcessState;
201  std::vector<G4shared_ptr<G4ProcessState_Lock> > fProcessState;
202 
203  //_______________________________________________________
205 
206  //_______________________________________________________
211 
217 };
218 
219 inline
221 {
222  fpStepProcessorState = state;
223 }
224 
226 {
227  return fpStepProcessorState;
228 }
229 /*
230  inline void G4TrackingInformation::RecordProcessState(G4ProcessState_Lock* state,
231  size_t index)
232  {
233  // G4cout << "G4TrackingInformation::RecordProcessState" << G4endl;
234  fProcessState[index] = state;
235  }*/
236 
237 inline
238 void G4TrackingInformation::RecordProcessState(G4shared_ptr<G4ProcessState_Lock> state,
239  size_t index)
240 {
241  // G4cout << "G4TrackingInformation::RecordProcessState" << G4endl;
242  fProcessState[index] = state;
243 }
244 
246 {
248 }
249 
251 {
253 }
254 
256 {
257  return fRecordedTrackPosition;
258 }
259 
260 inline void G4TrackingInformation::SetNavigatorState(G4ITNavigatorState_Lock* state)
261 {
262  // G4cout << "Set Navigator state : " << state << G4endl;
263  fNavigatorState = state;
264 }
265 
266 inline G4ITNavigatorState_Lock* G4TrackingInformation::GetNavigatorState() const
267 {
268  return fNavigatorState;
269 }
270 
271 #endif // G4TRACKINGINFORMATION_HH
bool IsLeadingStep()
If the track is the one having the minimum step time, then it "leads" the step.
Its role is the same as G4StepManager :
G4Trajectory_Lock * fpTrajectory_Lock
G4double GetPreStepGlobalTime() const
CLHEP::Hep3Vector G4ThreeVector
void SetStepProcessorState(G4ITStepProcessorState_Lock *)
std::vector< G4int > G4SelectedAtRestDoItVector
void RecordCurrentPositionNTime(G4Track *)
G4TrackingInformation & operator=(const G4TrackingInformation &other)
Assignment operator.
The class G4TrackingInformation (hold by G4IT) emcompasses processes informations computed at the PS/...
std::vector< G4int > G4SelectedPostStepDoItVector
void SetTrajectory_Lock(G4Trajectory_Lock *trajLock)
bool G4bool
Definition: G4Types.hh:79
G4ITStepProcessorState_Lock * GetStepProcessorState()
G4double GetPreStepLocalTime() const
const G4ThreeVector & GetPreStepPosition() const
G4ITNavigatorState_Lock * GetNavigatorState() const
G4TrackStateManager fTrackStateManager
std::vector< G4int > G4SelectedAlongStepDoItVector
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
Every process should store the information computed at the InteractionLegth stage in the track...
G4ITNavigatorState_Lock * fNavigatorState
std::vector< G4int > G4SelectedPostStepAtTimeDoItVector
void SetNavigatorState(G4ITNavigatorState_Lock *)
G4ITStepProcessorState_Lock * fpStepProcessorState
G4Trajectory_Lock * GetTrajectory_Lock()
double G4double
Definition: G4Types.hh:76
void RecordProcessState(G4shared_ptr< G4ProcessState_Lock >, size_t index)
G4TrackStateManager & GetTrackStateManager()
std::vector< G4shared_ptr< G4ProcessState_Lock > > fProcessState
Holds the information related to processes Indexed on GetPhysIntVector (cf.