Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UCons.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 G4UCons wrapper class
27//
28// 30.10.13 G.Cosmo, CERN
29// --------------------------------------------------------------------
30
31#include "G4Cons.hh"
32#include "G4UCons.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
48G4UCons::G4UCons( const G4String& pName,
49 G4double pRmin1, G4double pRmax1,
50 G4double pRmin2, G4double pRmax2,
51 G4double pDz,
52 G4double pSPhi, G4double pDPhi)
53 : Base_t(pName, pRmin1, pRmax1, pRmin2, pRmax2, pDz, pSPhi, pDPhi)
54{
55}
56
57//////////////////////////////////////////////////////////////////////////
58//
59// Copy constructor
60
61G4UCons::G4UCons(const G4UCons& rhs)
62 : Base_t(rhs)
63{
64}
65
66//////////////////////////////////////////////////////////////////////////
67//
68// Assignment operator
69
70G4UCons& G4UCons::operator = (const G4UCons& rhs)
71{
72 // Check assignment to self
73 //
74 if (this == &rhs) { return *this; }
75
76 // Copy base class data
77 //
78 Base_t::operator=(rhs);
79
80 return *this;
81}
82
83/////////////////////////////////////////////////////////////////////////
84//
85// Accessors and modifiers
86
87G4double G4UCons::GetInnerRadiusMinusZ() const
88{
89 return GetRmin1();
90}
91G4double G4UCons::GetOuterRadiusMinusZ() const
92{
93 return GetRmax1();
94}
95G4double G4UCons::GetInnerRadiusPlusZ() const
96{
97 return GetRmin2();
98}
99G4double G4UCons::GetOuterRadiusPlusZ() const
100{
101 return GetRmax2();
102}
103G4double G4UCons::GetZHalfLength() const
104{
105 return GetDz();
106}
107G4double G4UCons::GetStartPhiAngle() const
108{
109 return GetSPhi();
110}
111G4double G4UCons::GetDeltaPhiAngle() const
112{
113 return GetDPhi();
114}
115G4double G4UCons::GetSinStartPhi() const
116{
117 G4double phi = GetStartPhiAngle();
118 return std::sin(phi);
119}
120G4double G4UCons::GetCosStartPhi() const
121{
122 G4double phi = GetStartPhiAngle();
123 return std::cos(phi);
124}
125G4double G4UCons::GetSinEndPhi() const
126{
127 G4double phi = GetStartPhiAngle() + GetDeltaPhiAngle();
128 return std::sin(phi);
129}
130G4double G4UCons::GetCosEndPhi() const
131{
132 G4double phi = GetStartPhiAngle() + GetDeltaPhiAngle();
133 return std::cos(phi);
134}
135
136void G4UCons::SetInnerRadiusMinusZ(G4double Rmin1)
137{
138 SetRmin1(Rmin1);
139 fRebuildPolyhedron = true;
140}
141void G4UCons::SetOuterRadiusMinusZ(G4double Rmax1)
142{
143 SetRmax1(Rmax1);
144 fRebuildPolyhedron = true;
145}
146void G4UCons::SetInnerRadiusPlusZ(G4double Rmin2)
147{
148 SetRmin2(Rmin2);
149 fRebuildPolyhedron = true;
150}
151void G4UCons::SetOuterRadiusPlusZ(G4double Rmax2)
152{
153 SetRmax2(Rmax2);
154 fRebuildPolyhedron = true;
155}
156void G4UCons::SetZHalfLength(G4double newDz)
157{
158 SetDz(newDz);
159 fRebuildPolyhedron = true;
160}
161void G4UCons::SetStartPhiAngle(G4double newSPhi, G4bool)
162{
163 SetSPhi(newSPhi);
164 fRebuildPolyhedron = true;
165}
166void G4UCons::SetDeltaPhiAngle(G4double newDPhi)
167{
168 SetDPhi(newDPhi);
169 fRebuildPolyhedron = true;
170}
171
172/////////////////////////////////////////////////////////////////////////
173//
174// Dispatch to parameterisation for replication mechanism dimension
175// computation & modification.
176
177void G4UCons::ComputeDimensions( G4VPVParameterisation* p,
178 const G4int n,
179 const G4VPhysicalVolume* pRep )
180{
181 p->ComputeDimensions(*(G4Cons*)this,n,pRep);
182}
183
184//////////////////////////////////////////////////////////////////////////
185//
186// Make a clone of the object
187
188G4VSolid* G4UCons::Clone() const
189{
190 return new G4UCons(*this);
191}
192
193//////////////////////////////////////////////////////////////////////////
194//
195// Get bounding box
196
197void G4UCons::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
198{
199 static G4bool checkBBox = true;
200
201 G4double rmin = std::min(GetInnerRadiusMinusZ(),GetInnerRadiusPlusZ());
202 G4double rmax = std::max(GetOuterRadiusMinusZ(),GetOuterRadiusPlusZ());
203 G4double dz = GetZHalfLength();
204
205 // Find bounding box
206 //
207 if (GetDeltaPhiAngle() < twopi)
208 {
209 G4TwoVector vmin,vmax;
210 G4GeomTools::DiskExtent(rmin,rmax,
211 GetSinStartPhi(),GetCosStartPhi(),
212 GetSinEndPhi(),GetCosEndPhi(),
213 vmin,vmax);
214 pMin.set(vmin.x(),vmin.y(),-dz);
215 pMax.set(vmax.x(),vmax.y(), dz);
216 }
217 else
218 {
219 pMin.set(-rmax,-rmax,-dz);
220 pMax.set( rmax, rmax, dz);
221 }
222
223 // Check correctness of the bounding box
224 //
225 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
226 {
227 std::ostringstream message;
228 message << "Bad bounding box (min >= max) for solid: "
229 << GetName() << " !"
230 << "\npMin = " << pMin
231 << "\npMax = " << pMax;
232 G4Exception("G4UCons::BoundingLimits()", "GeomMgt0001",
233 JustWarning, message);
234 StreamInfo(G4cout);
235 }
236
237 // Check consistency of bounding boxes
238 //
239 if (checkBBox)
240 {
241 U3Vector vmin, vmax;
242 Extent(vmin,vmax);
243 if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
244 std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
245 std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
246 std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
247 std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
248 std::abs(pMax.z()-vmax.z()) > kCarTolerance)
249 {
250 std::ostringstream message;
251 message << "Inconsistency in bounding boxes for solid: "
252 << GetName() << " !"
253 << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
254 << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
255 G4Exception("G4UCons::BoundingLimits()", "GeomMgt0001",
256 JustWarning, message);
257 checkBBox = false;
258 }
259 }
260}
261
262/////////////////////////////////////////////////////////////////////////
263//
264// Calculate extent under transform and specified limit
265
266G4bool
267G4UCons::CalculateExtent(const EAxis pAxis,
268 const G4VoxelLimits& pVoxelLimit,
269 const G4AffineTransform& pTransform,
270 G4double& pMin, G4double& pMax) const
271{
272 G4ThreeVector bmin, bmax;
273 G4bool exist;
274
275 // Get bounding box
276 BoundingLimits(bmin,bmax);
277
278 // Check bounding box
279 G4BoundingEnvelope bbox(bmin,bmax);
280#ifdef G4BBOX_EXTENT
281 if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
282#endif
283 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
284 {
285 return exist = pMin < pMax;
286 }
287
288 // Get parameters of the solid
289 G4double rmin1 = GetInnerRadiusMinusZ();
290 G4double rmax1 = GetOuterRadiusMinusZ();
291 G4double rmin2 = GetInnerRadiusPlusZ();
292 G4double rmax2 = GetOuterRadiusPlusZ();
293 G4double dz = GetZHalfLength();
294 G4double dphi = GetDeltaPhiAngle();
295
296 // Find bounding envelope and calculate extent
297 //
298 const G4int NSTEPS = 24; // number of steps for whole circle
299 G4double astep = twopi/NSTEPS; // max angle for one step
300 G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
301 G4double ang = dphi/ksteps;
302
303 G4double sinHalf = std::sin(0.5*ang);
304 G4double cosHalf = std::cos(0.5*ang);
305 G4double sinStep = 2.*sinHalf*cosHalf;
306 G4double cosStep = 1. - 2.*sinHalf*sinHalf;
307 G4double rext1 = rmax1/cosHalf;
308 G4double rext2 = rmax2/cosHalf;
309
310 // bounding envelope for full cone without hole consists of two polygons,
311 // in other cases it is a sequence of quadrilaterals
312 if (rmin1 == 0 && rmin2 == 0 && dphi == twopi)
313 {
314 G4double sinCur = sinHalf;
315 G4double cosCur = cosHalf;
316
317 G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
318 for (G4int k=0; k<NSTEPS; ++k)
319 {
320 baseA[k].set(rext1*cosCur,rext1*sinCur,-dz);
321 baseB[k].set(rext2*cosCur,rext2*sinCur, dz);
322
323 G4double sinTmp = sinCur;
324 sinCur = sinCur*cosStep + cosCur*sinStep;
325 cosCur = cosCur*cosStep - sinTmp*sinStep;
326 }
327 std::vector<const G4ThreeVectorList *> polygons(2);
328 polygons[0] = &baseA;
329 polygons[1] = &baseB;
330 G4BoundingEnvelope benv(bmin,bmax,polygons);
331 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
332 }
333 else
334 {
335 G4double sinStart = GetSinStartPhi();
336 G4double cosStart = GetCosStartPhi();
337 G4double sinEnd = GetSinEndPhi();
338 G4double cosEnd = GetCosEndPhi();
339 G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
340 G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
341
342 // set quadrilaterals
343 G4ThreeVectorList pols[NSTEPS+2];
344 for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(4);
345 pols[0][0].set(rmin2*cosStart,rmin2*sinStart, dz);
346 pols[0][1].set(rmin1*cosStart,rmin1*sinStart,-dz);
347 pols[0][2].set(rmax1*cosStart,rmax1*sinStart,-dz);
348 pols[0][3].set(rmax2*cosStart,rmax2*sinStart, dz);
349 for (G4int k=1; k<ksteps+1; ++k)
350 {
351 pols[k][0].set(rmin2*cosCur,rmin2*sinCur, dz);
352 pols[k][1].set(rmin1*cosCur,rmin1*sinCur,-dz);
353 pols[k][2].set(rext1*cosCur,rext1*sinCur,-dz);
354 pols[k][3].set(rext2*cosCur,rext2*sinCur, dz);
355
356 G4double sinTmp = sinCur;
357 sinCur = sinCur*cosStep + cosCur*sinStep;
358 cosCur = cosCur*cosStep - sinTmp*sinStep;
359 }
360 pols[ksteps+1][0].set(rmin2*cosEnd,rmin2*sinEnd, dz);
361 pols[ksteps+1][1].set(rmin1*cosEnd,rmin1*sinEnd,-dz);
362 pols[ksteps+1][2].set(rmax1*cosEnd,rmax1*sinEnd,-dz);
363 pols[ksteps+1][3].set(rmax2*cosEnd,rmax2*sinEnd, dz);
364
365 // set envelope and calculate extent
366 std::vector<const G4ThreeVectorList *> polygons;
367 polygons.resize(ksteps+2);
368 for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
369 G4BoundingEnvelope benv(bmin,bmax,polygons);
370 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
371 }
372 return exist;
373}
374
375//////////////////////////////////////////////////////////////////////////
376//
377// Create polyhedron for visualization
378
379G4Polyhedron* G4UCons::CreatePolyhedron() const
380{
381 return new G4PolyhedronCons(GetInnerRadiusMinusZ(),
382 GetOuterRadiusMinusZ(),
383 GetInnerRadiusPlusZ(),
384 GetOuterRadiusPlusZ(),
385 GetZHalfLength(),
386 GetStartPhiAngle(),
387 GetDeltaPhiAngle());
388}
389
390#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...
G4Cons is, in the general case, a Phi segment of a cone, with half-length fDz, inner and outer radii ...
Definition G4Cons.hh:85
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
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