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