Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/specific/src/G4UGenericTrap.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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