Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4USphere.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 //
27 // $Id:$
28 //
29 //
30 // Implementation for G4USphere wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4Sphere.hh"
34 #include "G4USphere.hh"
35 
36 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
37 
38 #include "G4GeomTools.hh"
39 #include "G4AffineTransform.hh"
40 #include "G4VPVParameterisation.hh"
41 #include "G4BoundingEnvelope.hh"
42 
43 using namespace CLHEP;
44 
46 //
47 // constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
48 // - note if pDPhi>2PI then reset to 2PI
49 
50 G4USphere::G4USphere( const G4String& pName,
51  G4double pRmin, G4double pRmax,
52  G4double pSPhi, G4double pDPhi,
53  G4double pSTheta, G4double pDTheta )
54  : Base_t(pName, pRmin, pRmax, pSPhi, pDPhi, pSTheta, pDTheta)
55 {
56 }
57 
59 //
60 // Fake default constructor - sets only member data and allocates memory
61 // for usage restricted to object persistency.
62 //
63 G4USphere::G4USphere( __void__& a )
64  : Base_t(a)
65 {
66 }
67 
69 //
70 // Destructor
71 
72 G4USphere::~G4USphere()
73 {
74 }
75 
77 //
78 // Copy constructor
79 
80 G4USphere::G4USphere(const G4USphere& rhs)
81  : Base_t(rhs)
82 {
83 }
84 
86 //
87 // Assignment operator
88 
89 G4USphere& G4USphere::operator = (const G4USphere& rhs)
90 {
91  // Check assignment to self
92  //
93  if (this == &rhs) { return *this; }
94 
95  // Copy base class data
96  //
97  Base_t::operator=(rhs);
98 
99  return *this;
100 }
101 
103 //
104 // Accessors & modifiers
105 
106 G4double G4USphere::GetInnerRadius() const
107 {
108  return Base_t::GetInnerRadius();
109 }
110 G4double G4USphere::GetOuterRadius() const
111 {
112  return Base_t::GetOuterRadius();
113 }
114 G4double G4USphere::GetStartPhiAngle() const
115 {
116  return Base_t::GetStartPhiAngle();
117 }
118 G4double G4USphere::GetDeltaPhiAngle() const
119 {
120  return Base_t::GetDeltaPhiAngle();
121 }
122 G4double G4USphere::GetStartThetaAngle() const
123 {
124  return Base_t::GetStartThetaAngle();
125 }
126 G4double G4USphere::GetDeltaThetaAngle() const
127 {
128  return Base_t::GetDeltaThetaAngle();
129 }
130 G4double G4USphere::GetSinStartPhi() const
131 {
132  return Base_t::GetSinSPhi();
133 }
134 G4double G4USphere::GetCosStartPhi() const
135 {
136  return Base_t::GetCosSPhi();
137 }
138 G4double G4USphere::GetSinEndPhi() const
139 {
140  return Base_t::GetSinEPhi();
141 }
142 G4double G4USphere::GetCosEndPhi() const
143 {
144  return Base_t::GetCosEPhi();
145 }
146 G4double G4USphere::GetSinStartTheta() const
147 {
148  return Base_t::GetSinSTheta();
149 }
150 G4double G4USphere::GetCosStartTheta() const
151 {
152  return Base_t::GetCosSTheta();
153 }
154 G4double G4USphere::GetSinEndTheta() const
155 {
156  return Base_t::GetSinETheta();
157 }
158 G4double G4USphere::GetCosEndTheta() const
159 {
160  return Base_t::GetCosETheta();
161 }
162 
163 void G4USphere::SetInnerRadius(G4double newRMin)
164 {
165  Base_t::SetInnerRadius(newRMin);
166  fRebuildPolyhedron = true;
167 }
168 void G4USphere::SetOuterRadius(G4double newRmax)
169 {
170  Base_t::SetOuterRadius(newRmax);
171  fRebuildPolyhedron = true;
172 }
173 void G4USphere::SetStartPhiAngle(G4double newSphi, G4bool trig)
174 {
175  Base_t::SetStartPhiAngle(newSphi, trig);
176  fRebuildPolyhedron = true;
177 }
178 void G4USphere::SetDeltaPhiAngle(G4double newDphi)
179 {
180  Base_t::SetDeltaPhiAngle(newDphi);
181  fRebuildPolyhedron = true;
182 }
183 void G4USphere::SetStartThetaAngle(G4double newSTheta)
184 {
185  Base_t::SetStartThetaAngle(newSTheta);
186  fRebuildPolyhedron = true;
187 }
188 void G4USphere::SetDeltaThetaAngle(G4double newDTheta)
189 {
190  Base_t::SetDeltaThetaAngle(newDTheta);
191  fRebuildPolyhedron = true;
192 }
193 
195 //
196 // Dispatch to parameterisation for replication mechanism dimension
197 // computation & modification.
198 
199 void G4USphere::ComputeDimensions( G4VPVParameterisation* p,
200  const G4int n,
201  const G4VPhysicalVolume* pRep)
202 {
203  p->ComputeDimensions(*(G4Sphere*)this,n,pRep);
204 }
205 
207 //
208 // Make a clone of the object
209 
210 G4VSolid* G4USphere::Clone() const
211 {
212  return new G4USphere(*this);
213 }
214 
216 //
217 // Get bounding box
218 
219 void G4USphere::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
220 {
221  static G4bool checkBBox = true;
222 
223  G4double rmin = GetInnerRadius();
224  G4double rmax = GetOuterRadius();
225 
226  // Find bounding box
227  //
228  if (GetDeltaThetaAngle() >= pi && GetDeltaPhiAngle() >= twopi)
229  {
230  pMin.set(-rmax,-rmax,-rmax);
231  pMax.set( rmax, rmax, rmax);
232  }
233  else
234  {
235  G4double sinStart = GetSinStartTheta();
236  G4double cosStart = GetCosStartTheta();
237  G4double sinEnd = GetSinEndTheta();
238  G4double cosEnd = GetCosEndTheta();
239 
240  G4double stheta = GetStartThetaAngle();
241  G4double etheta = stheta + GetDeltaThetaAngle();
242  G4double rhomin = rmin*std::min(sinStart,sinEnd);
243  G4double rhomax = rmax;
244  if (stheta > halfpi) rhomax = rmax*sinStart;
245  if (etheta < halfpi) rhomax = rmax*sinEnd;
246 
247  G4TwoVector xymin,xymax;
248  G4GeomTools::DiskExtent(rhomin,rhomax,
249  GetSinStartPhi(),GetCosStartPhi(),
250  GetSinEndPhi(),GetCosEndPhi(),
251  xymin,xymax);
252 
253  G4double zmin = std::min(rmin*cosEnd,rmax*cosEnd);
254  G4double zmax = std::max(rmin*cosStart,rmax*cosStart);
255  pMin.set(xymin.x(),xymin.y(),zmin);
256  pMax.set(xymax.x(),xymax.y(),zmax);
257  }
258 
259  // Check correctness of the bounding box
260  //
261  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
262  {
263  std::ostringstream message;
264  message << "Bad bounding box (min >= max) for solid: "
265  << GetName() << " !"
266  << "\npMin = " << pMin
267  << "\npMax = " << pMax;
268  G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
269  JustWarning, message);
270  StreamInfo(G4cout);
271  }
272 
273  // Check consistency of bounding boxes
274  //
275  if (checkBBox)
276  {
277  U3Vector vmin, vmax;
278  Extent(vmin,vmax);
279  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
280  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
281  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
282  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
283  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
284  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
285  {
286  std::ostringstream message;
287  message << "Inconsistency in bounding boxes for solid: "
288  << GetName() << " !"
289  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
290  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
291  G4Exception("G4USphere::BoundingLimits()", "GeomMgt0001",
292  JustWarning, message);
293  checkBBox = false;
294  }
295  }
296 }
297 
299 //
300 // Calculate extent under transform and specified limit
301 
302 G4bool G4USphere::CalculateExtent(const EAxis pAxis,
303  const G4VoxelLimits& pVoxelLimit,
304  const G4AffineTransform& pTransform,
305  G4double& pMin, G4double& pMax) const
306 {
307  G4ThreeVector bmin, bmax;
308 
309  // Get bounding box
310  BoundingLimits(bmin,bmax);
311 
312  // Find extent
313  G4BoundingEnvelope bbox(bmin,bmax);
314  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
315 }
316 
318 //
319 // Create polyhedron for visualization
320 
321 G4Polyhedron* G4USphere::CreatePolyhedron() const
322 {
323  return new G4PolyhedronSphere(GetInnerRadius(),
324  GetOuterRadius(),
325  GetStartPhiAngle(),
326  GetDeltaPhiAngle(),
327  GetStartThetaAngle(),
328  GetDeltaThetaAngle());
329 }
330 
331 #endif // G4GEOM_USE_USOLIDS
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 ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
const char * p
Definition: xmltok.h:285
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
double z() const
const G4double kCarTolerance
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:396
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
double x() const
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double halfpi
Definition: G4SIunits.hh:77
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
double y() const
G4GLOB_DLL std::ostream G4cout
double x() const
Char_t n[5]
static constexpr double pi
Definition: G4SIunits.hh:75
double y() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments