Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27// Author: Mathieu Karamitros
28
29// The code is developed in the framework of the ESA AO7146
30//
31// We would be very happy hearing from you, send us your feedback! :)
32//
33// In order for Geant4-DNA to be maintained and still open-source,
34// article citations are crucial.
35// If you use Geant4-DNA chemistry and you publish papers about your software,
36// in addition to the general paper on Geant4-DNA:
37//
38// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39//
40// we would be very happy if you could please also cite the following
41// reference papers on chemistry:
42//
43// J. Comput. Phys. 274 (2014) 841-882
44// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45
46#ifndef G4ITSTEPPROCESSOR_H
47#define G4ITSTEPPROCESSOR_H
48
49#include "G4ios.hh" // Include from 'system'
50#include "globals.hh" // Include from 'global'
51#include "Randomize.hh" // Include from 'global'
52
53#include "G4LogicalVolume.hh" // Include from 'geometry'
54#include "G4VPhysicalVolume.hh" // Include from 'geometry'
55#include "G4ProcessManager.hh" // Include from 'piim'
56
57#include "G4Track.hh" // Include from 'track'
58#include "G4TrackVector.hh" // Include from 'track'
59#include "G4TrackStatus.hh" // Include from 'track'
60#include "G4StepStatus.hh" // Include from 'track'
61//#include "G4UserSteppingAction.hh" // Include from 'tracking'
62//#include "G4UserTrackingAction.hh" // Include from 'tracking'
63#include "G4Step.hh" // Include from 'track'
64#include "G4StepPoint.hh" // Include from 'track'
65#include "G4TouchableHandle.hh" // Include from 'geometry'
66#include "G4TouchableHistoryHandle.hh" // Include from 'geometry'
67
69#include "G4ITLeadingTracks.hh"
70
71class G4ITNavigator;
74class G4IT;
77class G4VITProcess;
79class G4ITTrackHolder;
80typedef class std::vector<int, std::allocator<int> > G4SelectedAtRestDoItVector;
81typedef class std::vector<int, std::allocator<int> > G4SelectedAlongStepDoItVector;
82typedef class std::vector<int, std::allocator<int> > G4SelectedPostStepDoItVector;
83
84//________________________________________________
85//
86// Members related to ParticleDefinition and not
87// proper to a track
88//________________________________________________
90{
94
98 //
99 // Note: DoItVector has inverse order against GetPhysIntVector
100 // and SelectedPostStepDoItVector.
101 //
102 // * Max number of processes
103 std::size_t MAXofAtRestLoops;
106 // Maximum number of processes for each type of process
107 // These depend on the G4ParticleDefinition, so on the track
108
109 // * Transportation process
111};
112
113//________________________________________________
114//
115// Members proper to a track
116//________________________________________________
118{
119public:
121 virtual ~G4ITStepProcessorState();
122
125
126 // * Max Number of Process
129
133
135
136 // * Safety
138 // This keeps the minimum safety value proposed by AlongStepGPILs.
141 // To get the true safety value at the PostStepPoint, you have
142 // to subtract the distance to 'endpointSafOrigin' from this value.
143
145};
146
147/**
148 * Its role is the same as G4StepManager :
149 * - Find the minimum physical length and corresponding time step
150 * - Step one track BUT on a given time step.
151 */
152
154{
155friend class G4Scheduler;
156public:
158 virtual ~G4ITStepProcessor();
159
160 inline void SetPreviousStepTime(G4double);
161
163 {
164 return fpTrack;
165 }
166 inline G4Step* GetStep()
167 {
168 return fpStep;
169 }
170 inline const G4Step* GetStep() const
171 {
172 return fpStep;
173 }
174 inline void SetStep(G4Step* val)
175 {
176 fpStep = val;
177 }
178
180 {
181 return fpSecondary;
182 }
184 {
185 fpTrackingManager = trackMan;
186 }
188 {
189 return fpTrackingManager;
190 }
191
192 //___________________________________
193
194 virtual void Initialize();
196
197 void ResetLeadingTracks();
199
200 //___________________________________
201 G4double ComputeInteractionLength(double previousTimeStep);
204 {
205 return fILTimeStep;
206 }
207
208 //___________________________________
209 // DoIt
210 void DoIt(double timeStep);
211 void ExtractDoItData();
212 void Stepping(G4Track*, const double&);
214 //___________________________________
215
216 inline double GetInteractionTime();
217 inline const G4Track* GetTrack() const;
218 inline void CleanProcessor();
219
220 std::size_t GetAtRestDoItProcTriggered() const
221 {
222 return fAtRestDoItProcTriggered;
223 }
224
226 {
227 return fGPILSelection;
228 }
229
231 {
232 return fN2ndariesAlongStepDoIt;
233 }
234
236 {
237 return fN2ndariesAtRestDoIt;
238 }
239
241 {
242 return fN2ndariesPostStepDoIt;
243 }
244
246 {
247 return fpCurrentProcess;
248 }
249
251 {
252 return fPhysIntLength;
253 }
254
256 {
257 return fPostStepAtTimeDoItProcTriggered;
258 }
259
261 {
262 return fPostStepDoItProcTriggered;
263 }
264
266 {
267 return fpProcessInfo;
268 }
269
271 {
272 return fpState;
273 }
274
276 {
277 return fpParticleChange;
278 }
279
281 {
282 return fpCurrentVolume;
283 }
284
286 {
287 return fCondition;
288 }
289
290protected:
291
292 void ExtractILData();
293
295 void ClearProcessInfo();
296 void SetTrack(G4Track*);
297
298 void GetProcessInfo();
299
300 void SetupMembers();
301 void ResetSecondaries();
302 void InitDefineStep();
303
304 void SetInitialStep();
305
306 void GetAtRestIL();
308 void DoStepping();
309 void PushSecondaries();
310
311 // void CalculateStep();
312 // void DoCalculateStep();
313
314 // void CloneProcesses();
315 void ActiveOnlyITProcess();
317
322 void InvokePSDIP(std::size_t); //
324 void SetNavigator(G4ITNavigator *value);
326
327 // Return the estimated safety value at the PostStepPoint
329
332
333private:
334 //________________________________________________
335 //
336 // General members
337 //________________________________________________
338
339 G4bool fInitialized;
340
341 G4ITTrackingManager* fpTrackingManager;
342
344 // Cached geometrical tolerance on surface
345
346 G4ITNavigator* fpNavigator;
347 G4int fStoreTrajectory;
348 G4VITSteppingVerbose* fpVerbose;
349
350 G4ITTrackHolder* fpTrackContainer;
351 G4ITLeadingTracks fLeadingTracks;
352
353 //________________________________________________
354 //
355 // Members used as temporaries (= not proper to a track)
356 //________________________________________________
357
358 G4double fTimeStep; // not proper to a track
359 G4double fILTimeStep; // proper to a track ensemble
360
361 G4double fPreviousTimeStep;
362 G4TrackVector* fpSecondary; // get from fpStep at every configuration setup
363 G4VParticleChange* fpParticleChange;
364
365 G4VITProcess* fpCurrentProcess;
366 // The pointer to the process of which DoIt or
367 // GetPhysicalInteractionLength has been just executed
368
369 // * Secondaries
370 G4int fN2ndariesAtRestDoIt;
371 G4int fN2ndariesAlongStepDoIt;
372 G4int fN2ndariesPostStepDoIt;
373 // These are the numbers of secondaries generated by the process
374 // just executed.
375
376 // * Process selection
377 std::size_t fAtRestDoItProcTriggered;
378 std::size_t fPostStepDoItProcTriggered;
379 std::size_t fPostStepAtTimeDoItProcTriggered;
380 // Record the selected process
381
382 G4ForceCondition fCondition;
383 G4GPILSelection fGPILSelection;
384 // Above three variables are for the method
385 // DefinePhysicalStepLength(). To pass these information to
386 // the method Verbose, they are kept at here. Need a more
387 // elegant mechanism.
388
389 G4double fPhysIntLength;
390 // The minimum physical interaction length over all possible processes
391
392 // * Sensitive detector
393 // G4SteppingControl StepControlFlag;
394 // G4VSensitiveDetector* fpSensitive;
395
396 G4VPhysicalVolume* fpCurrentVolume;
397 // Get from fpStep or touchable, keep as member for user interface
398
399 //________________________________________________
400 //
401 // Members related to ParticleDefinition and not
402 // proper to a track
403 //________________________________________________
404
405 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> fProcessGeneralInfoMap;
406 ProcessGeneralInfo* fpProcessInfo;
407 G4ITTransportation* fpTransportation;
408
409 //________________________________________________
410 //
411 // Members used for setting up the processor
412 //________________________________________________
413
414 G4Track* fpTrack; // Set track
415 G4IT* fpITrack; // Set track
416 G4TrackingInformation* fpTrackingInfo; // Set track
417
418 G4ITStepProcessorState* fpState; // SetupMembers or InitDefineStep
419 G4Step* fpStep; // Set track or InitDefineStep
420
421 G4StepPoint* fpPreStepPoint; // SetupMembers
422 G4StepPoint* fpPostStepPoint; // SetupMembers
423};
424
425//______________________________________________________________________________
426
428{
429 fPreviousTimeStep = previousTimeStep;
430}
431
432//______________________________________________________________________________
433
435{
436 return fpTrack;
437}
438
439//______________________________________________________________________________
440
442{
443 return std::max(fpState->fEndpointSafety - (fpState->fEndpointSafOrigin
444 - fpPostStepPoint->GetPosition()).mag(),
446}
447
448//______________________________________________________________________________
449
450inline void G4ITStepProcessor::SetNavigator(G4ITNavigator *value)
451{
452 fpNavigator = value;
453}
454
455//______________________________________________________________________________
456
458{
459 fTimeStep = DBL_MAX;
460 fPhysIntLength = DBL_MAX;
461
462 fpState = 0;
463 fpTrack = 0;
464 fpTrackingInfo = 0;
465 fpITrack = 0;
466 fpStep = 0;
467 fpPreStepPoint = 0;
468 fpPostStepPoint = 0;
469
470 fpParticleChange = 0;
471
472 fpCurrentVolume = 0;
473 // fpSensitive = 0;
474
475 fpSecondary = 0;
476
477 fpTransportation = 0;
478
479 fpCurrentProcess= 0;
480 fpProcessInfo = 0;
481
482 fAtRestDoItProcTriggered = INT_MAX;
483 fPostStepDoItProcTriggered = INT_MAX;
484 fPostStepAtTimeDoItProcTriggered = INT_MAX;
485 fGPILSelection = NotCandidateForSelection;
486 fCondition = NotForced;
487}
488
489//______________________________________________________________________________
490
492{
493 return fTimeStep;
494}
495
496#endif // G4ITSTEPPROCESSOR_H
const G4double kCarTolerance
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:40
std::vector< G4Track * > G4TrackVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4TouchableHandle fTouchableHandle
G4ThreeVector fEndpointSafOrigin
G4SelectedPostStepDoItVector fSelectedPostStepDoItVector
G4SelectedAtRestDoItVector fSelectedAtRestDoItVector
G4ITStepProcessorState & operator=(const G4ITStepProcessorState &)
G4int GetN2ndariesAtRestDoIt() const
G4ITTrackingManager * GetTrackingManager()
void SetTrackingManager(G4ITTrackingManager *trackMan)
void DealWithSecondaries(G4int &)
std::size_t GetAtRestDoItProcTriggered() const
void SetStep(G4Step *val)
void SetNavigator(G4ITNavigator *value)
G4int GetN2ndariesAlongStepDoIt() const
const ProcessGeneralInfo * GetCurrentProcessInfo() const
G4GPILSelection GetGPILSelection() const
void SetTrack(G4Track *)
G4ITStepProcessor & operator=(const G4ITStepProcessor &other)
G4ForceCondition GetCondition() const
void Stepping(G4Track *, const double &)
G4double ComputeInteractionLength(double previousTimeStep)
void SetPreviousStepTime(G4double)
G4int GetN2ndariesPostStepDoIt() const
G4TrackVector * GetSecondaries() const
void DefinePhysicalStepLength(G4Track *)
void DoIt(double timeStep)
const G4Step * GetStep() const
void SetupGeneralProcessInfo(G4ParticleDefinition *, G4ProcessManager *)
void ApplyProductionCut(G4Track *)
G4double GetPhysIntLength() const
const G4ITStepProcessorState * GetProcessorState() const
const G4VPhysicalVolume * GetCurrentVolume() const
std::size_t GetPostStepAtTimeDoItProcTriggered() const
const G4VParticleChange * GetParticleChange() const
virtual void Initialize()
const G4VITProcess * GetCurrentProcess() const
void InvokePSDIP(std::size_t)
std::size_t GetPostStepDoItProcTriggered() const
Definition: G4IT.hh:88
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:62
G4ProcessVector * fpPostStepDoItVector
G4ProcessVector * fpPostStepGetPhysIntVector
std::size_t MAXofAlongStepLoops
std::size_t MAXofAtRestLoops
G4ProcessVector * fpAlongStepGetPhysIntVector
G4ProcessVector * fpAlongStepDoItVector
G4ProcessVector * fpAtRestGetPhysIntVector
G4ITTransportation * fpTransportation
std::size_t MAXofPostStepLoops
G4ProcessVector * fpAtRestDoItVector
#define INT_MAX
Definition: templates.hh:90
#define DBL_MAX
Definition: templates.hh:62