Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UOrb.cc
이 파일의 문서화 페이지로 가기
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  : Base_t(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  : Base_t(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  : Base_t(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  Base_t::operator=(rhs);
97 
98  return *this;
99 }
100 
102 //
103 // Accessors & modifiers
104 
105 G4double G4UOrb::GetRadius() const
106 {
107  return Base_t::GetRadius();
108 }
109 
110 void G4UOrb::SetRadius(G4double newRmax)
111 {
112  Base_t::SetRadius(newRmax);
113  fRebuildPolyhedron = true;
114 }
115 
116 G4double G4UOrb::GetRadialTolerance() const
117 {
118  return Base_t::GetRadialTolerance();
119 }
120 
122 //
123 // Dispatch to parameterisation for replication mechanism dimension
124 // computation & modification.
125 
126 void G4UOrb::ComputeDimensions( G4VPVParameterisation* p,
127  const G4int n,
128  const G4VPhysicalVolume* pRep )
129 {
130  p->ComputeDimensions(*(G4Orb*)this,n,pRep);
131 }
132 
134 //
135 // Make a clone of the object
136 
137 G4VSolid* G4UOrb::Clone() const
138 {
139  return new G4UOrb(*this);
140 }
141 
143 //
144 // Get bounding box
145 
146 void G4UOrb::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
147 {
148  G4double radius = GetRadius();
149  pMin.set(-radius,-radius,-radius);
150  pMax.set( radius, radius, radius);
151 
152  // Check correctness of the bounding box
153  //
154  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
155  {
156  std::ostringstream message;
157  message << "Bad bounding box (min >= max) for solid: "
158  << GetName() << " !"
159  << "\npMin = " << pMin
160  << "\npMax = " << pMax;
161  G4Exception("G4UOrb::BoundingLimits()", "GeomMgt0001",
162  JustWarning, message);
163  StreamInfo(G4cout);
164  }
165 }
166 
168 //
169 // Calculate extent under transform and specified limit
170 
171 G4bool
172 G4UOrb::CalculateExtent(const EAxis pAxis,
173  const G4VoxelLimits& pVoxelLimit,
174  const G4AffineTransform& pTransform,
175  G4double& pMin, G4double& pMax) const
176 {
177  G4ThreeVector bmin, bmax;
178  G4bool exist;
179 
180  // Get bounding box
181  BoundingLimits(bmin,bmax);
182 
183  // Check bounding box
184  G4BoundingEnvelope bbox(bmin,bmax);
185 #ifdef G4BBOX_EXTENT
186  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
187 #endif
188  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
189  {
190  return exist = (pMin < pMax) ? true : false;
191  }
192 
193  // Find bounding envelope and calculate extent
194  //
195  static const G4int NTHETA = 8; // number of steps along Theta
196  static const G4int NPHI = 16; // number of steps along Phi
197  static const G4double sinHalfTheta = std::sin(halfpi/NTHETA);
198  static const G4double cosHalfTheta = std::cos(halfpi/NTHETA);
199  static const G4double sinHalfPhi = std::sin(pi/NPHI);
200  static const G4double cosHalfPhi = std::cos(pi/NPHI);
201  static const G4double sinStepTheta = 2.*sinHalfTheta*cosHalfTheta;
202  static const G4double cosStepTheta = 1. - 2.*sinHalfTheta*sinHalfTheta;
203  static const G4double sinStepPhi = 2.*sinHalfPhi*cosHalfPhi;
204  static const G4double cosStepPhi = 1. - 2.*sinHalfPhi*sinHalfPhi;
205 
206  G4double radius = GetRadius();
207  G4double rtheta = radius/cosHalfTheta;
208  G4double rphi = rtheta/cosHalfPhi;
209 
210  // set reference circle
211  G4TwoVector xy[NPHI];
212  G4double sinCurPhi = sinHalfPhi;
213  G4double cosCurPhi = cosHalfPhi;
214  for (G4int k=0; k<NPHI; ++k)
215  {
216  xy[k].set(cosCurPhi,sinCurPhi);
217  G4double sinTmpPhi = sinCurPhi;
218  sinCurPhi = sinCurPhi*cosStepPhi + cosCurPhi*sinStepPhi;
219  cosCurPhi = cosCurPhi*cosStepPhi - sinTmpPhi*sinStepPhi;
220  }
221 
222  // set bounding circles
223  G4ThreeVectorList circles[NTHETA];
224  for (G4int i=0; i<NTHETA; ++i) circles[i].resize(NPHI);
225 
226  G4double sinCurTheta = sinHalfTheta;
227  G4double cosCurTheta = cosHalfTheta;
228  for (G4int i=0; i<NTHETA; ++i)
229  {
230  G4double z = rtheta*cosCurTheta;
231  G4double rho = rphi*sinCurTheta;
232  for (G4int k=0; k<NPHI; ++k)
233  {
234  circles[i][k].set(rho*xy[k].x(),rho*xy[k].y(),z);
235  }
236  G4double sinTmpTheta = sinCurTheta;
237  sinCurTheta = sinCurTheta*cosStepTheta + cosCurTheta*sinStepTheta;
238  cosCurTheta = cosCurTheta*cosStepTheta - sinTmpTheta*sinStepTheta;
239  }
240 
241  // set envelope and calculate extent
242  std::vector<const G4ThreeVectorList *> polygons;
243  polygons.resize(NTHETA);
244  for (G4int i=0; i<NTHETA; ++i) polygons[i] = &circles[i];
245 
246  G4BoundingEnvelope benv(bmin,bmax,polygons);
247  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
248  return exist;
249 }
250 
252 //
253 // Create polyhedron for visualization
254 
255 G4Polyhedron* G4UOrb::CreatePolyhedron() const
256 {
257  return new G4PolyhedronSphere(0., GetRadius(), 0., twopi, 0., pi);
258 }
259 
260 #endif // G4GEOM_USE_USOLIDS
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
Double_t z
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
double z() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double halfpi
Definition: G4SIunits.hh:77
Definition: G4Orb.hh:62
Double_t radius
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
EAxis
Definition: geomdefs.hh:54
void set(double x, double y)
G4GLOB_DLL std::ostream G4cout
double x() const
Char_t n[5]
static constexpr double pi
Definition: G4SIunits.hh:75
double y() const