Geant4  10.03
G4UExtrudedSolid.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 G4UExtrudedSolid wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4ExtrudedSolid.hh"
34 #include "G4UExtrudedSolid.hh"
35 
36 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
37 
38 #include "G4GeomTools.hh"
39 #include "G4AffineTransform.hh"
40 #include "G4BoundingEnvelope.hh"
41 
42 #include "G4PolyhedronArbitrary.hh"
43 
45 //
46 // Constructors
47 //
48 G4UExtrudedSolid::G4UExtrudedSolid(const G4String& name,
49  std::vector<G4TwoVector> polygon,
50  std::vector<ZSection> zsections)
51  : G4USolid(name, new UExtrudedSolid())
52 {
53  GetShape()->SetName(name);
54  std::vector<UVector2> pvec;
55  for (unsigned int i=0; i<polygon.size(); ++i)
56  {
57  pvec.push_back(UVector2(polygon[i].x(), polygon[i].y()));
58  }
59  std::vector<UExtrudedSolid::ZSection> svec;
60  for (unsigned int i=0; i<zsections.size(); ++i)
61  {
62  ZSection sec = zsections[i];
63  svec.push_back(UExtrudedSolid::ZSection(sec.fZ,
64  UVector2(sec.fOffset.x(), sec.fOffset.y()), sec.fScale));
65  }
66  GetShape()->Initialise(pvec, svec);
67 }
68 
69 
70 G4UExtrudedSolid::G4UExtrudedSolid(const G4String& name,
71  std::vector<G4TwoVector> polygon,
72  G4double halfZ,
73  G4TwoVector off1, G4double scale1,
74  G4TwoVector off2, G4double scale2)
75  : G4USolid(name, new UExtrudedSolid())
76 {
77  GetShape()->SetName(name);
78  std::vector<UVector2> pvec;
79  for (unsigned int i=0; i<polygon.size(); ++i)
80  {
81  pvec.push_back(UVector2(polygon[i].x(), polygon[i].y()));
82  }
83  GetShape()->Initialise(pvec, halfZ, UVector2(off1.x(), off1.y()), scale1,
84  UVector2(off2.x(), off2.y()), scale2);
85 }
86 
88 //
89 // Fake default constructor - sets only member data and allocates memory
90 // for usage restricted to object persistency.
91 //
92 G4UExtrudedSolid::G4UExtrudedSolid(__void__& a)
93  : G4USolid(a)
94 {
95 }
96 
97 
99 //
100 // Destructor
101 //
102 G4UExtrudedSolid::~G4UExtrudedSolid()
103 {
104 }
105 
106 
108 //
109 // Copy constructor
110 //
111 G4UExtrudedSolid::G4UExtrudedSolid(const G4UExtrudedSolid &source)
112  : G4USolid(source)
113 {
114 }
115 
116 
118 //
119 // Assignment operator
120 //
121 G4UExtrudedSolid&
122 G4UExtrudedSolid::operator=(const G4UExtrudedSolid &source)
123 {
124  if (this == &source) return *this;
125 
126  G4USolid::operator=( source );
127 
128  return *this;
129 }
130 
131 
133 //
134 // Accessors
135 
136 G4int G4UExtrudedSolid::GetNofVertices() const
137 {
138  return GetShape()->GetNofVertices();
139 }
140 G4TwoVector G4UExtrudedSolid::GetVertex(G4int i) const
141 {
142  UVector2 v = GetShape()->GetVertex(i);
143  return G4TwoVector(v.x, v.y);
144 }
145 std::vector<G4TwoVector> G4UExtrudedSolid::GetPolygon() const
146 {
147  std::vector<UVector2> pol = GetShape()->GetPolygon();
148  std::vector<G4TwoVector> v;
149  for (unsigned int i=0; i<pol.size(); ++i)
150  {
151  v.push_back(G4TwoVector(pol[i].x, pol[i].y));
152  }
153  return v;
154 }
155 G4int G4UExtrudedSolid::GetNofZSections() const
156 {
157  return GetShape()->GetNofZSections();
158 }
159 G4UExtrudedSolid::ZSection G4UExtrudedSolid::GetZSection(G4int i) const
160 {
161  return ZSection(GetShape()->GetZSection(i));
162 }
163 std::vector<G4UExtrudedSolid::ZSection> G4UExtrudedSolid::GetZSections() const
164 {
165  std::vector<UExtrudedSolid::ZSection> sv = GetShape()->GetZSections();
166  std::vector<G4UExtrudedSolid::ZSection> vec;
167  for (unsigned int i=0; i<sv.size(); ++i)
168  {
169  vec.push_back(ZSection(sv[i]));
170  }
171  return vec;
172 }
173 
174 
176 //
177 // Get bounding box
178 
179 void G4UExtrudedSolid::Extent(G4ThreeVector& pMin, G4ThreeVector& pMax) const
180 {
181  static G4bool checkBBox = true;
182 
183  G4double xmin0 = kInfinity, xmax0 = -kInfinity;
184  G4double ymin0 = kInfinity, ymax0 = -kInfinity;
185 
186  for (G4int i=0; i<GetNofVertices(); ++i)
187  {
188  G4TwoVector vertex = GetVertex(i);
189  G4double x = vertex.x();
190  if (x < xmin0) xmin0 = x;
191  if (x > xmax0) xmax0 = x;
192  G4double y = vertex.y();
193  if (y < ymin0) ymin0 = y;
194  if (y > ymax0) ymax0 = y;
195  }
196 
197  G4double xmin = kInfinity, xmax = -kInfinity;
198  G4double ymin = kInfinity, ymax = -kInfinity;
199 
200  G4int nsect = GetNofZSections();
201  for (G4int i=0; i<nsect; ++i)
202  {
203  ZSection zsect = GetZSection(i);
204  G4double dx = zsect.fOffset.x();
205  G4double dy = zsect.fOffset.y();
206  G4double scale = zsect.fScale;
207  xmin = std::min(xmin,xmin0*scale+dx);
208  xmax = std::max(xmax,xmax0*scale+dx);
209  ymin = std::min(ymin,ymin0*scale+dy);
210  ymax = std::max(ymax,ymax0*scale+dy);
211  }
212 
213  G4double zmin = GetZSection(0).fZ;
214  G4double zmax = GetZSection(nsect-1).fZ;
215 
216  pMin.set(xmin,ymin,zmin);
217  pMax.set(xmax,ymax,zmax);
218 
219  // Check correctness of the bounding box
220  //
221  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
222  {
223  std::ostringstream message;
224  message << "Bad bounding box (min >= max) for solid: "
225  << GetName() << " !"
226  << "\npMin = " << pMin
227  << "\npMax = " << pMax;
228  G4Exception("G4UExtrudedSolid::Extent()", "GeomMgt0001",
229  JustWarning, message);
230  StreamInfo(G4cout);
231  }
232 
233  // Check consistency of bounding boxes
234  //
235  if (checkBBox)
236  {
237  UVector3 vmin, vmax;
238  GetShape()->Extent(vmin,vmax);
239  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
240  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
241  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
242  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
243  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
244  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
245  {
246  std::ostringstream message;
247  message << "Inconsistency in bounding boxes for solid: "
248  << GetName() << " !"
249  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
250  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
251  G4Exception("G4UExtrudedSolid::Extent()", "GeomMgt0001", JustWarning, message);
252  checkBBox = false;
253  }
254  }
255 }
256 
257 
259 //
260 // Calculate extent under transform and specified limit
261 
262 G4bool
264  const G4VoxelLimits& pVoxelLimit,
265  const G4AffineTransform& pTransform,
266  G4double& pMin, G4double& pMax) const
267 {
268  G4ThreeVector bmin, bmax;
269  G4bool exist;
270 
271  // Check bounding box (bbox)
272  //
273  Extent(bmin,bmax);
274  G4BoundingEnvelope bbox(bmin,bmax);
275 #ifdef G4BBOX_EXTENT
276  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
277 #endif
278  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
279  {
280  return exist = (pMin < pMax) ? true : false;
281  }
282 
283  // To find the extent, the base polygon is subdivided in triangles.
284  // The extent is calculated as cumulative extent of the parts
285  // formed by extrusion of the triangles
286  //
287  G4TwoVectorList basePolygon = GetPolygon();
288  G4TwoVectorList triangles;
289  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
290  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
291 
292  // triangulate the base polygon
293  if (!G4GeomTools::TriangulatePolygon(basePolygon,triangles))
294  {
295  std::ostringstream message;
296  message << "Triangulation of the base polygon has failed for solid: "
297  << GetName() << " !"
298  << "\nExtent has been calculated using boundary box";
299  G4Exception("G4UExtrudedSolid::CalculateExtent()",
300  "GeomMgt1002",JustWarning,message);
301  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
302  }
303 
304  // allocate vector lists
305  G4int nsect = GetNofZSections();
306  std::vector<const G4ThreeVectorList *> polygons;
307  polygons.resize(nsect);
308  for (G4int k=0; k<nsect; ++k) { polygons[k] = new G4ThreeVectorList(3); }
309 
310  // main loop along triangles
311  pMin = kInfinity;
312  pMax = -kInfinity;
313  G4int ntria = triangles.size()/3;
314  for (G4int i=0; i<ntria; ++i)
315  {
316  G4int i3 = i*3;
317  for (G4int k=0; k<nsect; ++k) // extrude triangle
318  {
319  ZSection zsect = GetZSection(k);
320  G4double z = zsect.fZ;
321  G4double dx = zsect.fOffset.x();
322  G4double dy = zsect.fOffset.y();
323  G4double scale = zsect.fScale;
324 
325  G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
326  G4ThreeVectorList::iterator iter = ptr->begin();
327  G4double x0 = triangles[i3+0].x()*scale+dx;
328  G4double y0 = triangles[i3+0].y()*scale+dy;
329  iter->set(x0,y0,z);
330  iter++;
331  G4double x1 = triangles[i3+1].x()*scale+dx;
332  G4double y1 = triangles[i3+1].y()*scale+dy;
333  iter->set(x1,y1,z);
334  iter++;
335  G4double x2 = triangles[i3+2].x()*scale+dx;
336  G4double y2 = triangles[i3+2].y()*scale+dy;
337  iter->set(x2,y2,z);
338  }
339 
340  // set sub-envelope and adjust extent
341  G4double emin,emax;
342  G4BoundingEnvelope benv(polygons);
343  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
344  if (emin < pMin) pMin = emin;
345  if (emax > pMax) pMax = emax;
346  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
347  }
348  // free memory
349  for (G4int k=0; k<nsect; ++k) { delete polygons[k]; polygons[k]=0;}
350  return (pMin < pMax);
351 }
352 
353 
355 //
356 // CreatePolyhedron()
357 //
358 G4Polyhedron* G4UExtrudedSolid::CreatePolyhedron () const
359 {
360  G4int nFacets = GetShape()->GetNumberOfFacets();
361  G4int nVertices = 0;
362  for (G4int l = 0; l<nFacets; ++l) // compute total number of vertices first
363  {
364  VUFacet* facet = GetShape()->GetFacet(l);
365  G4int n = facet->GetNumberOfVertices();
366  nVertices += n;
367  }
368 
369  G4PolyhedronArbitrary *polyhedron =
370  new G4PolyhedronArbitrary (nVertices,nFacets);
371 
372  for (G4int i = 0; i<nFacets; ++i)
373  {
374  VUFacet* facet = GetShape()->GetFacet(i);
375  G4int v[4];
376  G4int n = facet->GetNumberOfVertices();
377  for (G4int m = 0; m<n; ++m)
378  {
379  UVector3 vtx = facet->GetVertex(m);
380  polyhedron->AddVertex(G4ThreeVector(vtx.x(), vtx.y(), vtx.z()));
381  }
382  if (n > 4) n = 4;
383  else if (n == 3) v[3] = 0;
384  for (G4int j=0; j<n; ++j)
385  {
386  G4int k = facet->GetVertexIndex(j);
387  v[j] = k+1;
388  }
389  polyhedron->AddFacet(v[0],v[1],v[2],v[3]);
390  }
391  polyhedron->SetReferences();
392 
393  return (G4Polyhedron*) polyhedron;
394 }
395 
396 #endif // G4GEOM_USE_USOLIDS
G4String GetName() const
std::vector< G4TwoVector > GetPolygon() const
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:50
const char * name(G4int ptype)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
bool G4bool
Definition: G4Types.hh:79
const G4double kCarTolerance
const G4int n
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:178
std::vector< G4ThreeVector > G4ThreeVectorList
void AddVertex(const G4ThreeVector &v)
G4int GetNofZSections() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EAxis
Definition: geomdefs.hh:54
ZSection GetZSection(G4int index) const
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
double G4double
Definition: G4Types.hh:76
G4double GetMaxExtent(const EAxis pAxis) const
G4double GetMinExtent(const EAxis pAxis) const