Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolyhedraSide.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//
30// --------------------------------------------------------------------
31// GEANT 4 class header file
32//
33//
34// G4PolyhedraSide
35//
36// Class description:
37//
38// Class implementing a face that represents one segmented side
39// of a polyhedra:
40//
41// G4PolyhedraSide( const G4PolyhedraSideRZ *prevRZ,
42// const G4PolyhedraSideRZ *tail,
43// const G4PolyhedraSideRZ *head,
44// const G4PolyhedraSideRZ *nextRZ,
45// G4int numSide,
46// G4double phiStart, G4double phiTotal,
47// G4bool phiIsOpen, G4bool isAllBehind=false )
48//
49// Values for r1,z1 and r2,z2 should be specified in clockwise
50// order in (r,z).
51
52// Author:
53// David C. Williams ([email protected])
54// --------------------------------------------------------------------
55
56#ifndef G4PolyhedraSide_hh
57#define G4PolyhedraSide_hh
58
59#include "G4VCSGface.hh"
60
62
64{
65 G4double r, z; // start of vector
66};
67
69{
70
71 public: // with description
72
74 const G4PolyhedraSideRZ *tail,
75 const G4PolyhedraSideRZ *head,
76 const G4PolyhedraSideRZ *nextRZ,
78 G4double phiStart, G4double phiTotal,
79 G4bool phiIsOpen, G4bool isAllBehind=false );
80 virtual ~G4PolyhedraSide();
81
82 G4PolyhedraSide( const G4PolyhedraSide &source );
84
85 G4bool Intersect( const G4ThreeVector &p, const G4ThreeVector &v,
86 G4bool outgoing, G4double surfTolerance,
87 G4double &distance, G4double &distFromSurface,
88 G4ThreeVector &normal, G4bool &allBehind );
89
90 G4double Distance( const G4ThreeVector &p, G4bool outgoing );
91
92 EInside Inside( const G4ThreeVector &p, G4double tolerance,
93 G4double *bestDistance );
94
95 G4ThreeVector Normal( const G4ThreeVector &p, G4double *bestDistance );
96
97 G4double Extent( const G4ThreeVector axis );
98
99 void CalculateExtent( const EAxis axis,
100 const G4VoxelLimits &voxelLimit,
101 const G4AffineTransform &tranform,
102 G4SolidExtentList &extentList );
103
104 G4VCSGface *Clone() { return new G4PolyhedraSide( *this ); }
105
106 public: // without description
107
108 // Methods used for GetPointOnSurface()
109
111 G4ThreeVector p2,
112 G4ThreeVector p3,
113 G4ThreeVector *p4 );
116 G4double *Area );
119
120 public: // without description
121
122 G4PolyhedraSide(__void__&);
123 // Fake default constructor for usage restricted to direct object
124 // persistency for clients requiring preallocation of memory for
125 // persistifiable objects.
126
127 protected:
128
129 //
130 // A couple internal data structures
131 //
132 struct sG4PolyhedraSideVec; // Secret recipe for allowing
133 friend struct sG4PolyhedraSideVec; // protected nested structures
134
135 typedef struct sG4PolyhedraSideEdge
136 {
137 G4ThreeVector normal; // Unit normal to this edge
138 G4ThreeVector corner[2]; // The two corners of this phi edge
139 G4ThreeVector cornNorm[2]; // The normals of these corners
141
142 typedef struct sG4PolyhedraSideVec
143 {
144 G4ThreeVector normal, // Normal (point out of the shape)
145 center, // Point in center of side
146 surfPhi, // Unit vector on surface pointing along phi
147 surfRZ; // Unit vector on surface pointing along R/Z
148 G4PolyhedraSideEdge *edges[2]; // The phi boundary edges to this side
149 // [0]=low phi [1]=high phi
150 G4ThreeVector edgeNorm[2]; // RZ edge normals [i] at {r[i],z[i]}
152
154 const G4PolyhedraSideVec& vec,
155 G4double normSign,
156 G4double surfTolerance,
157 G4double &distance,
158 G4double &distFromSurface );
159
161 const G4ThreeVector &v,
162 G4int *i1, G4int *i2 );
163
165
167
168 G4double GetPhi( const G4ThreeVector& p );
169
171 const G4PolyhedraSideVec &vec,
172 G4double *normDist );
173
175 const G4PolyhedraSideVec &vec,
176 G4double *normDist );
177
178 void CopyStuff( const G4PolyhedraSide &source );
179
180 protected:
181
182 G4int numSide; // Number sides
183 G4double r[2], z[2]; // r, z parameters, in specified order
184 G4double startPhi, // Start phi (0 to 2pi), if phiIsOpen
185 deltaPhi, // Delta phi (0 to 2pi), if phiIsOpen
186 endPhi; // End phi (>startPhi), if phiIsOpen
187 G4bool phiIsOpen; // True if there is a phi slice
188 G4bool allBehind; // True if the entire solid is "behind" this face
189
190 G4IntersectingCone *cone; // Our intersecting cone
191
192 G4PolyhedraSideVec *vecs; // Vector set for each facet of our face
193 G4PolyhedraSideEdge *edges; // The edges belong to vecs
194 G4double lenRZ, // RZ length of each side
195 lenPhi[2]; // Phi dimensions of each side
196 G4double edgeNorm; // Normal in RZ/Phi space to each side
197
198 private:
199
200 std::pair<G4ThreeVector, G4double> fPhi; // Cached value for phi
201 G4double kCarTolerance; // Geometrical surface thickness
202 G4double fSurfaceArea; // Surface Area
203};
204
205#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)
G4PolyhedraSideVec * vecs
G4double Distance(const G4ThreeVector &p, G4bool outgoing)
G4PolyhedraSideEdge * edges
G4ThreeVector GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double *Area)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &allBehind)
G4VCSGface * Clone()
G4PolyhedraSide & operator=(const G4PolyhedraSide &source)
G4double SurfaceTriangle(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector *p4)
struct G4PolyhedraSide::sG4PolyhedraSideVec G4PolyhedraSideVec
G4double GetPhi(const G4ThreeVector &p)
G4bool IntersectSidePlane(const G4ThreeVector &p, const G4ThreeVector &v, const G4PolyhedraSideVec &vec, G4double normSign, G4double surfTolerance, G4double &distance, G4double &distFromSurface)
G4ThreeVector GetPointOnFace()
G4double Extent(const G4ThreeVector axis)
G4int ClosestPhiSegment(G4double phi)
G4double lenPhi[2]
struct G4PolyhedraSide::sG4PolyhedraSideEdge G4PolyhedraSideEdge
G4int PhiSegment(G4double phi)
G4int LineHitsSegments(const G4ThreeVector &p, const G4ThreeVector &v, G4int *i1, G4int *i2)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance)
G4double DistanceAway(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double DistanceToOneSide(const G4ThreeVector &p, const G4PolyhedraSideVec &vec, G4double *normDist)
G4double SurfaceArea()
virtual ~G4PolyhedraSide()
G4IntersectingCone * cone
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)
void CopyStuff(const G4PolyhedraSide &source)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58