Geant4  10.02.p02
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 94218 2015-11-09 08:24:48Z 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"
59 
61 
62 typedef std::vector<G4int> G4SelectedAtRestDoItVector;
63 typedef std::vector<G4int> G4SelectedAlongStepDoItVector;
64 typedef std::vector<G4int> G4SelectedPostStepDoItVector;
65 typedef std::vector<G4int> G4SelectedPostStepAtTimeDoItVector;
66 
67 class G4Trajectory_Lock;
68 class G4Track;
69 struct G4ProcessState_Lock;
71 class G4SaveNavigatorState_Lock;
72 struct G4ITNavigatorState_Lock;
73 
81 {
82 public:
85 
86  //________________________________________________
91  inline bool IsLeadingStep()
92  {
93  return fStepLeader;
94  }
95  inline void SetLeadingStep(bool value)
96  {
97  fStepLeader = value;
98  }
99 
100  //________________________________________________
105  G4shared_ptr<G4ProcessState_Lock> GetProcessState(size_t index);
106 
107  inline void RecordProcessState(G4shared_ptr<G4ProcessState_Lock>,
108  size_t index);
109 
110  //___________________________________________________
111 
114 
115  /*
116  std::map<int,G4VTrackStateHandle> fTrackStates;
117  std::map<void*,G4VTrackStateHandle> fMultipleTrackStates;
118 
119  void SetTrackState(void* adress, G4VTrackStateHandle state)
120  {
121  fMultipleTrackStates[adress] = state;
122  }
123  G4VTrackStateHandle GetTrackState(void* adress)
124  {
125  return fMultipleTrackStates[adress];
126  }
127 
128  void SetTrackState(G4VTrackStateHandle state)
129  {
130  fTrackStates[state->GetID()] = state;
131  }
132  template<typename T> G4VTrackStateHandle GetTrackState()
133  {
134  return fTrackStates[G4TrackStateID<T>::GetID()] ;
135  }
136  */
137 
139  {
140  return fTrackStateManager;
141  }
142  /*
143  G4TrackStateManager& GetTrackStateManager() const
144  {
145  return fTrackStateManager;
146  }
147  */
149  {
150  return fpTrajectory_Lock;
151  }
152 
153  inline void SetTrajectory_Lock(G4Trajectory_Lock* trajLock)
154  {
155  fpTrajectory_Lock = trajLock;
156  }
157 
159  inline const G4ThreeVector& GetPreStepPosition() const;
160  inline G4double GetPreStepLocalTime() const;
161  inline G4double GetPreStepGlobalTime() const;
162 
163  inline void SetNavigatorState(G4ITNavigatorState_Lock *);
164  inline G4ITNavigatorState_Lock* GetNavigatorState() const;
165 
166  //-------------
167 protected:
168  //-------------
169  friend class G4ITStepProcessor;
170  //_______________________________________________________
172  //_______________________________________________________
174 
176 
177  //_______________________________________________________
181 
182  //_______________________________________________________
183  G4ITNavigatorState_Lock* fNavigatorState;
184 // G4SaveNavigatorState_Lock* fNavigatorState;
185 
186 //_______________________________________________________
191 // std::vector<G4ProcessState_Lock*> fProcessState;
192  std::vector<G4shared_ptr<G4ProcessState_Lock> > fProcessState;
193 
194  //_______________________________________________________
196 
197  //_______________________________________________________
202 
208 };
209 
210 inline
212 {
213  fpStepProcessorState = state;
214 }
215 
217 {
218  return fpStepProcessorState;
219 }
220 /*
221  inline void G4TrackingInformation::RecordProcessState(G4ProcessState_Lock* state,
222  size_t index)
223  {
224  // G4cout << "G4TrackingInformation::RecordProcessState" << G4endl;
225  fProcessState[index] = state;
226  }*/
227 
228 inline
229 void G4TrackingInformation::RecordProcessState(G4shared_ptr<G4ProcessState_Lock> state,
230  size_t index)
231 {
232  // G4cout << "G4TrackingInformation::RecordProcessState" << G4endl;
233  fProcessState[index] = state;
234 }
235 
237 {
239 }
240 
242 {
244 }
245 
247 {
248  return fRecordedTrackPosition;
249 }
250 
251 inline void G4TrackingInformation::SetNavigatorState(G4ITNavigatorState_Lock* state)
252 {
253  // G4cout << "Set Navigator state : " << state << G4endl;
254  fNavigatorState = state;
255 }
256 
257 inline G4ITNavigatorState_Lock* G4TrackingInformation::GetNavigatorState() const
258 {
259  return fNavigatorState;
260 }
261 
262 #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
void SetLeadingStep(bool value)
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.