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
G4Polyhedron.hh
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
28#ifndef G4POLYHEDRON_HH
29#define G4POLYHEDRON_HH
30
31// Class Description:
32// G4Polyhedron is an intermediate class between G4 and visualization
33// systems. It is intended to provide some service like:
34// - polygonization of the G4 shapes with triangulization
35// (quadrilaterization) of complex polygons;
36// - calculation of normals for faces and vertices.
37//
38// Inherits from HepPolyhedron, to which reference should be made for
39// functionality.
40//
41// Public constructors:
42// G4PolyhedronBox(dx,dy,dz) - create G4Polyhedron for G4 Box;
43// G4PolyhedronTrd1(dx1,dx2,dy,dz) - create G4Polyhedron for G4 Trd1;
44// G4PolyhedronTrd2(dx1,dx2,dy1,dy2,dz) - create G4Polyhedron for G4 Trd2;
45// G4PolyhedronTrap(dz,theta,phi,
46// h1,bl1,tl1,alp1,
47// h2,bl2,tl2,alp2) - create G4Polyhedron for G4 Trap;
48// G4PolyhedronPara(dx,dy,dz,
49// alpha,theta,phi) - create G4Polyhedron for G4 Para;
50//
51// G4PolyhedronTube(rmin,rmax,dz) - create G4Polyhedron for G4 Tube;
52// G4PolyhedronTubs(rmin,rmax,dz,
53// phi1,dphi) - create G4Polyhedron for G4 Tubs;
54// G4PolyhedronCone(rmin1,rmax1,
55// rmin2,rmax2,dz) - create G4Polyhedron for G4 Cone;
56// G4PolyhedronCons(rmin1,rmax1,
57// rmin2,rmax2,dz,
58// phi1,dphi) - create G4Polyhedron for G4 Cons;
59//
60// G4PolyhedronPgon(phi,dphi,npdv,nz,
61// z(*),rmin(*),rmax(*)) - create G4Polyhedron for G4 Pgon;
62// G4PolyhedronPcon(phi,dphi,nz,
63// z(*),rmin(*),rmax(*)) - create G4Polyhedron for G4 Pcon;
64//
65// G4PolyhedronSphere(rmin,rmax,
66// phi,dphi,the,dthe) - create G4Polyhedron for Sphere;
67// G4PolyhedronTorus(rmin,rmax,rtor,
68// phi,dphi) - create G4Polyhedron for Torus;
69// G4PolyhedronTet(p0[3],p1[3],p2[3],p3[3]) - create polyhedron for Tet;
70//
71// G4PolyhedronEllipsoid(dx,dy,dz,
72// zcut1,zcut2) - create G4Polyhedron for Ellipsoid;
73// G4PolyhedronEllipticalCone(dx,dy,z,
74// zcut1) - create polyhedron for Elliptical cone;
75// G4PolyhedronParaboloid(r1,r2,dz,
76// phi,dphi) - create polyhedron for Paraboloid;
77// G4PolyhedronHype(r1,r2,
78// tan1,tan2,halfz) - create polyhedron for Hype;
79// G4PolyhedronHyperbolicMirror(a,h,r) - create polyhedron for Hyperbolic mirror;
80//
81// G4PolyhedronTetMesh(vector<p>) - create polyhedron for tetrahedron mesh;
82//
83// G4PolyhedronBoxMesh(sx,sy,sz,vector<p>) - create polyhedron for box mesh;
84//
85// Public functions inherited from HepPolyhedron (this list might be
86// incomplete):
87// GetNoVertices() - returns number of vertices
88// GetNoFacets() - returns number of faces
89// GetNextVertexIndex(index, edgeFlag) - get vertex indeces of the
90// quadrilaterals in order; returns false when
91// finished each face;
92// GetVertex(index) - returns vertex by index;
93// GetNextVertex(vertex, edgeFlag) - get vertices with edge visibility
94// of the quadrilaterals in order;
95// returns false when finished each face;
96// GetNextVertex(vertex, edgeFlag, normal) - get vertices with edge
97// visibility and normal of the quadrilaterals
98// in order; returns false when finished each face;
99// GetNextNormal(normal) - get normals of each face in order;
100// returns false when finished all faces;
101// GetNextUnitNormal(normal) - get normals of unit length of each face
102// in order; returns false when finished all faces;
103// GetNextEdgeIndeces(i1, i2, edgeFlag) - get indeces of the next edge;
104// returns false for the last edge;
105// GetNextEdge(p1, p2, edgeFlag) - get next edge;
106// returns false for the last edge;
107// SetVertex(index, v) - set vertex;
108// SetFacet(index,iv1,iv2,iv3,iv4) - set facet;
109// SetReferences() - set references to neighbouring facets;
110// JoinCoplanarFacets(tol) - join couples of triangles into quadrilaterals;
111// InvertFacets() - invert the order on nodes in facets;
112// SetNumberOfRotationSteps(G4int n) - set number of steps for whole circle;
113//
114// History:
115// 21st February 2000 Evgeni Chernaev, John Allison
116// - Re-written to inherit HepPolyhedron.
117//
118// 11.03.05 J.Allison
119// - Added fNumberOfRotationStepsAtTimeOfCreation and access method.
120// (NumberOfRotationSteps is also called number of sides per circle or
121// line segments per circle - see
122// /vis/viewer/set/lineSegmentsPerCircle.)
123// 20.06.05 G.Cosmo
124// - Added G4PolyhedronEllipsoid.
125// 09.03.06 J.Allison
126// - Added operator<<.
127
128#include "globals.hh"
129#include "HepPolyhedron.h"
130#include "G4Visible.hh"
131
132class G4Polyhedron : public HepPolyhedron, public G4Visible {
133public:
134 // Default constructor
136
137 // Constructors
138 G4Polyhedron (G4int Nvert, G4int Nface);
139 G4Polyhedron (const HepPolyhedron& from);
140
141 // Copy and move constructors
142 G4Polyhedron (const G4Polyhedron& from) = default;
143 G4Polyhedron (G4Polyhedron&& from) = default;
144
145 // Assignment and move assignment
146 G4Polyhedron & operator=(const G4Polyhedron & from) = default;
147 G4Polyhedron & operator=(G4Polyhedron && from) = default;
148
149 // Destructor
150 ~G4Polyhedron () override;
151
153 return fNumberOfRotationStepsAtTimeOfCreation;
154 }
155private:
156 G4int fNumberOfRotationStepsAtTimeOfCreation = fNumberOfRotationSteps;
157};
158
160public:
162 ~G4PolyhedronBox () override;
163};
164
166public:
168 G4double Rmn2, G4double Rmx2, G4double Dz);
169 ~G4PolyhedronCone () override;
170};
171
173public:
175 G4double Rmn2, G4double Rmx2, G4double Dz,
176 G4double Phi1, G4double Dphi);
177 ~G4PolyhedronCons () override;
178};
179
181public:
183 G4double Alpha, G4double Theta, G4double Phi);
184 ~G4PolyhedronPara () override;
185};
186
188public:
190 const G4double *z,
191 const G4double *rmin,
192 const G4double *rmax);
194 const std::vector<G4TwoVector> &rz);
195 ~G4PolyhedronPcon () override;
196};
197
199public:
200 G4PolyhedronPgon (G4double phi, G4double dphi, G4int npdv, G4int nz,
201 const G4double *z,
202 const G4double *rmin,
203 const G4double *rmax);
204 G4PolyhedronPgon (G4double phi, G4double dphi, G4int npdv,
205 const std::vector<G4TwoVector> &rz);
206
207 ~G4PolyhedronPgon () override;
208};
209
211public:
213 G4double phi, G4double dphi,
214 G4double the, G4double dthe);
216};
217
219public:
220 G4PolyhedronTet (const G4double p0[3],
221 const G4double p1[3],
222 const G4double p2[3],
223 const G4double p3[3]);
224 ~G4PolyhedronTet () override;
225};
226
228public:
230 G4double phi, G4double dphi);
232};
233
235public:
237 G4double Dy1,
238 G4double Dx1, G4double Dx2, G4double Alp1,
239 G4double Dy2,
240 G4double Dx3, G4double Dx4, G4double Alp2);
241 ~G4PolyhedronTrap () override;
242};
243
245public:
247 G4double Dy, G4double Dz);
248 ~G4PolyhedronTrd1 () override;
249};
250
252public:
254 G4double Dy1, G4double Dy2, G4double Dz);
255 ~G4PolyhedronTrd2 () override;
256};
257
259public:
261 ~G4PolyhedronTube () override;
262};
263
265public:
267 G4double Phi1, G4double Dphi);
268 ~G4PolyhedronTubs () override;
269};
270
272 public:
274 G4double sPhi, G4double dPhi);
276};
277
279 public:
281 G4double tan2, G4double halfZ);
282 ~G4PolyhedronHype () override;
283};
284
286 public:
288 G4double zcut1, G4double zcut2);
290};
291
293 public:
295 G4double zcut1);
297};
298
300 public:
303};
304
306 public:
307 G4PolyhedronTetMesh(const std::vector<G4ThreeVector>& tetrahedra);
309};
310
312 public:
313 G4PolyhedronBoxMesh(G4double sizeX, G4double sizeY, G4double sizeZ,
314 const std::vector<G4ThreeVector>& positions);
315
317};
318
319std::ostream& operator<<(std::ostream& os, const G4Polyhedron&);
320
321#endif /* G4POLYHEDRON_HH */
std::ostream & operator<<(std::ostream &os, const G4Polyhedron &)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
~G4PolyhedronBoxMesh() override
~G4PolyhedronBox() override
~G4PolyhedronCone() override
~G4PolyhedronCons() override
~G4PolyhedronEllipsoid() override
~G4PolyhedronEllipticalCone() override
~G4PolyhedronHype() override
~G4PolyhedronHyperbolicMirror() override
~G4PolyhedronPara() override
~G4PolyhedronParaboloid() override
~G4PolyhedronPcon() override
~G4PolyhedronPgon() override
~G4PolyhedronSphere() override
~G4PolyhedronTetMesh() override
~G4PolyhedronTet() override
~G4PolyhedronTorus() override
~G4PolyhedronTrap() override
~G4PolyhedronTrd1() override
~G4PolyhedronTrd2() override
~G4PolyhedronTube() override
~G4PolyhedronTubs() override
~G4Polyhedron() override
G4Polyhedron & operator=(G4Polyhedron &&from)=default
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4Polyhedron & operator=(const G4Polyhedron &from)=default
G4Polyhedron(G4Polyhedron &&from)=default
G4Polyhedron(const G4Polyhedron &from)=default
static G4ThreadLocal G4int fNumberOfRotationSteps