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