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
G4ITStepProcessor.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: G4ITStepProcessor.hh 64057 2012-10-30 15:04:49Z gcosmo $
27//
28// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
29//
30// WARNING : This class is released as a prototype.
31// It might strongly evolve or even disapear in the next releases.
32//
33// History:
34// -----------
35// 10 Oct 2011 M.Karamitros created
36//
37// -------------------------------------------------------------------
38
39#ifndef G4ITSTEPPROCESSOR_H
40#define G4ITSTEPPROCESSOR_H
41
42#include "G4ios.hh" // Include from 'system'
43#include "globals.hh" // Include from 'global'
44#include "Randomize.hh" // Include from 'global'
45
46#include "G4LogicalVolume.hh" // Include from 'geometry'
47#include "G4VPhysicalVolume.hh" // Include from 'geometry'
48#include "G4ProcessManager.hh" // Include from 'piim'
49
50#include "G4Track.hh" // Include from 'track'
51#include "G4TrackVector.hh" // Include from 'track'
52#include "G4TrackStatus.hh" // Include from 'track'
53#include "G4StepStatus.hh" // Include from 'track'
54//#include "G4UserSteppingAction.hh" // Include from 'tracking'
55//#include "G4UserTrackingAction.hh" // Include from 'tracking'
56#include "G4Step.hh" // Include from 'track'
57#include "G4StepPoint.hh" // Include from 'track'
58#include "G4TouchableHandle.hh" // Include from 'geometry'
59#include "G4TouchableHistoryHandle.hh" // Include from 'geometry'
60
61
63
64//class G4Navigator;
65class G4ITNavigator;
68class G4IT;
71class G4VITProcess;
72typedef class std::vector<int, std::allocator<int> > G4SelectedAtRestDoItVector;
73typedef class std::vector<int, std::allocator<int> > G4SelectedAlongStepDoItVector;
74typedef class std::vector<int, std::allocator<int> > G4SelectedPostStepDoItVector;
75
76
77/**
78 * Its role is the same as G4StepManager :
79 * - Find the minimum physical length and corresponding time step
80 * - Step one track BUT on a given time step.
81 */
82
84{
85
86public:
88 virtual ~G4ITStepProcessor();
89
90 inline void SetPreviousStepTime(G4double);
91
92 inline G4Track* GetTrack() {return fpTrack;}
93 inline G4Step* GetStep() {return fpStep;}
94 inline const G4Step* GetStep() const {return fpStep;}
95 inline void SetStep(G4Step* val) {fpStep = val;}
96
97 inline G4TrackVector* GetSecondaries() {return fpSecondary;}
98 inline void SetTrackingManager(G4ITTrackingManager* trackMan) {fpTrackingManager = trackMan;}
99 inline G4ITTrackingManager* GetTrackingManager() {return fpTrackingManager;}
100
101 virtual void Initialize();
103
105 void Stepping(G4Track*, const double&);
106 void CalculateStep(G4Track*, const double&);
108
109 void DoIt(G4Track*,double);
110
113
114 inline double GetInteractionTime();
115 inline const G4Track* GetTrack() const ;
116 inline void CleanProcessor();
117
118protected:
120 void ClearProcessInfo();
121 void SetTrack(G4Track*);
122
123 void GetProcessInfo();
124
125 void SetupMembers();
126 void ResetSecondaries();
127 void InitDefineStep();
128
129 void SetInitialStep();
130
131 void GetAtRestIL();
133 void DoStepping();
134
137
139 void ActiveOnlyITProcess();
141
146 void InvokePSDIP(size_t); //
148 void SetNavigator(G4ITNavigator *value);
150
151 // Return the estimated safety value at the PostStepPoint
153
154
157
158private:
159 //________________________________________________
160 //
161 // General members
162 //________________________________________________
163
164 G4bool fInitialized;
165
166 G4ITTrackingManager* fpTrackingManager;
167 // G4UserSteppingAction* fpUserSteppingAction;
168
169 G4double kCarTolerance;
170 // Cached geometrical tolerance on surface
171
172 G4ITNavigator* fpNavigator;
173// G4Navigator* fpNavigator;
174 G4int fStoreTrajectory;
175 G4int verboseLevel;
176
177 //________________________________________________
178 //
179 // Members used as temporaries (= not proper to a track)
180 //________________________________________________
181
182 G4double fTimeStep ; // not proper to a track
183 G4double fPreviousTimeStep;
184 G4TrackVector* fpSecondary ; // get from fpStep at every configuration setup
185 G4VParticleChange* fpParticleChange;
186
187 G4VITProcess* fpCurrentProcess;
188 // The pointer to the process of which DoIt or
189 // GetPhysicalInteractionLength has been just executed
190
191 // * Secondaries
192 G4int fN2ndariesAtRestDoIt;
193 G4int fN2ndariesAlongStepDoIt;
194 G4int fN2ndariesPostStepDoIt;
195 // These are the numbers of secondaries generated by the process
196 // just executed.
197
198 // * Process selection
199 size_t fAtRestDoItProcTriggered;
200 size_t fPostStepDoItProcTriggered;
201 size_t fPostStepAtTimeDoItProcTriggered;
202 // Record the selected process
203
204 G4ForceCondition fCondition;
205 G4GPILSelection fGPILSelection;
206 // Above three variables are for the method
207 // DefinePhysicalStepLength(). To pass these information to
208 // the method Verbose, they are kept at here. Need a more
209 // elegant mechanism.
210
211 G4double fPhysIntLength;
212 // The minimum physical interaction length over all possible processes
213
214 // * Sensitive detector
215// G4SteppingControl StepControlFlag;
216// G4VSensitiveDetector* fpSensitive;
217
218 G4VPhysicalVolume* fpCurrentVolume; // Get from fpStep or touchable, keep as member for user interface
219
220 //________________________________________________
221 //
222 // Members related to ParticleDefinition and not
223 // proper to a track
224 //________________________________________________
225 struct ProcessGeneralInfo
226 {
227 G4ProcessVector* fpAtRestDoItVector;
228 G4ProcessVector* fpAlongStepDoItVector;
229 G4ProcessVector* fpPostStepDoItVector;
230
231 G4ProcessVector* fpAtRestGetPhysIntVector;
232 G4ProcessVector* fpAlongStepGetPhysIntVector;
233 G4ProcessVector* fpPostStepGetPhysIntVector;
234 //
235 // Note: DoItVector has inverse order against GetPhysIntVector
236 // and SelectedPostStepDoItVector.
237 //
238 // * Max Number of Process
239 size_t MAXofAtRestLoops;
240 size_t MAXofAlongStepLoops;
241 size_t MAXofPostStepLoops;
242 // Maximum number of processes for each type of process
243 // These depend on the G4ParticleDefinition, so on the track
244
245 // * Transportation process
246 G4ITTransportation* fpTransportation ;
247 };
248
249 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> fProcessGeneralInfoMap;
250 ProcessGeneralInfo* fpProcessInfo;
251
252 G4ITTransportation* fpTransportation ;
253
254 //________________________________________________
255 //
256 // Members proper to a track
257 //________________________________________________
258 class G4ITStepProcessorState : public G4ITStepProcessorState_Lock
259 {
260 public:
261 G4ITStepProcessorState();
262 virtual ~G4ITStepProcessorState();
263
264 // * Max Number of Process
265 G4SelectedAtRestDoItVector fSelectedAtRestDoItVector;
266 G4SelectedPostStepDoItVector fSelectedPostStepDoItVector;
267
268 G4double fPhysicalStep;
269 G4double fPreviousStepSize;
270 G4double fSafety;
271
272 G4StepStatus fStepStatus;
273
274 // * Safety
275 G4double proposedSafety;
276 // This keeps the minimum safety value proposed by AlongStepGPILs.
277 G4ThreeVector endpointSafOrigin;
278 G4double endpointSafety;
279 // To get the true safety value at the PostStepPoint, you have
280 // to subtract the distance to 'endpointSafOrigin' from this value.
281
282 G4TouchableHandle fTouchableHandle;
283 private :
284 G4ITStepProcessorState(const G4ITStepProcessorState&);
285 G4ITStepProcessorState& operator=(const G4ITStepProcessorState&);
286 };
287
288 //________________________________________________
289 //
290 // Members used for configurating the processor
291 //________________________________________________
292
293 G4Track* fpTrack; // Set track
294 G4IT* fpITrack ; // Set track
295 G4TrackingInformation* fpTrackingInfo ; // Set track
296
297 G4ITStepProcessorState* fpState; // SetupMembers or InitDefineStep
298 G4Step* fpStep; // Set track or InitDefineStep
299
300 G4StepPoint* fpPreStepPoint; // SetupMembers
301 G4StepPoint* fpPostStepPoint; // SetupMembers
302};
303
305{
306 fPreviousTimeStep = previousTimeStep;
307}
308
310{
311 return fpTrack;
312}
313
315{
316 return std::max( fpState->endpointSafety -
317 (fpState->endpointSafOrigin - fpPostStepPoint->GetPosition()).mag(),
318 kCarTolerance );
319}
320
322{
323 fpNavigator = value;
324}
325
327{
328 fTimeStep = DBL_MAX ;
329 fPhysIntLength = DBL_MAX;
330
331 fpState = 0;
332 fpTrack = 0;
333 fpTrackingInfo = 0 ;
334 fpITrack = 0;
335 fpStep = 0;
336 fpPreStepPoint = 0;
337 fpPostStepPoint = 0;
338
339 fpParticleChange = 0;
340
341 fpCurrentVolume = 0;
342// fpSensitive = 0;
343
344 fpSecondary = 0 ;
345
346 fpTransportation = 0;
347
348 fpCurrentProcess= 0;
349 fpProcessInfo = 0;
350
351 fAtRestDoItProcTriggered = INT_MAX;
352 fPostStepDoItProcTriggered = INT_MAX;
353 fPostStepAtTimeDoItProcTriggered = INT_MAX;
354 fGPILSelection = NotCandidateForSelection ;
355 fCondition = NotForced;
356}
357
358//______________________________________________________________________________
360{
361 return fTimeStep ;
362}
363
364
365#endif // G4ITSTEPPROCESSOR_H
G4ForceCondition
@ NotForced
G4GPILSelection
@ NotCandidateForSelection
class std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
class std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
class std::vector< int, std::allocator< int > > G4SelectedAlongStepDoItVector
G4StepStatus
Definition: G4StepStatus.hh:51
std::vector< G4Track * > G4TrackVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void CalculateStep(G4Track *)
G4ITTrackingManager * GetTrackingManager()
void SetTrackingManager(G4ITTrackingManager *trackMan)
void DealWithSecondaries(G4int &)
void SetStep(G4Step *val)
void SetNavigator(G4ITNavigator *value)
void UpdateTrack(G4Track *)
void SetTrack(G4Track *)
G4ITStepProcessor & operator=(const G4ITStepProcessor &other)
void Stepping(G4Track *, const double &)
void SetPreviousStepTime(G4double)
void DoIt(G4Track *, double)
G4TrackVector * GetSecondaries()
void CalculateStep(G4Track *, const double &)
void DefinePhysicalStepLength(G4Track *)
const G4Step * GetStep() const
void SetupGeneralProcessInfo(G4ParticleDefinition *, G4ProcessManager *)
void ApplyProductionCut(G4Track *)
virtual void Initialize()
Definition: G4IT.hh:83
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:78
#define INT_MAX
Definition: templates.hh:111
#define DBL_MAX
Definition: templates.hh:83