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