Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// $Id$
28
29#ifndef G4POLYHEDRON_HH
30#define G4POLYHEDRON_HH
31
32// Class Description:
33// G4Polyhedron is an intermediate class between G4 and visualization
34// systems. It is intended to provide some service like:
35// - polygonization of the G4 shapes with triangulization
36// (quadrilaterization) of complex polygons;
37// - calculation of normals for faces and vertices.
38//
39// Inherits from HepPolyhedron, to which reference should be made for
40// functionality.
41//
42// Public constructors:
43// G4PolyhedronBox(dx,dy,dz) - create G4Polyhedron for G4 Box;
44// G4PolyhedronTrd1(dx1,dx2,dy,dz) - create G4Polyhedron for G4 Trd1;
45// G4PolyhedronTrd2(dx1,dx2,dy1,dy2,dz) - create G4Polyhedron for G4 Trd2;
46// G4PolyhedronTrap(dz,theta,phi,
47// h1,bl1,tl1,alp1,
48// h2,bl2,tl2,alp2) - create G4Polyhedron for G4 Trap;
49// G4PolyhedronPara(dx,dy,dz,
50// alpha,theta,phi) - create G4Polyhedron for G4 Para;
51//
52// G4PolyhedronTube(rmin,rmax,dz) - create G4Polyhedron for G4 Tube;
53// G4PolyhedronTubs(rmin,rmax,dz,
54// phi1,dphi) - create G4Polyhedron for G4 Tubs;
55// G4PolyhedronCone(rmin1,rmax1,
56// rmin2,rmax2,dz) - create G4Polyhedron for G4 Cone;
57// G4PolyhedronCons(rmin1,rmax1,
58// rmin2,rmax2,dz,
59// phi1,dphi) - create G4Polyhedron for G4 Cons;
60//
61// G4PolyhedronPgon(phi,dphi,npdv,nz,
62// z(*),rmin(*),rmax(*)) - create G4Polyhedron for G4 Pgon;
63// G4PolyhedronPcon(phi,dphi,nz,
64// z(*),rmin(*),rmax(*)) - create G4Polyhedron for G4 Pcon;
65//
66// G4PolyhedronSphere(rmin,rmax,
67// phi,dphi,the,dthe) - create G4Polyhedron for Sphere;
68// G4PolyhedronTorus(rmin,rmax,rtor,
69// phi,dphi) - create G4Polyhedron for Torus;
70// G4PolyhedronEllipsoid(dx,dy,dz,
71// zcut1,zcut2) - create G4Polyhedron for Ellipsoid;
72//
73// Public functions inherited from HepPolyhedron (this list might be
74// incomplete):
75// GetNoVertices() - returns number of vertices
76// GetNoFacets() - returns number of faces
77// GetNextVertexIndex(index, edgeFlag) - get vertex indeces of the
78// quadrilaterals in order; returns false when
79// finished each face;
80// GetVertex(index) - returns vertex by index;
81// GetNextVertex(vertex, edgeFlag) - get vertices with edge visibility
82// of the quadrilaterals in order;
83// returns false when finished each face;
84// GetNextVertex(vertex, edgeFlag, normal) - get vertices with edge
85// visibility and normal of the quadrilaterals
86// in order; returns false when finished each face;
87// GetNextNormal(normal) - get normals of each face in order;
88// returns false when finished all faces;
89// GetNextUnitNormal(normal) - get normals of unit length of each face
90// in order; returns false when finished all faces;
91// GetNextEdgeIndeces(i1, i2, edgeFlag) - get indeces of the next edge;
92// returns false for the last edge;
93// GetNextEdge(p1, p2, edgeFlag) - get next edge;
94// returns false for the last edge;
95// SetNumberOfRotationSteps(G4int n) - Set number of steps for whole circle;
96
97// History:
98// 21st February 2000 Evgeni Chernaev, John Allison
99// - Re-written to inherit HepPolyhedron.
100//
101// 11.03.05 J.Allison
102// - Added fNumberOfRotationStepsAtTimeOfCreation and access method.
103// (NumberOfRotationSteps is also called number of sides per circle or
104// line segments per circle - see
105// /vis/viewer/set/lineSegmentsPerCircle.)
106// 20.06.05 G.Cosmo
107// - Added G4PolyhedronEllipsoid.
108// 09.03.06 J.Allison
109// - Added operator<<.
110
111#include "globals.hh"
112#include "HepPolyhedron.h"
113#include "G4Visible.hh"
114
115class G4Polyhedron : public HepPolyhedron, public G4Visible {
116public:
117 G4Polyhedron ();
118 G4Polyhedron (const HepPolyhedron& from);
119 // Use compiler defaults for copy contructor and assignment. (They
120 // invoke their counterparts in HepPolyhedron and G4Visible.)
121 virtual ~G4Polyhedron ();
122
124 return fNumberOfRotationStepsAtTimeOfCreation;
125 }
126private:
127 G4int fNumberOfRotationStepsAtTimeOfCreation;
128};
129
131public:
133 virtual ~G4PolyhedronBox ();
134};
135
137public:
139 G4double Rmn2, G4double Rmx2, G4double Dz);
140 virtual ~G4PolyhedronCone ();
141};
142
144public:
146 G4double Rmn2, G4double Rmx2, G4double Dz,
147 G4double Phi1, G4double Dphi);
148 virtual ~G4PolyhedronCons ();
149};
150
152public:
154 G4double Alpha, G4double Theta, G4double Phi);
155 virtual ~G4PolyhedronPara ();
156};
157
159public:
161 const G4double *z,
162 const G4double *rmin,
163 const G4double *rmax);
164 virtual ~G4PolyhedronPcon ();
165};
166
168public:
169 G4PolyhedronPgon (G4double phi, G4double dphi, G4int npdv, G4int nz,
170 const G4double *z,
171 const G4double *rmin,
172 const G4double *rmax);
173 virtual ~G4PolyhedronPgon ();
174};
175
177public:
179 G4double phi, G4double dphi,
180 G4double the, G4double dthe);
181 virtual ~G4PolyhedronSphere ();
182};
183
185public:
187 G4double phi, G4double dphi);
188 virtual ~G4PolyhedronTorus ();
189};
190
192public:
194 G4double Dy1,
195 G4double Dx1, G4double Dx2, G4double Alp1,
196 G4double Dy2,
197 G4double Dx3, G4double Dx4, G4double Alp2);
198 virtual ~G4PolyhedronTrap ();
199};
200
202public:
204 G4double Dy, G4double Dz);
205 virtual ~G4PolyhedronTrd1 ();
206};
207
209public:
211 G4double Dy1, G4double Dy2, G4double Dz);
212 virtual ~G4PolyhedronTrd2 ();
213};
214
216public:
218 virtual ~G4PolyhedronTube ();
219};
220
222public:
224 G4double Phi1, G4double Dphi);
225 virtual ~G4PolyhedronTubs ();
226};
227
229 public:
231 G4double sPhi, G4double dPhi);
232 virtual ~G4PolyhedronParaboloid ();
233};
234
236 public:
238 G4double tan2, G4double halfZ);
239 virtual ~G4PolyhedronHype ();
240};
241
243 public:
245 G4double zcut1, G4double zcut2);
246 virtual ~G4PolyhedronEllipsoid ();
247};
248
250 public:
252 G4double zcut1);
254};
255
256std::ostream& operator<<(std::ostream& os, const G4Polyhedron&);
257
258#endif /* G4POLYHEDRON_HH */
std::ostream & operator<<(std::ostream &os, const G4Polyhedron &)
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
virtual ~G4PolyhedronBox()
Definition: G4Polyhedron.cc:47
virtual ~G4PolyhedronCone()
Definition: G4Polyhedron.cc:53
virtual ~G4PolyhedronCons()
Definition: G4Polyhedron.cc:60
virtual ~G4PolyhedronEllipsoid()
virtual ~G4PolyhedronHype()
virtual ~G4PolyhedronPara()
Definition: G4Polyhedron.cc:67
virtual ~G4PolyhedronParaboloid()
virtual ~G4PolyhedronPcon()
Definition: G4Polyhedron.cc:75
virtual ~G4PolyhedronPgon()
Definition: G4Polyhedron.cc:84
virtual ~G4PolyhedronSphere()
Definition: G4Polyhedron.cc:91
virtual ~G4PolyhedronTorus()
Definition: G4Polyhedron.cc:98
virtual ~G4PolyhedronTrap()
virtual ~G4PolyhedronTrd1()
virtual ~G4PolyhedronTrd2()
virtual ~G4PolyhedronTube()
virtual ~G4PolyhedronTubs()
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
virtual ~G4Polyhedron()
Definition: G4Polyhedron.cc:35