Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UTrap.cc
이 파일의 문서화 페이지로 가기
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  : Base_t(pName, pdz, pTheta, pPhi, pdy1, pdx1, pdx2,
56  pAlp1, pdy2, pdx3, pdx4, pAlp2)
57 {
58 }
59 
60 G4UTrap::G4UTrap( const G4String& pName,
61  const G4ThreeVector pt[8] )
62  : Base_t(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  : Base_t(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  : Base_t(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  : Base_t(pName, pdx, pdy, pdz, pAlpha, pTheta, pPhi)
87 {
88 }
89 
90 G4UTrap::G4UTrap( const G4String& pName )
91  : Base_t(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  : Base_t(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  : Base_t(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  Base_t::operator=(rhs);
135 
136  return *this;
137 }
138 
140 //
141 // Accessors & modifiers
142 
143 G4double G4UTrap::GetZHalfLength() const
144 {
145  return GetDz();
146 }
147 G4double G4UTrap::GetYHalfLength1() const
148 {
149  return GetDy1();
150 }
151 G4double G4UTrap::GetXHalfLength1() const
152 {
153  return GetDx1();
154 }
155 G4double G4UTrap::GetXHalfLength2() const
156 {
157  return GetDx2();
158 }
159 G4double G4UTrap::GetYHalfLength2() const
160 {
161  return GetDy2();
162 }
163 G4double G4UTrap::GetXHalfLength3() const
164 {
165  return GetDx3();
166 }
167 G4double G4UTrap::GetXHalfLength4() const
168 {
169  return GetDx4();
170 }
171 G4double G4UTrap::GetThetaCphi() const
172 {
173  return GetTanThetaCosPhi();
174 }
175 G4double G4UTrap::GetThetaSphi() const
176 {
177  return GetTanThetaSinPhi();
178 }
179 TrapSidePlane G4UTrap::GetSidePlane(G4int n) const
180 {
181  TrapSidePlane plane;
182  plane.a = GetStruct().GetPlane(n).fA;
183  plane.b = GetStruct().GetPlane(n).fB;
184  plane.c = GetStruct().GetPlane(n).fC;
185  plane.d = GetStruct().GetPlane(n).fD;
186  return plane;
187 }
188 G4ThreeVector G4UTrap::GetSymAxis() const
189 {
190  G4double tanThetaSphi = GetTanThetaSinPhi();
191  G4double tanThetaCphi = GetTanThetaCosPhi();
192  G4double tan2Theta = tanThetaSphi*tanThetaSphi + tanThetaCphi*tanThetaCphi;
193  G4double cosTheta = 1.0 / std::sqrt(1 + tan2Theta);
194  return G4ThreeVector(tanThetaCphi*cosTheta, tanThetaSphi*cosTheta, cosTheta);
195 }
196 
197 void G4UTrap::SetAllParameters(G4double pDz, G4double pTheta, G4double pPhi,
198  G4double pDy1, G4double pDx1, G4double pDx2,
199  G4double pAlp1,
200  G4double pDy2, G4double pDx3, G4double pDx4,
201  G4double pAlp2)
202 {
203  SetDz(pDz);
204  SetDy1(pDy1);
205  SetDy2(pDy2);
206  SetDx1(pDx1);
207  SetDx2(pDx2);
208  SetDx3(pDx3);
209  SetDx4(pDx4);
210  SetTanAlpha1(std::tan(pAlp1));
211  SetTanAlpha1(std::tan(pAlp2));
212  // last two will also reset cached variables
213  SetTheta(pTheta);
214  SetPhi(pPhi);
215  fRebuildPolyhedron = true;
216 }
217 
218 void G4UTrap::SetPlanes(const G4ThreeVector pt[8])
219 {
220  U3Vector upt[8];
221  for (unsigned int i=0; i<8; ++i)
222  {
223  upt[i] = U3Vector(pt[i].x(), pt[i].y(), pt[i].z());
224  }
225  fromCornersToParameters(upt);
226  fRebuildPolyhedron = true;
227 }
228 
230 //
231 // Dispatch to parameterisation for replication mechanism dimension
232 // computation & modification.
233 //
234 void G4UTrap::ComputeDimensions( G4VPVParameterisation* p,
235  const G4int n,
236  const G4VPhysicalVolume* pRep)
237 {
238  p->ComputeDimensions(*(G4Trap*)this,n,pRep);
239 }
240 
242 //
243 // Make a clone of the object
244 //
245 G4VSolid* G4UTrap::Clone() const
246 {
247  return new G4UTrap(*this);
248 }
249 
251 //
252 // Get bounding box
253 
254 void G4UTrap::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
255 {
256  static G4bool checkBBox = true;
257 
258  TrapSidePlane planes[4];
259  for (G4int i=0; i<4; ++i) { planes[i] = GetSidePlane(i); }
260 
261  G4double xmin = kInfinity, xmax = -kInfinity;
263  G4double dz = GetZHalfLength();
264  for (G4int i=0; i<8; ++i)
265  {
266  G4int iy = (i==0 || i==1 || i==4 || i==5) ? 0 : 1;
267  G4int ix = (i==0 || i==2 || i==4 || i==6) ? 2 : 3;
268  G4double z = (i < 4) ? -dz : dz;
269  G4double y = -(planes[iy].c*z + planes[iy].d)/planes[iy].b;
270  G4double x = -(planes[ix].b*y + planes[ix].c*z + planes[ix].d)/planes[ix].a;
271  if (x < xmin) xmin = x;
272  if (x > xmax) xmax = x;
273  if (y < ymin) ymin = y;
274  if (y > ymax) ymax = y;
275  }
276 
277  pMin.set(xmin,ymin,-dz);
278  pMax.set(xmax,ymax, dz);
279 
280  // Check correctness of the bounding box
281  //
282  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
283  {
284  std::ostringstream message;
285  message << "Bad bounding box (min >= max) for solid: "
286  << GetName() << " !"
287  << "\npMin = " << pMin
288  << "\npMax = " << pMax;
289  G4Exception("G4UTrap::BoundingLimits()", "GeomMgt0001",
290  JustWarning, message);
291  StreamInfo(G4cout);
292  }
293 
294  // Check consistency of bounding boxes
295  //
296  if (checkBBox)
297  {
298  G4double tolerance = kCarTolerance;
299  U3Vector vmin, vmax;
300  Extent(vmin,vmax);
301  if (std::abs(pMin.x()-vmin.x()) > tolerance ||
302  std::abs(pMin.y()-vmin.y()) > tolerance ||
303  std::abs(pMin.z()-vmin.z()) > tolerance ||
304  std::abs(pMax.x()-vmax.x()) > tolerance ||
305  std::abs(pMax.y()-vmax.y()) > tolerance ||
306  std::abs(pMax.z()-vmax.z()) > tolerance)
307  {
308  std::ostringstream message;
309  message << "Inconsistency in bounding boxes for solid: "
310  << GetName() << " !"
311  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
312  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
313  G4Exception("G4UTrap::BoundingLimits()", "GeomMgt0001",
314  JustWarning, message);
315  checkBBox = false;
316  }
317  }
318 }
319 
321 //
322 // Calculate extent under transform and specified limit
323 
324 G4bool
325 G4UTrap::CalculateExtent(const EAxis pAxis,
326  const G4VoxelLimits& pVoxelLimit,
327  const G4AffineTransform& pTransform,
328  G4double& pMin, G4double& pMax) const
329 {
330  G4ThreeVector bmin, bmax;
331  G4bool exist;
332 
333  // Check bounding box (bbox)
334  //
335  BoundingLimits(bmin,bmax);
336  G4BoundingEnvelope bbox(bmin,bmax);
337 #ifdef G4BBOX_EXTENT
338  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
339 #endif
340  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
341  {
342  return exist = (pMin < pMax) ? true : false;
343  }
344 
345  // Set bounding envelope (benv) and calculate extent
346  //
347  TrapSidePlane planes[4];
348  for (G4int i=0; i<4; ++i) { planes[i] = GetSidePlane(i); }
349 
350  G4ThreeVector pt[8];
351  G4double dz = GetZHalfLength();
352  for (G4int i=0; i<8; ++i)
353  {
354  G4int iy = (i==0 || i==1 || i==4 || i==5) ? 0 : 1;
355  G4int ix = (i==0 || i==2 || i==4 || i==6) ? 2 : 3;
356  G4double z = (i < 4) ? -dz : dz;
357  G4double y = -(planes[iy].c*z + planes[iy].d)/planes[iy].b;
358  G4double x = -(planes[ix].b*y + planes[ix].c*z + planes[ix].d)/planes[ix].a;
359  pt[i].set(x,y,z);
360  }
361 
362  G4ThreeVectorList baseA(4), baseB(4);
363  baseA[0] = pt[0];
364  baseA[1] = pt[1];
365  baseA[2] = pt[3];
366  baseA[3] = pt[2];
367 
368  baseB[0] = pt[4];
369  baseB[1] = pt[5];
370  baseB[2] = pt[7];
371  baseB[3] = pt[6];
372 
373  std::vector<const G4ThreeVectorList *> polygons(2);
374  polygons[0] = &baseA;
375  polygons[1] = &baseB;
376 
377  G4BoundingEnvelope benv(bmin,bmax,polygons);
378  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
379  return exist;
380 }
381 
383 //
384 // Create polyhedron for visualization
385 //
386 G4Polyhedron* G4UTrap::CreatePolyhedron() const
387 {
388  G4double fTthetaSphi = GetThetaSphi();
389  G4double fTthetaCphi = GetThetaCphi();
390  G4double phi = std::atan2(fTthetaSphi, fTthetaCphi);
391  G4double alpha1 = std::atan(GetTanAlpha1());
392  G4double alpha2 = std::atan(GetTanAlpha2());
393  G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi+fTthetaSphi*fTthetaSphi));
394 
395  return new G4PolyhedronTrap(GetZHalfLength(), theta, phi,
396  GetYHalfLength1(),
397  GetXHalfLength1(), GetXHalfLength2(), alpha1,
398  GetYHalfLength2(),
399  GetXHalfLength3(), GetXHalfLength4(), alpha2);
400 }
401 
402 #endif // G4GEOM_USE_USOLIDS
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static const G4double kInfinity
Definition: geomdefs.hh:42
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4double d
Definition: G4Trap.hh:103
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
Double_t z
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
G4double c
Definition: G4Trap.hh:103
double z() const
const G4double kCarTolerance
G4double a
Definition: G4Trap.hh:103
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
TMarker * pt
Definition: egs.C:25
G4double b
Definition: G4Trap.hh:103
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
EAxis
Definition: geomdefs.hh:54
G4GLOB_DLL std::ostream G4cout
double x() const
Char_t n[5]
double y() const