Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 "G4Step.hh" // Include from 'track'
62#include "G4StepPoint.hh" // Include from 'track'
63#include "G4TouchableHandle.hh" // Include from 'geometry'
64
66#include "G4ITLeadingTracks.hh"
67
68#include <vector>
69
70class G4ITNavigator;
73class G4IT;
76class G4VITProcess;
78class G4ITTrackHolder;
79using G4SelectedAtRestDoItVector = std::vector<int, std::allocator<int>>;
80using G4SelectedAlongStepDoItVector = std::vector<int, std::allocator<int>>;
81using G4SelectedPostStepDoItVector = std::vector<int, std::allocator<int>>;
82
83//________________________________________________
84//
85// Members related to ParticleDefinition and not
86// proper to a track
87//________________________________________________
89{
93
97 //
98 // Note: DoItVector has inverse order against GetPhysIntVector
99 // and SelectedPostStepDoItVector.
100 //
101 // * Max number of processes
102 std::size_t MAXofAtRestLoops;
105 // Maximum number of processes for each type of process
106 // These depend on the G4ParticleDefinition, so on the track
107
108 // * Transportation process
110};
111
112//________________________________________________
113//
114// Members proper to a track
115//________________________________________________
117{
118public:
120 ~G4ITStepProcessorState() override;
121
124
125 // * Max Number of Process
128
132
134
135 // * Safety
137 // This keeps the minimum safety value proposed by AlongStepGPILs.
140 // To get the true safety value at the PostStepPoint, you have
141 // to subtract the distance to 'endpointSafOrigin' from this value.
142
144};
145
146/**
147 * Its role is the same as G4StepManager :
148 * - Find the minimum physical length and corresponding time step
149 * - Step one track BUT on a given time step.
150 */
151
153{
154friend class G4Scheduler;
155public:
157 virtual ~G4ITStepProcessor();
158
159 inline void SetPreviousStepTime(G4double);
160
162 {
163 return fpTrack;
164 }
165 inline G4Step* GetStep()
166 {
167 return fpStep;
168 }
169 inline const G4Step* GetStep() const
170 {
171 return fpStep;
172 }
173 inline void SetStep(G4Step* val)
174 {
175 fpStep = val;
176 }
177
179 {
180 return fpSecondary;
181 }
183 {
184 fpTrackingManager = trackMan;
185 }
187 {
188 return fpTrackingManager;
189 }
190
191 //___________________________________
192
193 virtual void Initialize();
195
196 void ResetLeadingTracks();
198
199 //___________________________________
200 G4double ComputeInteractionLength(double previousTimeStep);
203 {
204 return fILTimeStep;
205 }
206
207 //___________________________________
208 // DoIt
209 void DoIt(double timeStep);
210 void ExtractDoItData();
211 void Stepping(G4Track*, const double&);
213 //___________________________________
214
215 inline double GetInteractionTime();
216 inline const G4Track* GetTrack() const;
217 inline void CleanProcessor();
218
219 std::size_t GetAtRestDoItProcTriggered() const
220 {
221 return fAtRestDoItProcTriggered;
222 }
223
225 {
226 return fGPILSelection;
227 }
228
230 {
231 return fN2ndariesAlongStepDoIt;
232 }
233
235 {
236 return fN2ndariesAtRestDoIt;
237 }
238
240 {
241 return fN2ndariesPostStepDoIt;
242 }
243
245 {
246 return fpCurrentProcess;
247 }
248
250 {
251 return fPhysIntLength;
252 }
253
255 {
256 return fPostStepAtTimeDoItProcTriggered;
257 }
258
260 {
261 return fPostStepDoItProcTriggered;
262 }
263
265 {
266 return fpProcessInfo;
267 }
268
270 {
271 return fpState;
272 }
273
275 {
276 return fpParticleChange;
277 }
278
280 {
281 return fpCurrentVolume;
282 }
283
285 {
286 return fCondition;
287 }
288
289protected:
290
291 void ExtractILData();
292
294 void ClearProcessInfo();
295 void SetTrack(G4Track*);
296
297 void GetProcessInfo();
298
299 void SetupMembers();
300 void ResetSecondaries();
301 void InitDefineStep();
302
303 void SetInitialStep();
304
305 void GetAtRestIL();
307 void DoStepping();
308 void PushSecondaries();
309
310 // void CalculateStep();
311 // void DoCalculateStep();
312
313 // void CloneProcesses();
314 void ActiveOnlyITProcess();
316
321 void InvokePSDIP(std::size_t); //
323 void SetNavigator(G4ITNavigator *value);
325
326 // Return the estimated safety value at the PostStepPoint
328
331
332private:
333 //________________________________________________
334 //
335 // General members
336 //________________________________________________
337
338 G4bool fInitialized;
339
340 G4ITTrackingManager* fpTrackingManager;
341
343 // Cached geometrical tolerance on surface
344
345 G4ITNavigator* fpNavigator;
346 G4int fStoreTrajectory;
347 G4VITSteppingVerbose* fpVerbose;
348
349 G4ITTrackHolder* fpTrackContainer;
350 G4ITLeadingTracks fLeadingTracks;
351
352 //________________________________________________
353 //
354 // Members used as temporaries (= not proper to a track)
355 //________________________________________________
356
357 G4double fTimeStep; // not proper to a track
358 G4double fILTimeStep; // proper to a track ensemble
359
360 G4double fPreviousTimeStep;
361 G4TrackVector* fpSecondary; // get from fpStep at every configuration setup
362 G4VParticleChange* fpParticleChange;
363
364 G4VITProcess* fpCurrentProcess;
365 // The pointer to the process of which DoIt or
366 // GetPhysicalInteractionLength has been just executed
367
368 // * Secondaries
369 G4int fN2ndariesAtRestDoIt;
370 G4int fN2ndariesAlongStepDoIt;
371 G4int fN2ndariesPostStepDoIt;
372 // These are the numbers of secondaries generated by the process
373 // just executed.
374
375 // * Process selection
376 std::size_t fAtRestDoItProcTriggered;
377 std::size_t fPostStepDoItProcTriggered;
378 std::size_t fPostStepAtTimeDoItProcTriggered;
379 // Record the selected process
380
381 G4ForceCondition fCondition;
382 G4GPILSelection fGPILSelection;
383 // Above three variables are for the method
384 // DefinePhysicalStepLength(). To pass these information to
385 // the method Verbose, they are kept at here. Need a more
386 // elegant mechanism.
387
388 G4double fPhysIntLength;
389 // The minimum physical interaction length over all possible processes
390
391 // * Sensitive detector
392 // G4SteppingControl StepControlFlag;
393 // G4VSensitiveDetector* fpSensitive;
394
395 G4VPhysicalVolume* fpCurrentVolume;
396 // Get from fpStep or touchable, keep as member for user interface
397
398 //________________________________________________
399 //
400 // Members related to ParticleDefinition and not
401 // proper to a track
402 //________________________________________________
403
404 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> fProcessGeneralInfoMap;
405 ProcessGeneralInfo* fpProcessInfo;
406 G4ITTransportation* fpTransportation;
407
408 //________________________________________________
409 //
410 // Members used for setting up the processor
411 //________________________________________________
412
413 G4Track* fpTrack; // Set track
414 G4IT* fpITrack; // Set track
415 G4TrackingInformation* fpTrackingInfo; // Set track
416
417 G4ITStepProcessorState* fpState; // SetupMembers or InitDefineStep
418 G4Step* fpStep; // Set track or InitDefineStep
419
420 G4StepPoint* fpPreStepPoint; // SetupMembers
421 G4StepPoint* fpPostStepPoint; // SetupMembers
422};
423
424//______________________________________________________________________________
425
427{
428 fPreviousTimeStep = previousTimeStep;
429}
430
431//______________________________________________________________________________
432
434{
435 return fpTrack;
436}
437
438//______________________________________________________________________________
439
441{
442 return std::max(fpState->fEndpointSafety - (fpState->fEndpointSafOrigin
443 - fpPostStepPoint->GetPosition()).mag(),
445}
446
447//______________________________________________________________________________
448
449inline void G4ITStepProcessor::SetNavigator(G4ITNavigator *value)
450{
451 fpNavigator = value;
452}
453
454//______________________________________________________________________________
455
457{
458 fTimeStep = DBL_MAX;
459 fPhysIntLength = DBL_MAX;
460
461 fpState = nullptr;
462 fpTrack = nullptr;
463 fpTrackingInfo = nullptr;
464 fpITrack = nullptr;
465 fpStep = nullptr;
466 fpPreStepPoint = nullptr;
467 fpPostStepPoint = nullptr;
468
469 fpParticleChange = nullptr;
470
471 fpCurrentVolume = nullptr;
472 // fpSensitive = 0;
473
474 fpSecondary = nullptr;
475
476 fpTransportation = nullptr;
477
478 fpCurrentProcess= nullptr;
479 fpProcessInfo = nullptr;
480
481 fAtRestDoItProcTriggered = INT_MAX;
482 fPostStepDoItProcTriggered = INT_MAX;
483 fPostStepAtTimeDoItProcTriggered = INT_MAX;
484 fGPILSelection = NotCandidateForSelection;
485 fCondition = NotForced;
486}
487
488//______________________________________________________________________________
489
491{
492 return fTimeStep;
493}
494
495#endif // G4ITSTEPPROCESSOR_H
const G4double kCarTolerance
G4ForceCondition
@ NotForced
G4GPILSelection
@ NotCandidateForSelection
std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
std::vector< int, std::allocator< int > > G4SelectedAlongStepDoItVector
G4StepStatus
std::vector< G4Track * > G4TrackVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4TouchableHandle fTouchableHandle
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
G4ProcessVector * fpPostStepDoItVector
G4ProcessVector * fpPostStepGetPhysIntVector
std::size_t MAXofAlongStepLoops
G4ProcessVector * fpAlongStepGetPhysIntVector
G4ProcessVector * fpAlongStepDoItVector
G4ProcessVector * fpAtRestGetPhysIntVector
G4ITTransportation * fpTransportation
G4ProcessVector * fpAtRestDoItVector
#define INT_MAX
Definition templates.hh:90
#define DBL_MAX
Definition templates.hh:62