Geant4  10.03
G4UTrap.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 G4UTrap wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4Trap.hh"
34 #include "G4UTrap.hh"
35 
36 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
37 
38 #include "G4AffineTransform.hh"
39 #include "G4VPVParameterisation.hh"
40 #include "G4BoundingEnvelope.hh"
41 
42 using namespace CLHEP;
43 
45 //
46 // Constructors
47 //
48 G4UTrap::G4UTrap( const G4String& pName,
49  G4double pdz,
50  G4double pTheta, G4double pPhi,
51  G4double pdy1, G4double pdx1, G4double pdx2,
52  G4double pAlp1,
53  G4double pdy2, G4double pdx3, G4double pdx4,
54  G4double pAlp2 )
55  : G4USolid(pName, new UTrap(pName, pdz, pTheta, pPhi,
56  pdy1, pdx1, pdx2, pAlp1, pdy2, pdx3, pdx4, pAlp2))
57 {
58 }
59 
60 G4UTrap::G4UTrap( const G4String& pName,
61  const G4ThreeVector pt[8] )
62  : G4USolid(pName, new UTrap(pName))
63 {
64  SetPlanes(pt);
65 }
66 
67 G4UTrap::G4UTrap( const G4String& pName,
68  G4double pZ,
69  G4double pY,
70  G4double pX, G4double pLTX )
71  : G4USolid(pName, new UTrap(pName, pZ, pY, pX, pLTX))
72 {
73 }
74 
75 G4UTrap::G4UTrap( const G4String& pName,
76  G4double pdx1, G4double pdx2,
77  G4double pdy1, G4double pdy2,
78  G4double pdz )
79  : G4USolid(pName, new UTrap(pName, pdx1, pdx2, pdy1, pdy2, pdz))
80 {
81 }
82 
83 G4UTrap::G4UTrap(const G4String& pName,
84  G4double pdx, G4double pdy, G4double pdz,
85  G4double pAlpha, G4double pTheta, G4double pPhi )
86  : G4USolid(pName, new UTrap(pName, pdx, pdy, pdz, pAlpha, pTheta, pPhi))
87 {
88 }
89 
90 G4UTrap::G4UTrap( const G4String& pName )
91  : G4USolid(pName, new UTrap(pName))
92 {
93 }
94 
96 //
97 // Fake default constructor - sets only member data and allocates memory
98 // for usage restricted to object persistency.
99 //
100 G4UTrap::G4UTrap( __void__& a )
101  : G4USolid(a)
102 {
103 }
104 
106 //
107 // Destructor
108 //
109 G4UTrap::~G4UTrap()
110 {
111 }
112 
114 //
115 // Copy constructor
116 //
117 G4UTrap::G4UTrap(const G4UTrap& rhs)
118  : G4USolid(rhs)
119 {
120 }
121 
123 //
124 // Assignment operator
125 //
126 G4UTrap& G4UTrap::operator = (const G4UTrap& rhs)
127 {
128  // Check assignment to self
129  //
130  if (this == &rhs) { return *this; }
131 
132  // Copy base class data
133  //
134  G4USolid::operator=(rhs);
135 
136  return *this;
137 }
138 
140 //
141 // Accessors & modifiers
142 
143 G4double G4UTrap::GetZHalfLength() const
144 {
145  return GetShape()->GetZHalfLength();
146 }
147 G4double G4UTrap::GetYHalfLength1() const
148 {
149  return GetShape()->GetYHalfLength1();
150 }
151 G4double G4UTrap::GetXHalfLength1() const
152 {
153  return GetShape()->GetXHalfLength1();
154 }
155 G4double G4UTrap::GetXHalfLength2() const
156 {
157  return GetShape()->GetXHalfLength2();
158 }
159 G4double G4UTrap::GetTanAlpha1() const
160 {
161  return GetShape()->GetTanAlpha1();
162 }
163 G4double G4UTrap::GetYHalfLength2() const
164 {
165  return GetShape()->GetYHalfLength2();
166 }
167 G4double G4UTrap::GetXHalfLength3() const
168 {
169  return GetShape()->GetXHalfLength3();
170 }
171 G4double G4UTrap::GetXHalfLength4() const
172 {
173  return GetShape()->GetXHalfLength4();
174 }
175 G4double G4UTrap::GetTanAlpha2() const
176 {
177  return GetShape()->GetTanAlpha2();
178 }
179 TrapSidePlane G4UTrap::GetSidePlane(G4int n) const
180 {
181  UTrapSidePlane iplane = GetShape()->GetSidePlane(n);
182  TrapSidePlane oplane = {iplane.a, iplane.b, iplane.c, iplane.d };
183  return oplane;
184 }
185 G4ThreeVector G4UTrap::GetSymAxis() const
186 {
187  UVector3 axis = GetShape()->GetSymAxis();
188  return G4ThreeVector(axis.x(), axis.y(), axis.z());
189 }
190 
191 void G4UTrap::SetAllParameters(G4double pDz, G4double pTheta, G4double pPhi,
192  G4double pDy1, G4double pDx1, G4double pDx2,
193  G4double pAlp1,
194  G4double pDy2, G4double pDx3, G4double pDx4,
195  G4double pAlp2)
196 {
197  GetShape()->SetAllParameters(pDz, pTheta, pPhi,
198  pDy1, pDx1, pDx2, pAlp1,
199  pDy2, pDx3, pDx4, pAlp2);
200  fRebuildPolyhedron = true;
201 }
202 
203 void G4UTrap::SetPlanes(const G4ThreeVector pt[8])
204 {
205  UVector3 upt[8];
206  for (unsigned int i=0; i<8; ++i)
207  {
208  upt[i] = UVector3(pt[i].x(), pt[i].y(), pt[i].z());
209  }
210  GetShape()->SetPlanes(upt);
211  fRebuildPolyhedron = true;
212 }
213 
215 //
216 // Dispatch to parameterisation for replication mechanism dimension
217 // computation & modification.
218 //
219 void G4UTrap::ComputeDimensions( G4VPVParameterisation* p,
220  const G4int n,
221  const G4VPhysicalVolume* pRep)
222 {
223  p->ComputeDimensions(*(G4Trap*)this,n,pRep);
224 }
225 
227 //
228 // Make a clone of the object
229 //
230 G4VSolid* G4UTrap::Clone() const
231 {
232  return new G4UTrap(*this);
233 }
234 
236 //
237 // Get bounding box
238 
239 void G4UTrap::Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const
240 {
241  static G4bool checkBBox = true;
242 
243  G4double dz = GetZHalfLength();
244  G4double dx1 = GetXHalfLength1();
245  G4double dx2 = GetXHalfLength2();
246  G4double dx3 = GetXHalfLength3();
247  G4double dx4 = GetXHalfLength4();
248  G4double dy1 = GetYHalfLength1();
249  G4double dy2 = GetYHalfLength2();
250  G4double fTthetaSphi = GetShape()->GetThetaSphi();
251  G4double fTthetaCphi = GetShape()->GetThetaCphi();
252 
253  G4double x0 = dz*fTthetaCphi;
254  G4double x1 = dy1*GetTanAlpha1();
255  G4double x2 = dy2*GetTanAlpha2();
256  G4double xmin =
257  std::min(
258  std::min(
259  std::min(-x0-x1-dx1,-x0+x1-dx2),x0-x2-dx3),x0+x2-dx4);
260  G4double xmax =
261  std::max(
262  std::max(
263  std::max(-x0-x1+dx1,-x0+x1+dx2),x0-x2+dx3),x0+x2+dx4);
264 
265  G4double y0 = dz*fTthetaSphi;
266  G4double ymin = std::min(-y0-dy1,y0-dy2);
267  G4double ymax = std::max(-y0+dy1,y0+dy2);
268 
269  pMin.set(xmin,ymin,-dz);
270  pMax.set(xmax,ymax, dz);
271 
272  // Check correctness of the bounding box
273  //
274  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
275  {
276  std::ostringstream message;
277  message << "Bad bounding box (min >= max) for solid: "
278  << GetName() << " !"
279  << "\npMin = " << pMin
280  << "\npMax = " << pMax;
281  G4Exception("G4UTrap::Extent()", "GeomMgt0001", JustWarning, message);
282  StreamInfo(G4cout);
283  }
284 
285  // Check consistency of bounding boxes
286  //
287  if (checkBBox)
288  {
289  UVector3 vmin, vmax;
290  GetShape()->Extent(vmin,vmax);
291  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
292  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
293  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
294  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
295  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
296  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
297  {
298  std::ostringstream message;
299  message << "Inconsistency in bounding boxes for solid: "
300  << GetName() << " !"
301  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
302  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
303  G4Exception("G4UTrap::Extent()", "GeomMgt0001", JustWarning, message);
304  checkBBox = false;
305  }
306  }
307 }
308 
310 //
311 // Calculate extent under transform and specified limit
312 
313 G4bool
314 G4UTrap::CalculateExtent(const EAxis pAxis,
315  const G4VoxelLimits& pVoxelLimit,
316  const G4AffineTransform& pTransform,
317  G4double& pMin, G4double& pMax) const
318 {
319  G4ThreeVector bmin, bmax;
320  G4bool exist;
321 
322  // Check bounding box (bbox)
323  //
324  Extent(bmin,bmax);
325  G4BoundingEnvelope bbox(bmin,bmax);
326 #ifdef G4BBOX_EXTENT
327  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
328 #endif
329  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
330  {
331  return exist = (pMin < pMax) ? true : false;
332  }
333 
334  // Set bounding envelope (benv) and calculate extent
335  //
336  G4double dz = GetZHalfLength();
337  G4double dx1 = GetXHalfLength1();
338  G4double dx2 = GetXHalfLength2();
339  G4double dx3 = GetXHalfLength3();
340  G4double dx4 = GetXHalfLength4();
341  G4double dy1 = GetYHalfLength1();
342  G4double dy2 = GetYHalfLength2();
343  G4double fTthetaSphi = GetShape()->GetThetaSphi();
344  G4double fTthetaCphi = GetShape()->GetThetaCphi();
345 
346  G4double x0 = dz*fTthetaCphi;
347  G4double x1 = dy1*GetTanAlpha1();
348  G4double x2 = dy2*GetTanAlpha2();
349  G4double y0 = dz*fTthetaSphi;
350 
351  G4ThreeVectorList baseA(4), baseB(4);
352  baseA[0].set(-x0-x1-dx1,-y0-dy1,-dz);
353  baseA[1].set(-x0-x1+dx1,-y0-dy1,-dz);
354  baseA[2].set(-x0+x1+dx2,-y0+dy1,-dz);
355  baseA[3].set(-x0+x1-dx2,-y0+dy1,-dz);
356 
357  baseB[0].set( x0-x2-dx3, y0-dy2, dz);
358  baseB[1].set( x0-x2+dx3, y0-dy2, dz);
359  baseB[2].set( x0+x2+dx4, y0+dy2, dz);
360  baseB[3].set( x0+x2-dx4, y0+dy2, dz);
361 
362  std::vector<const G4ThreeVectorList *> polygons(2);
363  polygons[0] = &baseA;
364  polygons[1] = &baseB;
365 
366  G4BoundingEnvelope benv(bmin,bmax,polygons);
367  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
368  return exist;
369 }
370 
372 //
373 // Create polyhedron for visualization
374 //
375 G4Polyhedron* G4UTrap::CreatePolyhedron() const
376 {
377  G4double fTthetaSphi = GetShape()->GetThetaSphi();
378  G4double fTthetaCphi = GetShape()->GetThetaCphi();
379  G4double phi = std::atan2(fTthetaSphi, fTthetaCphi);
380  G4double alpha1 = std::atan(GetTanAlpha1());
381  G4double alpha2 = std::atan(GetTanAlpha2());
382  G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi+fTthetaSphi*fTthetaSphi));
383 
384  return new G4PolyhedronTrap(GetZHalfLength(), theta, phi,
385  GetYHalfLength1(),
386  GetXHalfLength1(), GetXHalfLength2(), alpha1,
387  GetYHalfLength2(),
388  GetXHalfLength3(), GetXHalfLength4(), alpha2);
389 }
390 
391 #endif // G4GEOM_USE_USOLIDS
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
int G4int
Definition: G4Types.hh:78
G4double a
Definition: G4Trap.hh:103
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
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EAxis
Definition: geomdefs.hh:54
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76