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
G4ReplicaNavigation.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 G4ReplicaNavigation
27//
28// Class description:
29//
30// Utility for navigation in volumes containing a single G4PVParameterised
31// volume for which voxels for the replicated volumes have been constructed.
32// [Voxels MUST be along one axis only: NOT refined]
33
34// History:
35// - Created: Paul Kent, Aug 96
36// --------------------------------------------------------------------
37#ifndef G4REPLICANAVIGATION_HH
38#define G4REPLICANAVIGATION_HH
39
41
42#include "G4Types.hh"
44#include "G4LogicalVolume.hh"
45#include "G4VPhysicalVolume.hh"
46#include "G4ThreeVector.hh"
47#include "G4BlockingList.hh"
48
49// Required for voxel handling
50//
51#include "G4SmartVoxelHeader.hh"
52
53class G4VSolid;
54
56{
57 // Bucket to hold value of Normal (3-vector),
58 // bools for calculated and leave-behind or 'validConvex',
59 // and exiting side.
60
62 // Identity of 'Side' of Replicas. Used by DistanceToOut methods.
63
65 G4bool calculated; // Normal
66 G4bool validConvex; // Solid locally convex
68
69 public:
70
72 G4bool calc = false,
73 G4bool valid= false,
74 ESide side = kNull )
75 { exitNormal= norm; calculated= calc; validConvex=valid; exitSide=side;}
76};
77
79{
80 public: // with description
81
84
86 const G4VPhysicalVolume* blockedVol,
87 const G4int blockedNum,
88 const G4ThreeVector& globalPoint,
89 const G4ThreeVector* globalDirection,
90 const G4bool pLocatedOnEdge,
91 G4ThreeVector& localPoint );
92
93 G4double ComputeStep( const G4ThreeVector& globalPoint,
94 const G4ThreeVector& globalDirection,
95 const G4ThreeVector& localPoint,
96 const G4ThreeVector& localDirection,
97 const G4double currentProposedStepLength,
98 G4double& newSafety,
99 G4NavigationHistory& history,
100 G4bool& validExitNormal,
101 G4bool& calculatedExitNormal,
102 G4ThreeVector &exitNormal,
103 G4bool& exiting,
104 G4bool& entering,
105 G4VPhysicalVolume* (*pBlockedPhysical),
106 G4int &blockedReplicaNo );
107
108 G4double ComputeSafety( const G4ThreeVector& globalPoint,
109 const G4ThreeVector& localPoint,
110 G4NavigationHistory& history,
111 const G4double pProposedMaxLength = DBL_MAX );
112
114 const G4ThreeVector& globalPoint,
115 G4ThreeVector& localPoint,
116 const G4bool& exiting,
117 G4bool& notKnownInside ) const;
118
119 void ComputeTransformation( const G4int replicaNo,
120 G4VPhysicalVolume* pVol,
121 G4ThreeVector& point ) const;
122 void ComputeTransformation( const G4int replicaNo,
123 G4VPhysicalVolume* pVol ) const;
124
125 EInside Inside( const G4VPhysicalVolume* pVol,
126 const G4int replicaNo,
127 const G4ThreeVector& localPoint ) const;
129 const G4int replicaNo,
130 const G4ThreeVector &localPoint ) const;
132 const G4int replicaNo,
133 const G4ThreeVector& localPoint,
134 const G4ThreeVector& localDirection,
135 G4ExitNormal& candidateNormal ) const;
136
137 inline G4int GetVerboseLevel() const;
138 inline void SetVerboseLevel(G4int level);
139 // Get/Set Verbose(ness) level.
140 // [if level>0 && G4VERBOSE, printout can occur]
141
142 inline void CheckMode(G4bool mode);
143 // Run navigation in "check-mode", therefore using additional
144 // verifications and more strict correctness conditions.
145 // Is effective only with G4VERBOSE set.
146
147 private:
148
149 inline G4int VoxelLocate( const G4SmartVoxelHeader* pHead,
150 const G4ThreeVector& localPoint,
151 const G4int blocked=-1 ) const;
152
153 G4double DistanceToOutPhi( const G4ThreeVector& localPoint,
154 const G4ThreeVector& localDirection,
155 const G4double width,
156 G4ExitNormal& foundNormal ) const;
157
158 G4double DistanceToOutRad( const G4ThreeVector& localPoint,
159 const G4ThreeVector& localDirection,
160 const G4double width,
161 const G4double offset,
162 const G4int replicaNo,
163 G4ExitNormal& foundNormal ) const;
164 inline void SetPhiTransformation( const G4double ang,
165 G4VPhysicalVolume* pVol=0 ) const;
166 private:
167
168 // Invariants - unaltered during navigation
169 // **********
170
171 G4bool fCheck = false;
172 G4int fVerbose = 0;
173 // Configuration parameters
174
175 G4double kCarTolerance, kRadTolerance, kAngTolerance,
176 halfkCarTolerance, halfkRadTolerance, halfkAngTolerance,
177 fMinStep;
178 // Local copy of constants
179};
180
181#include "G4ReplicaNavigation.icc"
182
183#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double DistanceToOut(const G4VPhysicalVolume *pVol, const G4int replicaNo, const G4ThreeVector &localPoint) const
G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4bool &calculatedExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
G4double ComputeSafety(const G4ThreeVector &globalPoint, const G4ThreeVector &localPoint, G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
void ComputeTransformation(const G4int replicaNo, G4VPhysicalVolume *pVol, G4ThreeVector &point) const
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
void CheckMode(G4bool mode)
EInside Inside(const G4VPhysicalVolume *pVol, const G4int replicaNo, const G4ThreeVector &localPoint) const
G4int GetVerboseLevel() const
void SetVerboseLevel(G4int level)
EInside BackLocate(G4NavigationHistory &history, const G4ThreeVector &globalPoint, G4ThreeVector &localPoint, const G4bool &exiting, G4bool &notKnownInside) const
EInside
Definition: geomdefs.hh:67
G4ThreeVector exitNormal
G4ExitNormal(G4ThreeVector norm=G4ThreeVector(0., 0., 0.), G4bool calc=false, G4bool valid=false, ESide side=kNull)
#define DBL_MAX
Definition: templates.hh:62