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