Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UGenericTrap.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 of G4UGenericTrap wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4GenericTrap.hh"
34 #include "G4UGenericTrap.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 #include "G4Polyhedron.hh"
43 #include "G4PolyhedronArbitrary.hh"
44 
45 using namespace CLHEP;
46 
48 //
49 // Constructor (generic parameters)
50 //
51 G4UGenericTrap::G4UGenericTrap(const G4String& name, G4double halfZ,
52  const std::vector<G4TwoVector>& vertices)
53  : Base_t(name), fVisSubdivisions(0)
54 {
55  SetZHalfLength(halfZ);
56  Initialise(vertices);
57 }
58 
59 
61 //
62 // Fake default constructor - sets only member data and allocates memory
63 // for usage restricted to object persistency.
64 //
65 G4UGenericTrap::G4UGenericTrap(__void__& a)
66  : Base_t(a), fVisSubdivisions(0), fVertices()
67 {
68 }
69 
70 
72 //
73 // Destructor
74 //
75 G4UGenericTrap::~G4UGenericTrap()
76 {
77 }
78 
79 
81 //
82 // Copy constructor
83 //
84 G4UGenericTrap::G4UGenericTrap(const G4UGenericTrap &source)
85  : Base_t(source), fVisSubdivisions(source.fVisSubdivisions),
86  fVertices(source.fVertices)
87 
88 {
89 }
90 
91 
93 //
94 // Assignment operator
95 //
96 G4UGenericTrap&
97 G4UGenericTrap::operator=(const G4UGenericTrap &source)
98 {
99  if (this == &source) return *this;
100 
101  Base_t::operator=( source );
102  fVertices = source.fVertices;
103  fVisSubdivisions = source.fVisSubdivisions;
104 
105  return *this;
106 }
107 
109 //
110 // Accessors & modifiers
111 //
112 G4double G4UGenericTrap::GetZHalfLength() const
113 {
114  return GetDZ();
115 }
116 G4int G4UGenericTrap::GetNofVertices() const
117 {
118  return fVertices.size();
119 }
120 G4TwoVector G4UGenericTrap::GetVertex(G4int index) const
121 {
122  return G4TwoVector(GetVerticesX()[index], GetVerticesY()[index]);
123 }
124 const std::vector<G4TwoVector>& G4UGenericTrap::GetVertices() const
125 {
126  return fVertices;
127 }
128 G4double G4UGenericTrap::GetTwistAngle(G4int index) const
129 {
130  return GetTwist(index);
131 }
132 G4bool G4UGenericTrap::IsTwisted() const
133 {
134  return !IsPlanar();
135 }
136 G4int G4UGenericTrap::GetVisSubdivisions() const
137 {
138  return fVisSubdivisions;
139 }
140 
141 void G4UGenericTrap::SetVisSubdivisions(G4int subdiv)
142 {
143  fVisSubdivisions = subdiv;
144 }
145 
146 void G4UGenericTrap::SetZHalfLength(G4double halfZ)
147 {
148  SetDZ(halfZ);
149 }
150 
151 void G4UGenericTrap::Initialise(const std::vector<G4TwoVector>& v)
152 {
153  G4double verticesx[8], verticesy[8];
154  for (G4int i=0; i<8; ++i)
155  {
156  fVertices.push_back(v[i]);
157  verticesx[i] = v[i].x();
158  verticesy[i] = v[i].y();
159  }
160  Initialize(verticesx, verticesy, GetZHalfLength());
161 }
162 
164 //
165 // Get bounding box
166 
167 void G4UGenericTrap::BoundingLimits(G4ThreeVector& pMin,
168  G4ThreeVector& pMax) const
169 {
170  U3Vector vmin, vmax;
171  Extent(vmin,vmax);
172  pMin.set(vmin.x(),vmin.y(),vmin.z());
173  pMax.set(vmax.x(),vmax.y(),vmax.z());
174 
175  // Check correctness of the bounding box
176  //
177  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
178  {
179  std::ostringstream message;
180  message << "Bad bounding box (min >= max) for solid: "
181  << GetName() << " !"
182  << "\npMin = " << pMin
183  << "\npMax = " << pMax;
184  G4Exception("G4UGenericTrap::BoundingLimits()", "GeomMgt0001",
185  JustWarning, message);
186  StreamInfo(G4cout);
187  }
188 }
189 
191 //
192 // Calculate extent under transform and specified limit
193 
194 G4bool
195 G4UGenericTrap::CalculateExtent(const EAxis pAxis,
196  const G4VoxelLimits& pVoxelLimit,
197  const G4AffineTransform& pTransform,
198  G4double& pMin, G4double& pMax) const
199 {
200  G4ThreeVector bmin, bmax;
201  G4bool exist;
202 
203  // Check bounding box (bbox)
204  //
205  BoundingLimits(bmin,bmax);
206  G4BoundingEnvelope bbox(bmin,bmax);
207 #ifdef G4BBOX_EXTENT
208  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
209 #endif
210  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
211  {
212  return exist = (pMin < pMax) ? true : false;
213  }
214 
215  // Set bounding envelope (benv) and calculate extent
216  //
217  // To build the bounding envelope with plane faces each side face of
218  // the trapezoid is subdivided in triangles. Subdivision is done by
219  // duplication of vertices in the bases in a way that the envelope be
220  // a convex polyhedron (some faces of the envelope can be degenerate)
221  //
222  G4double dz = GetZHalfLength();
223  G4ThreeVectorList baseA(8), baseB(8);
224  for (G4int i=0; i<4; ++i)
225  {
226  G4TwoVector va = GetVertex(i);
227  G4TwoVector vb = GetVertex(i+4);
228  baseA[2*i].set(va.x(),va.y(),-dz);
229  baseB[2*i].set(vb.x(),vb.y(), dz);
230  }
231  for (G4int i=0; i<4; ++i)
232  {
233  G4int k1=2*i, k2=(2*i+2)%8;
234  G4double ax = (baseA[k2].x()-baseA[k1].x());
235  G4double ay = (baseA[k2].y()-baseA[k1].y());
236  G4double bx = (baseB[k2].x()-baseB[k1].x());
237  G4double by = (baseB[k2].y()-baseB[k1].y());
238  G4double znorm = ax*by - ay*bx;
239  baseA[k1+1] = (znorm < 0.0) ? baseA[k2] : baseA[k1];
240  baseB[k1+1] = (znorm < 0.0) ? baseB[k1] : baseB[k2];
241  }
242 
243  std::vector<const G4ThreeVectorList *> polygons(2);
244  polygons[0] = &baseA;
245  polygons[1] = &baseB;
246 
247  G4BoundingEnvelope benv(bmin,bmax,polygons);
248  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
249  return exist;
250 }
251 
253 //
254 // CreatePolyhedron()
255 //
256 G4Polyhedron* G4UGenericTrap::CreatePolyhedron() const
257 {
258  // Approximation of Twisted Side
259  // Construct extra Points, if Twisted Side
260  //
261  G4PolyhedronArbitrary* polyhedron;
262  size_t nVertices, nFacets;
263  G4double fDz = GetZHalfLength();
264 
265  G4int subdivisions=0;
266  G4int i;
267  if(IsTwisted())
268  {
269  if ( GetVisSubdivisions()!= 0 )
270  {
271  subdivisions=GetVisSubdivisions();
272  }
273  else
274  {
275  // Estimation of Number of Subdivisions for smooth visualisation
276  //
277  G4double maxTwist=0.;
278  for(i=0; i<4; i++)
279  {
280  if(GetTwistAngle(i)>maxTwist) { maxTwist=GetTwistAngle(i); }
281  }
282 
283  // Computes bounding vectors for the shape
284  //
285  G4double Dx,Dy;
286  G4ThreeVector minVec, maxVec;
287  BoundingLimits(minVec,maxVec);
288  Dx = 0.5*(maxVec.x()- minVec.y());
289  Dy = 0.5*(maxVec.y()- minVec.y());
290  if (Dy > Dx) { Dx=Dy; }
291 
292  subdivisions=8*G4int(maxTwist/(Dx*Dx*Dx)*fDz);
293  if (subdivisions<4) { subdivisions=4; }
294  if (subdivisions>30) { subdivisions=30; }
295  }
296  }
297  G4int sub4=4*subdivisions;
298  nVertices = 8+subdivisions*4;
299  nFacets = 6+subdivisions*4;
300  G4double cf=1./(subdivisions+1);
301  polyhedron = new G4PolyhedronArbitrary (nVertices, nFacets);
302 
303  // Add Vertex
304  //
305  for (i=0;i<4;i++)
306  {
307  polyhedron->AddVertex(G4ThreeVector(GetVertex(i).x(),
308  GetVertex(i).y(),-fDz));
309  }
310  for( i=0;i<subdivisions;i++)
311  {
312  for(G4int j=0;j<4;j++)
313  {
314  G4TwoVector u=GetVertex(j)+cf*(i+1)*( GetVertex(j+4)-GetVertex(j));
315  polyhedron->AddVertex(G4ThreeVector(u.x(),u.y(),-fDz+cf*2*fDz*(i+1)));
316  }
317  }
318  for (i=4;i<8;i++)
319  {
320  polyhedron->AddVertex(G4ThreeVector(GetVertex(i).x(),
321  GetVertex(i).y(),fDz));
322  }
323 
324  // Add Facets
325  //
326  polyhedron->AddFacet(1,4,3,2); //Z-plane
327  for (i=0;i<subdivisions+1;i++)
328  {
329  G4int is=i*4;
330  polyhedron->AddFacet(5+is,8+is,4+is,1+is);
331  polyhedron->AddFacet(8+is,7+is,3+is,4+is);
332  polyhedron->AddFacet(7+is,6+is,2+is,3+is);
333  polyhedron->AddFacet(6+is,5+is,1+is,2+is);
334  }
335  polyhedron->AddFacet(5+sub4,6+sub4,7+sub4,8+sub4); //Z-plane
336 
337  polyhedron->SetReferences();
338  polyhedron->InvertFacets();
339 
340  return (G4Polyhedron*) polyhedron;
341 }
342 
343 #endif // G4GEOM_USE_USOLIDS
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
const XML_Char * name
Definition: expat.h:151
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
void AddVertex(const G4ThreeVector &v)
Float_t y
Definition: compare.C:6
void Initialize()
Definition: errprop.cc:101
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
double z() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
double x() const
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
double y() const
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const