Geant4-11
Public Member Functions | Private Attributes
G4RegularNavigation Class Reference

#include <G4RegularNavigation.hh>

Public Member Functions

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

Private Attributes

G4int fAbandonThreshold_NoZeroSteps = 25
 
G4int fActionThreshold_NoZeroSteps = 2
 
G4bool fcheck = false
 
G4bool fLastStepWasZero = false
 
G4double fMinStep
 
G4NormalNavigationfnormalNav = nullptr
 
G4int fNoStepsAllowed = 10000
 
G4int fNumberZeroSteps = 0
 
G4int fverbose = false
 
G4double kCarTolerance
 

Detailed Description

Definition at line 50 of file G4RegularNavigation.hh.

Constructor & Destructor Documentation

◆ G4RegularNavigation()

G4RegularNavigation::G4RegularNavigation ( )

◆ ~G4RegularNavigation()

G4RegularNavigation::~G4RegularNavigation ( )

Definition at line 50 of file G4RegularNavigation.cc.

51{
52}

Member Function Documentation

◆ CheckMode()

void G4RegularNavigation::CheckMode ( G4bool  mode)
inline

Definition at line 114 of file G4RegularNavigation.hh.

114{ fcheck = mode; }

References fcheck.

◆ ComputeSafety()

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

Definition at line 351 of file G4RegularNavigation.cc.

354{
355 // This method is never called because to be called the daughter has to be a
356 // regular structure. This would only happen if the track is in the mother of
357 // voxels volume. But the voxels fill completely their mother, so when a
358 // track enters the mother it automatically enters a voxel. Only precision
359 // problems would make this method to be called
360
361 // Compute step in voxel
362 //
363 return fnormalNav->ComputeSafety(localPoint,
364 history,
365 pMaxLength );
366}
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
G4NormalNavigation * fnormalNav
def history()
Definition: g4zmq.py:84

References G4NormalNavigation::ComputeSafety(), fnormalNav, and g4zmq::history().

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

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

References G4NormalNavigation::ComputeStep(), fnormalNav, G4LogicalVolume::GetDaughter(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetParameterisation(), G4PhantomParameterisation::GetReplicaNo(), G4PhantomParameterisation::GetTranslation(), g4zmq::history(), and LevelLocate().

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

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

References G4RegularNavigationHelper::AddStepLength(), G4RegularNavigationHelper::ClearStepLengths(), G4PhantomParameterisation::ComputeMaterial(), G4NormalNavigation::ComputeStep(), G4VSolid::DistanceToOut(), EventMustBeAborted, fAbandonThreshold_NoZeroSteps, fActionThreshold_NoZeroSteps, fLastStepWasZero, fMinStep, fnormalNav, fNoStepsAllowed, fNumberZeroSteps, G4endl, G4Exception(), G4PhantomParameterisation::GetContainerSolid(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetParameterisation(), G4PhantomParameterisation::GetReplicaNo(), G4LogicalVolume::GetSolid(), G4PhantomParameterisation::GetTranslation(), g4zmq::history(), G4VSolid::Inside(), G4RegularNavigationHelper::Instance(), JustWarning, kCarTolerance, kInside, G4INCL::Math::min(), and G4PhantomParameterisation::SkipEqualMaterials().

◆ 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 371 of file G4RegularNavigation.cc.

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

References G4PhantomParameterisation::ComputeMaterial(), G4PhantomParameterisation::ComputeTransformation(), G4LogicalVolume::GetDaughter(), G4VPhysicalVolume::GetLogicalVolume(), G4PhantomParameterisation::GetNoVoxels(), G4VPhysicalVolume::GetParameterisation(), G4PhantomParameterisation::GetReplicaNo(), g4zmq::history(), kParameterised, G4VPhysicalVolume::SetCopyNo(), and G4LogicalVolume::UpdateMaterial().

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

◆ SetNormalNavigation()

void G4RegularNavigation::SetNormalNavigation ( G4NormalNavigation fnormnav)
inline

◆ SetVerboseLevel()

void G4RegularNavigation::SetVerboseLevel ( G4int  level)
inline

Definition at line 113 of file G4RegularNavigation.hh.

113{ fverbose = level; }

References fverbose.

Field Documentation

◆ fAbandonThreshold_NoZeroSteps

G4int G4RegularNavigation::fAbandonThreshold_NoZeroSteps = 25
private

Definition at line 133 of file G4RegularNavigation.hh.

Referenced by ComputeStepSkippingEqualMaterials().

◆ fActionThreshold_NoZeroSteps

G4int G4RegularNavigation::fActionThreshold_NoZeroSteps = 2
private

Definition at line 131 of file G4RegularNavigation.hh.

Referenced by ComputeStepSkippingEqualMaterials().

◆ fcheck

G4bool G4RegularNavigation::fcheck = false
private

Definition at line 121 of file G4RegularNavigation.hh.

Referenced by CheckMode().

◆ fLastStepWasZero

G4bool G4RegularNavigation::fLastStepWasZero = false
private

Definition at line 127 of file G4RegularNavigation.hh.

Referenced by ComputeStepSkippingEqualMaterials().

◆ fMinStep

G4double G4RegularNavigation::fMinStep
private

◆ fnormalNav

G4NormalNavigation* G4RegularNavigation::fnormalNav = nullptr
private

◆ fNoStepsAllowed

G4int G4RegularNavigation::fNoStepsAllowed = 10000
private

Definition at line 135 of file G4RegularNavigation.hh.

Referenced by ComputeStepSkippingEqualMaterials().

◆ fNumberZeroSteps

G4int G4RegularNavigation::fNumberZeroSteps = 0
private

Definition at line 129 of file G4RegularNavigation.hh.

Referenced by ComputeStepSkippingEqualMaterials().

◆ fverbose

G4int G4RegularNavigation::fverbose = false
private

Definition at line 120 of file G4RegularNavigation.hh.

Referenced by SetVerboseLevel().

◆ kCarTolerance

G4double G4RegularNavigation::kCarTolerance
private

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