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

#include <G4RegularNavigation.hh>

+ Inheritance diagram for G4RegularNavigation:

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) final
 
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) final
 
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) final
 
void SetNormalNavigation (G4NormalNavigation *fnormnav)
 
- Public Member Functions inherited from G4VNavigation
virtual ~G4VNavigation ()
 
virtual void RelocateWithinVolume (G4VPhysicalVolume *, const G4ThreeVector &)
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int level)
 
void CheckMode (G4bool mode)
 

Additional Inherited Members

- Protected Attributes inherited from G4VNavigation
G4int fVerbose = 0
 
G4bool fCheck = false
 

Detailed Description

Definition at line 51 of file G4RegularNavigation.hh.

Constructor & Destructor Documentation

◆ G4RegularNavigation()

G4RegularNavigation::G4RegularNavigation ( )

Definition at line 42 of file G4RegularNavigation.cc.

43{
45 fMinStep = 101*kCarTolerance;
46}
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()

◆ ~G4RegularNavigation()

G4RegularNavigation::~G4RegularNavigation ( )
default

Member Function Documentation

◆ ComputeSafety()

G4double G4RegularNavigation::ComputeSafety ( const G4ThreeVector & globalpoint,
const G4NavigationHistory & history,
const G4double pMaxLength = DBL_MAX )
finalvirtual

Compute the distance to the closest surface.

Parameters
[in]globalPointGlobal point.
[in]historyNavigation history.
[in]pMaxLengthMaximum step length beyond which volumes need not be checked.
Returns
Length from current point to closest surface.

Implements G4VNavigation.

Definition at line 345 of file G4RegularNavigation.cc.

348{
349 // This method is never called because to be called the daughter has to be a
350 // regular structure. This would only happen if the track is in the mother of
351 // voxels volume. But the voxels fill completely their mother, so when a
352 // track enters the mother it automatically enters a voxel. Only precision
353 // problems would make this method to be called
354
355 // Compute step in voxel
356 //
357 return fnormalNav->ComputeSafety(localPoint,
358 history,
359 pMaxLength );
360}
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX) final

Referenced by G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), and G4SafetyCalculator::SafetyInCurrentVolume().

◆ ComputeStep()

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

Compute the length of a step to the next boundary. Do not test against pBlockedPhysical. Identify the next candidate volume (if a daughter of current volume), and return it in pBlockedPhysical, blockedReplicaNo.

Parameters
[in]localPointLocal point
[in]localDirectionPointer to local direction or null pointer.
[in]currentProposedStepLengthCurrent proposed step length.
[in,out]newSafetyNew safety.
[in,out]historyNavigation history.
[in,out]validExitNormalFlag to indicate whether exit normal is valid or not.
[in,out]exitNormalExit normal.
[in,out]enteringFlag to indicate whether we are entering a volume.
[in,out]exitingFlag to indicate whether we are exiting a volume.
[in,out]pBlockedPhysicalBlocked physical volume that should be ignored in queries.
[in,out]blockedReplicaNoCopy number for blocked replica volumes.
Returns
Length from current point to next boundary surface along localDirection.

Implements G4VNavigation.

Definition at line 54 of file G4RegularNavigation.cc.

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

Referenced by G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), and G4Navigator::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 117 of file G4RegularNavigation.cc.

