Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ParallelGeometriesLimiterProcess.hh
이 파일의 문서화 페이지로 가기
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 // $Id: G4ParallelGeometriesLimiterProcess.hh $
28 //
29 //
30 //---------------------------------------------------------------
31 //
32 // G4ParallelGeometriesLimiterProcess.hh
33 //
34 // Description:
35 // Process dedicated to limiting the step on boudaries of
36 // parallel geometries used for the generic biasing.
37 //
38 // History:
39 // Sep 16: Creation.
40 //
41 //---------------------------------------------------------------
42 
43 
44 #ifndef G4ParallelGeometriesLimiterProcess_hh
45 #define G4ParallelGeometriesLimiterProcess_hh
46 
47 #include "globals.hh"
48 #include "G4VProcess.hh"
49 #include "G4Step.hh"
50 #include "G4Navigator.hh"
51 #include "G4VPhysicalVolume.hh"
53 #include "G4FieldTrack.hh"
54 class G4PathFinder;
56 
57 
58 // Class Description:
59 // -- G4VProcess limiting the step on parallel geometries used by generic biasing.
60 
61 
63 {
64 public:
65  // --------------------------
66  // -- Constructor/Destructor:
67  // --------------------------
68  G4ParallelGeometriesLimiterProcess(const G4String& processName = "biasLimiter");
69 
70 public:
72  {}
73 
74  // ----------------------------------------------------
75  // -- Registration / deregistration of parallel worlds.
76  // -- Access to the list of registered parallel worlds.
77  // ----------------------------------------------------
78  void AddParallelWorld(const G4String& parallelWorldName);
79  void RemoveParallelWorld(const G4String& parallelWorldName);
80  // -- The list of registered parallel worlds:
81  const std::vector< G4VPhysicalVolume* >& GetParallelWorlds() const { return fParallelWorlds; }
82 
83  // -- Get the parallel world volume index in the list of world volumes handled
84  // -- by the process. This index can then be used to access current volume and step
85  // -- limitation status below.
86  // -- If the world passed is unknown to the process, -1 is returned.
87  G4int GetParallelWorldIndex( const G4VPhysicalVolume* parallelWorld ) const;
88  G4int GetParallelWorldIndex( G4String parallelWorldName ) const;
89 
90 
91  // ---------------------
92  // -- Active navigators:
93  // ---------------------
94  // -- The list of navigators handled by the process:
95  const std::vector< G4Navigator* >& GetActiveNavigators() const { return fParallelWorldNavigators; }
96  // -- The navigator used for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
97  // -- Note that no boundary checks are done on the index passed.
98  const G4Navigator* GetNavigator( G4int worldIndex ) const { return fParallelWorldNavigators[size_t(worldIndex)]; }
99  // ---------------------------------------------------
100  // -- Previous and current steps geometry information:
101  // ---------------------------------------------------
102  // -- The "switch" between the previous and current step is done in the PostStepGPIL
103  // -- The update on the current step is done:
104  // -- - in the PostStepGPIL for the volumes
105  // -- - in the AlongStepGPIL for the step limitations
106  // --
107  // -- The list of previous step and current step volumes:
108  const std::vector< const G4VPhysicalVolume* >& GetCurrentVolumes() const { return fCurrentVolumes; }
109  const std::vector< const G4VPhysicalVolume* >& GetPreviousVolumes() const { return fPreviousVolumes; }
110  // -- The current and previous volume for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
111  // -- Note that no boundary checks are done on the index passed.
112  const G4VPhysicalVolume* GetCurrentVolume( G4int worldIndex ) const { return fCurrentVolumes[size_t(worldIndex)]; }
113  const G4VPhysicalVolume* GetPreviousVolume( G4int worldIndex ) const { return fPreviousVolumes[size_t(worldIndex)]; }
114  // -- Flags telling about step limitation in previous and current step:
115  const std::vector< G4bool >& GetIsLimiting() const { return fParallelWorldIsLimiting; }
116  const std::vector< G4bool >& GetWasLimiting() const { return fParallelWorldWasLimiting; }
117  // -- The current and previous step limitation status for the passed parallel world index (obtained with GetParallelWorldIndex(...) above)
118  // -- Note that no boundary checks are done on the index passed.
119  G4bool GetIsLimiting( G4int worldIndex ) const { return fParallelWorldIsLimiting[size_t(worldIndex)]; }
120  G4bool GetWasLimiting( G4int worldIndex ) const { return fParallelWorldWasLimiting[size_t(worldIndex)]; }
121 
122 
123  // --------------------------------------------------------------
124  // From process interface
125  // --------------------------------------------------------------
126  // -- Start/End tracking:
127  void StartTracking(G4Track*);
128  void EndTracking();
129 
130  // --------------------
131  // -- PostStep methods:
132  // --------------------
133  // -- PostStepGPIL is used to collect up to date volumes in the parallel geometries:
135  // -- PostStepDoIt is not used (never called):
137  { return nullptr; }
138 
139  // ---------------------------------------------------------------------------
140  // -- Along step used for limiting the step on parallel geometries boundaries:
141  // ---------------------------------------------------------------------------
143  G4double previousStepSize,
144  G4double currentMinimumStep,
145  G4double& proposedSafety,
146  G4GPILSelection* selection);
148  const G4Step& step);
149 
150  // ---------------------------
151  // -- AtRest methods not used:
152  // ---------------------------
155  { return DBL_MAX; }
156 
158  { return nullptr; }
159 
160 
161  // --
162  virtual void SetProcessManager(const G4ProcessManager*);
163 
164 
165 
166 private:
167  std::vector< G4VPhysicalVolume* > fParallelWorlds;
168  std::vector< G4Navigator* > fParallelWorldNavigators;
169  std::vector< G4int > fParallelWorldNavigatorIndeces;
170  std::vector< G4double > fParallelWorldSafeties;
171  std::vector< G4bool > fParallelWorldIsLimiting;
172  std::vector< G4bool > fParallelWorldWasLimiting;
173  std::vector< const G4VPhysicalVolume* > fCurrentVolumes;
174  std::vector< const G4VPhysicalVolume* > fPreviousVolumes;
181 };
182 
183 #endif
G4ParallelGeometriesLimiterProcess(const G4String &processName="biasLimiter")
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *)
const std::vector< G4VPhysicalVolume * > & GetParallelWorlds() const
G4int GetParallelWorldIndex(const G4VPhysicalVolume *parallelWorld) const
const G4VPhysicalVolume * GetPreviousVolume(G4int worldIndex) const
const std::vector< G4Navigator * > & GetActiveNavigators() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
const std::vector< G4bool > & GetIsLimiting() const
void AddParallelWorld(const G4String &parallelWorldName)
const G4VPhysicalVolume * GetCurrentVolume(G4int worldIndex) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
const std::vector< const G4VPhysicalVolume * > & GetCurrentVolumes() const
virtual void SetProcessManager(const G4ProcessManager *)
std::vector< G4VPhysicalVolume * > fParallelWorlds
G4GPILSelection
const std::vector< G4bool > & GetWasLimiting() const
Definition: G4Step.hh:76
int G4int
Definition: G4Types.hh:78
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4ForceCondition
std::vector< const G4VPhysicalVolume * > fPreviousVolumes
std::vector< const G4VPhysicalVolume * > fCurrentVolumes
const std::vector< const G4VPhysicalVolume * > & GetPreviousVolumes() const
#define DBL_MAX
Definition: templates.hh:83
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
const G4Navigator * GetNavigator(G4int worldIndex) const
void RemoveParallelWorld(const G4String &parallelWorldName)