Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RegularNavigation Class Reference

#include <G4RegularNavigation.hh>

Public Member Functions

 G4RegularNavigation ()
 
 ~G4RegularNavigation ()
 
G4bool LevelLocate (G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
 
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)
 
G4double ComputeStepSkippingEqualMaterials (G4ThreeVector &localPoint, const G4ThreeVector &globalDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo, G4VPhysicalVolume *pCurrentPhysical)
 
G4double ComputeSafety (const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
 
void SetVerboseLevel (G4int level)
 
void CheckMode (G4bool mode)
 
void SetNormalNavigation (G4NormalNavigation *fnormnav)
 

Detailed Description

Definition at line 54 of file G4RegularNavigation.hh.

Constructor & Destructor Documentation

◆ G4RegularNavigation()

G4RegularNavigation::G4RegularNavigation ( )

Definition at line 46 of file G4RegularNavigation.cc.

47 : fverbose(false), fcheck(false), fnormalNav(0)
48{
50}
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()

◆ ~G4RegularNavigation()

G4RegularNavigation::~G4RegularNavigation ( )

Definition at line 54 of file G4RegularNavigation.cc.

55{
56}

Member Function Documentation

◆ CheckMode()

void G4RegularNavigation::CheckMode ( G4bool  mode)
inline

Definition at line 118 of file G4RegularNavigation.hh.

118{ fcheck = mode; }

◆ ComputeSafety()

G4double G4RegularNavigation::ComputeSafety ( const G4ThreeVector localPoint,
const G4NavigationHistory history,
const G4double  pProposedMaxLength = DBL_MAX 
)

Definition at line 266 of file G4RegularNavigation.cc.

269{
270 // This method is never called because to be called the daughter has to be a
271 // regular structure. This would only happen if the track is in the mother of
272 // voxels volume. But the voxels fill completely their mother, so when a
273 // track enters the mother it automatically enters a voxel. Only precision
274 // problems would make this method to be called
275
276 // Compute step in voxel
277 //
278 return fnormalNav->ComputeSafety(localPoint,
279 history,
280 pMaxLength );
281}
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)

Referenced by G4ITNavigator::ComputeSafety(), and G4Navigator::ComputeSafety().

◆ ComputeStep()

G4double G4RegularNavigation::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 
)

Definition at line 60 of file G4RegularNavigation.cc.

