Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VoxelSafety.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// G4VoxelSafety
27//
28// Class description:
29//
30// Utility for isotropic safety in volumes containing only G4PVPlacement
31// daughter volumes for which voxels have been constructed.
32// Implementation extracted and modified/adapted from G4VoxelNavigation class.
33
34// Author: John Apostolakis (CERN), 30 April 2010
35// --------------------------------------------------------------------
36#ifndef G4VOXELSAFETY_HH
37#define G4VOXELSAFETY_HH 1
38
39#include "geomdefs.hh"
41#include "G4AffineTransform.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4LogicalVolume.hh"
44#include "G4VSolid.hh"
45#include "G4ThreeVector.hh"
46
47#include "G4BlockingList.hh"
48
49#include <vector> // Required for voxel handling & voxel stack
50
53
54/**
55 * @brief G4VoxelSafety is an utility class for the handling isotropic safety
56 * in volumes containing only G4PVPlacement daughter volumes for which voxels
57 * have been constructed.
58 */
59
61{
62 public:
63
64 /**
65 * Constructor and Destructor.
66 */
69
70 /**
71 * Calculates the isotropic distance to the nearest boundary from the
72 * specified point in the local coordinate system.
73 * The localpoint utilised must be within the current volume.
74 * @param[in] localPoint Local point.
75 * @param[in] currentPhysical Current physical volume.
76 * @param[in] maxLength Maximum length beyond which volumes are not checked.
77 * @returns Isotropic distance of given point to closest surface.
78 */
79 G4double ComputeSafety( const G4ThreeVector& localPoint,
80 const G4VPhysicalVolume& currentPhysical,
81 G4double maxLength = DBL_MAX );
82
83 /**
84 * Verbosity control.
85 * @note If level>0 && G4VERBOSE, printout can occur.
86 */
87 inline G4int GetVerboseLevel() const { return fVerbose; }
88 inline void SetVerboseLevel(G4int level) { fVerbose = level; }
89
90 protected:
91
92 /**
93 * Cycles through levels of headers to process each node level.
94 * @param[in] pHead Voxel header.
95 * @param[in] localPoint Local point.
96 * @param[in] maxLength Maximum length beyond which volumes are not checked.
97 * @param[in] currentPhysical Current volume (used for debug printout).
98 * @param[in] distUpperDepth Upper square distance from voxel.
99 * @param[in] previousMinSafety Minimum distance beyond which not to look.
100 * @returns Isotropic distance of the point to closest volume in all nodes.
101 */
103 const G4ThreeVector& localPoint,
104 G4double maxLength,
105 const G4VPhysicalVolume& currentPhysical,
106 G4double distUpperDepth = 0.0,
107 G4double previousMinSafety = DBL_MAX );
108
109 /**
110 * Calculates the safety for volumes included in current Voxel Node.
111 * @param[in] curVoxelNode Voxel node.
112 * @param[in] localPoint Local point.
113 * @returns Isotropic distance of given point to closest volume in node.
114 */
115 G4double SafetyForVoxelNode( const G4SmartVoxelNode* curVoxelNode,
116 const G4ThreeVector& localPoint );
117
118 private:
119
120 // ---- BEGIN State - values used during computation of Safety ------------
121
122 /** Blocked volumes */
123 G4BlockingList fBlockList;
124
125 /** Cached pointer to mother logical volume */
126 G4LogicalVolume* fpMotherLogical = nullptr;
127
128 // ---- BEGIN Voxel Stack information -------------------------------------
129
130 /** Voxel depth
131 * @note fVoxelDepth==0+ => fVoxelAxisStack(0+) contains axes of voxel
132 * fVoxelDepth==-1 -> not in voxel.
133 */
134 G4int fVoxelDepth = -1;
135
136 /** Voxel axes */
137 std::vector<EAxis> fVoxelAxisStack;
138
139 /** No slices per voxel at each level */
140 std::vector<G4int> fVoxelNoSlicesStack;
141
142 /** Width of voxels at each level */
143 std::vector<G4double> fVoxelSliceWidthStack;
144
145 /** Node no point is inside at each level */
146 std::vector<G4int> fVoxelNodeNoStack;
147
148 /** Voxel headers at each level */
149 std::vector<const G4SmartVoxelHeader*> fVoxelHeaderStack;
150
151 // ----- END Voxel Stack information --------------------------------------
152
153 G4bool fCheck = false;
154 G4int fVerbose = 0;
156};
157
158#endif
const G4double kCarTolerance
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...
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
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...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4int GetVerboseLevel() const
void SetVerboseLevel(G4int level)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4VPhysicalVolume &currentPhysical, G4double maxLength=DBL_MAX)
G4double SafetyForVoxelHeader(const G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint, G4double maxLength, const G4VPhysicalVolume &currentPhysical, G4double distUpperDepth=0.0, G4double previousMinSafety=DBL_MAX)
G4double SafetyForVoxelNode(const G4SmartVoxelNode *curVoxelNode, const G4ThreeVector &localPoint)
#define DBL_MAX
Definition templates.hh:62