Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// Implementation for G4UTubs wrapper class
27//
28// 30.10.13 G.Cosmo, CERN
29// --------------------------------------------------------------------
30
31#include "G4Tubs.hh"
32#include "G4UTubs.hh"
33
34#if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35
36#include "G4GeomTools.hh"
37#include "G4AffineTransform.hh"
39#include "G4BoundingEnvelope.hh"
40
41using namespace CLHEP;
42
43/////////////////////////////////////////////////////////////////////////
44//
45// Constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
46// - note if pdphi>2PI then reset to 2PI
47
48G4UTubs::G4UTubs( const G4String& pName,
49 G4double pRMin, G4double pRMax,
50 G4double pDz,
51 G4double pSPhi, G4double pDPhi )
52 : Base_t(pName, pRMin, pRMax, pDz, pSPhi, pDPhi)
53{
54}
55
56//////////////////////////////////////////////////////////////////////////
57//
58// Copy constructor
59
60G4UTubs::G4UTubs(const G4UTubs& rhs)
61 : Base_t(rhs)
62{
63}
64
65//////////////////////////////////////////////////////////////////////////
66//
67// Assignment operator
68
69G4UTubs& G4UTubs::operator = (const G4UTubs& rhs)
70{
71 // Check assignment to self
72 //
73 if (this == &rhs) { return *this; }
74
75 // Copy base class data
76 //
77 Base_t::operator=(rhs);
78
79 return *this;
80}
81
82/////////////////////////////////////////////////////////////////////////
83//
84// Accessors and modifiers
85
86G4double G4UTubs::GetInnerRadius() const
87{
88 return rmin();
89}
90G4double G4UTubs::GetOuterRadius() const
91{
92 return rmax();
93}
94G4double G4UTubs::GetZHalfLength() const
95{
96 return z();
97}
98G4double G4UTubs::GetStartPhiAngle() const
99{
100 return sphi();
101}
102G4double G4UTubs::GetDeltaPhiAngle() const
103{
104 return dphi();
105}
106G4double G4UTubs::GetSinStartPhi() const
107{
108 return std::sin(GetStartPhiAngle());
109}
110G4double G4UTubs::GetCosStartPhi() const
111{
112 return std::cos(GetStartPhiAngle());
113}
114G4double G4UTubs::GetSinEndPhi() const
115{
116 return std::sin(GetStartPhiAngle()+GetDeltaPhiAngle());
117}
118G4double G4UTubs::GetCosEndPhi() const
119{
120 return std::cos(GetStartPhiAngle()+GetDeltaPhiAngle());
121}
122
123void G4UTubs::SetInnerRadius(G4double newRMin)
124{
125 SetRMin(newRMin);
126 fRebuildPolyhedron = true;
127}
128void G4UTubs::SetOuterRadius(G4double newRMax)
129{
130 SetRMax(newRMax);
131 fRebuildPolyhedron = true;
132}
133void G4UTubs::SetZHalfLength(G4double newDz)
134{
135 SetDz(newDz);
136 fRebuildPolyhedron = true;
137}
138void G4UTubs::SetStartPhiAngle(G4double newSPhi, G4bool)
139{
140 SetSPhi(newSPhi);
141 fRebuildPolyhedron = true;
142}
143void G4UTubs::SetDeltaPhiAngle(G4double newDPhi)
144{
145 SetDPhi(newDPhi);
146 fRebuildPolyhedron = true;
147}
148
149/////////////////////////////////////////////////////////////////////////
150//
151// Dispatch to parameterisation for replication mechanism dimension
152// computation & modification.
153
154void G4UTubs::ComputeDimensions( G4VPVParameterisation* p,
155 const G4int n,
156 const G4VPhysicalVolume* pRep )
157{
158 p->ComputeDimensions(*(G4Tubs*)this,n,pRep) ;
159}
160
161/////////////////////////////////////////////////////////////////////////
162//
163// Make a clone of the object
164
165G4VSolid* G4UTubs::Clone() const
166{
167 return new G4UTubs(*this);
168}
169
170//////////////////////////////////////////////////////////////////////////
171//
172// Get bounding box
173
174void G4UTubs::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
175{
176 static G4bool checkBBox = true;
177
178 G4double rmin = GetInnerRadius();
179 G4double rmax = GetOuterRadius();
180 G4double dz = GetZHalfLength();
181
182 // Find bounding box
183 //
184 if (GetDeltaPhiAngle() < twopi)
185 {
186 G4TwoVector vmin,vmax;
187 G4GeomTools::DiskExtent(rmin,rmax,
188 GetSinStartPhi(),GetCosStartPhi(),
189 GetSinEndPhi(),GetCosEndPhi(),
190 vmin,vmax);
191 pMin.set(vmin.x(),vmin.y(),-dz);
192 pMax.set(vmax.x(),vmax.y(), dz);
193 }
194 else
195 {
196 pMin.set(-rmax,-rmax,-dz);
197 pMax.set( rmax, rmax, dz);
198 }
199
200 // Check correctness of the bounding box
201 //
202 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
203 {
204 std::ostringstream message;
205 message << "Bad bounding box (min >= max) for solid: "
206 << GetName() << " !"
207 << "\npMin = " << pMin
208 << "\npMax = " << pMax;
209 G4Exception("G4UTubs::BoundingLimits()", "GeomMgt0001",
210 JustWarning, message);
211 StreamInfo(G4cout);
212 }
213
214 // Check consistency of bounding boxes
215 //
216 if (checkBBox)
217 {
218 U3Vector vmin, vmax;
219 Extent(vmin,vmax);
220 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
221 std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
222 std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
223 std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
224 std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
225 std::abs(pMax.z()-vmax.z()) > kCarTolerance)
226 {
227 std::ostringstream message;
228 message << "Inconsistency in bounding boxes for solid: "
229 << GetName() << " !"
230 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
231 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
232 G4Exception("G4UTubs::BoundingLimits()", "GeomMgt0001",
233 JustWarning, message);
234 checkBBox = false;
235 }
236 }
237}
238
239//////////////////////////////////////////////////////////////////////////
240//
241// Calculate extent under transform and specified limit
242
243G4bool
244G4UTubs::CalculateExtent(const EAxis pAxis,
245 const G4VoxelLimits& pVoxelLimit,
246 const G4AffineTransform& pTransform,
247 G4double& pMin, G4double& pMax) const
248{
249 G4ThreeVector bmin, bmax;
250 G4bool exist;
251
252 // Get bounding box
253 BoundingLimits(bmin,bmax);
254
255 // Check bounding box
256 G4BoundingEnvelope bbox(bmin,bmax);
257#ifdef G4BBOX_EXTENT
258 if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
259#endif
260 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
261 {
262 return exist = pMin < pMax;
263 }
264
265 // Get parameters of the solid
266 G4double rmin = GetInnerRadius();
267 G4double rmax = GetOuterRadius();
268 G4double dz = GetZHalfLength();
269 G4double dphi = GetDeltaPhiAngle();
270
271 // Find bounding envelope and calculate extent
272 //
273 const G4int NSTEPS = 24; // number of steps for whole circle
274 G4double astep = twopi/NSTEPS; // max angle for one step
275 G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
276 G4double ang = dphi/ksteps;
277
278 G4double sinHalf = std::sin(0.5*ang);
279 G4double cosHalf = std::cos(0.5*ang);
280 G4double sinStep = 2.*sinHalf*cosHalf;
281 G4double cosStep = 1. - 2.*sinHalf*sinHalf;
282 G4double rext = rmax/cosHalf;
283
284 // bounding envelope for full cylinder consists of two polygons,
285 // in other cases it is a sequence of quadrilaterals
286 if (rmin == 0 && dphi == twopi)
287 {
288 G4double sinCur = sinHalf;
289 G4double cosCur = cosHalf;
290
291 G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
292 for (G4int k=0; k<NSTEPS; ++k)
293 {
294 baseA[k].set(rext*cosCur,rext*sinCur,-dz);
295 baseB[k].set(rext*cosCur,rext*sinCur, dz);
296
297 G4double sinTmp = sinCur;
298 sinCur = sinCur*cosStep + cosCur*sinStep;
299 cosCur = cosCur*cosStep - sinTmp*sinStep;
300 }
301 std::vector<const G4ThreeVectorList *> polygons(2);
302 polygons[0] = &baseA;
303 polygons[1] = &baseB;
304 G4BoundingEnvelope benv(bmin,bmax,polygons);
305 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
306 }
307 else
308 {
309 G4double sinStart = GetSinStartPhi();
310 G4double cosStart = GetCosStartPhi();
311 G4double sinEnd = GetSinEndPhi();
312 G4double cosEnd = GetCosEndPhi();
313 G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
314 G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
315
316 // set quadrilaterals
317 G4ThreeVectorList pols[NSTEPS+2];
318 for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(4);
319 pols[0][0].set(rmin*cosStart,rmin*sinStart, dz);
320 pols[0][1].set(rmin*cosStart,rmin*sinStart,-dz);
321 pols[0][2].set(rmax*cosStart,rmax*sinStart,-dz);
322 pols[0][3].set(rmax*cosStart,rmax*sinStart, dz);
323 for (G4int k=1; k<ksteps+1; ++k)
324 {
325 pols[k][0].set(rmin*cosCur,rmin*sinCur, dz);
326 pols[k][1].set(rmin*cosCur,rmin*sinCur,-dz);
327 pols[k][2].set(rext*cosCur,rext*sinCur,-dz);
328 pols[k][3].set(rext*cosCur,rext*sinCur, dz);
329
330 G4double sinTmp = sinCur;
331 sinCur = sinCur*cosStep + cosCur*sinStep;
332 cosCur = cosCur*cosStep - sinTmp*sinStep;
333 }
334 pols[ksteps+1][0].set(rmin*cosEnd,rmin*sinEnd, dz);
335 pols[ksteps+1][1].set(rmin*cosEnd,rmin*sinEnd,-dz);
336 pols[ksteps+1][2].set(rmax*cosEnd,rmax*sinEnd,-dz);
337 pols[ksteps+1][3].set(rmax*cosEnd,rmax*sinEnd, dz);
338
339 // set envelope and calculate extent
340 std::vector<const G4ThreeVectorList *> polygons;
341 polygons.resize(ksteps+2);
342 for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
343 G4BoundingEnvelope benv(bmin,bmax,polygons);
344 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
345 }
346 return exist;
347}
348
349//////////////////////////////////////////////////////////////////////////
350//
351// Create polyhedron for visualization
352//
353G4Polyhedron* G4UTubs::CreatePolyhedron() const
354{
355 return new G4PolyhedronTubs(GetInnerRadius(),
356 GetOuterRadius(),
357 GetZHalfLength(),
358 GetStartPhiAngle(),
359 GetDeltaPhiAngle());
360}
361
362#endif // G4GEOM_USE_USOLIDS
const G4double kCarTolerance
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
G4Tubs is a tube or tube segment with curved sides parallel to the Z-axis. The tube has a specified h...
Definition G4Tubs.hh:85
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54