72{
73 // This method is never called because to be called the daughter has to be
74 // a regular structure. This would only happen if the track is in the mother
75 // of voxels volume. But the voxels fill completely their mother, so when a
76 // track enters the mother it automatically enters a voxel. Only precision
77 // problems would make this method to be called
78
79 G4ThreeVector globalPoint =
80 history.GetTopTransform().Inverse().TransformPoint(localPoint);
81 G4ThreeVector globalDirection =
82 history.GetTopTransform().Inverse().TransformAxis(localDirection);
83
84 G4ThreeVector localPoint2 = localPoint; // take away constantness
85
86 LevelLocate( history, *pBlockedPhysical, blockedReplicaNo,
87 globalPoint, &globalDirection, true, localPoint2 );
88
89
90 // Get in which voxel it is
91 //
92 G4VPhysicalVolume *motherPhysical, *daughterPhysical;
93 G4LogicalVolume *motherLogical;
94 motherPhysical = history.GetTopVolume();
95 motherLogical = motherPhysical->GetLogicalVolume();
96 daughterPhysical = motherLogical->GetDaughter(0);
97
98 G4PhantomParameterisation * daughterParam =
99 (G4PhantomParameterisation*)(daughterPhysical->GetParameterisation());
100 G4int copyNo = daughterParam ->GetReplicaNo(localPoint,localDirection);
101
102 G4ThreeVector voxelTranslation = daughterParam->GetTranslation( copyNo );
103 G4ThreeVector daughterPoint = localPoint - voxelTranslation;
104
105
106 // Compute step in voxel
107 //
108 return fnormalNav->ComputeStep(daughterPoint,
109 localDirection,
110 currentProposedStepLength,
111 newSafety,
112 history,
113 validExitNormal,
114 exitNormal,
115 exiting,
116 entering,
117 pBlockedPhysical,
118 blockedReplicaNo);
119}
int G4int
Definition: G4Types.hh:66
G4AffineTransform Inverse() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4VPhysicalVolume * GetDaughter(const G4int i) const
const G4AffineTransform & GetTopTransform() const
G4VPhysicalVolume * GetTopVolume() const
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)
virtual G4int GetReplicaNo(const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
G4ThreeVector GetTranslation(const G4int copyNo) const
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0

Referenced by G4Navigator::ComputeStep(), and G4ITNavigator::ComputeStep().

◆ ComputeStepSkippingEqualMaterials()

G4double G4RegularNavigation::ComputeStepSkippingEqualMaterials ( G4ThreeVector localPoint,
const G4ThreeVector globalDirection,
const G4double  currentProposedStepLength,
G4double newSafety,
G4NavigationHistory history,
G4bool validExitNormal,
G4ThreeVector exitNormal,
G4bool exiting,
G4bool entering,
G4VPhysicalVolume **  pBlockedPhysical,
G4int blockedReplicaNo,
G4VPhysicalVolume pCurrentPhysical 
)

Definition at line 123 of file G4RegularNavigation.cc.

136{
138
140 (G4PhantomParameterisation*)(pCurrentPhysical->GetParameterisation());
141
142 if( !param->SkipEqualMaterials() )
143 {
144 return fnormalNav->ComputeStep(localPoint,
145 localDirection,
146 currentProposedStepLength,
147 newSafety,
148 history,
149 validExitNormal,
150 exitNormal,
151 exiting,
152 entering,
153 pBlockedPhysical,
154 blockedReplicaNo);
155 }
156
157
158 G4double ourStep = 0.;
159
160 // To get replica No: transform local point to the reference system of the
161 // param container volume
162 //
163 G4int ide = history.GetDepth();
164 G4ThreeVector containerPoint = history.GetTransform(ide).Inverse().TransformPoint(localPoint);
165
166 // Point in global frame
167 //
168 containerPoint = history.GetTransform(ide).Inverse().TransformPoint(localPoint);
169
170 // Point in voxel parent volume frame
171 //
172 containerPoint = history.GetTransform(ide-1).TransformPoint(containerPoint);
173
174 // Store previous voxel translation to move localPoint by the difference
175 // with the new one
176 //
177 G4ThreeVector prevVoxelTranslation = containerPoint - localPoint;
178
179 // Do not use the expression below: There are cases where the
180 // fLastLocatedPointLocal does not give the correct answer
181 // (particle reaching a wall and bounced back, particle travelling through
182 // the wall that is deviated in an step, ...; these are pathological cases
183 // that give wrong answers in G4PhantomParameterisation::GetReplicaNo()
184 //
185 // G4ThreeVector prevVoxelTranslation = param->GetTranslation( copyNo );
186
187 G4int copyNo = param->GetReplicaNo(containerPoint,localDirection);
188
189 G4Material* currentMate = param->ComputeMaterial( copyNo, 0, 0 );
190 G4VSolid* voxelBox = pCurrentPhysical->GetLogicalVolume()->GetSolid();
191
192 G4VSolid* containerSolid = param->GetContainerSolid();
193 G4Material* nextMate;
194 G4bool bFirstStep = true;
195 G4double newStep;
196 G4double totalNewStep = 0.;
197
198 // Loop while same material is found
199 //
200 for( ;; )
201 {
202 newStep = voxelBox->DistanceToOut( localPoint, localDirection );
203
204 if( (bFirstStep) && (newStep < currentProposedStepLength) )
205 {
206 exiting = true;
207 }
208 bFirstStep = false;
209
210 newStep += kCarTolerance; // Avoid precision problems
211 ourStep += newStep;
212 totalNewStep += newStep;
213
214 // Physical process is limiting the step, don't continue
215 //
216 if(std::fabs(totalNewStep-currentProposedStepLength) < kCarTolerance)
217 {
218 return currentProposedStepLength;
219 }
220 if(totalNewStep > currentProposedStepLength)
221 {
223 AddStepLength(copyNo, newStep-totalNewStep+currentProposedStepLength);
224 return currentProposedStepLength;
225 }
226 else
227 {
229 }
230
231
232 // Move container point until wall of voxel
233 //
234 containerPoint += newStep*localDirection;
235 if( containerSolid->Inside( containerPoint ) != kInside )
236 {
237 break;
238 }
239
240 // Get copyNo and translation of new voxel
241 //
242 copyNo = param->GetReplicaNo(containerPoint,localDirection);
243 G4ThreeVector voxelTranslation = param->GetTranslation( copyNo );
244
245 // G4cout << " copyNo " << copyNo << " = " << pCurrentPhysical->GetCopyNo() << G4endl;
246 // Move local point until wall of voxel and then put it in the new voxel
247 // local coordinates
248 //
249 localPoint += newStep*localDirection;
250 localPoint += prevVoxelTranslation - voxelTranslation;
251
252 prevVoxelTranslation = voxelTranslation;
253
254 // Check if material of next voxel is the same as that of the current voxel
255 nextMate = param->ComputeMaterial( copyNo, 0, 0 );
256
257 if( currentMate != nextMate ) { break; }
258 }
259
260 return ourStep;
261}
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
G4VSolid * GetSolid() const
G4int GetDepth() const
const G4AffineTransform & GetTransform(G4int n) const
G4VSolid * GetContainerSolid() const
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
G4bool SkipEqualMaterials() const
static G4RegularNavigationHelper * Instance()
void AddStepLength(G4int copyNo, G4double slen)
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
@ kInside
Definition: geomdefs.hh:58

◆ LevelLocate()

G4bool G4RegularNavigation::LevelLocate ( G4NavigationHistory history,
const G4VPhysicalVolume blockedVol,
const G4int  blockedNum,
const G4ThreeVector globalPoint,
const G4ThreeVector globalDirection,
const G4bool  pLocatedOnEdge,
G4ThreeVector localPoint 
)

Definition at line 286 of file G4RegularNavigation.cc.

293{
294 G4VPhysicalVolume *motherPhysical, *pPhysical;
296 G4LogicalVolume *motherLogical;
297 G4ThreeVector localDir;
298 G4int replicaNo;
299
300 motherPhysical = history.GetTopVolume();
301 motherLogical = motherPhysical->GetLogicalVolume();
302
303 pPhysical = motherLogical->GetDaughter(0);
304 pParam = (G4PhantomParameterisation*)(pPhysical->GetParameterisation());
305
306 // Save parent history in touchable history
307 // ... for use as parent t-h in ComputeMaterial method of param
308 //
309 G4TouchableHistory parentTouchable( history );
310
311 // Get local direction
312 //
313 if( globalDirection )
314 {
315 localDir = history.GetTopTransform().TransformAxis(*globalDirection);
316 }
317 else
318 {
319 localDir = G4ThreeVector(0.,0.,0.);
320 }
321
322 // Enter this daughter
323 //
324 replicaNo = pParam->GetReplicaNo( localPoint, localDir );
325
326 if( replicaNo < 0 || replicaNo >= G4int(pParam->GetNoVoxel()) )
327 {
328 return false;
329 }
330
331 // Set the correct copy number in physical
332 //
333 pPhysical->SetCopyNo(replicaNo);
334 pParam->ComputeTransformation(replicaNo,pPhysical);
335
336 history.NewLevel(pPhysical, kParameterised, replicaNo );
337 localPoint = history.GetTopTransform().TransformPoint(globalPoint);
338
339 // Set the correct solid and material in Logical Volume
340 //
341 G4LogicalVolume *pLogical = pPhysical->GetLogicalVolume();
342
343 pLogical->UpdateMaterial(pParam->ComputeMaterial(replicaNo,
344 pPhysical, &parentTouchable) );
345 return true;
346}
CLHEP::Hep3Vector G4ThreeVector
void UpdateMaterial(G4Material *pMaterial)
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const
size_t GetNoVoxel() const
virtual void SetCopyNo(G4int CopyNo)=0
@ kParameterised
Definition: geomdefs.hh:68

Referenced by ComputeStep(), G4Navigator::LocateGlobalPointAndSetup(), and G4ITNavigator::LocateGlobalPointAndSetup().

◆ SetNormalNavigation()

void G4RegularNavigation::SetNormalNavigation ( G4NormalNavigation fnormnav)
inline

Definition at line 119 of file G4RegularNavigation.hh.

120 { fnormalNav = fnormnav; }

Referenced by G4ITNavigator::G4ITNavigator(), and G4Navigator::G4Navigator().

◆ SetVerboseLevel()

void G4RegularNavigation::SetVerboseLevel ( G4int  level)
inline

Definition at line 117 of file G4RegularNavigation.hh.

117{ fverbose = level; }

The documentation for this class was generated from the following files: