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
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// class 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// History:
34// - Created: Paul Kent, Aug 96
35// --------------------------------------------------------------------
36#ifndef G4VOXELNAVIGATION_HH
37#define G4VOXELNAVIGATION_HH
38
39#include "geomdefs.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
64{
65 public: // with description
66
68 virtual ~G4VoxelNavigation();
69
71 const G4ThreeVector& localPoint );
72
74 const G4VPhysicalVolume* blockedVol,
75 const G4int blockedNum,
76 const G4ThreeVector& globalPoint,
77 const G4ThreeVector* globalDirection,
78 const G4bool pLocatedOnEdge,
79 G4ThreeVector& localPoint );
80
81 virtual G4double ComputeStep( const G4ThreeVector& globalPoint,
82 const G4ThreeVector& globalDirection,
83 const G4double currentProposedStepLength,
84 G4double& newSafety,
85 G4NavigationHistory& history,
86 G4bool& validExitNormal,
87 G4ThreeVector& exitNormal,
88 G4bool& exiting,
89 G4bool& entering,
90 G4VPhysicalVolume* (*pBlockedPhysical),
91 G4int& blockedReplicaNo );
92
93 virtual G4double ComputeSafety( const G4ThreeVector& globalpoint,
94 const G4NavigationHistory& history,
95 const G4double pMaxLength = DBL_MAX );
96
97 inline G4int GetVerboseLevel() const;
98 void SetVerboseLevel(G4int level);
99 // Get/Set Verbose(ness) level.
100 // [if level>0 && G4VERBOSE, printout can occur]
101
102 inline void CheckMode(G4bool mode);
103 // Run navigation in "check-mode", therefore using additional
104 // verifications and more strict correctness conditions.
105 // Is effective only with G4VERBOSE set.
106
107 inline void EnableBestSafety( G4bool flag = false );
108 // Enable best-possible evaluation of isotropic safety
109
110 protected:
111
112 G4double ComputeVoxelSafety( const G4ThreeVector& localPoint ) const;
113 G4bool LocateNextVoxel( const G4ThreeVector& localPoint,
114 const G4ThreeVector& localDirection,
115 const G4double currentStep );
116
118 const G4ThreeVector& localPoint ) const;
119
120 private: // Logging functions
121
122 void PreComputeStepLog (const G4VPhysicalVolume* motherPhysical,
123 G4double motherSafety,
124 const G4ThreeVector& localPoint);
125 void AlongComputeStepLog(const G4VSolid* sampleSolid,
126 const G4ThreeVector& samplePoint,
127 const G4ThreeVector& sampleDirection,
128 const G4ThreeVector& localDirection,
129 G4double sampleSafety,
130 G4double sampleStep);
131 void PostComputeStepLog (const G4VSolid* motherSolid,
132 const G4ThreeVector& localPoint,
133 const G4ThreeVector& localDirection,
134 G4double motherStep,
135 G4double motherSafety);
136 void ComputeSafetyLog (const G4VSolid* solid,
137 const G4ThreeVector& point,
138 G4double safety,
139 G4bool banner);
140 inline void PrintDaughterLog (const G4VSolid* sampleSolid,
141 const G4ThreeVector& samplePoint,
142 G4double sampleSafety,
143 G4double sampleStep);
144 protected:
145
147 // Blocked volumes
148
149 //
150 // BEGIN Voxel Stack information
151 //
152
154 // Note: fVoxelDepth==0+ => fVoxelAxisStack(0+) contains axes of voxel
155 // fVoxelDepth==-1 -> not in voxel
156
157 std::vector<EAxis> fVoxelAxisStack;
158 // Voxel axes
159
160 std::vector<G4int> fVoxelNoSlicesStack;
161 // No slices per voxel at each level
162
163 std::vector<G4double> fVoxelSliceWidthStack;
164 // Width of voxels at each level
165
166 std::vector<G4int> fVoxelNodeNoStack;
167 // Node no point is inside at each level
168
169 std::vector<G4SmartVoxelHeader*> fVoxelHeaderStack;
170 // Voxel headers at each level
171
173 // Node containing last located point
174
175 //
176 // END Voxel Stack information
177 //
178
180 // Helper object for Voxel Safety
181
183 // Surface tolerance
184
185 G4bool fCheck = false;
187
189 // Verbosity logger
190};
191
192#include "G4VoxelNavigation.icc"
193
194#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4NavigationLogger * fLogger
void EnableBestSafety(G4bool flag=false)
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
virtual G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
std::vector< G4double > fVoxelSliceWidthStack
std::vector< EAxis > fVoxelAxisStack
void CheckMode(G4bool mode)
G4SmartVoxelNode * VoxelLocateLight(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint) const
G4VoxelSafety * fpVoxelSafety
G4SmartVoxelNode * fVoxelNode
G4bool LocateNextVoxel(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentStep)
G4SmartVoxelNode * VoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
virtual ~G4VoxelNavigation()
G4int GetVerboseLevel() const
virtual G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
std::vector< G4int > fVoxelNodeNoStack
G4BlockingList fBList
std::vector< G4SmartVoxelHeader * > fVoxelHeaderStack
std::vector< G4int > fVoxelNoSlicesStack
G4double ComputeVoxelSafety(const G4ThreeVector &localPoint) const
void SetVerboseLevel(G4int level)
#define DBL_MAX
Definition: templates.hh:62