Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4UOrb.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 // $Id:$
27 //
28 //
29 // Implementation for G4UOrb wrapper class
30 // --------------------------------------------------------------------
31 
32 #include "G4Orb.hh"
33 #include "G4UOrb.hh"
34 
35 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
36 
37 #include "G4TwoVector.hh"
38 #include "G4AffineTransform.hh"
39 #include "G4GeometryTolerance.hh"
40 #include "G4BoundingEnvelope.hh"
41 
42 #include "G4VPVParameterisation.hh"
43 #include "G4PhysicalConstants.hh"
44 
45 using namespace CLHEP;
46 
48 //
49 // constructor - check positive radius
50 //
51 
52 G4UOrb::G4UOrb( const G4String& pName, G4double pRmax )
53  : G4USolid(pName, new UOrb(pName, pRmax))
54 {
55 }
56 
58 //
59 // Fake default constructor - sets only member data and allocates memory
60 // for usage restricted to object persistency.
61 //
62 G4UOrb::G4UOrb( __void__& a )
63  : G4USolid(a)
64 {
65 }
66 
68 //
69 // Destructor
70 
71 G4UOrb::~G4UOrb()
72 {
73 }
74 
76 //
77 // Copy constructor
78 
79 G4UOrb::G4UOrb(const G4UOrb& rhs)
80  : G4USolid(rhs)
81 {
82 }
83 
85 //
86 // Assignment operator
87 
88 G4UOrb& G4UOrb::operator = (const G4UOrb& rhs)
89 {
90  // Check assignment to self
91  //
92  if (this == &rhs) { return *this; }
93 
94  // Copy base class data
95  //
96  G4USolid::operator=(rhs);
97 
98  return *this;
99 }
100 
102 //
103 // Accessors & modifiers
104 
105 G4double G4UOrb::GetRadius() const
106 {
107  return GetShape()->GetRadius();
108 }
109 
110 void G4UOrb::SetRadius(G4double newRmax)
111 {
112  GetShape()->SetRadius(newRmax);
113  fRebuildPolyhedron = true;
114 }
115 
117 //
118 // Dispatch to parameterisation for replication mechanism dimension
119 // computation & modification.
120 
121 void G4UOrb::ComputeDimensions( G4VPVParameterisation* p,
122  const G4int n,
123  const G4VPhysicalVolume* pRep )
124 {
125  p->ComputeDimensions(*(G4Orb*)this,n,pRep);
126 }
127 
129 //
130 // Make a clone of the object
131 
132 G4VSolid* G4UOrb::Clone() const
133 {
134  return new G4UOrb(*this);
135 }
136 
138 //
139 // Get bounding box
140 
141 void G4UOrb::Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const
142 {
143  G4double radius = GetRadius();
144  pMin.set(-radius,-radius,-radius);
145  pMax.set( radius, radius, radius);
146 
147  // Check correctness of the bounding box
148  //
149  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
150  {
151  std::ostringstream message;
152  message << "Bad bounding box (min >= max) for solid: "
153  << GetName() << " !"
154  << "\npMin = " << pMin
155  << "\npMax = " << pMax;
156  G4Exception("G4UOrb::Extent()", "GeomMgt0001", JustWarning, message);
157  StreamInfo(G4cout);
158  }
159 }
160 
162 //
163 // Calculate extent under transform and specified limit
164 
165 G4bool
166 G4UOrb::CalculateExtent(const EAxis pAxis,
167  const G4VoxelLimits& pVoxelLimit,
168  const G4AffineTransform& pTransform,
169  G4double& pMin, G4double& pMax) const
170 {
171  G4ThreeVector bmin, bmax;
172  G4bool exist;
173 
174  // Get bounding box
175  Extent(bmin,bmax);
176 
177  // Check bounding box
178  G4BoundingEnvelope bbox(bmin,bmax);
179 #ifdef G4BBOX_EXTENT
180  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
181 #endif
182  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
183  {
184  return exist = (pMin < pMax) ? true : false;
185  }
186 
187  // Find bounding envelope and calculate extent
188  //
189  static const G4int NTHETA = 8; // number of steps along Theta
190  static const G4int NPHI = 16; // number of steps along Phi
191  static const G4double sinHalfTheta = std::sin(halfpi/NTHETA);
192  static const G4double cosHalfTheta = std::cos(halfpi/NTHETA);
193  static const G4double sinHalfPhi = std::sin(pi/NPHI);
194  static const G4double cosHalfPhi = std::cos(pi/NPHI);
195  static const G4double sinStepTheta = 2.*sinHalfTheta*cosHalfTheta;
196  static const G4double cosStepTheta = 1. - 2.*sinHalfTheta*sinHalfTheta;
197  static const G4double sinStepPhi = 2.*sinHalfPhi*cosHalfPhi;
198  static const G4double cosStepPhi = 1. - 2.*sinHalfPhi*sinHalfPhi;
199 
200  G4double radius = GetRadius();
201  G4double rtheta = radius/cosHalfTheta;
202  G4double rphi = rtheta/cosHalfPhi;
203 
204  // set reference circle
205  G4TwoVector xy[NPHI];
206  G4double sinCurPhi = sinHalfPhi;
207  G4double cosCurPhi = cosHalfPhi;
208  for (G4int k=0; k<NPHI; ++k)
209  {
210  xy[k].set(cosCurPhi,sinCurPhi);
211  G4double sinTmpPhi = sinCurPhi;
212  sinCurPhi = sinCurPhi*cosStepPhi + cosCurPhi*sinStepPhi;
213  cosCurPhi = cosCurPhi*cosStepPhi - sinTmpPhi*sinStepPhi;
214  }
215 
216  // set bounding circles
217  G4ThreeVectorList circles[NTHETA];
218  for (G4int i=0; i<NTHETA; ++i) circles[i].resize(NPHI);
219 
220  G4double sinCurTheta = sinHalfTheta;
221  G4double cosCurTheta = cosHalfTheta;
222  for (G4int i=0; i<NTHETA; ++i)
223  {
224  G4double z = rtheta*cosCurTheta;
225  G4double rho = rphi*sinCurTheta;
226  for (G4int k=0; k<NPHI; ++k)
227  {
228  circles[i][k].set(rho*xy[k].x(),rho*xy[k].y(),z);
229  }
230  G4double sinTmpTheta = sinCurTheta;
231  sinCurTheta = sinCurTheta*cosStepTheta + cosCurTheta*sinStepTheta;
232  cosCurTheta = cosCurTheta*cosStepTheta - sinTmpTheta*sinStepTheta;
233  }
234 
235  // set envelope and calculate extent
236  std::vector<const G4ThreeVectorList *> polygons;
237  polygons.resize(NTHETA);
238  for (G4int i=0; i<NTHETA; ++i) polygons[i] = &circles[i];
239 
240  G4BoundingEnvelope benv(bmin,bmax,polygons);
241  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
242  return exist;
243 }
244 
246 //
247 // Create polyhedron for visualization
248 
249 G4Polyhedron* G4UOrb::CreatePolyhedron() const
250 {
251  return new G4PolyhedronSphere(0., GetRadius(), 0., twopi, 0., pi);
252 }
253 
254 #endif // G4GEOM_USE_USOLIDS
void set(double x, double y, double z)
double x() const
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
double z() const
static constexpr double twopi
Definition: G4SIunits.hh:76
G4GLOB_DLL std::ostream G4cout
void set(double x, double y)
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4Orb.hh:61
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
EAxis
Definition: geomdefs.hh:54
double y() const
static constexpr double pi
Definition: G4SIunits.hh:75
static constexpr double halfpi
Definition: G4SIunits.hh:77
double G4double
Definition: G4Types.hh:76