Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VoxelNavigation.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// G4VoxelNavigation
27//
28// Class description:
29//
30// Utility for navigation in volumes containing only G4PVPlacement
31// daughter volumes for which voxels have been constructed.
32
33// Author: Paul Kent (CERN), August 1996
34// --------------------------------------------------------------------
35#ifndef G4VOXELNAVIGATION_HH
36#define G4VOXELNAVIGATION_HH 1
37
38#include "geomdefs.hh"
39#include "G4VNavigation.hh"
41#include "G4NavigationLogger.hh"
42#include "G4AffineTransform.hh"
43#include "G4VPhysicalVolume.hh"
44#include "G4LogicalVolume.hh"
45#include "G4VSolid.hh"
46#include "G4ThreeVector.hh"
47
48#include "G4BlockingList.hh"
49
50class G4VoxelSafety;
51
52// Required for inline implementation
53//
55
56// Required for voxel handling & voxel stack
57//
58#include <vector>
59#include "G4SmartVoxelProxy.hh"
60#include "G4SmartVoxelNode.hh"
61#include "G4SmartVoxelHeader.hh"
62
63/**
64 * @brief G4VoxelNavigation is a concrete utility class for navigation in
65 * volumes containing only G4PVPlacement daughter volumes for which voxels
66 * have been constructed.
67 */
68
70{
71 public:
72
73 /**
74 * Constructor and Destructor.
75 */
77 ~G4VoxelNavigation() override;
78
79 /**
80 * Locates voxel node based on given point.
81 * @param[in] pHead Pointer to header of nodes to look through.
82 * @param[in] localPoint Local point
83 * @returns Pointer to the node where the given point is located.
84 */
86 const G4ThreeVector& localPoint );
87
88 /**
89 * Searches positioned volumes in mother at current top level of @p history
90 * for volume containing @p globalPoint. Do not test against @p blockedVol.
91 * If a containing volume is found, push it onto navigation history state.
92 * @param[in,out] history Navigation history.
93 * @param[in,out] blockedVol Blocked volume to be ignored in queries.
94 * @param[in,out] blockedNum Copy number for blocked replica volumes.
95 * @param[in,out] globalPoint Point in global coordinates system.
96 * @param[in,out] globalDirection Pointer to global direction or null.
97 * @param[in] pLocatedOnEdge Flag specifying if point is located on edge.
98 * @param[in,out] localPoint Point in local coordinates system.
99 * @returns Whether a containing volume has been found.
100 */
102 const G4VPhysicalVolume* blockedVol,
103 const G4int blockedNum,
104 const G4ThreeVector& globalPoint,
105 const G4ThreeVector* globalDirection,
106 const G4bool pLocatedOnEdge,
107 G4ThreeVector& localPoint ) override;
108
109 /**
110 * Computes the length of a step to the next boundary.
111 * Does not test against @p pBlockedPhysical. Identifies the next candidate
112 * volume (if a daughter of the current volume), and returns it in:
113 * pBlockedPhysical, blockedReplicaNo.
114 * @param[in] localPoint Local point.
115 * @param[in] localDirection Local direction vector.
116 * @param[in] currentProposedStepLength Current proposed step length.
117 * @param[in,out] newSafety New safety.
118 * @param[in,out] history Navigation history.
119 * @param[in,out] validExitNormal Flag to indicate whether exit normal is
120 * valid or not.
121 * @param[in,out] exitNormal Exit normal.
122 * @param[in,out] exiting Flag to indicate whether exiting a volume.
123 * @param[in,out] entering Flag to indicate whether entering a volume.
124 * @param[in,out] pBlockedPhysical Blocked physical volume that should be
125 * ignored in queries.
126 * @param[in,out] blockedReplicaNo Copy number for blocked replica volumes.
127 * @returns Length from current point to next boundary surface along
128 * @p localDirection.
129 */
130 G4double ComputeStep( const G4ThreeVector& localPoint,
131 const G4ThreeVector& localDirection,
132 const G4double currentProposedStepLength,
133 G4double& newSafety,
134 G4NavigationHistory& history,
135 G4bool& validExitNormal,
136 G4ThreeVector& exitNormal,
137 G4bool& exiting,
138 G4bool& entering,
139 G4VPhysicalVolume* (*pBlockedPhysical),
140 G4int& blockedReplicaNo ) override;
141
142 /**
143 * Calculates the isotropic distance to the nearest boundary from the
144 * specified point in the local coordinate system.
145 * The localpoint utilised must be within the current volume.
146 * @param[in] localPoint Local point.
147 * @param[in] history Navigation history.
148 * @param[in] pMaxLength Maximum step length beyond which volumes
149 * need not be checked.
150 * @returns Length from current point to closest surface.
151 */
152 G4double ComputeSafety( const G4ThreeVector& localpoint,
153 const G4NavigationHistory& history,
154 const G4double pMaxLength = DBL_MAX ) override;
155
156 /**
157 * Updates internal navigation state to take into account that location
158 * has been moved, but remains within the @p motherPhysical volume.
159 * @param[in] motherPhysical Current physical volume.
160 * @param[in] localPoint Local point.
161 */
162 void RelocateWithinVolume( G4VPhysicalVolume* motherPhysical,
163 const G4ThreeVector& localPoint ) override;
164
165 /**
166 * Verbosity control.
167 * @note If level>0 && G4VERBOSE, printout can occur.
168 */
169 inline G4int GetVerboseLevel() const override;
170 void SetVerboseLevel(G4int level) override;
171
172 /**
173 * Enables best-possible evaluation of isotropic safety.
174 */
175 inline void EnableBestSafety( G4bool flag = false );
176
177 protected:
178
179 /**
180 * Computes safety from specified point to voxel boundaries using already
181 * located point.
182 * @param[in] localPoint Local point.
183 * @returns Safety length from current point to voxel boundary.
184 */
185 G4double ComputeVoxelSafety( const G4ThreeVector& localPoint ) const;
186
187 /**
188 * Finds the next voxel from the current voxel and point in the specified
189 * direction.
190 * @param[in] localPoint Local point.
191 * @param[in] localDirection Direction along which compute the distance.
192 * @param[in] currentStep Current step size.
193 * @returns false if all voxels considered
194 * [current Step ends inside same voxel or leaves all voxels]
195 * true otherwise
196 * [the information on the next voxel is saved].
197 */
198 G4bool LocateNextVoxel( const G4ThreeVector& localPoint,
199 const G4ThreeVector& localDirection,
200 const G4double currentStep );
201
202 private: // Logging functions
203
204 void PreComputeStepLog (const G4VPhysicalVolume* motherPhysical,
205 G4double motherSafety,
206 const G4ThreeVector& localPoint);
207 void AlongComputeStepLog(const G4VSolid* sampleSolid,
208 const G4ThreeVector& samplePoint,
209 const G4ThreeVector& sampleDirection,
210 const G4ThreeVector& localDirection,
211 G4double sampleSafety,
212 G4double sampleStep);
213 void PostComputeStepLog (const G4VSolid* motherSolid,
214 const G4ThreeVector& localPoint,
215 const G4ThreeVector& localDirection,
216 G4double motherStep,
217 G4double motherSafety);
218 void ComputeSafetyLog (const G4VSolid* solid,
219 const G4ThreeVector& point,
220 G4double safety,
221 G4bool banner);
222 inline void PrintDaughterLog (const G4VSolid* sampleSolid,
223 const G4ThreeVector& samplePoint,
224 G4double sampleSafety,
225 G4double sampleStep);
226 protected:
227
228 /** Blocked volumes. */
230
231 // -----------------------------------------------------------------------
232 // BEGIN Voxel Stack information
233
234 /** Voxels depth.
235 * @note fVoxelDepth==0+ => fVoxelAxisStack(0+) contains axes of voxel
236 * fVoxelDepth==-1 -> not in voxel.
237 */
239
240 /** Voxel axes. */
241 std::vector<EAxis> fVoxelAxisStack;
242
243 /** No slices per voxel at each level. */
244 std::vector<G4int> fVoxelNoSlicesStack;
245
246 /** Width of voxels at each level. */
247 std::vector<G4double> fVoxelSliceWidthStack;
248
249 /** Node no point is inside at each level. */
250 std::vector<G4int> fVoxelNodeNoStack;
251
252 /** Voxel headers at each level. */
253 std::vector<G4SmartVoxelHeader*> fVoxelHeaderStack;
254
255 /** Node containing last located point. */
257
258 // END Voxel Stack information
259 // -----------------------------------------------------------------------
260
261 /** Helper object for Voxel Safety. */
263
264 /** Surface tolerance. */
266
267 /** Flag for best safety. */
269
270 /** Verbosity logger. */
272};
273
274#include "G4VoxelNavigation.icc"
275
276#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4BlockingList is an utility class responsible for (efficiently) maintaining a list of blocked volume...
G4NavigationHistory is a class responsible for the maintenance of the history of the path taken throu...
G4NavigationLogger is a simple utility class for use by the navigation systems for verbosity and chec...
G4SmartVoxelHeader represents a set of voxels, created by a single axis of virtual division....
G4SmartVoxelNode defines a node in the smart voxel hierarchy, i.e. a 'slice' of space along a given a...
G4VNavigation class holds the common navigation interface for all geometry navigator types.
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
G4double ComputeStep(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo) override
void SetVerboseLevel(G4int level) override
G4NavigationLogger * fLogger
void EnableBestSafety(G4bool flag=false)
std::vector< G4double > fVoxelSliceWidthStack
std::vector< EAxis > fVoxelAxisStack
G4int GetVerboseLevel() const override
G4VoxelSafety * fpVoxelSafety
G4double ComputeSafety(const G4ThreeVector &localpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX) override
G4SmartVoxelNode * fVoxelNode
G4bool LocateNextVoxel(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentStep)
G4SmartVoxelNode * VoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
void RelocateWithinVolume(G4VPhysicalVolume *motherPhysical, const G4ThreeVector &localPoint) override
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint) override
std::vector< G4int > fVoxelNodeNoStack
std::vector< G4SmartVoxelHeader * > fVoxelHeaderStack
std::vector< G4int > fVoxelNoSlicesStack
~G4VoxelNavigation() override
G4double ComputeVoxelSafety(const G4ThreeVector &localPoint) const
G4VoxelSafety is an utility class for the handling isotropic safety in volumes containing only G4PVPl...
#define DBL_MAX
Definition templates.hh:62