Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UTubs.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 //
27 // $Id:$
28 //
29 //
30 // Implementation for G4UTubs wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4Tubs.hh"
34 #include "G4UTubs.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 G4UTubs::G4UTubs( const G4String& pName,
51  G4double pRMin, G4double pRMax,
52  G4double pDz,
53  G4double pSPhi, G4double pDPhi )
54  : G4USolid(pName, new UTubs(pName, pRMin, pRMax, pDz, pSPhi, pDPhi))
55 {
56 }
57 
59 //
60 // Fake default constructor - sets only member data and allocates memory
61 // for usage restricted to object persistency.
62 //
63 G4UTubs::G4UTubs( __void__& a )
64  : G4USolid(a)
65 {
66 }
67 
69 //
70 // Destructor
71 
72 G4UTubs::~G4UTubs()
73 {
74 }
75 
77 //
78 // Copy constructor
79 
80 G4UTubs::G4UTubs(const G4UTubs& rhs)
81  : G4USolid(rhs)
82 {
83 }
84 
86 //
87 // Assignment operator
88 
89 G4UTubs& G4UTubs::operator = (const G4UTubs& rhs)
90 {
91  // Check assignment to self
92  //
93  if (this == &rhs) { return *this; }
94 
95  // Copy base class data
96  //
97  G4USolid::operator=(rhs);
98 
99  return *this;
100 }
101 
103 //
104 // Accessors and modifiers
105 
106 G4double G4UTubs::GetInnerRadius() const
107 {
108  return GetShape()->GetInnerRadius();
109 }
110 G4double G4UTubs::GetOuterRadius() const
111 {
112  return GetShape()->GetOuterRadius();
113 }
114 G4double G4UTubs::GetZHalfLength() const
115 {
116  return GetShape()->GetZHalfLength();
117 }
118 G4double G4UTubs::GetStartPhiAngle() const
119 {
120  return GetShape()->GetStartPhiAngle();
121 }
122 G4double G4UTubs::GetDeltaPhiAngle() const
123 {
124  return GetShape()->GetDeltaPhiAngle();
125 }
126 G4double G4UTubs::GetSinStartPhi() const
127 {
128  G4double phi = GetShape()->GetStartPhiAngle();
129  return std::sin(phi);
130 }
131 G4double G4UTubs::GetCosStartPhi() const
132 {
133  G4double phi = GetShape()->GetStartPhiAngle();
134  return std::cos(phi);
135 }
136 G4double G4UTubs::GetSinEndPhi() const
137 {
138  G4double phi = GetShape()->GetStartPhiAngle() +
139  GetShape()->GetDeltaPhiAngle();
140  return std::sin(phi);
141 }
142 G4double G4UTubs::GetCosEndPhi() const
143 {
144  G4double phi = GetShape()->GetStartPhiAngle() +
145  GetShape()->GetDeltaPhiAngle();
146  return std::cos(phi);
147 }
148 
149 void G4UTubs::SetInnerRadius(G4double newRMin)
150 {
151  GetShape()->SetInnerRadius(newRMin);
152  fRebuildPolyhedron = true;
153 }
154 void G4UTubs::SetOuterRadius(G4double newRMax)
155 {
156  GetShape()->SetOuterRadius(newRMax);
157  fRebuildPolyhedron = true;
158 }
159 void G4UTubs::SetZHalfLength(G4double newDz)
160 {
161  GetShape()->SetZHalfLength(newDz);
162  fRebuildPolyhedron = true;
163 }
164 void G4UTubs::SetStartPhiAngle(G4double newSPhi, G4bool trig)
165 {
166  GetShape()->SetStartPhiAngle(newSPhi, trig);
167  fRebuildPolyhedron = true;
168 }
169 void G4UTubs::SetDeltaPhiAngle(G4double newDPhi)
170 {
171  GetShape()->SetDeltaPhiAngle(newDPhi);
172  fRebuildPolyhedron = true;
173 }
174 
176 //
177 // Dispatch to parameterisation for replication mechanism dimension
178 // computation & modification.
179 
180 void G4UTubs::ComputeDimensions( G4VPVParameterisation* p,
181  const G4int n,
182  const G4VPhysicalVolume* pRep )
183 {
184  p->ComputeDimensions(*(G4Tubs*)this,n,pRep) ;
185 }
186 
188 //
189 // Make a clone of the object
190 
191 G4VSolid* G4UTubs::Clone() const
192 {
193  return new G4UTubs(*this);
194 }
195 
197 //
198 // Get bounding box
199 
200 void G4UTubs::Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const
201 {
202  static G4bool checkBBox = true;
203 
204  G4double rmin = GetInnerRadius();
205  G4double rmax = GetOuterRadius();
206  G4double dz = GetZHalfLength();
207 
208  // Find bounding box
209  //
210  if (GetDeltaPhiAngle() < twopi)
211  {
212  G4TwoVector vmin,vmax;
213  G4GeomTools::DiskExtent(rmin,rmax,
214  GetSinStartPhi(),GetCosStartPhi(),
215  GetSinEndPhi(),GetCosEndPhi(),
216  vmin,vmax);
217  pMin.set(vmin.x(),vmin.y(),-dz);
218  pMax.set(vmax.x(),vmax.y(), dz);
219  }
220  else
221  {
222  pMin.set(-rmax,-rmax,-dz);
223  pMax.set( rmax, rmax, dz);
224  }
225 
226  // Check correctness of the bounding box
227  //
228  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
229  {
230  std::ostringstream message;
231  message << "Bad bounding box (min >= max) for solid: "
232  << GetName() << " !"
233  << "\npMin = " << pMin
234  << "\npMax = " << pMax;
235  G4Exception("G4UTubs::Extent()", "GeomMgt0001", JustWarning, message);
236  StreamInfo(G4cout);
237  }
238 
239  // Check consistency of bounding boxes
240  //
241  if (checkBBox)
242  {
243  UVector3 vmin, vmax;
244  GetShape()->Extent(vmin,vmax);
245  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
246  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
247  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
248  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
249  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
250  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
251  {
252  std::ostringstream message;
253  message << "Inconsistency in bounding boxes for solid: "
254  << GetName() << " !"
255  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
256  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
257  G4Exception("G4UTubs::Extent()", "GeomMgt0001", JustWarning, message);
258  checkBBox = false;
259  }
260  }
261 }
262 
264 //
265 // Calculate extent under transform and specified limit
266 
267 G4bool
268 G4UTubs::CalculateExtent(const EAxis pAxis,
269  const G4VoxelLimits& pVoxelLimit,
270  const G4AffineTransform& pTransform,
271  G4double& pMin, G4double& pMax) const
272 {
273  G4ThreeVector bmin, bmax;
274  G4bool exist;
275 
276  // Get bounding box
277  Extent(bmin,bmax);
278 
279  // Check bounding box
280  G4BoundingEnvelope bbox(bmin,bmax);
281 #ifdef G4BBOX_EXTENT
282  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
283 #endif
284  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
285  {
286  return exist = (pMin < pMax) ? true : false;
287  }
288 
289  // Get parameters of the solid
290  G4double rmin = GetInnerRadius();
291  G4double rmax = GetOuterRadius();
292  G4double dz = GetZHalfLength();
293  G4double dphi = GetDeltaPhiAngle();
294 
295  // Find bounding envelope and calculate extent
296  //
297  const G4int NSTEPS = 24; // number of steps for whole circle
298  G4double astep = twopi/NSTEPS; // max angle for one step
299  G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
300  G4double ang = dphi/ksteps;
301 
302  G4double sinHalf = std::sin(0.5*ang);
303  G4double cosHalf = std::cos(0.5*ang);
304  G4double sinStep = 2.*sinHalf*cosHalf;
305  G4double cosStep = 1. - 2.*sinHalf*sinHalf;
306  G4double rext = rmax/cosHalf;
307 
308  // bounding envelope for full cylinder consists of two polygons,
309  // in other cases it is a sequence of quadrilaterals
310  if (rmin == 0 && dphi == twopi)
311  {
312  G4double sinCur = sinHalf;
313  G4double cosCur = cosHalf;
314 
315  G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
316  for (G4int k=0; k<NSTEPS; ++k)
317  {
318  baseA[k].set(rext*cosCur,rext*sinCur,-dz);
319  baseB[k].set(rext*cosCur,rext*sinCur, dz);
320 
321  G4double sinTmp = sinCur;
322  sinCur = sinCur*cosStep + cosCur*sinStep;
323  cosCur = cosCur*cosStep - sinTmp*sinStep;
324  }
325  std::vector<const G4ThreeVectorList *> polygons(2);
326  polygons[0] = &baseA;
327  polygons[1] = &baseB;
328  G4BoundingEnvelope benv(bmin,bmax,polygons);
329  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
330  }
331  else
332  {
333  G4double sinStart = GetSinStartPhi();
334  G4double cosStart = GetCosStartPhi();
335  G4double sinEnd = GetSinEndPhi();
336  G4double cosEnd = GetCosEndPhi();
337  G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
338  G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
339 
340  // set quadrilaterals
341  G4ThreeVectorList pols[NSTEPS+2];
342  for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(4);
343  pols[0][0].set(rmin*cosStart,rmin*sinStart, dz);
344  pols[0][1].set(rmin*cosStart,rmin*sinStart,-dz);
345  pols[0][2].set(rmax*cosStart,rmax*sinStart,-dz);
346  pols[0][3].set(rmax*cosStart,rmax*sinStart, dz);
347  for (G4int k=1; k<ksteps+1; ++k)
348  {
349  pols[k][0].set(rmin*cosCur,rmin*sinCur, dz);
350  pols[k][1].set(rmin*cosCur,rmin*sinCur,-dz);
351  pols[k][2].set(rext*cosCur,rext*sinCur,-dz);
352  pols[k][3].set(rext*cosCur,rext*sinCur, dz);
353 
354  G4double sinTmp = sinCur;
355  sinCur = sinCur*cosStep + cosCur*sinStep;
356  cosCur = cosCur*cosStep - sinTmp*sinStep;
357  }
358  pols[ksteps+1][0].set(rmin*cosEnd,rmin*sinEnd, dz);
359  pols[ksteps+1][1].set(rmin*cosEnd,rmin*sinEnd,-dz);
360  pols[ksteps+1][2].set(rmax*cosEnd,rmax*sinEnd,-dz);
361  pols[ksteps+1][3].set(rmax*cosEnd,rmax*sinEnd, dz);
362 
363  // set envelope and calculate extent
364  std::vector<const G4ThreeVectorList *> polygons;
365  polygons.resize(ksteps+2);
366  for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
367  G4BoundingEnvelope benv(bmin,bmax,polygons);
368  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
369  }
370  return exist;
371 }
372 
374 //
375 // Create polyhedron for visualization
376 //
377 G4Polyhedron* G4UTubs::CreatePolyhedron() const
378 {
379  return new G4PolyhedronTubs(GetInnerRadius(),
380  GetOuterRadius(),
381  GetZHalfLength(),
382  GetStartPhiAngle(),
383  GetDeltaPhiAngle());
384 }
385 
386 #endif // G4GEOM_USE_USOLIDS
void set(double x, double y, double z)
double y() const
double x() const
double x() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
const char * p
Definition: xmltok.h:285
Definition: G4Tubs.hh:85
int G4int
Definition: G4Types.hh:78
double z() const
static constexpr double twopi
Definition: G4SIunits.hh:76
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4double kCarTolerance
const G4int n
std::vector< G4ThreeVector > G4ThreeVectorList
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
double G4double
Definition: G4Types.hh:76
static constexpr double deg
Definition: G4SIunits.hh:152
tuple astep
Definition: test1.py:13
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:378