130{
132
133 auto param =
134 (G4PhantomParameterisation*)(pCurrentPhysical->GetParameterisation());
135
136 if( !param->SkipEqualMaterials() )
137 {
138 return fnormalNav->ComputeStep(localPoint,
139 localDirection,
140 currentProposedStepLength,
141 newSafety,
142 history,
143 validExitNormal,
144 exitNormal,
145 exiting,
146 entering,
147 pBlockedPhysical,
148 blockedReplicaNo);
149 }
150
151
152 G4double ourStep = 0.;
153
154 // To get replica No: transform local point to the reference system of the
155 // param container volume
156 //
157 auto ide = (G4int)history.GetDepth();
158 G4ThreeVector containerPoint = history.GetTransform(ide)
159 .InverseTransformPoint(localPoint);
160
161 // Point in global frame
162 //
163 containerPoint = history.GetTransform(ide).InverseTransformPoint(localPoint);
164
165 // Point in voxel parent volume frame
166 //
167 containerPoint = history.GetTransform(ide-1).TransformPoint(containerPoint);
168
169 // Store previous voxel translation to move localPoint by the difference
170 // with the new one
171 //
172 G4ThreeVector prevVoxelTranslation = containerPoint - localPoint;
173
174 // Do not use the expression below: There are cases where the
175 // fLastLocatedPointLocal does not give the correct answer
176 // (particle reaching a wall and bounced back, particle travelling through
177 // the wall that is deviated in an step, ...; these are pathological cases
178 // that give wrong answers in G4PhantomParameterisation::GetReplicaNo()
179 //
180 // G4ThreeVector prevVoxelTranslation = param->GetTranslation( copyNo );
181
182 G4int copyNo = param->GetReplicaNo(containerPoint,localDirection);
183
184 G4Material* currentMate = param->ComputeMaterial( copyNo, nullptr, nullptr );
185 G4VSolid* voxelBox = pCurrentPhysical->GetLogicalVolume()->GetSolid();
186
187 G4VSolid* containerSolid = param->GetContainerSolid();
188 G4Material* nextMate;
189 G4bool bFirstStep = true;
190 G4double newStep;
191 G4double totalNewStep = 0.;
192
193 // Loop while same material is found
194 //
195 //
196 fNumberZeroSteps = 0;
197 for( G4int ii = 0; ii < fNoStepsAllowed+1; ++ii )
198 {
199 if( ii == fNoStepsAllowed ) {
200 // Must kill this stuck track
201 //
202 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
203 .InverseTransformPoint(localPoint);
204 std::ostringstream message;
205 message << "G4RegularNavigation::ComputeStepSkippingEqualMaterials()"
206 << "Stuck Track: potential geometry or navigation problem."
207 << G4endl
208 << " Track stuck, moving for more than "
209 << ii << " steps" << G4endl
210 << "- at point " << pGlobalpoint << G4endl
211 << " local direction: " << localDirection << G4endl;
212 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
213 "GeomRegNav1001",
215 message);
216 }
217 newStep = voxelBox->DistanceToOut( localPoint, localDirection );
218 fLastStepWasZero = (newStep<fMinStep);
219 if( fLastStepWasZero )
220 {
221 ++fNumberZeroSteps;
222#ifdef G4DEBUG_NAVIGATION
223 if( fNumberZeroSteps > 1 )
224 {
225 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
226 .InverseTransformPoint(localPoint);
227 std::ostringstream message;
228 message.precision(16);
229 message << "G4RegularNavigation::ComputeStepSkippingEqualMaterials(): another 'zero' step, # "
230 << fNumberZeroSteps
231 << ", at " << pGlobalpoint
232 << ", nav-comp-step calls # " << ii
233 << ", Step= " << newStep;
234 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
235 "GeomRegNav1002", JustWarning, message,
236 "Potential overlap in geometry!");
237 }
238#endif
239 if( fNumberZeroSteps > fActionThreshold_NoZeroSteps-1 )
240 {
241 // Act to recover this stuck track. Pushing it along direction
242 //
243 newStep = std::min(101*kCarTolerance*std::pow(10,fNumberZeroSteps-2),0.1);
244#ifdef G4DEBUG_NAVIGATION
245 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
246 .InverseTransformPoint(localPoint);
247 std::ostringstream message;
248 message.precision(16);
249 message << "Track stuck or not moving." << G4endl
250 << " Track stuck, not moving for "
251 << fNumberZeroSteps << " steps" << G4endl
252 << "- at point " << pGlobalpoint
253 << " (local point " << localPoint << ")" << G4endl
254 << " local direction: " << localDirection
255 << " Potential geometry or navigation problem !"
256 << G4endl
257 << " Trying pushing it of " << newStep << " mm ...";
258 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
259 "GeomRegNav1003", JustWarning, message,
260 "Potential overlap in geometry!");
261#endif
262 }
263 if( fNumberZeroSteps > fAbandonThreshold_NoZeroSteps-1 )
264 {
265 // Must kill this stuck track
266 //
267 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
268 .InverseTransformPoint(localPoint);
269 std::ostringstream message;
270 message << "G4RegularNavigation::ComputeStepSkippingEqualMaterials()"
271 << "Stuck Track: potential geometry or navigation problem."
272 << G4endl
273 << " Track stuck, not moving for "
274 << fNumberZeroSteps << " steps" << G4endl
275 << "- at point " << pGlobalpoint << G4endl
276 << " local direction: " << localDirection << G4endl;
277 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
278 "GeomRegNav1004",
280 message);
281 }
282 }
283 else
284 {
285 // reset the zero step counter when a non-zero step was performed
286 fNumberZeroSteps = 0;
287 }
288 if( (bFirstStep) && (newStep < currentProposedStepLength) )
289 {
290 exiting = true;
291 }
292 bFirstStep = false;
293
294 newStep += kCarTolerance; // Avoid precision problems
295 ourStep += newStep;
296 totalNewStep += newStep;
297
298 // Physical process is limiting the step, don't continue
299 //
300 if(std::fabs(totalNewStep-currentProposedStepLength) < kCarTolerance)
301 {
302 return currentProposedStepLength;
303 }
304 if(totalNewStep > currentProposedStepLength)
305 {
307 AddStepLength(copyNo, newStep-totalNewStep+currentProposedStepLength);
308 return currentProposedStepLength;
309 }
311
312 // Move container point until wall of voxel
313 //
314 containerPoint += newStep*localDirection;
315 if( containerSolid->Inside( containerPoint ) != kInside )
316 {
317 break;
318 }
319
320 // Get copyNo and translation of new voxel
321 //
322 copyNo = param->GetReplicaNo(containerPoint, localDirection);
323 G4ThreeVector voxelTranslation = param->GetTranslation( copyNo );
324
325 // Move local point until wall of voxel and then put it in the new voxel
326 // local coordinates
327 //
328 localPoint += newStep*localDirection;
329 localPoint += prevVoxelTranslation - voxelTranslation;
330
331 prevVoxelTranslation = voxelTranslation;
332
333 // Check if material of next voxel is the same as that of the current voxel
334 nextMate = param->ComputeMaterial( copyNo, nullptr, nullptr );
335
336 if( currentMate != nextMate ) { break; }
337 }
338
339 return ourStep;
340}
@ JustWarning
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4VSolid * GetSolid() const
std::size_t GetDepth() const
const G4AffineTransform & GetTransform(G4int n) 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=nullptr, G4ThreeVector *n=nullptr) const =0
@ kInside
Definition geomdefs.hh:70

