Geant4  10.00.p02
G4Box.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: G4Box.cc 81636 2014-06-04 09:06:08Z gcosmo $
28 //
29 //
30 //
31 // Implementation for G4Box class
32 //
33 // 24.06.98 - V.Grichine: insideEdge in DistanceToIn(p,v)
34 // 20.09.98 - V.Grichine: new algorithm of DistanceToIn(p,v)
35 // 07.05.00 - V.Grichine: d= DistanceToIn(p,v), if d<e/2, d=0
36 // 09.06.00 - V.Grichine: safety in DistanceToIn(p) against Inside(p)=kOutside
37 // and information before exception in DistanceToOut(p,v,...)
38 // 15.11.00 - D.Williams, V.Grichine: bug fixed in CalculateExtent - change
39 // algorithm for rotated vertices
40 // --------------------------------------------------------------------
41 
42 #include "G4Box.hh"
43 
44 #if !defined(G4GEOM_USE_UBOX)
45 
46 #include "G4SystemOfUnits.hh"
47 #include "G4VoxelLimits.hh"
48 #include "G4AffineTransform.hh"
49 #include "Randomize.hh"
50 
51 #include "G4VPVParameterisation.hh"
52 
53 #include "G4VGraphicsScene.hh"
54 #include "G4VisExtent.hh"
55 
57 //
58 // Constructor - check & set half widths
59 
60 G4Box::G4Box(const G4String& pName,
61  G4double pX,
62  G4double pY,
63  G4double pZ)
64  : G4CSGSolid(pName), fDx(pX), fDy(pY), fDz(pZ)
65 {
66  delta = 0.5*kCarTolerance;
67  if ( (pX < 2*kCarTolerance)
68  || (pY < 2*kCarTolerance)
69  || (pZ < 2*kCarTolerance) ) // limit to thickness of surfaces
70  {
71  std::ostringstream message;
72  message << "Dimensions too small for Solid: " << GetName() << "!" << G4endl
73  << " hX, hY, hZ = " << pX << ", " << pY << ", " << pZ;
74  G4Exception("G4Box::G4Box()", "GeomSolids0002", FatalException, message);
75  }
76 }
77 
79 //
80 // Fake default constructor - sets only member data and allocates memory
81 // for usage restricted to object persistency.
82 
83 G4Box::G4Box( __void__& a )
84  : G4CSGSolid(a), fDx(0.), fDy(0.), fDz(0.), delta(0.)
85 {
86 }
87 
89 //
90 // Destructor
91 
93 {
94 }
95 
97 //
98 // Copy constructor
99 
100 G4Box::G4Box(const G4Box& rhs)
101  : G4CSGSolid(rhs), fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), delta(rhs.delta)
102 {
104 }
105 
107 //
108 // Assignment operator
109 
111 {
112  // Check assignment to self
113  //
114  if (this == &rhs) { return *this; }
115 
116  // Copy base class data
117  //
119 
120  // Copy data
121  //
122  fDx = rhs.fDx;
123  fDy = rhs.fDy;
124  fDz = rhs.fDz;
125  delta = rhs.delta;
127 
128  return *this;
129 }
130 
132 
134 {
135  if(dx > 2*kCarTolerance) // limit to thickness of surfaces
136  {
137  fDx = dx;
138  }
139  else
140  {
141  std::ostringstream message;
142  message << "Dimension X too small for solid: " << GetName() << "!"
143  << G4endl
144  << " hX = " << dx;
145  G4Exception("G4Box::SetXHalfLength()", "GeomSolids0002",
146  FatalException, message);
147  }
148  fCubicVolume= 0.;
149  fSurfaceArea= 0.;
150  delete fpPolyhedron; fpPolyhedron = 0;
151 }
152 
154 {
155  if(dy > 2*kCarTolerance) // limit to thickness of surfaces
156  {
157  fDy = dy;
158  }
159  else
160  {
161  std::ostringstream message;
162  message << "Dimension Y too small for solid: " << GetName() << "!"
163  << G4endl
164  << " hY = " << dy;
165  G4Exception("G4Box::SetYHalfLength()", "GeomSolids0002",
166  FatalException, message);
167  }
168  fCubicVolume= 0.;
169  fSurfaceArea= 0.;
170  delete fpPolyhedron; fpPolyhedron = 0;
171 }
172 
174 {
175  if(dz > 2*kCarTolerance) // limit to thickness of surfaces
176  {
177  fDz = dz;
178  }
179  else
180  {
181  std::ostringstream message;
182  message << "Dimension Z too small for solid: " << GetName() << "!"
183  << G4endl
184  << " hZ = " << dz;
185  G4Exception("G4Box::SetZHalfLength()", "GeomSolids0002",
186  FatalException, message);
187  }
188  fCubicVolume= 0.;
189  fSurfaceArea= 0.;
190  delete fpPolyhedron; fpPolyhedron = 0;
191 }
192 
194 //
195 // Dispatch to parameterisation for replication mechanism dimension
196 // computation & modification.
197 
199  const G4int n,
200  const G4VPhysicalVolume* pRep)
201 {
202  p->ComputeDimensions(*this,n,pRep);
203 }
204 
206 //
207 // Calculate extent under transform and specified limit
208 
210  const G4VoxelLimits& pVoxelLimit,
211  const G4AffineTransform& pTransform,
212  G4double& pMin, G4double& pMax) const
213 {
214  if (!pTransform.IsRotated())
215  {
216  // Special case handling for unrotated boxes
217  // Compute x/y/z mins and maxs respecting limits, with early returns
218  // if outside limits. Then switch() on pAxis
219 
220  G4double xoffset,xMin,xMax;
221  G4double yoffset,yMin,yMax;
222  G4double zoffset,zMin,zMax;
223 
224  xoffset = pTransform.NetTranslation().x() ;
225  xMin = xoffset - fDx ;
226  xMax = xoffset + fDx ;
227 
228  if (pVoxelLimit.IsXLimited())
229  {
230  if ((xMin > pVoxelLimit.GetMaxXExtent()+kCarTolerance) ||
231  (xMax < pVoxelLimit.GetMinXExtent()-kCarTolerance)) { return false ; }
232  else
233  {
234  xMin = std::max(xMin, pVoxelLimit.GetMinXExtent());
235  xMax = std::min(xMax, pVoxelLimit.GetMaxXExtent());
236  }
237  }
238  yoffset = pTransform.NetTranslation().y() ;
239  yMin = yoffset - fDy ;
240  yMax = yoffset + fDy ;
241 
242  if (pVoxelLimit.IsYLimited())
243  {
244  if ((yMin > pVoxelLimit.GetMaxYExtent()+kCarTolerance) ||
245  (yMax < pVoxelLimit.GetMinYExtent()-kCarTolerance)) { return false ; }
246  else
247  {
248  yMin = std::max(yMin, pVoxelLimit.GetMinYExtent());
249  yMax = std::min(yMax, pVoxelLimit.GetMaxYExtent());
250  }
251  }
252  zoffset = pTransform.NetTranslation().z() ;
253  zMin = zoffset - fDz ;
254  zMax = zoffset + fDz ;
255 
256  if (pVoxelLimit.IsZLimited())
257  {
258  if ((zMin > pVoxelLimit.GetMaxZExtent()+kCarTolerance) ||
259  (zMax < pVoxelLimit.GetMinZExtent()-kCarTolerance)) { return false ; }
260  else
261  {
262  zMin = std::max(zMin, pVoxelLimit.GetMinZExtent());
263  zMax = std::min(zMax, pVoxelLimit.GetMaxZExtent());
264  }
265  }
266  switch (pAxis)
267  {
268  case kXAxis:
269  pMin = xMin ;
270  pMax = xMax ;
271  break ;
272  case kYAxis:
273  pMin=yMin;
274  pMax=yMax;
275  break;
276  case kZAxis:
277  pMin=zMin;
278  pMax=zMax;
279  break;
280  default:
281  break;
282  }
283  pMin -= kCarTolerance ;
284  pMax += kCarTolerance ;
285 
286  return true;
287  }
288  else // General rotated case - create and clip mesh to boundaries
289  {
290  G4bool existsAfterClip = false ;
291  G4ThreeVectorList* vertices ;
292 
293  pMin = +kInfinity ;
294  pMax = -kInfinity ;
295 
296  // Calculate rotated vertex coordinates
297 
298  vertices = CreateRotatedVertices(pTransform) ;
299  ClipCrossSection(vertices,0,pVoxelLimit,pAxis,pMin,pMax) ;
300  ClipCrossSection(vertices,4,pVoxelLimit,pAxis,pMin,pMax) ;
301  ClipBetweenSections(vertices,0,pVoxelLimit,pAxis,pMin,pMax) ;
302 
303  if (pVoxelLimit.IsLimited(pAxis) == false)
304  {
305  if ( (pMin != kInfinity) || (pMax != -kInfinity) )
306  {
307  existsAfterClip = true ;
308 
309  // Add 2*tolerance to avoid precision troubles
310 
311  pMin -= kCarTolerance;
312  pMax += kCarTolerance;
313  }
314  }
315  else
316  {
317  G4ThreeVector clipCentre(
318  ( pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
319  ( pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
320  ( pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
321 
322  if ( (pMin != kInfinity) || (pMax != -kInfinity) )
323  {
324  existsAfterClip = true ;
325 
326 
327  // Check to see if endpoints are in the solid
328 
329  clipCentre(pAxis) = pVoxelLimit.GetMinExtent(pAxis);
330 
331  if (Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside)
332  {
333  pMin = pVoxelLimit.GetMinExtent(pAxis);
334  }
335  else
336  {
337  pMin -= kCarTolerance;
338  }
339  clipCentre(pAxis) = pVoxelLimit.GetMaxExtent(pAxis);
340 
341  if (Inside(pTransform.Inverse().TransformPoint(clipCentre)) != kOutside)
342  {
343  pMax = pVoxelLimit.GetMaxExtent(pAxis);
344  }
345  else
346  {
347  pMax += kCarTolerance;
348  }
349  }
350 
351  // Check for case where completely enveloping clipping volume
352  // If point inside then we are confident that the solid completely
353  // envelopes the clipping volume. Hence set min/max extents according
354  // to clipping volume extents along the specified axis.
355 
356  else if (Inside(pTransform.Inverse().TransformPoint(clipCentre))
357  != kOutside)
358  {
359  existsAfterClip = true ;
360  pMin = pVoxelLimit.GetMinExtent(pAxis) ;
361  pMax = pVoxelLimit.GetMaxExtent(pAxis) ;
362  }
363  }
364  delete vertices;
365  return existsAfterClip;
366  }
367 }
368 
370 //
371 // Return whether point inside/outside/on surface, using tolerance
372 
374 {
375  EInside in = kOutside ;
376  G4ThreeVector q(std::fabs(p.x()), std::fabs(p.y()), std::fabs(p.z()));
377 
378  if ( q.x() <= (fDx - delta) )
379  {
380  if (q.y() <= (fDy - delta) )
381  {
382  if ( q.z() <= (fDz - delta) ) { in = kInside ; }
383  else if ( q.z() <= (fDz + delta) ) { in = kSurface ; }
384  }
385  else if ( q.y() <= (fDy + delta) )
386  {
387  if ( q.z() <= (fDz + delta) ) { in = kSurface ; }
388  }
389  }
390  else if ( q.x() <= (fDx + delta) )
391  {
392  if ( q.y() <= (fDy + delta) )
393  {
394  if ( q.z() <= (fDz + delta) ) { in = kSurface ; }
395  }
396  }
397  return in ;
398 }
399 
401 //
402 // Calculate side nearest to p, and return normal
403 // If two sides are equidistant, normal of first side (x/y/z)
404 // encountered returned
405 
407 {
408  G4double distx, disty, distz ;
409  G4ThreeVector norm(0.,0.,0.);
410 
411  // Calculate distances as if in 1st octant
412 
413  distx = std::fabs(std::fabs(p.x()) - fDx) ;
414  disty = std::fabs(std::fabs(p.y()) - fDy) ;
415  distz = std::fabs(std::fabs(p.z()) - fDz) ;
416 
417  // New code for particle on surface including edges and corners with specific
418  // normals
419 
420  const G4ThreeVector nX = G4ThreeVector( 1.0, 0,0 );
421  const G4ThreeVector nmX = G4ThreeVector(-1.0, 0,0 );
422  const G4ThreeVector nY = G4ThreeVector( 0, 1.0,0 );
423  const G4ThreeVector nmY = G4ThreeVector( 0,-1.0,0 );
424  const G4ThreeVector nZ = G4ThreeVector( 0, 0, 1.0);
425  const G4ThreeVector nmZ = G4ThreeVector( 0, 0,- 1.0);
426 
427  G4ThreeVector normX(0.,0.,0.), normY(0.,0.,0.), normZ(0.,0.,0.);
428  G4ThreeVector sumnorm(0., 0., 0.);
429  G4int noSurfaces=0;
430 
431  if (distx <= delta) // on X/mX surface and around
432  {
433  noSurfaces ++;
434  if ( p.x() >= 0. ) { normX= nX ; } // on +X surface : (1,0,0)
435  else { normX= nmX; } // (-1,0,0)
436  sumnorm= normX;
437  }
438 
439  if (disty <= delta) // on one of the +Y or -Y surfaces
440  {
441  noSurfaces ++;
442  if ( p.y() >= 0. ) { normY= nY; } // on +Y surface
443  else { normY= nmY; }
444  sumnorm += normY;
445  }
446 
447  if (distz <= delta) // on one of the +Z or -Z surfaces
448  {
449  noSurfaces ++;
450  if ( p.z() >= 0. ) { normZ= nZ; } // on +Z surface
451  else { normZ= nmZ; }
452  sumnorm += normZ;
453  }
454 
455  static const G4double invSqrt2 = 1.0 / std::sqrt(2.0);
456  static const G4double invSqrt3 = 1.0 / std::sqrt(3.0);
457 
458  if( noSurfaces > 0 )
459  {
460  if( noSurfaces == 1 )
461  {
462  norm= sumnorm;
463  }
464  else
465  {
466  // norm = sumnorm . unit();
467  if( noSurfaces == 2 )
468  {
469  // 2 surfaces -> on edge
470  norm = invSqrt2 * sumnorm;
471  }
472  else
473  {
474  // 3 surfaces (on corner)
475  norm = invSqrt3 * sumnorm;
476  }
477  }
478  }
479  else
480  {
481 #ifdef G4CSGDEBUG
482  G4Exception("G4Box::SurfaceNormal(p)", "Notification", JustWarning,
483  "Point p is not on surface !?" );
484 #endif
485  norm = ApproxSurfaceNormal(p);
486  }
487 
488  return norm;
489 }
490 
492 //
493 // Algorithm for SurfaceNormal() following the original specification
494 // for points not on the surface
495 
497 {
498  G4double distx, disty, distz ;
499  G4ThreeVector norm(0.,0.,0.);
500 
501  // Calculate distances as if in 1st octant
502 
503  distx = std::fabs(std::fabs(p.x()) - fDx) ;
504  disty = std::fabs(std::fabs(p.y()) - fDy) ;
505  distz = std::fabs(std::fabs(p.z()) - fDz) ;
506 
507  if ( distx <= disty )
508  {
509  if ( distx <= distz ) // Closest to X
510  {
511  if ( p.x() < 0 ) { norm = G4ThreeVector(-1.0,0,0) ; }
512  else { norm = G4ThreeVector( 1.0,0,0) ; }
513  }
514  else // Closest to Z
515  {
516  if ( p.z() < 0 ) { norm = G4ThreeVector(0,0,-1.0) ; }
517  else { norm = G4ThreeVector(0,0, 1.0) ; }
518  }
519  }
520  else
521  {
522  if ( disty <= distz ) // Closest to Y
523  {
524  if ( p.y() < 0 ) { norm = G4ThreeVector(0,-1.0,0) ; }
525  else { norm = G4ThreeVector(0, 1.0,0) ; }
526  }
527  else // Closest to Z
528  {
529  if ( p.z() < 0 ) { norm = G4ThreeVector(0,0,-1.0) ; }
530  else { norm = G4ThreeVector(0,0, 1.0) ; }
531  }
532  }
533  return norm;
534 }
535 
537 //
538 // Calculate distance to box from an outside point
539 // - return kInfinity if no intersection.
540 //
541 // ALGORITHM:
542 //
543 // Check that if point lies outside x/y/z extent of box, travel is towards
544 // the box (ie. there is a possibility of an intersection)
545 //
546 // Calculate pairs of minimum and maximum distances for x/y/z travel for
547 // intersection with the box's x/y/z extent.
548 // If there is a valid intersection, it is given by the maximum min distance
549 // (ie. distance to satisfy x/y/z intersections) *if* <= minimum max distance
550 // (ie. distance after which 1+ of x/y/z intersections not satisfied)
551 //
552 // NOTE:
553 //
554 // `Inside' safe - meaningful answers given if point is inside the exact
555 // shape.
556 
558  const G4ThreeVector& v) const
559 {
560  G4double safx, safy, safz ;
561  G4double smin=0.0, sminy, sminz ; // , sminx ;
562  G4double smax=kInfinity, smaxy, smaxz ; // , smaxx ; // they always > 0
563  G4double stmp ;
564  G4double sOut=kInfinity, sOuty=kInfinity, sOutz=kInfinity ;
565 
566  safx = std::fabs(p.x()) - fDx ; // minimum distance to x surface of shape
567  safy = std::fabs(p.y()) - fDy ;
568  safz = std::fabs(p.z()) - fDz ;
569 
570  // Will we intersect?
571  // If safx/y/z is >-tol/2 the point is outside/on the box's x/y/z extent.
572  // If both p.x/y/z and v.x/y/z repectively are both positive/negative,
573  // travel is in a direction away from the shape.
574 
575  if ( ((p.x()*v.x() >= 0.0) && (safx > -delta))
576  || ((p.y()*v.y() >= 0.0) && (safy > -delta))
577  || ((p.z()*v.z() >= 0.0) && (safz > -delta)) )
578  {
579  return kInfinity ; // travel away or parallel within tolerance
580  }
581 
582  // Compute min / max distances for x/y/z travel:
583  // X Planes
584 
585  if ( v.x() ) // != 0
586  {
587  stmp = 1.0/std::fabs(v.x()) ;
588 
589  if (safx >= 0.0)
590  {
591  smin = safx*stmp ;
592  smax = (fDx+std::fabs(p.x()))*stmp ;
593  }
594  else
595  {
596  if (v.x() < 0) { sOut = (fDx + p.x())*stmp ; }
597  else { sOut = (fDx - p.x())*stmp ; }
598  }
599  }
600 
601  // Y Planes
602 
603  if ( v.y() ) // != 0
604  {
605  stmp = 1.0/std::fabs(v.y()) ;
606 
607  if (safy >= 0.0)
608  {
609  sminy = safy*stmp ;
610  smaxy = (fDy+std::fabs(p.y()))*stmp ;
611 
612  if (sminy > smin) { smin=sminy ; }
613  if (smaxy < smax) { smax=smaxy ; }
614 
615  if (smin >= (smax-delta))
616  {
617  return kInfinity ; // touch XY corner
618  }
619  }
620  else
621  {
622  if (v.y() < 0) { sOuty = (fDy + p.y())*stmp ; }
623  else { sOuty = (fDy - p.y())*stmp ; }
624  if( sOuty < sOut ) { sOut = sOuty ; }
625  }
626  }
627 
628  // Z planes
629 
630  if ( v.z() ) // != 0
631  {
632  stmp = 1.0/std::fabs(v.z()) ;
633 
634  if ( safz >= 0.0 )
635  {
636  sminz = safz*stmp ;
637  smaxz = (fDz+std::fabs(p.z()))*stmp ;
638 
639  if (sminz > smin) { smin = sminz ; }
640  if (smaxz < smax) { smax = smaxz ; }
641 
642  if (smin >= (smax-delta))
643  {
644  return kInfinity ; // touch ZX or ZY corners
645  }
646  }
647  else
648  {
649  if (v.z() < 0) { sOutz = (fDz + p.z())*stmp ; }
650  else { sOutz = (fDz - p.z())*stmp ; }
651  if( sOutz < sOut ) { sOut = sOutz ; }
652  }
653  }
654 
655  if (sOut <= (smin + delta)) // travel over edge
656  {
657  return kInfinity ;
658  }
659  if (smin < delta) { smin = 0.0 ; }
660 
661  return smin ;
662 }
663 
665 //
666 // Appoximate distance to box.
667 // Returns largest perpendicular distance to the closest x/y/z sides of
668 // the box, which is the most fast estimation of the shortest distance to box
669 // - If inside return 0
670 
672 {
673  G4double safex, safey, safez, safe = 0.0 ;
674 
675  safex = std::fabs(p.x()) - fDx ;
676  safey = std::fabs(p.y()) - fDy ;
677  safez = std::fabs(p.z()) - fDz ;
678 
679  if (safex > safe) { safe = safex ; }
680  if (safey > safe) { safe = safey ; }
681  if (safez > safe) { safe = safez ; }
682 
683  return safe ;
684 }
685 
687 //
688 // Calculate distance to surface of box from inside
689 // by calculating distances to box's x/y/z planes.
690 // Smallest distance is exact distance to exiting.
691 // - Eliminate one side of each pair by considering direction of v
692 // - when leaving a surface & v.close, return 0
693 
695  const G4bool calcNorm,
696  G4bool *validNorm,G4ThreeVector *n) const
697 {
698  ESide side = kUndefined ;
699  G4double pdist,stmp,snxt=kInfinity;
700 
701  if (calcNorm) { *validNorm = true ; } // All normals are valid
702 
703  if (v.x() > 0) // X planes
704  {
705  pdist = fDx - p.x() ;
706 
707  if (pdist > delta)
708  {
709  snxt = pdist/v.x() ;
710  side = kPX ;
711  }
712  else
713  {
714  if (calcNorm) { *n = G4ThreeVector(1,0,0) ; }
715  return snxt = 0 ;
716  }
717  }
718  else if (v.x() < 0)
719  {
720  pdist = fDx + p.x() ;
721 
722  if (pdist > delta)
723  {
724  snxt = -pdist/v.x() ;
725  side = kMX ;
726  }
727  else
728  {
729  if (calcNorm) { *n = G4ThreeVector(-1,0,0) ; }
730  return snxt = 0 ;
731  }
732  }
733 
734  if (v.y() > 0) // Y planes
735  {
736  pdist = fDy-p.y();
737 
738  if (pdist > delta)
739  {
740  stmp = pdist/v.y();
741 
742  if (stmp < snxt)
743  {
744  snxt = stmp;
745  side = kPY;
746  }
747  }
748  else
749  {
750  if (calcNorm) { *n = G4ThreeVector(0,1,0) ; }
751  return snxt = 0 ;
752  }
753  }
754  else if (v.y() < 0)
755  {
756  pdist = fDy + p.y() ;
757 
758  if (pdist > delta)
759  {
760  stmp = -pdist/v.y();
761 
762  if ( stmp < snxt )
763  {
764  snxt = stmp;
765  side = kMY;
766  }
767  }
768  else
769  {
770  if (calcNorm) { *n = G4ThreeVector(0,-1,0) ; }
771  return snxt = 0 ;
772  }
773  }
774 
775  if (v.z() > 0) // Z planes
776  {
777  pdist = fDz-p.z();
778 
779  if ( pdist > delta )
780  {
781  stmp = pdist/v.z();
782 
783  if ( stmp < snxt )
784  {
785  snxt = stmp;
786  side = kPZ;
787  }
788  }
789  else
790  {
791  if (calcNorm) { *n = G4ThreeVector(0,0,1) ; }
792  return snxt = 0 ;
793  }
794  }
795  else if (v.z() < 0)
796  {
797  pdist = fDz + p.z();
798 
799  if ( pdist > delta )
800  {
801  stmp = -pdist/v.z();
802 
803  if ( stmp < snxt )
804  {
805  snxt = stmp;
806  side = kMZ;
807  }
808  }
809  else
810  {
811  if (calcNorm) { *n = G4ThreeVector(0,0,-1) ; }
812  return snxt = 0 ;
813  }
814  }
815 
816  if (calcNorm)
817  {
818  switch (side)
819  {
820  case kPX:
821  *n=G4ThreeVector(1,0,0);
822  break;
823  case kMX:
824  *n=G4ThreeVector(-1,0,0);
825  break;
826  case kPY:
827  *n=G4ThreeVector(0,1,0);
828  break;
829  case kMY:
830  *n=G4ThreeVector(0,-1,0);
831  break;
832  case kPZ:
833  *n=G4ThreeVector(0,0,1);
834  break;
835  case kMZ:
836  *n=G4ThreeVector(0,0,-1);
837  break;
838  default:
839  G4cout << G4endl;
840  DumpInfo();
841  std::ostringstream message;
842  G4int oldprc = message.precision(16);
843  message << "Undefined side for valid surface normal to solid."
844  << G4endl
845  << "Position:" << G4endl << G4endl
846  << "p.x() = " << p.x()/mm << " mm" << G4endl
847  << "p.y() = " << p.y()/mm << " mm" << G4endl
848  << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl
849  << "Direction:" << G4endl << G4endl
850  << "v.x() = " << v.x() << G4endl
851  << "v.y() = " << v.y() << G4endl
852  << "v.z() = " << v.z() << G4endl << G4endl
853  << "Proposed distance :" << G4endl << G4endl
854  << "snxt = " << snxt/mm << " mm" << G4endl;
855  message.precision(oldprc);
856  G4Exception("G4Box::DistanceToOut(p,v,..)", "GeomSolids1002",
857  JustWarning, message);
858  break;
859  }
860  }
861  return snxt;
862 }
863 
865 //
866 // Calculate exact shortest distance to any boundary from inside
867 // - If outside return 0
868 
870 {
871  G4double safx1,safx2,safy1,safy2,safz1,safz2,safe=0.0;
872 
873 #ifdef G4CSGDEBUG
874  if( Inside(p) == kOutside )
875  {
876  G4int oldprc = G4cout.precision(16) ;
877  G4cout << G4endl ;
878  DumpInfo();
879  G4cout << "Position:" << G4endl << G4endl ;
880  G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
881  G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
882  G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
883  G4cout.precision(oldprc) ;
884  G4Exception("G4Box::DistanceToOut(p)", "GeomSolids1002",
885  JustWarning, "Point p is outside !?" );
886  }
887 #endif
888 
889  safx1 = fDx - p.x() ;
890  safx2 = fDx + p.x() ;
891  safy1 = fDy - p.y() ;
892  safy2 = fDy + p.y() ;
893  safz1 = fDz - p.z() ;
894  safz2 = fDz + p.z() ;
895 
896  // shortest Dist to any boundary now MIN(safx1,safx2,safy1..)
897 
898  if (safx2 < safx1) { safe = safx2; }
899  else { safe = safx1; }
900  if (safy1 < safe) { safe = safy1; }
901  if (safy2 < safe) { safe = safy2; }
902  if (safz1 < safe) { safe = safz1; }
903  if (safz2 < safe) { safe = safz2; }
904 
905  if (safe < 0) { safe = 0 ; }
906  return safe ;
907 }
908 
910 //
911 // Create a List containing the transformed vertices
912 // Ordering [0-3] -fDz cross section
913 // [4-7] +fDz cross section such that [0] is below [4],
914 // [1] below [5] etc.
915 // Note:
916 // Caller has deletion resposibility
917 
920 {
921  G4ThreeVectorList* vertices = new G4ThreeVectorList();
922 
923  if (vertices)
924  {
925  vertices->reserve(8);
926  G4ThreeVector vertex0(-fDx,-fDy,-fDz) ;
927  G4ThreeVector vertex1(fDx,-fDy,-fDz) ;
928  G4ThreeVector vertex2(fDx,fDy,-fDz) ;
929  G4ThreeVector vertex3(-fDx,fDy,-fDz) ;
930  G4ThreeVector vertex4(-fDx,-fDy,fDz) ;
931  G4ThreeVector vertex5(fDx,-fDy,fDz) ;
932  G4ThreeVector vertex6(fDx,fDy,fDz) ;
933  G4ThreeVector vertex7(-fDx,fDy,fDz) ;
934 
935  vertices->push_back(pTransform.TransformPoint(vertex0));
936  vertices->push_back(pTransform.TransformPoint(vertex1));
937  vertices->push_back(pTransform.TransformPoint(vertex2));
938  vertices->push_back(pTransform.TransformPoint(vertex3));
939  vertices->push_back(pTransform.TransformPoint(vertex4));
940  vertices->push_back(pTransform.TransformPoint(vertex5));
941  vertices->push_back(pTransform.TransformPoint(vertex6));
942  vertices->push_back(pTransform.TransformPoint(vertex7));
943  }
944  else
945  {
946  DumpInfo();
947  G4Exception("G4Box::CreateRotatedVertices()",
948  "GeomSolids0003", FatalException,
949  "Error in allocation of vertices. Out of memory !");
950  }
951  return vertices;
952 }
953 
955 //
956 // GetEntityType
957 
959 {
960  return G4String("G4Box");
961 }
962 
964 //
965 // Stream object contents to an output stream
966 
967 std::ostream& G4Box::StreamInfo(std::ostream& os) const
968 {
969  G4int oldprc = os.precision(16);
970  os << "-----------------------------------------------------------\n"
971  << " *** Dump for solid - " << GetName() << " ***\n"
972  << " ===================================================\n"
973  << " Solid type: G4Box\n"
974  << " Parameters: \n"
975  << " half length X: " << fDx/mm << " mm \n"
976  << " half length Y: " << fDy/mm << " mm \n"
977  << " half length Z: " << fDz/mm << " mm \n"
978  << "-----------------------------------------------------------\n";
979  os.precision(oldprc);
980 
981  return os;
982 }
983 
985 //
986 // GetPointOnSurface
987 //
988 // Return a point (G4ThreeVector) randomly and uniformly selected
989 // on the solid surface
990 
992 {
993  G4double px, py, pz, select, sumS;
994  G4double Sxy = fDx*fDy, Sxz = fDx*fDz, Syz = fDy*fDz;
995 
996  sumS = Sxy + Sxz + Syz;
997  select = sumS*G4UniformRand();
998 
999  if( select < Sxy )
1000  {
1001  px = -fDx +2*fDx*G4UniformRand();
1002  py = -fDy +2*fDy*G4UniformRand();
1003 
1004  if(G4UniformRand() > 0.5) { pz = fDz; }
1005  else { pz = -fDz; }
1006  }
1007  else if ( ( select - Sxy ) < Sxz )
1008  {
1009  px = -fDx +2*fDx*G4UniformRand();
1010  pz = -fDz +2*fDz*G4UniformRand();
1011 
1012  if(G4UniformRand() > 0.5) { py = fDy; }
1013  else { py = -fDy; }
1014  }
1015  else
1016  {
1017  py = -fDy +2*fDy*G4UniformRand();
1018  pz = -fDz +2*fDz*G4UniformRand();
1019 
1020  if(G4UniformRand() > 0.5) { px = fDx; }
1021  else { px = -fDx; }
1022  }
1023  return G4ThreeVector(px,py,pz);
1024 }
1025 
1027 //
1028 // Make a clone of the object
1029 //
1031 {
1032  return new G4Box(*this);
1033 }
1034 
1036 //
1037 // Methods for visualisation
1038 
1040 {
1041  scene.AddSolid (*this);
1042 }
1043 
1045 {
1046  return G4VisExtent (-fDx, fDx, -fDy, fDy, -fDz, fDz);
1047 }
1048 
1050 {
1051  return new G4PolyhedronBox (fDx, fDy, fDz);
1052 }
1053 #endif
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Box.cc:198
G4ThreeVector GetPointOnSurface() const
Definition: G4Box.cc:991
G4String GetName() const
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:347
void SetZHalfLength(G4double dz)
Definition: G4Box.cc:173
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Box.cc:694
G4double fDy
Definition: G4Box.hh:151
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Box.cc:209
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
CLHEP::Hep3Vector G4ThreeVector
G4Box & operator=(const G4Box &rhs)
Definition: G4Box.cc:110
G4double delta
Definition: G4Box.hh:152
G4AffineTransform Inverse() const
Definition: G4Box.hh:64
G4bool IsYLimited() const
G4GeometryType GetEntityType() const
Definition: G4Box.cc:958
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Box.cc:1039
G4bool IsRotated() const
G4double fCubicVolume
Definition: G4CSGSolid.hh:78
G4ThreeVector NetTranslation() const
G4Box(const G4String &pName, G4double pX, G4double pY, G4double pZ)
Definition: G4Box.cc:60
G4bool IsXLimited() const
G4double a
Definition: TRTMaterials.hh:39
G4VisExtent GetExtent() const
Definition: G4Box.cc:1044
virtual void AddSolid(const G4Box &)=0
G4VSolid * Clone() const
Definition: G4Box.cc:1030
int G4int
Definition: G4Types.hh:78
G4Polyhedron * CreatePolyhedron() const
Definition: G4Box.cc:1049
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Box.cc:967
const G4int smax
void DumpInfo() const
G4double fDz
Definition: G4Box.hh:151
G4double GetMaxXExtent() const
G4double GetMinZExtent() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4double fSurfaceArea
Definition: G4CSGSolid.hh:79
G4GLOB_DLL std::ostream G4cout
G4bool IsLimited() const
G4Polyhedron * fpPolyhedron
Definition: G4CSGSolid.hh:80
bool G4bool
Definition: G4Types.hh:79
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4CSGSolid.cc:124
G4double fDx
Definition: G4Box.hh:151
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
const G4int n
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Box.cc:919
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
ESide
Definition: G4Box.hh:140
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Box.cc:406
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4double GetMinXExtent() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Box.cc:557
EInside Inside(const G4ThreeVector &p) const
Definition: G4Box.cc:373
G4double GetMaxZExtent() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
virtual ~G4Box()
Definition: G4Box.cc:92
EAxis
Definition: geomdefs.hh:54
void SetYHalfLength(G4double dy)
Definition: G4Box.cc:153
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetXHalfLength(G4double dx)
Definition: G4Box.cc:133
#define G4endl
Definition: G4ios.hh:61
G4double GetMaxYExtent() const
G4double kCarTolerance
Definition: G4VSolid.hh:305
double G4double
Definition: G4Types.hh:76
G4double GetMaxExtent(const EAxis pAxis) const
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:82
G4bool IsZLimited() const
static const double mm
Definition: G4SIunits.hh:102
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:378
G4double GetMinExtent(const EAxis pAxis) const
G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector &p) const
Definition: G4Box.cc:496