Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
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 // $Id: G4Orb.cc 105834 2017-08-23 08:14:34Z gcosmo $
26 //
27 // class G4Orb
28 //
29 // Implementation for G4Orb class
30 //
31 // History:
32 //
33 // 08.08.17 E.Tcherniaev - complete revision, speed-up
34 // 27.10.16 E.Tcherniaev - reimplemented CalculateExtent()
35 // 05.04.12 M.Kelsey - GetPointOnSurface() throw flat in cos(theta)
36 // 30.06.04 V.Grichine - bug fixed in DistanceToIn(p,v) on Rmax surface
37 // 20.08.03 V.Grichine - created
38 // --------------------------------------------------------------------
39 
40 #include "G4Orb.hh"
41 
42 #if !defined(G4GEOM_USE_UORB)
43 
44 #include "G4TwoVector.hh"
45 #include "G4VoxelLimits.hh"
46 #include "G4AffineTransform.hh"
47 #include "G4GeometryTolerance.hh"
48 #include "G4BoundingEnvelope.hh"
49 
50 #include "G4VPVParameterisation.hh"
51 
52 #include "G4RandomDirection.hh"
53 #include "Randomize.hh"
54 
55 #include "G4VGraphicsScene.hh"
56 #include "G4VisExtent.hh"
57 
58 using namespace CLHEP;
59 
61 //
62 // Constructor
63 
64 G4Orb::G4Orb( const G4String& pName, G4double pRmax )
65  : G4CSGSolid(pName), fRmax(pRmax)
66 {
67  Initialize();
68 }
69 
71 //
72 // Fake default constructor - sets only member data and allocates memory
73 // for usage restricted to object persistency
74 
75 G4Orb::G4Orb( __void__& a )
76  : G4CSGSolid(a), fRmax(0.), halfRmaxTol(0.),
77  sqrRmaxPlusTol(0.), sqrRmaxMinusTol(0.)
78 {
79 }
80 
82 //
83 // Destructor
84 
86 {
87 }
88 
90 //
91 // Copy constructor
92 
93 G4Orb::G4Orb(const G4Orb& rhs)
94  : G4CSGSolid(rhs), fRmax(rhs.fRmax), halfRmaxTol(rhs.halfRmaxTol),
95  sqrRmaxPlusTol(rhs.sqrRmaxPlusTol), sqrRmaxMinusTol(rhs.sqrRmaxMinusTol)
96 {
97 }
98 
100 //
101 // Assignment operator
102 
104 {
105  // Check assignment to self
106  //
107  if (this == &rhs) { return *this; }
108 
109  // Copy base class data
110  //
112 
113  // Copy data
114  //
115  fRmax = rhs.fRmax;
116  halfRmaxTol = rhs.halfRmaxTol;
119 
120  return *this;
121 }
122 
124 //
125 // Check radius and initialize dada members
126 
128 {
129  const G4double fEpsilon = 2.e-11; // relative tolerance of fRmax
130 
131  // Check radius
132  //
133  if ( fRmax < 10*kCarTolerance )
134  {
135  G4Exception("G4Orb::Initialize()", "GeomSolids0002", FatalException,
136  "Invalid radius < 10*kCarTolerance.");
137  }
138  halfRmaxTol = 0.5 * std::max(kCarTolerance, fEpsilon*fRmax);
139  G4double rmaxPlusTol = fRmax + halfRmaxTol;
140  G4double rmaxMinusTol = fRmax - halfRmaxTol;
141  sqrRmaxPlusTol = rmaxPlusTol*rmaxPlusTol;
142  sqrRmaxMinusTol = rmaxMinusTol*rmaxMinusTol;
143 }
144 
146 //
147 // Dispatch to parameterisation for replication mechanism dimension
148 // computation & modification
149 
151  const G4int n,
152  const G4VPhysicalVolume* pRep )
153 {
154  p->ComputeDimensions(*this,n,pRep);
155 }
156 
158 //
159 // Get bounding box
160 
162 {
164  pMin.set(-radius,-radius,-radius);
165  pMax.set( radius, radius, radius);
166 
167  // Check correctness of the bounding box
168  //
169  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
170  {
171  std::ostringstream message;
172  message << "Bad bounding box (min >= max) for solid: "
173  << GetName() << " !"
174  << "\npMin = " << pMin
175  << "\npMax = " << pMax;
176  G4Exception("G4Orb::BoundingLimits()", "GeomMgt0001",
177  JustWarning, message);
178  DumpInfo();
179  }
180 }
181 
183 //
184 // Calculate extent under transform and specified limit
185 
187  const G4VoxelLimits& pVoxelLimit,
188  const G4AffineTransform& pTransform,
189  G4double& pMin, G4double& pMax) const
190 {
191  G4ThreeVector bmin, bmax;
192  G4bool exist;
193 
194  // Get bounding box
195  BoundingLimits(bmin,bmax);
196 
197  // Check bounding box
198  G4BoundingEnvelope bbox(bmin,bmax);
199 #ifdef G4BBOX_EXTENT
200  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
201 #endif
202  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
203  {
204  return exist = (pMin < pMax) ? true : false;
205  }
206 
207  // Find bounding envelope and calculate extent
208  //
209  static const G4int NTHETA = 8; // number of steps along Theta
210  static const G4int NPHI = 16; // number of steps along Phi
211  static const G4double sinHalfTheta = std::sin(halfpi/NTHETA);
212  static const G4double cosHalfTheta = std::cos(halfpi/NTHETA);
213  static const G4double sinHalfPhi = std::sin(pi/NPHI);
214  static const G4double cosHalfPhi = std::cos(pi/NPHI);
215  static const G4double sinStepTheta = 2.*sinHalfTheta*cosHalfTheta;
216  static const G4double cosStepTheta = 1. - 2.*sinHalfTheta*sinHalfTheta;
217  static const G4double sinStepPhi = 2.*sinHalfPhi*cosHalfPhi;
218  static const G4double cosStepPhi = 1. - 2.*sinHalfPhi*sinHalfPhi;
219 
221  G4double rtheta = radius/cosHalfTheta;
222  G4double rphi = rtheta/cosHalfPhi;
223 
224  // set reference circle
225  G4TwoVector xy[NPHI];
226  G4double sinCurPhi = sinHalfPhi;
227  G4double cosCurPhi = cosHalfPhi;
228  for (G4int k=0; k<NPHI; ++k)
229  {
230  xy[k].set(cosCurPhi,sinCurPhi);
231  G4double sinTmpPhi = sinCurPhi;
232  sinCurPhi = sinCurPhi*cosStepPhi + cosCurPhi*sinStepPhi;
233  cosCurPhi = cosCurPhi*cosStepPhi - sinTmpPhi*sinStepPhi;
234  }
235 
236  // set bounding circles
237  G4ThreeVectorList circles[NTHETA];
238  for (G4int i=0; i<NTHETA; ++i) circles[i].resize(NPHI);
239 
240  G4double sinCurTheta = sinHalfTheta;
241  G4double cosCurTheta = cosHalfTheta;
242  for (G4int i=0; i<NTHETA; ++i)
243  {
244  G4double z = rtheta*cosCurTheta;
245  G4double rho = rphi*sinCurTheta;
246  for (G4int k=0; k<NPHI; ++k)
247  {
248  circles[i][k].set(rho*xy[k].x(),rho*xy[k].y(),z);
249  }
250  G4double sinTmpTheta = sinCurTheta;
251  sinCurTheta = sinCurTheta*cosStepTheta + cosCurTheta*sinStepTheta;
252  cosCurTheta = cosCurTheta*cosStepTheta - sinTmpTheta*sinStepTheta;
253  }
254 
255  // set envelope and calculate extent
256  std::vector<const G4ThreeVectorList *> polygons;
257  polygons.resize(NTHETA);
258  for (G4int i=0; i<NTHETA; ++i) polygons[i] = &circles[i];
259 
260  G4BoundingEnvelope benv(bmin,bmax,polygons);
261  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
262  return exist;
263 }
264 
266 //
267 // Return whether point is inside/outside/on surface
268 
270 {
271  G4double rr = p.mag2();
272  if (rr > sqrRmaxPlusTol) return kOutside;
273  return (rr > sqrRmaxMinusTol) ? kSurface : kInside;
274 }
275 
277 //
278 // Return unit normal of surface closest to p
279 
281 {
282  return (1/p.mag())*p;
283 }
284 
286 //
287 // Calculate distance to the surface of the orb from outside
288 // - return kInfinity if no intersection or
289 // intersection distance <= tolerance
290 
292  const G4ThreeVector& v ) const
293 {
294  // Check if point is on the surface and traveling away
295  //
296  G4double rr = p.mag2();
297  G4double pv = p.dot(v);
298  if (rr >= sqrRmaxMinusTol && pv >= 0) return kInfinity;
299 
300  // Find intersection
301  //
302  // Sphere eqn: x^2 + y^2 + z^2 = R^2
303  //
304  // => (px + t*vx)^2 + (py + t*vy)^2 + (pz + t*vz)^2 = R^2
305  // => r^2 + 2t(p.v) + t^2 = R^2
306  // => tmin = -(p.v) - Sqrt((p.v)^2 - (r^2 - R^2))
307  //
308  G4double D = pv*pv - rr + fRmax*fRmax;
309  if (D < 0) return kInfinity; // no intersection
310 
311  G4double sqrtD = std::sqrt(D);
312  G4double tmin = -pv - sqrtD;
313 
314  // Avoid rounding errors due to precision issues seen on 64 bits systems.
315  // Split long distances and recompute
316  //
317  G4double Dmax = 32*fRmax;
318  if (tmin > Dmax)
319  {
320  G4double tadd = tmin - Dmax + fRmax;
321  tadd += DistanceToIn(p + tadd*v, v);
322  return (tadd >= kInfinity) ? kInfinity : tadd;
323  }
324 
325  if (sqrtD*2 <= halfRmaxTol) return kInfinity; // touch
326  return (tmin < halfRmaxTol) ? 0. : tmin;
327 }
328 
330 //
331 // Calculate shortest distance to the boundary from outside
332 // - Return 0 if point is inside
333 
335 {
336  G4double dist = p.mag() - fRmax;
337  return (dist > 0) ? dist : 0.;
338 }
339 
341 //
342 // Calculate distance to the surface of the orb from inside and
343 // find normal at exit point, if required
344 // - when leaving the surface, return 0
345 
347  const G4ThreeVector& v,
348  const G4bool calcNorm,
349  G4bool *validNorm,
350  G4ThreeVector *n ) const
351 {
352  // Check if point is on the surface and traveling away
353  //
354  G4double rr = p.mag2();
355  G4double pv = p.dot(v);
356  if (rr >= sqrRmaxMinusTol && pv > 0)
357  {
358  if (calcNorm)
359  {
360  *validNorm = true;
361  *n = p*(1./std::sqrt(rr));
362  }
363  return 0.;
364  }
365 
366  // Find intersection
367  //
368  // Sphere eqn: x^2 + y^2 + z^2 = R^2
369  //
370  // => (px + t*vx)^2 + (py + t*vy)^2 + (pz + t*vz)^2 = R^2
371  // => r^2 + 2t(p.v) + t^2 = R^2
372  // => tmax = -(p.v) + Sqrt((p.v)^2 - (r^2 - R^2))
373  //
374  G4double D = pv*pv - rr + fRmax*fRmax;
375  G4double tmax = (D <= 0) ? 0. : std::sqrt(D) - pv;
376  if (tmax < halfRmaxTol) tmax = 0.;
377  if (calcNorm)
378  {
379  *validNorm = true;
380  G4ThreeVector ptmax = p + tmax*v;
381  *n = ptmax*(1./ptmax.mag());
382  }
383  return tmax;
384 }
385 
387 //
388 // Calculate distance (<=actual) to closest surface of shape from inside
389 
391 {
392 #ifdef G4CSGDEBUG
393  if( Inside(p) == kOutside )
394  {
395  std::ostringstream message;
396  G4int oldprc = message.precision(16);
397  message << "Point p is outside (!?) of solid: " << GetName() << "\n";
398  message << "Position:\n";
399  message << " p.x() = " << p.x()/mm << " mm\n";
400  message << " p.y() = " << p.y()/mm << " mm\n";
401  message << " p.z() = " << p.z()/mm << " mm";
402  G4cout.precision(oldprc);
403  G4Exception("G4Trap::DistanceToOut(p)", "GeomSolids1002",
404  JustWarning, message );
405  DumpInfo();
406  }
407 #endif
408  G4double dist = fRmax - p.mag();
409  return (dist > 0) ? dist : 0.;
410 }
411 
413 //
414 // G4EntityType
415 
417 {
418  return G4String("G4Orb");
419 }
420 
422 //
423 // Make a clone of the object
424 
426 {
427  return new G4Orb(*this);
428 }
429 
431 //
432 // Stream object contents to an output stream
433 
434 std::ostream& G4Orb::StreamInfo( std::ostream& os ) const
435 {
436  G4int oldprc = os.precision(16);
437  os << "-----------------------------------------------------------\n"
438  << " *** Dump for solid - " << GetName() << " ***\n"
439  << " ===================================================\n"
440  << " Solid type: G4Orb\n"
441  << " Parameters: \n"
442  << " outer radius: " << fRmax/mm << " mm \n"
443  << "-----------------------------------------------------------\n";
444  os.precision(oldprc);
445  return os;
446 }
447 
449 //
450 // GetPointOnSurface
451 
453 {
454  return fRmax * G4RandomDirection();
455 }
456 
458 //
459 // Methods for visualisation
460 
462 {
463  scene.AddSolid (*this);
464 }
465 
467 {
468  return G4VisExtent (-fRmax, fRmax, -fRmax, fRmax, -fRmax, fRmax);
469 }
470 
472 {
473  return new G4PolyhedronSphere (0., fRmax, 0., 2*pi, 0., pi);
474 }
475 
476 #endif
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
virtual void AddSolid(const G4Box &)=0
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetRadius() const
EInside Inside(const G4ThreeVector &p) const
Definition: G4Orb.cc:269
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Orb.cc:434
static constexpr double mm
Definition: G4SIunits.hh:115
G4double halfRmaxTol
Definition: G4Orb.hh:142
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition: G4CSGSolid.cc:92
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
Float_t y
Definition: compare.C:6
double D(double temp)
const char * p
Definition: xmltok.h:285
Double_t z
G4double sqrRmaxMinusTol
Definition: G4Orb.hh:143
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
G4ThreeVector GetPointOnSurface() const
Definition: G4Orb.cc:452
G4Polyhedron * CreatePolyhedron() const
Definition: G4Orb.cc:471
double z() const
double dot(const Hep3Vector &) const
void Initialize()
Definition: G4Orb.cc:127
G4ThreeVector G4RandomDirection()
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=G4bool(false), G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Orb.cc:346
G4GeometryType GetEntityType() const
Definition: G4Orb.cc:416
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4VSolid * Clone() const
Definition: G4Orb.cc:425
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Orb.cc:280
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Orb.cc:150
G4double kCarTolerance
Definition: G4VSolid.hh:307
static constexpr double halfpi
Definition: G4SIunits.hh:77
double mag2() const
Definition: G4Orb.hh:62
Double_t radius
G4Orb & operator=(const G4Orb &rhs)
Definition: G4Orb.cc:103
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
double mag() const
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
void set(double x, double y)
G4String GetName() const
void DumpInfo() const
~G4Orb()
Definition: G4Orb.cc:85
G4GLOB_DLL std::ostream G4cout
double x() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Orb.cc:186
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Orb.cc:161
G4double fRmax
Definition: G4Orb.hh:141
Char_t n[5]
static constexpr double pi
Definition: G4SIunits.hh:75
double y() const
G4double sqrRmaxPlusTol
Definition: G4Orb.hh:143
G4VisExtent GetExtent() const
Definition: G4Orb.cc:466
G4Orb(const G4String &pName, G4double pRmax)
Definition: G4Orb.cc:64
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Orb.cc:291
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Orb.cc:461