◆ LevelLocate()

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

Search positioned volumes in mother at current top level of history for volume containing globalPoint. Do not test against blockedVol. If a containing volume is found, push it onto navigation history state.

Parameters
[in,out]historyNavigation history.
[in,out]blockedVolBlocked volume that should be ignored in queries.
[in,out]blockedNumCopy number for blocked replica volumes.
[in,out]globalPointGlobal point
[in,out]globalDirectionPointer to global direction or null pointer.
[in,out]localPoint= global point in local system on entry, point in new system on exit.
Returns
Whether a containing volume has been found.

Implements G4VNavigation.

Definition at line 365 of file G4RegularNavigation.cc.

372{
373 G4VPhysicalVolume *motherPhysical, *pPhysical;
375 G4LogicalVolume *motherLogical;
376 G4ThreeVector localDir;
377 G4int replicaNo;
378
379 motherPhysical = history.GetTopVolume();
380 motherLogical = motherPhysical->GetLogicalVolume();
381
382 pPhysical = motherLogical->GetDaughter(0);
383 pParam = (G4PhantomParameterisation*)(pPhysical->GetParameterisation());
384
385 // Save parent history in touchable history
386 // ... for use as parent t-h in ComputeMaterial method of param
387 //
388 G4TouchableHistory parentTouchable( history );
389
390 // Get local direction
391 //
392 if( globalDirection != nullptr )
393 {
394 localDir = history.GetTopTransform().TransformAxis(*globalDirection);
395 }
396 else
397 {
398 localDir = G4ThreeVector(0.,0.,0.);
399 }
400
401 // Enter this daughter
402 //
403 replicaNo = pParam->GetReplicaNo( localPoint, localDir );
404
405 if( replicaNo < 0 || replicaNo >= G4int(pParam->GetNoVoxels()) )
406 {
407 return false;
408 }
409
410 // Set the correct copy number in physical
411 //
412 pPhysical->SetCopyNo(replicaNo);
413 pParam->ComputeTransformation(replicaNo,pPhysical);
414
415 history.NewLevel(pPhysical, kParameterised, replicaNo );
416 localPoint = history.GetTopTransform().TransformPoint(globalPoint);
417
418 // Set the correct solid and material in Logical Volume
419 //
420 G4LogicalVolume *pLogical = pPhysical->GetLogicalVolume();
421
422 pLogical->UpdateMaterial(pParam->ComputeMaterial(replicaNo,
423 pPhysical, &parentTouchable) );
424 return true;
425}
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
void UpdateMaterial(G4Material *pMaterial)
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
virtual G4int GetReplicaNo(const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr) override
void ComputeTransformation(const G4int, G4VPhysicalVolume *) const override
std::size_t GetNoVoxels() const
virtual void SetCopyNo(G4int CopyNo)=0
@ kParameterised
Definition geomdefs.hh:86

Referenced by ComputeStep(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), and G4Navigator::LocateGlobalPointAndSetup().

◆ SetNormalNavigation()

void G4RegularNavigation::SetNormalNavigation ( G4NormalNavigation * fnormnav)
inline

Definition at line 112 of file G4RegularNavigation.hh.

113 { fnormalNav = fnormnav; }

Referenced by G4ITNavigator1::G4ITNavigator1(), G4ITNavigator2::G4ITNavigator2(), and G4Navigator::G4Navigator().


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