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