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

#include <G4ParameterisedNavigation.hh>

+ Inheritance diagram for G4ParameterisedNavigation:

Public Member Functions

 G4ParameterisedNavigation ()
 
 ~G4ParameterisedNavigation () override
 
G4SmartVoxelNodeParamVoxelLocate (G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
 
G4bool LevelLocate (G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint) override
 
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) override
 
G4double ComputeSafety (const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX) override
 
void RelocateWithinVolume (G4VPhysicalVolume *motherPhysical, const G4ThreeVector &localPoint) override
 
- Public Member Functions inherited from G4VoxelNavigation
 G4VoxelNavigation ()
 
virtual ~G4VoxelNavigation ()
 
G4SmartVoxelNodeVoxelLocate (G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
 
virtual G4int GetVerboseLevel () const override
 
virtual void SetVerboseLevel (G4int level) override
 
void EnableBestSafety (G4bool flag=false)
 
- Public Member Functions inherited from G4VNavigation
virtual ~G4VNavigation ()
 
void CheckMode (G4bool mode)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VoxelNavigation
G4double ComputeVoxelSafety (const G4ThreeVector &localPoint) const
 
G4bool LocateNextVoxel (const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentStep)
 
G4SmartVoxelNodeVoxelLocateLight (G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint) const
 
- Protected Attributes inherited from G4VoxelNavigation
G4BlockingList fBList
 
G4int fVoxelDepth = -1
 
std::vector< EAxisfVoxelAxisStack
 
std::vector< G4intfVoxelNoSlicesStack
 
std::vector< G4doublefVoxelSliceWidthStack
 
std::vector< G4intfVoxelNodeNoStack
 
std::vector< G4SmartVoxelHeader * > fVoxelHeaderStack
 
G4SmartVoxelNodefVoxelNode = nullptr
 
G4VoxelSafetyfpVoxelSafety = nullptr
 
G4double fHalfTolerance
 
G4bool fBestSafety = false
 
G4NavigationLoggerfLogger
 
- Protected Attributes inherited from G4VNavigation
G4int fVerbose = 0
 
G4bool fCheck = false
 

Detailed Description

Definition at line 54 of file G4ParameterisedNavigation.hh.

Constructor & Destructor Documentation

◆ G4ParameterisedNavigation()

G4ParameterisedNavigation::G4ParameterisedNavigation ( )
default

◆ ~G4ParameterisedNavigation()

G4ParameterisedNavigation::~G4ParameterisedNavigation ( )
overridedefault

Member Function Documentation

◆ ComputeSafety()

G4double G4ParameterisedNavigation::ComputeSafety ( const G4ThreeVector & globalpoint,
const G4NavigationHistory & history,
const G4double pMaxLength = DBL_MAX )
overridevirtual

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.

Reimplemented from G4VoxelNavigation.

Definition at line 404 of file G4ParameterisedNavigation.cc.

407{
408 G4VPhysicalVolume *motherPhysical, *samplePhysical;
409 G4VPVParameterisation *sampleParam;
410 G4LogicalVolume *motherLogical;
411 G4VSolid *motherSolid, *sampleSolid;
412 G4double motherSafety, ourSafety;
413 G4int sampleNo, curVoxelNodeNo;
414
415 G4SmartVoxelNode *curVoxelNode;
416 G4long curNoVolumes, contentNo;
417 G4double voxelSafety;
418
419 // Replication data
420 //
421 EAxis axis;
422 G4int nReplicas;
423 G4double width, offset;
424 G4bool consuming;
425
426 motherPhysical = history.GetTopVolume();
427 motherLogical = motherPhysical->GetLogicalVolume();
428 motherSolid = motherLogical->GetSolid();
429
430 //
431 // Compute mother safety
432 //
433
434 motherSafety = motherSolid->DistanceToOut(localPoint);
435 ourSafety = motherSafety; // Working isotropic safety
436
437 //
438 // Compute daughter safeties
439 //
440
441 // By definition, parameterised volumes exist as first
442 // daughter of the mother volume
443 //
444 samplePhysical = motherLogical->GetDaughter(0);
445 samplePhysical->GetReplicationData(axis, nReplicas,
446 width, offset, consuming);
447 sampleParam = samplePhysical->GetParameterisation();
448
449 // Look inside the current Voxel only at the current point
450 //
451 if ( axis==kUndefined ) // 3D case: current voxel node is retrieved
452 { // from G4VoxelNavigation.
453 curVoxelNode = fVoxelNode;
454 }
455 else // 1D case: current voxel node is computed here.
456 {
457 curVoxelNodeNo = G4int((localPoint(fVoxelAxis)
458 -fVoxelHeader->GetMinExtent()) / fVoxelSliceWidth );
459 curVoxelNode = fVoxelHeader->GetSlice(curVoxelNodeNo)->GetNode();
460 fVoxelNodeNo = curVoxelNodeNo;
461 fVoxelNode = curVoxelNode;
462 }
463 curNoVolumes = curVoxelNode->GetNoContained();
464
465 for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
466 {
467 sampleNo = curVoxelNode->GetVolume((G4int)contentNo);
468
469 // Call virtual methods, and copy information if needed
470 //
471 sampleSolid= IdentifyAndPlaceSolid( sampleNo,samplePhysical,sampleParam );
472
473 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
474 samplePhysical->GetTranslation());
475 sampleTf.Invert();
476 const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
477 G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
478 if ( sampleSafety<ourSafety )
479 {
480 ourSafety = sampleSafety;
481 }
482 }
483
484 voxelSafety = ComputeVoxelSafety(localPoint,axis);
485 if ( voxelSafety<ourSafety )
486 {
487 ourSafety=voxelSafety;
488 }
489
490 return ourSafety;
491}
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
G4VPhysicalVolume * GetTopVolume() const
G4double GetMinExtent() const
G4SmartVoxelProxy * GetSlice(std::size_t n) const
G4int GetVolume(G4int pVolumeNo) const
std::size_t GetNoContained() const
G4SmartVoxelNode * GetNode() const
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const
G4LogicalVolume * GetLogicalVolume() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
virtual G4VPVParameterisation * GetParameterisation() const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4SmartVoxelNode * fVoxelNode
EAxis
Definition geomdefs.hh:54
@ kUndefined
Definition geomdefs.hh:61

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

