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