Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParallelGeometriesLimiterProcess.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// G4ParallelGeometriesLimiterProcess
27//
28// Description:
29//
30// Process dedicated to limiting the step on boudaries of
31// parallel geometries used for the generic biasing.
32//
33// Author: Marc Verderi, September 2016.
34// --------------------------------------------------------------------
35#ifndef G4ParallelGeometriesLimiterProcess_hh
36#define G4ParallelGeometriesLimiterProcess_hh 1
37
38#include "globals.hh"
39#include "G4VProcess.hh"
40#include "G4Step.hh"
41#include "G4Navigator.hh"
42#include "G4VPhysicalVolume.hh"
44#include "G4FieldTrack.hh"
45
46class G4PathFinder;
48
50{
51 public:
52
53 // --------------------------
54 // -- Constructor/Destructor:
55 // --------------------------
56 G4ParallelGeometriesLimiterProcess(const G4String& processName = "biasLimiter");
57
59
60 // ----------------------------------------------------
61 // -- Registration / deregistration of parallel worlds.
62 // -- Access to the list of registered parallel worlds.
63 // ----------------------------------------------------
64 void AddParallelWorld(const G4String& parallelWorldName);
65 void RemoveParallelWorld(const G4String& parallelWorldName);
66 // -- The list of registered parallel worlds:
67 const std::vector< G4VPhysicalVolume* >& GetParallelWorlds() const
68 { return fParallelWorlds; }
69
70 // -- Get the parallel world volume index in the list of world volumes handled
71 // -- by the process. This index can then be used to access current volume and step
72 // -- limitation status below.
73 // -- If the world passed is unknown to the process, -1 is returned.
74 G4int GetParallelWorldIndex( const G4VPhysicalVolume* parallelWorld ) const;
75 G4int GetParallelWorldIndex( const G4String& parallelWorldName ) const;
76
77 // ---------------------
78 // -- Active navigators:
79 // ---------------------
80 // -- The list of navigators handled by the process:
81 const std::vector< G4Navigator* >& GetActiveNavigators() const
82 { return fParallelWorldNavigators; }
83 // -- The navigator used for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
84 // -- Note that no boundary checks are done on the index passed.
85 const G4Navigator* GetNavigator( G4int worldIndex ) const
86 { return fParallelWorldNavigators[std::size_t(worldIndex)]; }
87
88 // ---------------------------------------------------
89 // -- Previous and current steps geometry information:
90 // ---------------------------------------------------
91 // -- The "switch" between the previous and current step is done in the PostStepGPIL
92 // -- The update on the current step is done:
93 // -- - in the PostStepGPIL for the volumes
94 // -- - in the AlongStepGPIL for the step limitations
95 // --
96 // -- The list of previous step and current step volumes:
97 const std::vector< const G4VPhysicalVolume* >& GetCurrentVolumes() const
98 { return fCurrentVolumes; }
99 const std::vector< const G4VPhysicalVolume* >& GetPreviousVolumes() const
100 { return fPreviousVolumes; }
101 // -- The current and previous volume for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
102 // -- Note that no boundary checks are done on the index passed.
103 const G4VPhysicalVolume* GetCurrentVolume( G4int worldIndex ) const
104 { return fCurrentVolumes[std::size_t(worldIndex)]; }
105 const G4VPhysicalVolume* GetPreviousVolume( G4int worldIndex ) const
106 { return fPreviousVolumes[std::size_t(worldIndex)]; }
107 // -- Flags telling about step limitation in previous and current step:
108 const std::vector< G4bool >& GetIsLimiting() const
109 { return fParallelWorldIsLimiting; }
110 const std::vector< G4bool >& GetWasLimiting() const
111 { return fParallelWorldWasLimiting; }
112 // -- The current and previous step limitation status for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
113 // -- Note that no boundary checks are done on the index passed.
114 G4bool GetIsLimiting( G4int worldIndex ) const
115 { return fParallelWorldIsLimiting[std::size_t(worldIndex)]; }
116 G4bool GetWasLimiting( G4int worldIndex ) const
117 { return fParallelWorldWasLimiting[std::size_t(worldIndex)]; }
118
119 // --------------------------------------------------------------
120 // From process interface
121 // --------------------------------------------------------------
122 // -- Start/End tracking:
123 void StartTracking(G4Track*);
124 void EndTracking();
125
126 // --------------------
127 // -- PostStep methods:
128 // --------------------
129 // -- PostStepGPIL is used to collect up to date volumes in the parallel geometries:
131 // -- PostStepDoIt is not used (never called):
132 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step& ) { return nullptr; }
133
134 // ---------------------------------------------------------------------------
135 // -- Along step used for limiting the step on parallel geometries boundaries:
136 // ---------------------------------------------------------------------------
138 G4double previousStepSize, G4double currentMinimumStep,
139 G4double& proposedSafety, G4GPILSelection* selection);
140 G4VParticleChange* AlongStepDoIt(const G4Track& track, const G4Step& step);
141
142 // ---------------------------
143 // -- AtRest methods not used:
144 // ---------------------------
148 { return nullptr; }
149
150 // --
151 virtual void SetProcessManager(const G4ProcessManager*);
152
153 private:
154
155 std::vector< G4VPhysicalVolume* > fParallelWorlds;
156 std::vector< G4Navigator* > fParallelWorldNavigators;
157 std::vector< G4int > fParallelWorldNavigatorIndeces;
158 std::vector< G4double > fParallelWorldSafeties;
159 std::vector< G4bool > fParallelWorldIsLimiting;
160 std::vector< G4bool > fParallelWorldWasLimiting;
161 std::vector< const G4VPhysicalVolume* > fCurrentVolumes;
162 std::vector< const G4VPhysicalVolume* > fPreviousVolumes;
163 G4double fParallelWorldSafety = 0.0;
164 G4bool fIsTrackingTime = false;
165 G4FieldTrack fFieldTrack = '0';
166 G4ParticleChangeForNothing fDummyParticleChange;
167 G4PathFinder* fPathFinder = nullptr;
168 G4TransportationManager* fTransportationManager = nullptr;
169};
170
171#endif
G4ForceCondition
G4GPILSelection
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const std::vector< const G4VPhysicalVolume * > & GetPreviousVolumes() const
const std::vector< G4Navigator * > & GetActiveNavigators() const
G4ParallelGeometriesLimiterProcess(const G4String &processName="biasLimiter")
void AddParallelWorld(const G4String &parallelWorldName)
G4int GetParallelWorldIndex(const G4VPhysicalVolume *parallelWorld) const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual ~G4ParallelGeometriesLimiterProcess()=default
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *)
const G4Navigator * GetNavigator(G4int worldIndex) const
const std::vector< G4VPhysicalVolume * > & GetParallelWorlds() const
void RemoveParallelWorld(const G4String &parallelWorldName)
const G4VPhysicalVolume * GetCurrentVolume(G4int worldIndex) const
const G4VPhysicalVolume * GetPreviousVolume(G4int worldIndex) const
const std::vector< const G4VPhysicalVolume * > & GetCurrentVolumes() const
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
virtual void SetProcessManager(const G4ProcessManager *)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
const std::vector< G4bool > & GetWasLimiting() const
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
const std::vector< G4bool > & GetIsLimiting() const
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46
#define DBL_MAX
Definition templates.hh:62