◆ ComputeStep()

G4double G4ParameterisedNavigation::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 )
overridevirtual

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.

Reimplemented from G4VoxelNavigation.

Definition at line 70 of file G4ParameterisedNavigation.cc.

82{
83 G4VPhysicalVolume *motherPhysical, *samplePhysical;
84 G4VPVParameterisation *sampleParam;
85 G4LogicalVolume *motherLogical;
86 G4VSolid *motherSolid, *sampleSolid;
87 G4ThreeVector sampleDirection;
88 G4double ourStep=currentProposedStepLength, ourSafety;
89 G4double motherSafety, motherStep = DBL_MAX;
90 G4bool motherValidExitNormal = false;
91 G4ThreeVector motherExitNormal;
92
93 G4int sampleNo;
94
95 G4bool initialNode, noStep;
96 G4SmartVoxelNode *curVoxelNode;
97 G4long curNoVolumes, contentNo;
98 G4double voxelSafety;
99
100 // Replication data
101 //
102 EAxis axis;
103 G4int nReplicas;
104 G4double width, offset;
105 G4bool consuming;
106
107 motherPhysical = history.GetTopVolume();
108 motherLogical = motherPhysical->GetLogicalVolume();
109 motherSolid = motherLogical->GetSolid();
110
111 //
112 // Compute mother safety
113 //
114
115 motherSafety = motherSolid->DistanceToOut(localPoint);
116 ourSafety = motherSafety; // Working isotropic safety
117
118#ifdef G4VERBOSE
119 if ( fCheck )
120 {
121 if( motherSafety < 0.0 )
122 {
123 motherSolid->DumpInfo();
124 std::ostringstream message;
125 message << "Negative Safety In Voxel Navigation !" << G4endl
126 << " Current solid " << motherSolid->GetName()
127 << " gave negative safety: " << motherSafety << G4endl
128 << " for the current (local) point " << localPoint;
129 G4Exception("G4ParameterisedNavigation::ComputeStep()",
130 "GeomNav0003", FatalException, message);
131 }
132 if( motherSolid->Inside(localPoint) == kOutside )
133 {
134 std::ostringstream message;
135 message << "Point is outside Current Volume !" << G4endl
136 << " Point " << localPoint
137 << " is outside current volume " << motherPhysical->GetName()
138 << G4endl;
139 G4double estDistToSolid = motherSolid->DistanceToIn(localPoint);
140 G4cout << " Estimated isotropic distance to solid (distToIn)= "
141 << estDistToSolid;
142 if( estDistToSolid > 100.0 * motherSolid->GetTolerance() )
143 {
144 motherSolid->DumpInfo();
145 G4Exception("G4ParameterisedNavigation::ComputeStep()",
146 "GeomNav0003", FatalException, message,
147 "Point is far outside Current Volume !");
148 }
149 else
150 {
151 G4Exception("G4ParameterisedNavigation::ComputeStep()",
152 "GeomNav1002", JustWarning, message,
153 "Point is a little outside Current Volume.");
154 }
155 }
156
157 // Compute early:
158 // a) to check whether point is (wrongly) outside
159 // (signaled if step < 0 or step == kInfinity )
160 // b) to check value against answer of daughters!
161 //
162 motherStep = motherSolid->DistanceToOut(localPoint,
163 localDirection,
164 true,
165 &motherValidExitNormal,
166 &motherExitNormal);
167
168 if( (motherStep >= kInfinity) || (motherStep < 0.0) )
169 {
170 // Error - indication of being outside solid !!
171 //
172 fLogger->ReportOutsideMother(localPoint, localDirection, motherPhysical);
173
174 ourStep = motherStep = 0.0;
175 exiting = true;
176 entering = false;
177
178 // If we are outside the solid does the normal make sense?
179 validExitNormal = motherValidExitNormal;
180 exitNormal = motherExitNormal;
181
182 *pBlockedPhysical = nullptr; // or motherPhysical ?
183 blockedReplicaNo = 0; // or motherReplicaNumber ?
184
185 newSafety = 0.0;
186 return ourStep;
187 }
188 }
189#endif
190
191 initialNode = true;
192 noStep = true;
193
194 // By definition, the parameterised volume is the first
195 // (and only) daughter of the mother volume
196 //
197 samplePhysical = motherLogical->GetDaughter(0);
198 samplePhysical->GetReplicationData(axis,nReplicas,width,offset,consuming);
199 fBList.Enlarge(nReplicas);
200 fBList.Reset();
201
202 // Exiting normal optimisation
203 //
204 if (exiting && (*pBlockedPhysical==samplePhysical) && validExitNormal)
205 {
206 if (localDirection.dot(exitNormal)>=kMinExitingNormalCosine)
207 {
208 // Block exited daughter replica; Must be on boundary => zero safety
209 //
210 fBList.BlockVolume(blockedReplicaNo);
211 ourSafety = 0;
212 }
213 }
214 exiting = false;
215 entering = false;
216
217 sampleParam = samplePhysical->GetParameterisation();
218
219 // Loop over voxels & compute daughter safeties & intersections
220
221 do
222 {
223 curVoxelNode = fVoxelNode;
224 curNoVolumes = curVoxelNode->GetNoContained();
225
226 for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
227 {
228 sampleNo = curVoxelNode->GetVolume((G4int)contentNo);
229 if ( !fBList.IsBlocked(sampleNo) )
230 {
231 fBList.BlockVolume(sampleNo);
232
233 // Call virtual methods, and copy information if needed
234 //
235 sampleSolid = IdentifyAndPlaceSolid( sampleNo, samplePhysical,
236 sampleParam );
237
238 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
239 samplePhysical->GetTranslation());
240 sampleTf.Invert();
241 const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
242 const G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
243 if ( sampleSafety<ourSafety )
244 {
245 ourSafety = sampleSafety;
246 }
247 if ( sampleSafety<=ourStep )
248 {
249 sampleDirection = sampleTf.TransformAxis(localDirection);
250 G4double sampleStep =
251 sampleSolid->DistanceToIn(samplePoint, sampleDirection);
252 if ( sampleStep<=ourStep )
253 {
254 ourStep = sampleStep;
255 entering = true;
256 exiting = false;
257 *pBlockedPhysical = samplePhysical;
258 blockedReplicaNo = sampleNo;
259#ifdef G4VERBOSE
260 // Check to see that the resulting point is indeed in/on volume.
261 // This check could eventually be made only for successful
262 // candidate.
263
264 if ( ( fCheck ) && ( sampleStep < kInfinity ) )
265 {
266 G4ThreeVector intersectionPoint;
267 intersectionPoint = samplePoint + sampleStep * sampleDirection;
268 EInside insideIntPt = sampleSolid->Inside(intersectionPoint);
269 if( insideIntPt != kSurface )
270 {
271 G4long oldcoutPrec = G4cout.precision(16);
272 std::ostringstream message;
273 message << "Navigator gets conflicting response from Solid."
274 << G4endl
275 << " Inaccurate solid DistanceToIn"
276 << " for solid " << sampleSolid->GetName() << G4endl
277 << " Solid gave DistanceToIn = "
278 << sampleStep << " yet returns " ;
279 if( insideIntPt == kInside )
280 {
281 message << "-kInside-";
282 }
283 else if( insideIntPt == kOutside )
284 {
285 message << "-kOutside-";
286 }
287 else
288 {
289 message << "-kSurface-";
290 }
291 message << " for this point !" << G4endl
292 << " Point = " << intersectionPoint
293 << G4endl;
294 if ( insideIntPt != kInside )
295 {
296 message << " DistanceToIn(p) = "
297 << sampleSolid->DistanceToIn(intersectionPoint);
298 }
299 if ( insideIntPt != kOutside )
300 {
301 message << " DistanceToOut(p) = "
302 << sampleSolid->DistanceToOut(intersectionPoint);
303 }
304 G4Exception("G4ParameterisedNavigation::ComputeStep()",
305 "GeomNav1002", JustWarning, message);
306 G4cout.precision(oldcoutPrec);
307 }
308 }
309#endif
310 }
311 }
312 }
313 }
314
315 if ( initialNode )
316 {
317 initialNode = false;
318 voxelSafety = ComputeVoxelSafety(localPoint,axis);
319 if ( voxelSafety<ourSafety )
320 {
321 ourSafety = voxelSafety;
322 }
323 if ( currentProposedStepLength<ourSafety )
324 {
325 // Guaranteed physics limited
326 //
327 noStep = false;
328 entering = false;
329 exiting = false;
330 *pBlockedPhysical = nullptr;
331 ourStep = kInfinity;
332 }
333 else
334 {
335 // Consider intersection with mother solid
336 //
337 if ( motherSafety<=ourStep )
338 {
339 if ( !fCheck )
340 {
341 motherStep = motherSolid->DistanceToOut(localPoint,
342 localDirection,
343 true,
344 &motherValidExitNormal,
345 &motherExitNormal);
346 }
347
348 if( ( motherStep < 0.0 ) || ( motherStep >= kInfinity) )
349 {
350#ifdef G4VERBOSE
351 fLogger->ReportOutsideMother(localPoint, localDirection,
352 motherPhysical);
353#endif
354 ourStep = motherStep = 0.0;
355 // Rely on the code below to set the remaining state, i.e.
356 // exiting, entering, exitNormal & validExitNormal,
357 // pBlockedPhysical etc.
358 }
359#ifdef G4VERBOSE
360 if( motherValidExitNormal && ( fCheck || (motherStep<=ourStep)) )
361 {
362 fLogger->CheckAndReportBadNormal(motherExitNormal,
363 localPoint, localDirection,
364 motherStep, motherSolid,
365 "From motherSolid::DistanceToOut");
366 }
367#endif
368 if ( motherStep<=ourStep )
369 {
370 ourStep = motherStep;
371 exiting = true;
372 entering = false;
373 if ( validExitNormal )
374 {
375 const G4RotationMatrix* rot = motherPhysical->GetRotation();
376 if (rot != nullptr)
377 {
378 exitNormal *= rot->inverse();
379 }
380 }
381 }
382 else
383 {
384 validExitNormal = false;
385 }
386 }
387 }
388 newSafety = ourSafety;
389 }
390 if (noStep)
391 {
392 noStep = LocateNextVoxel(localPoint, localDirection, ourStep, axis);
393 }
394 } while (noStep);
395
396 return ourStep;
397}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double dot(const Hep3Vector &) const
HepRotation inverse() const
void BlockVolume(const G4int v)
void Enlarge(const G4int nv)
G4bool IsBlocked(const G4int v) const
G4bool CheckAndReportBadNormal(const G4ThreeVector &unitNormal, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double step, const G4VSolid *solid, const char *msg) const
void ReportOutsideMother(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4VPhysicalVolume *motherPV, G4double tDist=30.0 *CLHEP::cm) const
const G4String & GetName() const
G4String GetName() const
G4double GetTolerance() const
virtual EInside Inside(const G4ThreeVector &p) const =0
void DumpInfo() const
G4NavigationLogger * fLogger
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
#define DBL_MAX
Definition templates.hh:62

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

