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