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