Geant4-11
G4GenericTrap.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// G4GenericTrap
27//
28// Class description:
29//
30// G4GenericTrap is a solid which represents an arbitrary trapezoid with
31// up to 8 vertices standing on two parallel planes perpendicular to Z axis.
32//
33// Parameters in the constructor:
34// - name - solid name
35// - halfZ - the solid half length in Z
36// - vertices - the (x,y) coordinates of vertices:
37// o first four points: vertices[i], i<4
38// are the vertices sitting on the -halfZ plane;
39// o last four points: vertices[i], i>=4
40// are the vertices sitting on the +halfZ plane.
41//
42// The order of defining the vertices of the solid is the following:
43// - point 0 is connected with points 1,3,4
44// - point 1 is connected with points 0,2,5
45// - point 2 is connected with points 1,3,6
46// - point 3 is connected with points 0,2,7
47// - point 4 is connected with points 0,5,7
48// - point 5 is connected with points 1,4,6
49// - point 6 is connected with points 2,5,7
50// - point 7 is connected with points 3,4,6
51// Points can be identical in order to create shapes with less than
52// 8 vertices.
53
54// Authors:
55// Tatiana Nikitina, CERN; Ivana Hrivnacova, IPN Orsay
56// Adapted from Root Arb8 implementation, author Andrei Gheata, CERN
57// -------------------------------------------------------------------
58#ifndef G4GENERICTRAP_HH
59#define G4GENERICTRAP_HH
60
61#include "G4GeomTypes.hh"
62
63#if defined(G4GEOM_USE_USOLIDS)
64#define G4GEOM_USE_UGENERICTRAP 1
65#endif
66
67#if defined(G4GEOM_USE_UGENERICTRAP)
68 #define G4UGenericTrap G4GenericTrap
69 #include "G4UGenericTrap.hh"
70#else
71
72#include <vector>
73
74#include "globals.hh"
75#include "G4TwoVector.hh"
76#include "G4VSolid.hh"
77#include "G4TessellatedSolid.hh"
78
79class G4GenericTrap : public G4VSolid
80{
81 public: // with description
82
83 G4GenericTrap( const G4String& name, G4double halfZ,
84 const std::vector<G4TwoVector>& vertices );
85 // Constructor
86
88 // Destructor
89
90 // Accessors
91
92 inline G4double GetZHalfLength() const;
93 inline G4int GetNofVertices() const;
94 inline G4TwoVector GetVertex(G4int index) const;
95 inline const std::vector<G4TwoVector>& GetVertices() const;
96 inline G4double GetTwistAngle(G4int index) const;
97 inline G4bool IsTwisted() const;
98 inline G4int GetVisSubdivisions() const;
99 inline void SetVisSubdivisions(G4int subdiv);
100
101 // Solid methods
102
103 EInside Inside(const G4ThreeVector& p) const;
106 const G4ThreeVector& v) const;
107 G4double DistanceToIn(const G4ThreeVector& p) const;
109 const G4ThreeVector& v,
110 const G4bool calcNorm = false,
111 G4bool* validNorm = nullptr,
112 G4ThreeVector* n = nullptr) const;
113 G4double DistanceToOut(const G4ThreeVector& p) const;
115 G4bool CalculateExtent(const EAxis pAxis,
116 const G4VoxelLimits& pVoxelLimit,
117 const G4AffineTransform& pTransform,
118 G4double& pmin, G4double& pmax) const;
119
121
122 G4VSolid* Clone() const;
123
124 std::ostream& StreamInfo(std::ostream& os) const;
125
127
130
131 // Visualisation functions
132
133 G4Polyhedron* GetPolyhedron () const;
134 void DescribeYourselfTo(G4VGraphicsScene& scene) const;
135 G4VisExtent GetExtent() const;
137
138 public: // without description
139
140 G4GenericTrap(__void__&);
141 // Fake default constructor for usage restricted to direct object
142 // persistency for clients requiring preallocation of memory for
143 // persistifiable objects.
144
145 G4GenericTrap(const G4GenericTrap& rhs);
147 // Copy constructor and assignment operator.
148
149 private:
150
151 // Internal methods
152
153 inline void SetTwistAngle(G4int index, G4double twist);
155 G4bool CheckOrder(const std::vector<G4TwoVector>& vertices) const;
156 G4bool IsSegCrossing(const G4TwoVector& a, const G4TwoVector& b,
157 const G4TwoVector& c, const G4TwoVector& d) const;
158 G4bool IsSegCrossingZ(const G4TwoVector& a, const G4TwoVector& b,
159 const G4TwoVector& c, const G4TwoVector& d) const;
160 void ReorderVertices(std::vector<G4ThreeVector>& vertices) const;
161 void ComputeBBox();
164
165 G4VFacet* MakeDownFacet(const std::vector<G4ThreeVector>& fromVertices,
166 G4int ind1, G4int ind2, G4int ind3) const;
167 G4VFacet* MakeUpFacet(const std::vector<G4ThreeVector>& fromVertices,
168 G4int ind1, G4int ind2, G4int ind3) const;
169 G4VFacet* MakeSideFacet(const G4ThreeVector& downVertex0,
170 const G4ThreeVector& downVertex1,
171 const G4ThreeVector& upVertex1,
172 const G4ThreeVector& upVertex0) const;
174
176 const std::vector<G4TwoVector>& poly) const;
178 const G4ThreeVector& v, const G4int ipl) const ;
180 const G4ThreeVector& v, const G4int ipl) const;
182 const G4int ipl) const;
183 G4double SafetyToFace(const G4ThreeVector& p, const G4int iseg) const;
185 const G4ThreeVector& p1,
186 const G4ThreeVector& p2,
187 const G4ThreeVector& p3) const;
189 const G4ThreeVector& p1,
190 const G4ThreeVector& p2,
191 const G4ThreeVector& p3) const;
193 const G4ThreeVector& p1,
194 const G4ThreeVector& p2,
195 const G4ThreeVector& p3) const;
196 protected:
197
198 mutable G4bool fRebuildPolyhedron = false;
199 mutable G4Polyhedron* fpPolyhedron = nullptr;
200
201 private:
202
203 // static data members
204
205 static const G4int fgkNofVertices;
206 static const G4double fgkTolerance;
207
209
210 // data members
211
213 std::vector<G4TwoVector> fVertices;
220
222 // Codes for faces (kXY[num]=num of lateral face,kMZ= minus z face etc)
223
226 // Surface and Volume
227};
228
229#include "G4GenericTrap.icc"
230
231#endif
232
233#endif
static const G4double pMax
static const G4double pMin
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4VFacet * MakeUpFacet(const std::vector< G4ThreeVector > &fromVertices, G4int ind1, G4int ind2, G4int ind3) const
G4VSolid * Clone() const
std::vector< G4TwoVector > fVertices
G4int GetNofVertices() const
G4Polyhedron * CreatePolyhedron() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4GenericTrap & operator=(const G4GenericTrap &rhs)
G4double SafetyToFace(const G4ThreeVector &p, const G4int iseg) const
G4double GetFaceCubicVolume(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
G4double GetZHalfLength() const
G4double GetCubicVolume()
G4TwoVector GetVertex(G4int index) const
G4double GetSurfaceArea()
G4int fVisSubdivisions
EInside InsidePolygone(const G4ThreeVector &p, const std::vector< G4TwoVector > &poly) const
static const G4int fgkNofVertices
G4bool fRebuildPolyhedron
G4ThreeVector GetPointOnSurface() const
G4double GetTwistAngle(G4int index) const
G4double halfCarTolerance
G4double GetTwistedFaceSurfaceArea(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
G4ThreeVector GetMinimumBBox() const
G4double fTwist[4]
G4double GetFaceSurfaceArea(const G4ThreeVector &p0, const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3) const
G4Polyhedron * GetPolyhedron() const
std::ostream & StreamInfo(std::ostream &os) const
G4ThreeVector fMinBBoxVector
G4ThreeVector fMaxBBoxVector
G4bool ComputeIsTwisted()
G4ThreeVector GetMaximumBBox() const
G4VisExtent GetExtent() const
G4VFacet * MakeDownFacet(const std::vector< G4ThreeVector > &fromVertices, G4int ind1, G4int ind2, G4int ind3) const
EInside Inside(const G4ThreeVector &p) const
void ReorderVertices(std::vector< G4ThreeVector > &vertices) const
G4double DistToTriangle(const G4ThreeVector &p, const G4ThreeVector &v, const G4int ipl) const
G4bool IsSegCrossingZ(const G4TwoVector &a, const G4TwoVector &b, const G4TwoVector &c, const G4TwoVector &d) const
G4VFacet * MakeSideFacet(const G4ThreeVector &downVertex0, const G4ThreeVector &downVertex1, const G4ThreeVector &upVertex1, const G4ThreeVector &upVertex0) const
G4GeometryType GetEntityType() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4TessellatedSolid * fTessellatedSolid
const std::vector< G4TwoVector > & GetVertices() const
static const G4double fgkTolerance
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4int GetVisSubdivisions() const
G4bool IsTwisted() const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double fCubicVolume
void SetVisSubdivisions(G4int subdiv)
G4double DistToPlane(const G4ThreeVector &p, const G4ThreeVector &v, const G4int ipl) const
G4TessellatedSolid * CreateTessellatedSolid() const
G4bool CheckOrder(const std::vector< G4TwoVector > &vertices) const
void SetTwistAngle(G4int index, G4double twist)
G4double fSurfaceArea
G4GenericTrap(const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
G4Polyhedron * fpPolyhedron
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
G4ThreeVector NormalToPlane(const G4ThreeVector &p, const G4int ipl) const
G4bool IsSegCrossing(const G4TwoVector &a, const G4TwoVector &b, const G4TwoVector &c, const G4TwoVector &d) const
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
const char * name(G4int ptype)