Geant4  10.03
G4UGenericTrap.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 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 #include "G4AutoLock.hh"
46 namespace { G4Mutex UGenericTrapMutex = G4MUTEX_INITIALIZER; }
47 using namespace CLHEP;
48 
50 //
51 // Constructor (generic parameters)
52 //
53 G4UGenericTrap::G4UGenericTrap(const G4String& name, G4double halfZ,
54  const std::vector<G4TwoVector>& vertices)
55  : G4USolid(name, new UGenericTrap())
56 {
57  SetZHalfLength(halfZ);
58  std::vector<UVector2> v;
59  for (size_t n=0; n<vertices.size(); ++n)
60  {
61  v.push_back(UVector2(vertices[n].x(),vertices[n].y()));
62  }
63  GetShape()->SetName(name);
64  GetShape()->Initialise(v);
65 }
66 
67 
69 //
70 // Fake default constructor - sets only member data and allocates memory
71 // for usage restricted to object persistency.
72 //
73 G4UGenericTrap::G4UGenericTrap(__void__& a)
74  : G4USolid(a)
75 {
76 }
77 
78 
80 //
81 // Destructor
82 //
83 G4UGenericTrap::~G4UGenericTrap()
84 {
85 }
86 
87 
89 //
90 // Copy constructor
91 //
92 G4UGenericTrap::G4UGenericTrap(const G4UGenericTrap &source)
93  : G4USolid(source)
94 {
95 }
96 
97 
99 //
100 // Assignment operator
101 //
102 G4UGenericTrap&
103 G4UGenericTrap::operator=(const G4UGenericTrap &source)
104 {
105  if (this == &source) return *this;
106 
107  G4USolid::operator=( source );
108 
109  return *this;
110 }
111 
113 //
114 // Accessors & modifiers
115 //
116 G4double G4UGenericTrap::GetZHalfLength() const
117 {
118  return GetShape()->GetZHalfLength();
119 }
120 G4int G4UGenericTrap::GetNofVertices() const
121 {
122  return GetShape()->GetNofVertices();
123 }
124 G4TwoVector G4UGenericTrap::GetVertex(G4int index) const
125 {
126  UVector2 v = GetShape()->GetVertex(index);
127  return G4TwoVector(v.x,v.y);
128 }
129 const std::vector<G4TwoVector>& G4UGenericTrap::GetVertices() const
130 {
131  G4AutoLock l(&UGenericTrapMutex);
132  std::vector<UVector2> v = GetShape()->GetVertices();
133  static std::vector<G4TwoVector> vertices; vertices.clear();
134  for (size_t n=0; n<v.size(); ++n)
135  {
136  vertices.push_back(G4TwoVector(v[n].x,v[n].y));
137  }
138  return vertices;
139 }
140 G4double G4UGenericTrap::GetTwistAngle(G4int index) const
141 {
142  return GetShape()->GetTwistAngle(index);
143 }
144 G4bool G4UGenericTrap::IsTwisted() const
145 {
146  return GetShape()->IsTwisted();
147 }
148 G4int G4UGenericTrap::GetVisSubdivisions() const
149 {
150  return GetShape()->GetVisSubdivisions();
151 }
152 
153 void G4UGenericTrap::SetVisSubdivisions(G4int subdiv)
154 {
155  GetShape()->SetVisSubdivisions(subdiv);
156 }
157 
158 void G4UGenericTrap::SetZHalfLength(G4double halfZ)
159 {
160  GetShape()->SetZHalfLength(halfZ);
161 }
162 
164 //
165 // Get bounding box
166 
167 void G4UGenericTrap::Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const
168 {
169  UVector3 vmin, vmax;
170  GetShape()->Extent(vmin,vmax);
171  pMin.set(vmin.x(),vmin.y(),vmin.z());
172  pMax.set(vmax.x(),vmax.y(),vmax.z());
173 
174  // Check correctness of the bounding box
175  //
176  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
177  {
178  std::ostringstream message;
179  message << "Bad bounding box (min >= max) for solid: "
180  << GetName() << " !"
181  << "\npMin = " << pMin
182  << "\npMax = " << pMax;
183  G4Exception("G4UGenericTrap::Extent()", "GeomMgt0001", JustWarning, message);
184  StreamInfo(G4cout);
185  }
186 }
187 
189 //
190 // Calculate extent under transform and specified limit
191 
192 G4bool
193 G4UGenericTrap::CalculateExtent(const EAxis pAxis,
194  const G4VoxelLimits& pVoxelLimit,
195  const G4AffineTransform& pTransform,
196  G4double& pMin, G4double& pMax) const
197 {
198  G4ThreeVector bmin, bmax;
199  G4bool exist;
200 
201  // Check bounding box (bbox)
202  //
203  Extent(bmin,bmax);
204  G4BoundingEnvelope bbox(bmin,bmax);
205 #ifdef G4BBOX_EXTENT
206  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
207 #endif
208  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
209  {
210  return exist = (pMin < pMax) ? true : false;
211  }
212 
213  // Set bounding envelope (benv) and calculate extent
214  //
215  // To build the bounding envelope with plane faces each side face of
216  // the trapezoid is subdivided in triangles. Subdivision is done by
217  // duplication of vertices in the bases in a way that the envelope be
218  // a convex polyhedron (some faces of the envelope can be degenerate)
219  //
220  G4double dz = GetZHalfLength();
221  G4ThreeVectorList baseA(8), baseB(8);
222  for (G4int i=0; i<4; ++i)
223  {
224  G4TwoVector va = GetVertex(i);
225  G4TwoVector vb = GetVertex(i+4);
226  baseA[2*i].set(va.x(),va.y(),-dz);
227  baseB[2*i].set(vb.x(),vb.y(), dz);
228  }
229  for (G4int i=0; i<4; ++i)
230  {
231  G4int k1=2*i, k2=(2*i+2)%8;
232  G4double ax = (baseA[k2].x()-baseA[k1].x());
233  G4double ay = (baseA[k2].y()-baseA[k1].y());
234  G4double bx = (baseB[k2].x()-baseB[k1].x());
235  G4double by = (baseB[k2].y()-baseB[k1].y());
236  G4double znorm = ax*by - ay*bx;
237  baseA[k1+1] = (znorm < 0.0) ? baseA[k2] : baseA[k1];
238  baseB[k1+1] = (znorm < 0.0) ? baseB[k1] : baseB[k2];
239  }
240 
241  std::vector<const G4ThreeVectorList *> polygons(2);
242  polygons[0] = &baseA;
243  polygons[1] = &baseB;
244 
245  G4BoundingEnvelope benv(bmin,bmax,polygons);
246  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
247  return exist;
248 }
249 
251 //
252 // CreatePolyhedron()
253 //
254 G4Polyhedron* G4UGenericTrap::CreatePolyhedron() const
255 {
256  // Approximation of Twisted Side
257  // Construct extra Points, if Twisted Side
258  //
259  G4PolyhedronArbitrary* polyhedron;
260  size_t nVertices, nFacets;
261  G4double fDz = GetZHalfLength();
262 
263  G4int subdivisions=0;
264  G4int i;
265  if(IsTwisted())
266  {
267  if ( GetVisSubdivisions()!= 0 )
268  {
269  subdivisions=GetVisSubdivisions();
270  }
271  else
272  {
273  // Estimation of Number of Subdivisions for smooth visualisation
274  //
275  G4double maxTwist=0.;
276  for(i=0; i<4; i++)
277  {
278  if(GetTwistAngle(i)>maxTwist) { maxTwist=GetTwistAngle(i); }
279  }
280 
281  // Computes bounding vectors for the shape
282  //
283  G4double Dx,Dy;
284  UVector3 minBox = GetShape()->GetMinimumBBox();
285  UVector3 maxBox = GetShape()->GetMaximumBBox();
286  G4ThreeVector minVec(minBox.x(), minBox.y(), minBox.z());
287  G4ThreeVector maxVec(maxBox.x(), maxBox.y(), maxBox.z());
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
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
const char * name(G4int ptype)
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4int n
std::vector< G4ThreeVector > G4ThreeVectorList
void AddVertex(const G4ThreeVector &v)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int G4Mutex
Definition: G4Threading.hh:173
EAxis
Definition: geomdefs.hh:54
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
double G4double
Definition: G4Types.hh:76