Geant4  10.03
G4UTorus.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 // Implementation for G4UTorus wrapper class
30 //
31 // 19-08-2015 Guilherme Lima, FNAL
32 //
33 // --------------------------------------------------------------------
34 
35 #include "G4Torus.hh"
36 #include "G4UTorus.hh"
37 
38 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
39 
40 #include "G4TwoVector.hh"
41 #include "G4GeomTools.hh"
42 #include "G4AffineTransform.hh"
43 #include "G4BoundingEnvelope.hh"
44 
45 #include "G4VPVParameterisation.hh"
46 
47 using namespace CLHEP;
48 
50 //
51 // Constructor - check & set half widths
52 
53 
54 G4UTorus::G4UTorus(const G4String& pName,
55  G4double rmin, G4double rmax, G4double rtor,
56  G4double sphi, G4double dphi)
57  : G4USolid(pName, new UTorus(pName, rmin, rmax, rtor, sphi, dphi))
58 { }
59 
61 //
62 // Fake default constructor - sets only member data and allocates memory
63 // for usage restricted to object persistency.
64 
65 G4UTorus::G4UTorus( __void__& a )
66  : G4USolid(a)
67 { }
68 
70 //
71 // Destructor
72 
73 G4UTorus::~G4UTorus() { }
74 
76 //
77 // Copy constructor
78 
79 G4UTorus::G4UTorus(const G4UTorus& rhs)
80  : G4USolid(rhs)
81 { }
82 
84 //
85 // Assignment operator
86 
87 G4UTorus& G4UTorus::operator = (const G4UTorus& rhs)
88 {
89  // Check assignment to self
90  //
91  if (this == &rhs) { return *this; }
92 
93  // Copy base class data
94  //
95  G4USolid::operator=(rhs);
96 
97  return *this;
98 }
99 
101 //
102 // Accessors & modifiers
103 
104 G4double G4UTorus::GetRmin() const
105 {
106  return GetShape()->GetRmin();
107 }
108 
109 G4double G4UTorus::GetRmax() const
110 {
111  return GetShape()->GetRmax();
112 }
113 
114 G4double G4UTorus::GetRtor() const
115 {
116  return GetShape()->GetRtor();
117 }
118 
119 G4double G4UTorus::GetSPhi() const
120 {
121  return GetShape()->GetSPhi();
122 }
123 
124 G4double G4UTorus::GetDPhi() const
125 {
126  return GetShape()->GetDPhi();
127 }
128 
129 G4double G4UTorus::GetSinStartPhi() const
130 {
131  G4double phi = GetShape()->GetSPhi();
132  return std::sin(phi);
133 }
134 
135 G4double G4UTorus::GetCosStartPhi() const
136 {
137  G4double phi = GetShape()->GetSPhi();
138  return std::cos(phi);
139 }
140 
141 G4double G4UTorus::GetSinEndPhi() const
142 {
143  G4double phi = GetShape()->GetSPhi() +
144  GetShape()->GetDPhi();
145  return std::sin(phi);
146 }
147 
148 G4double G4UTorus::GetCosEndPhi() const
149 {
150  G4double phi = GetShape()->GetSPhi() +
151  GetShape()->GetDPhi();
152  return std::cos(phi);
153 }
154 
155 void G4UTorus::SetRmin(G4double arg)
156 {
157  GetShape()->SetRmin(arg);
158  fRebuildPolyhedron = true;
159 }
160 
161 void G4UTorus::SetRmax(G4double arg)
162 {
163  GetShape()->SetRmax(arg);
164  fRebuildPolyhedron = true;
165 }
166 
167 void G4UTorus::SetRtor(G4double arg)
168 {
169  GetShape()->SetRtor(arg);
170  fRebuildPolyhedron = true;
171 }
172 
173 void G4UTorus::SetSPhi(G4double arg)
174 {
175  GetShape()->SetSPhi(arg);
176  fRebuildPolyhedron = true;
177 }
178 
179 void G4UTorus::SetDPhi(G4double arg)
180 {
181  GetShape()->SetDPhi(arg);
182  fRebuildPolyhedron = true;
183 }
184 
185 void G4UTorus::SetAllParameters(G4double arg1, G4double arg2,
186  G4double arg3, G4double arg4, G4double arg5)
187 {
188  GetShape()->SetRmin(arg1);
189  GetShape()->SetRmax(arg2);
190  GetShape()->SetRtor(arg3);
191  GetShape()->SetSPhi(arg4);
192  GetShape()->SetDPhi(arg5);
193  fRebuildPolyhedron = true;
194 }
195 
197 //
198 // Dispatch to parameterisation for replication mechanism dimension
199 // computation & modification.
200 
201 void G4UTorus::ComputeDimensions(G4VPVParameterisation* p,
202  const G4int n,
203  const G4VPhysicalVolume* pRep)
204 {
205  p->ComputeDimensions(*(G4Torus*)this,n,pRep);
206 }
207 
209 //
210 // Make a clone of the object
211 
212 G4VSolid* G4UTorus::Clone() const
213 {
214  return new G4UTorus(*this);
215 }
216 
218 //
219 // Get bounding box
220 
221 void G4UTorus::Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const
222 {
223  static G4bool checkBBox = true;
224 
225  G4double rmax = GetRmax();
226  G4double rtor = GetRtor();
227  G4double rint = rtor - rmax;
228  G4double rext = rtor + rmax;
229  G4double dz = rmax;
230 
231  // Find bounding box
232  //
233  if (GetDPhi() >= twopi)
234  {
235  pMin.set(-rext,-rext,-dz);
236  pMax.set( rext, rext, dz);
237  }
238  else
239  {
240  G4TwoVector vmin,vmax;
241  G4GeomTools::DiskExtent(rint,rext,
242  GetSinStartPhi(),GetCosStartPhi(),
243  GetSinEndPhi(),GetCosEndPhi(),
244  vmin,vmax);
245  pMin.set(vmin.x(),vmin.y(),-dz);
246  pMax.set(vmax.x(),vmax.y(), dz);
247  }
248 
249  // Check correctness of the bounding box
250  //
251  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
252  {
253  std::ostringstream message;
254  message << "Bad bounding box (min >= max) for solid: "
255  << GetName() << " !"
256  << "\npMin = " << pMin
257  << "\npMax = " << pMax;
258  G4Exception("G4UTorus::Extent()", "GeomMgt0001", JustWarning, message);
259  StreamInfo(G4cout);
260  }
261 
262  // Check consistency of bounding boxes
263  //
264  if (checkBBox)
265  {
266  UVector3 vmin, vmax;
267  GetShape()->Extent(vmin,vmax);
268  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
269  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
270  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
271  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
272  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
273  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
274  {
275  std::ostringstream message;
276  message << "Inconsistency in bounding boxes for solid: "
277  << GetName() << " !"
278  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
279  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
280  G4Exception("G4UTorus::Extent()", "GeomMgt0001", JustWarning, message);
281  checkBBox = false;
282  }
283  }
284 }
285 
287 //
288 // Calculate extent under transform and specified limit
289 
290 G4bool
291 G4UTorus::CalculateExtent(const EAxis pAxis,
292  const G4VoxelLimits& pVoxelLimit,
293  const G4AffineTransform& pTransform,
294  G4double& pMin, G4double& pMax) const
295 {
296  G4ThreeVector bmin, bmax;
297  G4bool exist;
298 
299  // Get bounding box
300  Extent(bmin,bmax);
301 
302  // Check bounding box
303  G4BoundingEnvelope bbox(bmin,bmax);
304 #ifdef G4BBOX_EXTENT
305  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
306 #endif
307  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
308  {
309  return exist = (pMin < pMax) ? true : false;
310  }
311 
312  // Get parameters of the solid
313  G4double rmin = GetRmin();
314  G4double rmax = GetRmax();
315  G4double rtor = GetRtor();
316  G4double dphi = GetDPhi();
317  G4double sinStart = GetSinStartPhi();
318  G4double cosStart = GetCosStartPhi();
319  G4double sinEnd = GetSinEndPhi();
320  G4double cosEnd = GetCosEndPhi();
321  G4double rint = rtor - rmax;
322  G4double rext = rtor + rmax;
323 
324  // Find bounding envelope and calculate extent
325  //
326  static const G4int NPHI = 24; // number of steps for whole torus
327  static const G4int NDISK = 16; // number of steps for disk
328  static const G4double sinHalfDisk = std::sin(pi/NDISK);
329  static const G4double cosHalfDisk = std::cos(pi/NDISK);
330  static const G4double sinStepDisk = 2.*sinHalfDisk*cosHalfDisk;
331  static const G4double cosStepDisk = 1. - 2.*sinHalfDisk*sinHalfDisk;
332 
333  G4double astep = (360/NPHI)*deg; // max angle for one slice in phi
334  G4int kphi = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
335  G4double ang = dphi/kphi;
336 
337  G4double sinHalf = std::sin(0.5*ang);
338  G4double cosHalf = std::cos(0.5*ang);
339  G4double sinStep = 2.*sinHalf*cosHalf;
340  G4double cosStep = 1. - 2.*sinHalf*sinHalf;
341 
342  // define vectors for bounding envelope
343  G4ThreeVectorList pols[NDISK+1];
344  for (G4int k=0; k<NDISK+1; ++k) pols[k].resize(4);
345 
346  std::vector<const G4ThreeVectorList *> polygons;
347  polygons.resize(NDISK+1);
348  for (G4int k=0; k<NDISK+1; ++k) polygons[k] = &pols[k];
349 
350  // set internal and external reference circles
351  G4TwoVector rzmin[NDISK];
352  G4TwoVector rzmax[NDISK];
353 
354  if ((rtor-rmin*sinHalfDisk)/cosHalf > (rtor+rmin*sinHalfDisk)) rmin = 0;
355  rmax /= cosHalfDisk;
356  G4double sinCurDisk = sinHalfDisk;
357  G4double cosCurDisk = cosHalfDisk;
358  for (G4int k=0; k<NDISK; ++k)
359  {
360  G4double rmincur = rtor + rmin*cosCurDisk;
361  if (cosCurDisk < 0 && rmin > 0) rmincur /= cosHalf;
362  rzmin[k].set(rmincur,rmin*sinCurDisk);
363 
364  G4double rmaxcur = rtor + rmax*cosCurDisk;
365  if (cosCurDisk > 0) rmaxcur /= cosHalf;
366  rzmax[k].set(rmaxcur,rmax*sinCurDisk);
367 
368  G4double sinTmpDisk = sinCurDisk;
369  sinCurDisk = sinCurDisk*cosStepDisk + cosCurDisk*sinStepDisk;
370  cosCurDisk = cosCurDisk*cosStepDisk - sinTmpDisk*sinStepDisk;
371  }
372 
373  // Loop along slices in Phi. The extent is calculated as cumulative
374  // extent of the slices
375  pMin = kInfinity;
376  pMax = -kInfinity;
377  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
378  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
379  G4double sinCur1 = 0, cosCur1 = 0, sinCur2 = 0, cosCur2 = 0;
380  for (G4int i=0; i<kphi+1; ++i)
381  {
382  if (i == 0)
383  {
384  sinCur1 = sinStart;
385  cosCur1 = cosStart;
386  sinCur2 = sinCur1*cosHalf + cosCur1*sinHalf;
387  cosCur2 = cosCur1*cosHalf - sinCur1*sinHalf;
388  }
389  else
390  {
391  sinCur1 = sinCur2;
392  cosCur1 = cosCur2;
393  sinCur2 = (i == kphi) ? sinEnd : sinCur1*cosStep + cosCur1*sinStep;
394  cosCur2 = (i == kphi) ? cosEnd : cosCur1*cosStep - sinCur1*sinStep;
395  }
396  for (G4int k=0; k<NDISK; ++k)
397  {
398  G4double r1 = rzmin[k].x(), r2 = rzmax[k].x();
399  G4double z1 = rzmin[k].y(), z2 = rzmax[k].y();
400  pols[k][0].set(r1*cosCur1,r1*sinCur1,z1);
401  pols[k][1].set(r2*cosCur1,r2*sinCur1,z2);
402  pols[k][2].set(r2*cosCur2,r2*sinCur2,z2);
403  pols[k][3].set(r1*cosCur2,r1*sinCur2,z1);
404  }
405  pols[NDISK] = pols[0];
406 
407  // get bounding box of current slice
408  G4TwoVector vmin,vmax;
410  DiskExtent(rint,rext,sinCur1,cosCur1,sinCur2,cosCur2,vmin,vmax);
411  bmin.setX(vmin.x()); bmin.setY(vmin.y());
412  bmax.setX(vmax.x()); bmax.setY(vmax.y());
413 
414  // set bounding envelope for current slice and adjust extent
415  G4double emin,emax;
416  G4BoundingEnvelope benv(bmin,bmax,polygons);
417  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
418  if (emin < pMin) pMin = emin;
419  if (emax > pMax) pMax = emax;
420  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
421  }
422  return (pMin < pMax);
423 }
424 
426 //
427 // Create polyhedron for visualization
428 
429 G4Polyhedron* G4UTorus::CreatePolyhedron() const
430 {
431  return new G4PolyhedronTorus(GetRmin(),
432  GetRmax(),
433  GetRtor(),
434  GetSPhi(),
435  GetDPhi());
436 }
437 
438 #endif // G4GEOM_USE_USOLIDS
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
int G4int
Definition: G4Types.hh:78
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
static const G4double emax
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
EAxis
Definition: geomdefs.hh:54
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
static constexpr double deg
Definition: G4SIunits.hh:152
G4double GetMaxExtent(const EAxis pAxis) const
G4double GetMinExtent(const EAxis pAxis) const
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:378