ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Orb.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4Orb.cc
1 // ********************************************************************
2 // * License and Disclaimer *
3 // * *
4 // * The Geant4 software is copyright of the Copyright Holders of *
5 // * the Geant4 Collaboration. It is provided under the terms and *
6 // * conditions of the Geant4 Software License, included in the file *
7 // * LICENSE and available at http://cern.ch/geant4/license . These *
8 // * include a list of copyright holders. *
9 // * *
10 // * Neither the authors of this software system, nor their employing *
11 // * institutes,nor the agencies providing financial support for this *
12 // * work make any representation or warranty, express or implied, *
13 // * regarding this software system or assume any liability for its *
14 // * use. Please see the license in the file LICENSE and URL above *
15 // * for the full disclaimer and the limitation of liability. *
16 // * *
17 // * This code implementation is the result of the scientific and *
18 // * technical work of the GEANT4 collaboration. *
19 // * By using, copying, modifying or distributing the software (or *
20 // * any work based on the software) you agree to acknowledge its *
21 // * use in resulting scientific publications, and indicate your *
22 // * acceptance of all terms of the Geant4 Software license. *
23 // ********************************************************************
24 //
25 // Implementation for G4Orb class
26 //
27 // 20.08.03 V.Grichine - created
28 // 08.08.17 E.Tcherniaev - complete revision, speed-up
29 // --------------------------------------------------------------------
30 
31 #include "G4Orb.hh"
32 
33 #if !defined(G4GEOM_USE_UORB)
34 
35 #include "G4TwoVector.hh"
36 #include "G4VoxelLimits.hh"
37 #include "G4AffineTransform.hh"
38 #include "G4GeometryTolerance.hh"
39 #include "G4BoundingEnvelope.hh"
40 
41 #include "G4VPVParameterisation.hh"
42 
43 #include "G4RandomDirection.hh"
44 #include "Randomize.hh"
45 
46 #include "G4VGraphicsScene.hh"
47 #include "G4VisExtent.hh"
48 
49 using namespace CLHEP;
50 
52 //
53 // Constructor
54 
55 G4Orb::G4Orb( const G4String& pName, G4double pRmax )
56  : G4CSGSolid(pName), fRmax(pRmax)
57 {
58  Initialize();
59 }
60 
62 //
63 // Fake default constructor - sets only member data and allocates memory
64 // for usage restricted to object persistency
65 
66 G4Orb::G4Orb( __void__& a )
67  : G4CSGSolid(a)
68 {
69 }
70 
72 //
73 // Destructor
74 
76 {
77 }
78 
80 //
81 // Copy constructor
82 
83 G4Orb::G4Orb(const G4Orb& rhs)
84  : G4CSGSolid(rhs), fRmax(rhs.fRmax), halfRmaxTol(rhs.halfRmaxTol),
85  sqrRmaxPlusTol(rhs.sqrRmaxPlusTol), sqrRmaxMinusTol(rhs.sqrRmaxMinusTol)
86 {
87 }
88 
90 //
91 // Assignment operator
92 
94 {
95  // Check assignment to self
96  //
97  if (this == &rhs) { return *this; }
98 
99  // Copy base class data
100  //
102 
103  // Copy data
104  //
105  fRmax = rhs.fRmax;
106  halfRmaxTol = rhs.halfRmaxTol;
109 
110  return *this;
111 }
112 
114 //
115 // Check radius and initialize dada members
116 
118 {
119  const G4double fEpsilon = 2.e-11; // relative tolerance of fRmax
120 
121  // Check radius
122  //
123  if ( fRmax < 10*kCarTolerance )
124  {
125  G4Exception("G4Orb::Initialize()", "GeomSolids0002", FatalException,
126  "Invalid radius < 10*kCarTolerance.");
127  }
128  halfRmaxTol = 0.5 * std::max(kCarTolerance, fEpsilon*fRmax);
129  G4double rmaxPlusTol = fRmax + halfRmaxTol;
130  G4double rmaxMinusTol = fRmax - halfRmaxTol;
131  sqrRmaxPlusTol = rmaxPlusTol*rmaxPlusTol;
132  sqrRmaxMinusTol = rmaxMinusTol*rmaxMinusTol;
133 }
134 
136 //
137 // Dispatch to parameterisation for replication mechanism dimension
138 // computation & modification
139 
141  const G4int n,
142  const G4VPhysicalVolume* pRep )
143 {
144  p->ComputeDimensions(*this,n,pRep);
145 }
146 
148 //
149 // Get bounding box
150 
152 {
154  pMin.set(-radius,-radius,-radius);
155  pMax.set( radius, radius, radius);
156 
157  // Check correctness of the bounding box
158  //
159  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
160  {
161  std::ostringstream message;
162  message << "Bad bounding box (min >= max) for solid: "
163  << GetName() << " !"
164  << "\npMin = " << pMin
165  << "\npMax = " << pMax;
166  G4Exception("G4Orb::BoundingLimits()", "GeomMgt0001",
167  JustWarning, message);
168  DumpInfo();
169  }
170 }
171 
173 //
174 // Calculate extent under transform and specified limit
175 
177  const G4VoxelLimits& pVoxelLimit,
178  const G4AffineTransform& pTransform,
179  G4double& pMin, G4double& pMax) const
180 {
181  G4ThreeVector bmin, bmax;
182  G4bool exist;
183 
184  // Get bounding box
185  BoundingLimits(bmin,bmax);
186 
187  // Check bounding box
188  G4BoundingEnvelope bbox(bmin,bmax);
189 #ifdef G4BBOX_EXTENT
190  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
191 #endif
192  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
193  {
194  return exist = (pMin < pMax) ? true : false;
195  }
196 
197  // Find bounding envelope and calculate extent
198  //
199  static const G4int NTHETA = 8; // number of steps along Theta
200  static const G4int NPHI = 16; // number of steps along Phi
201  static const G4double sinHalfTheta = std::sin(halfpi/NTHETA);
202  static const G4double cosHalfTheta = std::cos(halfpi/NTHETA);
203  static const G4double sinHalfPhi = std::sin(pi/NPHI);
204  static const G4double cosHalfPhi = std::cos(pi/NPHI);
205  static const G4double sinStepTheta = 2.*sinHalfTheta*cosHalfTheta;
206  static const G4double cosStepTheta = 1. - 2.*sinHalfTheta*sinHalfTheta;
207  static const G4double sinStepPhi = 2.*sinHalfPhi*cosHalfPhi;
208  static const G4double cosStepPhi = 1. - 2.*sinHalfPhi*sinHalfPhi;
209 
211  G4double rtheta = radius/cosHalfTheta;
212  G4double rphi = rtheta/cosHalfPhi;
213 
214  // set reference circle
215  G4TwoVector xy[NPHI];
216  G4double sinCurPhi = sinHalfPhi;
217  G4double cosCurPhi = cosHalfPhi;
218  for (G4int k=0; k<NPHI; ++k)
219  {
220  xy[k].set(cosCurPhi,sinCurPhi);
221  G4double sinTmpPhi = sinCurPhi;
222  sinCurPhi = sinCurPhi*cosStepPhi + cosCurPhi*sinStepPhi;
223  cosCurPhi = cosCurPhi*cosStepPhi - sinTmpPhi*sinStepPhi;
224  }
225 
226  // set bounding circles
227  G4ThreeVectorList circles[NTHETA];
228  for (G4int i=0; i<NTHETA; ++i) { circles[i].resize(NPHI); }
229 
230  G4double sinCurTheta = sinHalfTheta;
231  G4double cosCurTheta = cosHalfTheta;
232  for (G4int i=0; i<NTHETA; ++i)
233  {
234  G4double z = rtheta*cosCurTheta;
235  G4double rho = rphi*sinCurTheta;
236  for (G4int k=0; k<NPHI; ++k)
237  {
238  circles[i][k].set(rho*xy[k].x(),rho*xy[k].y(),z);
239  }
240  G4double sinTmpTheta = sinCurTheta;
241  sinCurTheta = sinCurTheta*cosStepTheta + cosCurTheta*sinStepTheta;
242  cosCurTheta = cosCurTheta*cosStepTheta - sinTmpTheta*sinStepTheta;
243  }
244 
245  // set envelope and calculate extent
246  std::vector<const G4ThreeVectorList *> polygons;
247  polygons.resize(NTHETA);
248  for (G4int i=0; i<NTHETA; ++i) { polygons[i] = &circles[i]; }
249 
250  G4BoundingEnvelope benv(bmin,bmax,polygons);
251  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
252  return exist;
253 }
254 
256 //
257 // Return whether point is inside/outside/on surface
258 
260 {
261  G4double rr = p.mag2();
262  if (rr > sqrRmaxPlusTol) return kOutside;
263  return (rr > sqrRmaxMinusTol) ? kSurface : kInside;
264 }
265 
267 //
268 // Return unit normal of surface closest to p
269 
271 {
272  return (1/p.mag())*p;
273 }
274 
276 //
277 // Calculate distance to the surface of the orb from outside
278 // - return kInfinity if no intersection or
279 // intersection distance <= tolerance
280 
282  const G4ThreeVector& v ) const
283 {
284  // Check if point is on the surface and traveling away
285  //
286  G4double rr = p.mag2();
287  G4double pv = p.dot(v);
288  if (rr >= sqrRmaxMinusTol && pv >= 0) return kInfinity;
289 
290  // Find intersection
291  //
292  // Sphere eqn: x^2 + y^2 + z^2 = R^2
293  //
294  // => (px + t*vx)^2 + (py + t*vy)^2 + (pz + t*vz)^2 = R^2
295  // => r^2 + 2t(p.v) + t^2 = R^2
296  // => tmin = -(p.v) - Sqrt((p.v)^2 - (r^2 - R^2))
297  //
298  G4double D = pv*pv - rr + fRmax*fRmax;
299  if (D < 0) return kInfinity; // no intersection
300 
301  G4double sqrtD = std::sqrt(D);
302  G4double dist = -pv - sqrtD;
303 
304  // Avoid rounding errors due to precision issues seen on 64 bits systems.
305  // Split long distances and recompute
306  //
307  G4double Dmax = 32*fRmax;
308  if (dist > Dmax)
309  {
310  dist = dist - 1.e-8*dist - fRmax; // to stay outside after the move
311  dist += DistanceToIn(p + dist*v, v);
312  return (dist >= kInfinity) ? kInfinity : dist;
313  }
314 
315  if (sqrtD*2 <= halfRmaxTol) return kInfinity; // touch
316  return (dist < halfRmaxTol) ? 0. : dist;
317 }
318 
320 //
321 // Calculate shortest distance to the boundary from outside
322 // - Return 0 if point is inside
323 
325 {
326  G4double dist = p.mag() - fRmax;
327  return (dist > 0) ? dist : 0.;
328 }
329 
331 //
332 // Calculate distance to the surface of the orb from inside and
333 // find normal at exit point, if required
334 // - when leaving the surface, return 0
335 
337  const G4ThreeVector& v,
338  const G4bool calcNorm,
339  G4bool* validNorm,
340  G4ThreeVector* n ) const
341 {
342  // Check if point is on the surface and traveling away
343  //
344  G4double rr = p.mag2();
345  G4double pv = p.dot(v);
346  if (rr >= sqrRmaxMinusTol && pv > 0)
347  {
348  if (calcNorm)
349  {
350  *validNorm = true;
351  *n = p*(1./std::sqrt(rr));
352  }
353  return 0.;
354  }
355 
356  // Find intersection
357  //
358  // Sphere eqn: x^2 + y^2 + z^2 = R^2
359  //
360  // => (px + t*vx)^2 + (py + t*vy)^2 + (pz + t*vz)^2 = R^2
361  // => r^2 + 2t(p.v) + t^2 = R^2
362  // => tmax = -(p.v) + Sqrt((p.v)^2 - (r^2 - R^2))
363  //
364  G4double D = pv*pv - rr + fRmax*fRmax;
365  G4double tmax = (D <= 0) ? 0. : std::sqrt(D) - pv;
366  if (tmax < halfRmaxTol) tmax = 0.;
367  if (calcNorm)
368  {
369  *validNorm = true;
370  G4ThreeVector ptmax = p + tmax*v;
371  *n = ptmax*(1./ptmax.mag());
372  }
373  return tmax;
374 }
375 
377 //
378 // Calculate distance (<=actual) to closest surface of shape from inside
379 
381 {
382 #ifdef G4CSGDEBUG
383  if( Inside(p) == kOutside )
384  {
385  std::ostringstream message;
386  G4int oldprc = message.precision(16);
387  message << "Point p is outside (!?) of solid: " << GetName() << "\n";
388  message << "Position:\n";
389  message << " p.x() = " << p.x()/mm << " mm\n";
390  message << " p.y() = " << p.y()/mm << " mm\n";
391  message << " p.z() = " << p.z()/mm << " mm";
392  G4cout.precision(oldprc);
393  G4Exception("G4Trap::DistanceToOut(p)", "GeomSolids1002",
394  JustWarning, message );
395  DumpInfo();
396  }
397 #endif
398  G4double dist = fRmax - p.mag();
399  return (dist > 0) ? dist : 0.;
400 }
401 
403 //
404 // G4EntityType
405 
407 {
408  return G4String("G4Orb");
409 }
410 
412 //
413 // Make a clone of the object
414 
416 {
417  return new G4Orb(*this);
418 }
419 
421 //
422 // Stream object contents to an output stream
423 
424 std::ostream& G4Orb::StreamInfo( std::ostream& os ) const
425 {
426  G4int oldprc = os.precision(16);
427  os << "-----------------------------------------------------------\n"
428  << " *** Dump for solid - " << GetName() << " ***\n"
429  << " ===================================================\n"
430  << " Solid type: G4Orb\n"
431  << " Parameters: \n"
432  << " outer radius: " << fRmax/mm << " mm \n"
433  << "-----------------------------------------------------------\n";
434  os.precision(oldprc);
435  return os;
436 }
437 
439 //
440 // GetPointOnSurface
441 
443 {
444  return fRmax * G4RandomDirection();
445 }
446 
448 //
449 // Methods for visualisation
450 
452 {
453  scene.AddSolid (*this);
454 }
455 
457 {
458  return G4VisExtent (-fRmax, fRmax, -fRmax, fRmax, -fRmax, fRmax);
459 }
460 
462 {
463  return new G4PolyhedronSphere (0., fRmax, 0., 2*pi, 0., pi);
464 }
465 
466 #endif