◆ LevelLocate()

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

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.

Reimplemented from G4VoxelNavigation.

Definition at line 610 of file G4ParameterisedNavigation.cc.

617{
618 G4SmartVoxelHeader *motherVoxelHeader;
619 G4SmartVoxelNode *motherVoxelNode;
620 G4VPhysicalVolume *motherPhysical, *pPhysical;
621 G4VPVParameterisation *pParam;
622 G4LogicalVolume *motherLogical;
623 G4VSolid *pSolid;
624 G4ThreeVector samplePoint;
625 G4int voxelNoDaughters, replicaNo;
626
627 motherPhysical = history.GetTopVolume();
628 motherLogical = motherPhysical->GetLogicalVolume();
629 motherVoxelHeader = motherLogical->GetVoxelHeader();
630
631 // Find the voxel containing the point
632 //
633 motherVoxelNode = ParamVoxelLocate(motherVoxelHeader,localPoint);
634
635 voxelNoDaughters = (G4int)motherVoxelNode->GetNoContained();
636 if ( voxelNoDaughters==0 ) { return false; }
637
638 pPhysical = motherLogical->GetDaughter(0);
639 pParam = pPhysical->GetParameterisation();
640
641 // Save parent history in touchable history
642 // ... for use as parent t-h in ComputeMaterial method of param
643 //
644 G4TouchableHistory parentTouchable( history );
645
646 // Search replicated daughter volume
647 //
648 for ( auto sampleNo=voxelNoDaughters-1; sampleNo>=0; sampleNo-- )
649 {
650 replicaNo = motherVoxelNode->GetVolume(sampleNo);
651 if ( (replicaNo!=blockedNum) || (pPhysical!=blockedVol) )
652 {
653 // Obtain solid (as it can vary) and obtain its parameters
654 //
655 pSolid = IdentifyAndPlaceSolid( replicaNo, pPhysical, pParam );
656
657 // Setup history
658 //
659 history.NewLevel(pPhysical, kParameterised, replicaNo);
660 samplePoint = history.GetTopTransform().TransformPoint(globalPoint);
661 if ( !G4AuxiliaryNavServices::CheckPointOnSurface( pSolid, samplePoint,
662 globalDirection, history.GetTopTransform(), pLocatedOnEdge) )
663 {
664 history.BackLevel();
665 }
666 else
667 {
668 // Enter this daughter
669 //
670 localPoint = samplePoint;
671
672 // Set the correct copy number in physical
673 //
674 pPhysical->SetCopyNo(replicaNo);
675
676 // Set the correct solid and material in Logical Volume
677 //
678 G4LogicalVolume *pLogical = pPhysical->GetLogicalVolume();
679 pLogical->SetSolid(pSolid);
680 pLogical->UpdateMaterial(pParam->ComputeMaterial(replicaNo,
681 pPhysical, &parentTouchable) );
682 return true;
683 }
684 }
685 }
686 return false;
687}
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
static G4bool CheckPointOnSurface(const G4VSolid *sampleSolid, const G4ThreeVector &localPoint, const G4ThreeVector *globalDirection, const G4AffineTransform &sampleTransform, const G4bool locatedOnEdge)
void SetSolid(G4VSolid *pSolid)
G4SmartVoxelHeader * GetVoxelHeader() const
void UpdateMaterial(G4Material *pMaterial)
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
const G4AffineTransform & GetTopTransform() const
G4SmartVoxelNode * ParamVoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
virtual void SetCopyNo(G4int CopyNo)=0
@ kParameterised
Definition geomdefs.hh:86

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

◆ ParamVoxelLocate()

◆ RelocateWithinVolume()

void G4ParameterisedNavigation::RelocateWithinVolume ( G4VPhysicalVolume * ,
const G4ThreeVector &  )
overridevirtual

Update internal navigation state to take into account that location has been moved, but remains within the motherPhysical volume.

Parameters
[in]motherPhysicalCurrent physical volume.
[in]localPointLocal point.

Reimplemented from G4VoxelNavigation.

Definition at line 689 of file G4ParameterisedNavigation.cc.

691{
692 auto motherLogical = motherPhysical->GetLogicalVolume();
693
694 /* this should only be called on parameterized volumes, which always satisfy the conditions below */
695 assert(motherPhysical->GetRegularStructureId() != 1);
696 assert(motherLogical->GetNoDaughters() == 1);
697
698 if ( auto pVoxelHeader = motherLogical->GetVoxelHeader() )
699 ParamVoxelLocate( pVoxelHeader, localPoint );
700}

Referenced by G4Navigator::LocateGlobalPointWithinVolume().


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