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