Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParameterisedNavigation.cc
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 G4ParameterisedNavigation Implementation
27//
28// Original author: Paul Kent (CERN), August 1996
29//
30// Revisions:
31// - J. Apostolakis 5 Mar 1998, Enabled parameterisation of mat & solid type
32// - G. Cosmo 15 May 2002, Extended to 3-d voxelisation, made subclass
33// - G. Cosmo 11 Mar 2004, Added Check mode
34// - J. Apostolakis 24 Nov 2005, Revised/fixed treatment of nested params
35// --------------------------------------------------------------------
36
37// Note: We cannot make the solid, dimensions and transformation dependent on
38// parent because the voxelisation will not have access to this.
39// So the following can NOT be done:
40// sampleSolid = curParam->ComputeSolid(num, curPhysical, pParentTouch);
41// sampleSolid->ComputeDimensions(curParam, num, curPhysical, pParentTouch);
42// curParam->ComputeTransformation(num, curPhysical, pParentTouch);
43
45#include "G4TouchableHistory.hh"
47
49
50// #include <cassert>
51
52// ********************************************************************
53// Constructor
54// ********************************************************************
55//
57
58// ***************************************************************************
59// Destructor
60// ***************************************************************************
61//
63
64// ***************************************************************************
65// ComputeStep
66// ***************************************************************************
67//
69 ComputeStep(const G4ThreeVector& localPoint,
70 const G4ThreeVector& localDirection,
71 const G4double currentProposedStepLength,
72 G4double& newSafety,
73 G4NavigationHistory& history,
74 G4bool& validExitNormal,
75 G4ThreeVector& exitNormal,
76 G4bool& exiting,
77 G4bool& entering,
78 G4VPhysicalVolume *(*pBlockedPhysical),
79 G4int& blockedReplicaNo)
80{
81 G4VPhysicalVolume *motherPhysical, *samplePhysical;
82 G4VPVParameterisation *sampleParam;
83 G4LogicalVolume *motherLogical;
84 G4VSolid *motherSolid, *sampleSolid;
85 G4ThreeVector sampleDirection;
86 G4double ourStep=currentProposedStepLength, ourSafety;
87 G4double motherSafety, motherStep = DBL_MAX;
88 G4bool motherValidExitNormal = false;
89 G4ThreeVector motherExitNormal;
90
91 G4int sampleNo;
92
93 G4bool initialNode, noStep;
94 G4SmartVoxelNode *curVoxelNode;
95 G4long curNoVolumes, contentNo;
96 G4double voxelSafety;
97
98 // Replication data
99 //
100 EAxis axis;
101 G4int nReplicas;
102 G4double width, offset;
103 G4bool consuming;
104
105 motherPhysical = history.GetTopVolume();
106 motherLogical = motherPhysical->GetLogicalVolume();
107 motherSolid = motherLogical->GetSolid();
108
109 //
110 // Compute mother safety
111 //
112
113 motherSafety = motherSolid->DistanceToOut(localPoint);
114 ourSafety = motherSafety; // Working isotropic safety
115
116#ifdef G4VERBOSE
117 if ( fCheck )
118 {
119 if( motherSafety < 0.0 )
120 {
121 motherSolid->DumpInfo();
122 std::ostringstream message;
123 message << "Negative Safety In Voxel Navigation !" << G4endl
124 << " Current solid " << motherSolid->GetName()
125 << " gave negative safety: " << motherSafety << G4endl
126 << " for the current (local) point " << localPoint;
127 G4Exception("G4ParameterisedNavigation::ComputeStep()",
128 "GeomNav0003", FatalException, message);
129 }
130 if( motherSolid->Inside(localPoint) == kOutside )
131 {
132 std::ostringstream message;
133 message << "Point is outside Current Volume !" << G4endl
134 << " Point " << localPoint
135 << " is outside current volume " << motherPhysical->GetName()
136 << G4endl;
137 G4double estDistToSolid = motherSolid->DistanceToIn(localPoint);
138 G4cout << " Estimated isotropic distance to solid (distToIn)= "
139 << estDistToSolid;
140 if( estDistToSolid > 100.0 * motherSolid->GetTolerance() )
141 {
142 motherSolid->DumpInfo();
143 G4Exception("G4ParameterisedNavigation::ComputeStep()",
144 "GeomNav0003", FatalException, message,
145 "Point is far outside Current Volume !");
146 }
147 else
148 {
149 G4Exception("G4ParameterisedNavigation::ComputeStep()",
150 "GeomNav1002", JustWarning, message,
151 "Point is a little outside Current Volume.");
152 }
153 }
154
155 // Compute early:
156 // a) to check whether point is (wrongly) outside
157 // (signaled if step < 0 or step == kInfinity )
158 // b) to check value against answer of daughters!
159 //
160 motherStep = motherSolid->DistanceToOut(localPoint,
161 localDirection,
162 true,
163 &motherValidExitNormal,
164 &motherExitNormal);
165
166 if( (motherStep >= kInfinity) || (motherStep < 0.0) )
167 {
168 // Error - indication of being outside solid !!
169 //
170 fLogger->ReportOutsideMother(localPoint, localDirection, motherPhysical);
171
172 ourStep = motherStep = 0.0;
173 exiting = true;
174 entering = false;
175
176 // If we are outside the solid does the normal make sense?
177 validExitNormal = motherValidExitNormal;
178 exitNormal = motherExitNormal;
179
180 *pBlockedPhysical = nullptr; // or motherPhysical ?
181 blockedReplicaNo = 0; // or motherReplicaNumber ?
182
183 newSafety = 0.0;
184 return ourStep;
185 }
186 }
187#endif
188
189 initialNode = true;
190 noStep = true;
191
192 // By definition, the parameterised volume is the first
193 // (and only) daughter of the mother volume
194 //
195 samplePhysical = motherLogical->GetDaughter(0);
196 samplePhysical->GetReplicationData(axis,nReplicas,width,offset,consuming);
197 fBList.Enlarge(nReplicas);
198 fBList.Reset();
199
200 // Exiting normal optimisation
201 //
202 if (exiting && (*pBlockedPhysical==samplePhysical) && validExitNormal)
203 {
204 if (localDirection.dot(exitNormal)>=kMinExitingNormalCosine)
205 {
206 // Block exited daughter replica; Must be on boundary => zero safety
207 //
208 fBList.BlockVolume(blockedReplicaNo);
209 ourSafety = 0;
210 }
211 }
212 exiting = false;
213 entering = false;
214
215 sampleParam = samplePhysical->GetParameterisation();
216
217 // Loop over voxels & compute daughter safeties & intersections
218
219 do
220 {
221 curVoxelNode = fVoxelNode;
222 curNoVolumes = curVoxelNode->GetNoContained();
223
224 for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
225 {
226 sampleNo = curVoxelNode->GetVolume((G4int)contentNo);
227 if ( !fBList.IsBlocked(sampleNo) )
228 {
229 fBList.BlockVolume(sampleNo);
230
231 // Call virtual methods, and copy information if needed
232 //
233 sampleSolid = IdentifyAndPlaceSolid( sampleNo, samplePhysical,
234 sampleParam );
235
236 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
237 samplePhysical->GetTranslation());
238 sampleTf.Invert();
239 const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
240 const G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
241 if ( sampleSafety<ourSafety )
242 {
243 ourSafety = sampleSafety;
244 }
245 if ( sampleSafety<=ourStep )
246 {
247 sampleDirection = sampleTf.TransformAxis(localDirection);
248 G4double sampleStep =
249 sampleSolid->DistanceToIn(samplePoint, sampleDirection);
250 if ( sampleStep<=ourStep )
251 {
252 ourStep = sampleStep;
253 entering = true;
254 exiting = false;
255 *pBlockedPhysical = samplePhysical;
256 blockedReplicaNo = sampleNo;
257#ifdef G4VERBOSE
258 // Check to see that the resulting point is indeed in/on volume.
259 // This check could eventually be made only for successful
260 // candidate.
261
262 if ( ( fCheck ) && ( sampleStep < kInfinity ) )
263 {
264 G4ThreeVector intersectionPoint;
265 intersectionPoint = samplePoint + sampleStep * sampleDirection;
266 EInside insideIntPt = sampleSolid->Inside(intersectionPoint);
267 if( insideIntPt != kSurface )
268 {
269 G4long oldcoutPrec = G4cout.precision(16);
270 std::ostringstream message;
271 message << "Navigator gets conflicting response from Solid."
272 << G4endl
273 << " Inaccurate solid DistanceToIn"
274 << " for solid " << sampleSolid->GetName() << G4endl
275 << " Solid gave DistanceToIn = "
276 << sampleStep << " yet returns " ;
277 if( insideIntPt == kInside )
278 {
279 message << "-kInside-";
280 }
281 else if( insideIntPt == kOutside )
282 {
283 message << "-kOutside-";
284 }
285 else
286 {
287 message << "-kSurface-";
288 }
289 message << " for this point !" << G4endl
290 << " Point = " << intersectionPoint
291 << G4endl;
292 if ( insideIntPt != kInside )
293 {
294 message << " DistanceToIn(p) = "
295 << sampleSolid->DistanceToIn(intersectionPoint);
296 }
297 if ( insideIntPt != kOutside )
298 {
299 message << " DistanceToOut(p) = "
300 << sampleSolid->DistanceToOut(intersectionPoint);
301 }
302 G4Exception("G4ParameterisedNavigation::ComputeStep()",
303 "GeomNav1002", JustWarning, message);
304 G4cout.precision(oldcoutPrec);
305 }
306 }
307#endif
308 }
309 }
310 }
311 }
312
313 if ( initialNode )
314 {
315 initialNode = false;
316 voxelSafety = ComputeVoxelSafety(localPoint,axis);
317 if ( voxelSafety<ourSafety )
318 {
319 ourSafety = voxelSafety;
320 }
321 if ( currentProposedStepLength<ourSafety )
322 {
323 // Guaranteed physics limited
324 //
325 noStep = false;
326 entering = false;
327 exiting = false;
328 *pBlockedPhysical = nullptr;
329 ourStep = kInfinity;
330 }
331 else
332 {
333 // Consider intersection with mother solid
334 //
335 if ( motherSafety<=ourStep )
336 {
337 if ( !fCheck )
338 {
339 motherStep = motherSolid->DistanceToOut(localPoint,
340 localDirection,
341 true,
342 &motherValidExitNormal,
343 &motherExitNormal);
344 }
345
346 if( ( motherStep < 0.0 ) || ( motherStep >= kInfinity) )
347 {
348#ifdef G4VERBOSE
349 fLogger->ReportOutsideMother(localPoint, localDirection,
350 motherPhysical);
351#endif
352 ourStep = motherStep = 0.0;
353 // Rely on the code below to set the remaining state, i.e.
354 // exiting, entering, exitNormal & validExitNormal,
355 // pBlockedPhysical etc.
356 }
357#ifdef G4VERBOSE
358 if( motherValidExitNormal && ( fCheck || (motherStep<=ourStep)) )
359 {
360 fLogger->CheckAndReportBadNormal(motherExitNormal,
361 localPoint, localDirection,
362 motherStep, motherSolid,
363 "From motherSolid::DistanceToOut");
364 }
365#endif
366 if ( motherStep<=ourStep )
367 {
368 ourStep = motherStep;
369 exiting = true;
370 entering = false;
371 if ( validExitNormal )
372 {
373 const G4RotationMatrix* rot = motherPhysical->GetRotation();
374 if (rot != nullptr)
375 {
376 exitNormal *= rot->inverse();
377 }
378 }
379 }
380 else
381 {
382 validExitNormal = false;
383 }
384 }
385 }
386 newSafety = ourSafety;
387 }
388 if (noStep)
389 {
390 noStep = LocateNextVoxel(localPoint, localDirection, ourStep, axis);
391 }
392 } while (noStep);
393
394 return ourStep;
395}
396
397// ***************************************************************************
398// ComputeSafety
399// ***************************************************************************
400//
403 const G4NavigationHistory& history,
404 const G4double )
405{
406 G4VPhysicalVolume *motherPhysical, *samplePhysical;
407 G4VPVParameterisation *sampleParam;
408 G4LogicalVolume *motherLogical;
409 G4VSolid *motherSolid, *sampleSolid;
410 G4double motherSafety, ourSafety;
411 G4int sampleNo, curVoxelNodeNo;
412
413 G4SmartVoxelNode *curVoxelNode;
414 G4long curNoVolumes, contentNo;
415 G4double voxelSafety;
416
417 // Replication data
418 //
419 EAxis axis;
420 G4int nReplicas;
421 G4double width, offset;
422 G4bool consuming;
423
424 motherPhysical = history.GetTopVolume();
425 motherLogical = motherPhysical->GetLogicalVolume();
426 motherSolid = motherLogical->GetSolid();
427
428 //
429 // Compute mother safety
430 //
431
432 motherSafety = motherSolid->DistanceToOut(localPoint);
433 ourSafety = motherSafety; // Working isotropic safety
434
435 //
436 // Compute daughter safeties
437 //
438
439 // By definition, parameterised volumes exist as first
440 // daughter of the mother volume
441 //
442 samplePhysical = motherLogical->GetDaughter(0);
443 samplePhysical->GetReplicationData(axis, nReplicas,
444 width, offset, consuming);
445 sampleParam = samplePhysical->GetParameterisation();
446
447 // Look inside the current Voxel only at the current point
448 //
449 if ( axis==kUndefined ) // 3D case: current voxel node is retrieved
450 { // from G4VoxelNavigation.
451 curVoxelNode = fVoxelNode;
452 }
453 else // 1D case: current voxel node is computed here.
454 {
455 curVoxelNodeNo = G4int((localPoint(fVoxelAxis)
456 -fVoxelHeader->GetMinExtent()) / fVoxelSliceWidth );
457 curVoxelNode = fVoxelHeader->GetSlice(curVoxelNodeNo)->GetNode();
458 fVoxelNodeNo = curVoxelNodeNo;
459 fVoxelNode = curVoxelNode;
460 }
461 curNoVolumes = curVoxelNode->GetNoContained();
462
463 for ( contentNo=curNoVolumes-1; contentNo>=0; contentNo-- )
464 {
465 sampleNo = curVoxelNode->GetVolume((G4int)contentNo);
466
467 // Call virtual methods, and copy information if needed
468 //
469 sampleSolid= IdentifyAndPlaceSolid( sampleNo,samplePhysical,sampleParam );
470
471 G4AffineTransform sampleTf(samplePhysical->GetRotation(),
472 samplePhysical->GetTranslation());
473 sampleTf.Invert();
474 const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
475 G4double sampleSafety = sampleSolid->DistanceToIn(samplePoint);
476 if ( sampleSafety<ourSafety )
477 {
478 ourSafety = sampleSafety;
479 }
480 }
481
482 voxelSafety = ComputeVoxelSafety(localPoint,axis);
483 if ( voxelSafety<ourSafety )
484 {
485 ourSafety=voxelSafety;
486 }
487
488 return ourSafety;
489}
490
491// ********************************************************************
492// ComputeVoxelSafety
493//
494// Computes safety from specified point to collected voxel boundaries
495// using already located point.
496// ********************************************************************
497//
498G4double G4ParameterisedNavigation::
499ComputeVoxelSafety(const G4ThreeVector& localPoint,
500 const EAxis pAxis) const
501{
502 // If no best axis is specified, adopt default
503 // strategy as for placements
504 //
505 if ( pAxis==kUndefined )
506 {
507 return G4VoxelNavigation::ComputeVoxelSafety(localPoint);
508 }
509
510 G4double voxelSafety, plusVoxelSafety, minusVoxelSafety;
511 G4double curNodeOffset, minCurCommonDelta, maxCurCommonDelta;
512 G4long minCurNodeNoDelta, maxCurNodeNoDelta;
513
514 // Compute linear intersection distance to boundaries of max/min
515 // to collected nodes at current level
516 //
517 curNodeOffset = fVoxelNodeNo*fVoxelSliceWidth;
518 minCurCommonDelta = localPoint(fVoxelAxis)
519 - fVoxelHeader->GetMinExtent()-curNodeOffset;
520 maxCurNodeNoDelta = fVoxelNode->GetMaxEquivalentSliceNo()-fVoxelNodeNo;
521 minCurNodeNoDelta = fVoxelNodeNo-fVoxelNode->GetMinEquivalentSliceNo();
522 maxCurCommonDelta = fVoxelSliceWidth-minCurCommonDelta;
523 plusVoxelSafety = minCurNodeNoDelta*fVoxelSliceWidth+minCurCommonDelta;
524 minusVoxelSafety = maxCurNodeNoDelta*fVoxelSliceWidth+maxCurCommonDelta;
525 voxelSafety = std::min(plusVoxelSafety,minusVoxelSafety);
526
527 if ( voxelSafety<0 )
528 {
529 voxelSafety = 0;
530 }
531
532 return voxelSafety;
533}
534
535// ********************************************************************
536// LocateNextVoxel
537//
538// Finds the next voxel from the current voxel and point
539// in the specified direction.
540//
541// Returns false if all voxels considered
542// true otherwise
543// [current Step ends inside same voxel or leaves all voxels]
544// ********************************************************************
545//
546G4bool G4ParameterisedNavigation::
547LocateNextVoxel( const G4ThreeVector& localPoint,
548 const G4ThreeVector& localDirection,
549 const G4double currentStep,
550 const EAxis pAxis)
551{
552 // If no best axis is specified, adopt default
553 // location strategy as for placements
554 //
555 if ( pAxis==kUndefined )
556 {
557 return G4VoxelNavigation::LocateNextVoxel(localPoint,
558 localDirection,
559 currentStep);
560 }
561
562 G4bool isNewVoxel;
563 G4int newNodeNo;
564 G4double minVal, maxVal, curMinExtent, curCoord;
565
566 curMinExtent = fVoxelHeader->GetMinExtent();
567 curCoord = localPoint(fVoxelAxis)+currentStep*localDirection(fVoxelAxis);
568 minVal = curMinExtent+fVoxelNode->GetMinEquivalentSliceNo()*fVoxelSliceWidth;
569 isNewVoxel = false;
570
571 if ( minVal<=curCoord )
572 {
573 maxVal = curMinExtent
574 + (fVoxelNode->GetMaxEquivalentSliceNo()+1)*fVoxelSliceWidth;
575 if ( maxVal<curCoord )
576 {
577 newNodeNo = fVoxelNode->GetMaxEquivalentSliceNo()+1;
578 if ( newNodeNo<G4int(fVoxelHeader->GetNoSlices()) )
579 {
580 fVoxelNodeNo = newNodeNo;
581 fVoxelNode = fVoxelHeader->GetSlice(newNodeNo)->GetNode();
582 isNewVoxel = true;
583 }
584 }
585 }
586 else
587 {
588 newNodeNo = fVoxelNode->GetMinEquivalentSliceNo()-1;
589
590 // Must locate from newNodeNo no and down to setup stack and fVoxelNode
591 // Repeat or earlier code...
592 //
593 if ( newNodeNo>=0 )
594 {
595 fVoxelNodeNo = newNodeNo;
596 fVoxelNode = fVoxelHeader->GetSlice(newNodeNo)->GetNode();
597 isNewVoxel = true;
598 }
599 }
600 return isNewVoxel;
601}
602
603// ********************************************************************
604// LevelLocate
605// ********************************************************************
606//
607G4bool
609 const G4VPhysicalVolume* blockedVol,
610 const G4int blockedNum,
611 const G4ThreeVector& globalPoint,
612 const G4ThreeVector* globalDirection,
613 const G4bool pLocatedOnEdge,
614 G4ThreeVector& localPoint )
615{
616 G4SmartVoxelHeader *motherVoxelHeader;
617 G4SmartVoxelNode *motherVoxelNode;
618 G4VPhysicalVolume *motherPhysical, *pPhysical;
619 G4VPVParameterisation *pParam;
620 G4LogicalVolume *motherLogical;
621 G4VSolid *pSolid;
622 G4ThreeVector samplePoint;
623 G4int voxelNoDaughters, replicaNo;
624
625 motherPhysical = history.GetTopVolume();
626 motherLogical = motherPhysical->GetLogicalVolume();
627 motherVoxelHeader = motherLogical->GetVoxelHeader();
628
629 // Find the voxel containing the point
630 //
631 motherVoxelNode = ParamVoxelLocate(motherVoxelHeader,localPoint);
632
633 voxelNoDaughters = (G4int)motherVoxelNode->GetNoContained();
634 if ( voxelNoDaughters==0 ) { return false; }
635
636 pPhysical = motherLogical->GetDaughter(0);
637 pParam = pPhysical->GetParameterisation();
638
639 // Save parent history in touchable history
640 // ... for use as parent t-h in ComputeMaterial method of param
641 //
642 G4TouchableHistory parentTouchable( history );
643
644 // Search replicated daughter volume
645 //
646 for ( auto sampleNo=voxelNoDaughters-1; sampleNo>=0; sampleNo-- )
647 {
648 replicaNo = motherVoxelNode->GetVolume(sampleNo);
649 if ( (replicaNo!=blockedNum) || (pPhysical!=blockedVol) )
650 {
651 // Obtain solid (as it can vary) and obtain its parameters
652 //
653 pSolid = IdentifyAndPlaceSolid( replicaNo, pPhysical, pParam );
654
655 // Setup history
656 //
657 history.NewLevel(pPhysical, kParameterised, replicaNo);
658 samplePoint = history.GetTopTransform().TransformPoint(globalPoint);
659 if ( !G4AuxiliaryNavServices::CheckPointOnSurface( pSolid, samplePoint,
660 globalDirection, history.GetTopTransform(), pLocatedOnEdge) )
661 {
662 history.BackLevel();
663 }
664 else
665 {
666 // Enter this daughter
667 //
668 localPoint = samplePoint;
669
670 // Set the correct copy number in physical
671 //
672 pPhysical->SetCopyNo(replicaNo);
673
674 // Set the correct solid and material in Logical Volume
675 //
676 G4LogicalVolume *pLogical = pPhysical->GetLogicalVolume();
677 pLogical->SetSolid(pSolid);
678 pLogical->UpdateMaterial(pParam->ComputeMaterial(replicaNo,
679 pPhysical, &parentTouchable) );
680 return true;
681 }
682 }
683 }
684 return false;
685}
686
688 const G4ThreeVector& localPoint )
689{
690 auto motherLogical = motherPhysical->GetLogicalVolume();
691
692 // this should only be called on parameterized volumes,
693 // which always satisfy the conditions below
694 // assert(motherPhysical->GetRegularStructureId() != 1);
695 // assert(motherLogical->GetNoDaughters() == 1);
696
697 if ( auto pVoxelHeader = motherLogical->GetVoxelHeader() )
698 {
699 ParamVoxelLocate( pVoxelHeader, localPoint );
700 }
701}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ThreadLocal T * G4GeomSplitter< T >::offset
CLHEP::HepRotation G4RotationMatrix
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double dot(const Hep3Vector &) const
HepRotation inverse() const
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
static G4bool CheckPointOnSurface(const G4VSolid *sampleSolid, const G4ThreeVector &localPoint, const G4ThreeVector *globalDirection, const G4AffineTransform &sampleTransform, const G4bool locatedOnEdge)
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
void SetSolid(G4VSolid *pSolid)
G4SmartVoxelHeader * GetVoxelHeader() const
void UpdateMaterial(G4Material *pMaterial)
G4NavigationHistory is a class responsible for the maintenance of the history of the path taken throu...
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
const G4AffineTransform & GetTopTransform() const
G4VPhysicalVolume * GetTopVolume() const
G4SmartVoxelNode * ParamVoxelLocate(G4SmartVoxelHeader *pHead, const G4ThreeVector &localPoint)
~G4ParameterisedNavigation() override
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX) override
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) override
void RelocateWithinVolume(G4VPhysicalVolume *motherPhysical, const G4ThreeVector &localPoint) override
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint) override
G4SmartVoxelHeader represents a set of voxels, created by a single axis of virtual division....
G4SmartVoxelNode defines a node in the smart voxel hierarchy, i.e. a 'slice' of space along a given a...
G4int GetVolume(G4int pVolumeNo) const
std::size_t GetNoContained() const
G4TouchableHistory is an object representing a touchable detector element, and its history in the geo...
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
virtual G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr)
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
const G4RotationMatrix * GetRotation() const
virtual void SetCopyNo(G4int CopyNo)=0
const G4ThreeVector GetTranslation() const
G4LogicalVolume * GetLogicalVolume() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
G4String GetName() const
G4double GetTolerance() const
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
void DumpInfo() const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4NavigationLogger * fLogger
G4SmartVoxelNode * fVoxelNode
G4bool LocateNextVoxel(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentStep)
G4double ComputeVoxelSafety(const G4ThreeVector &localPoint) const
EAxis
Definition geomdefs.hh:54
@ kUndefined
Definition geomdefs.hh:61
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
@ kParameterised
Definition geomdefs.hh:86
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
#define DBL_MAX
Definition templates.hh:62