Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolyconeSide.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// G4PolyconeSide
27//
28// Class description:
29//
30// Class implmenting a face that represents one conical side
31// of a polycone:
32//
33// G4PolyconeSide( const G4PolyconeSideRZ *prevRZ,
34// const G4PolyconeSideRZ *tail,
35// const G4PolyconeSideRZ *head,
36// const G4PolyconeSideRZ *nextRZ,
37// G4double phiStart, G4double deltaPhi,
38// G4bool phiIsOpen, G4bool isAllBehind=false )
39//
40// Values for r1,z1 and r2,z2 should be specified in clockwise
41// order in (r,z).
42
43// Author: David C. Williams ([email protected])
44// --------------------------------------------------------------------
45#ifndef G4POLYCONESIDE_HH
46#define G4POLYCONESIDE_HH
47
48#include "G4VCSGface.hh"
49
51
53{
54 G4double r, z; // start of vector
55};
56
57// ----------------------------------------------------------------------------
58// MT-specific utility code
59
60#include "G4GeomSplitter.hh"
61
62// The class G4PlSideData is introduced to encapsulate the
63// fields of the class G4PolyconeSide that may not be read-only.
64//
66{
67 public:
68
70 {
71 fPhix = 0.; fPhiy = 0.; fPhiz = 0.; fPhik = 0.;
72 }
73
74 G4double fPhix=0., fPhiy=0., fPhiz=0., fPhik=0.; // Cached values for phi
75};
76
77// The type G4PlSideManager is introduced to
78// encapsulate the methods used by both the master thread and
79// worker threads to allocate memory space for the fields encapsulated
80// by the class G4PlSideData.
81//
83
84//
85// ----------------------------------------------------------------------------
86
88{
89 public:
90
91 G4PolyconeSide( const G4PolyconeSideRZ* prevRZ,
92 const G4PolyconeSideRZ* tail,
93 const G4PolyconeSideRZ* head,
94 const G4PolyconeSideRZ* nextRZ,
95 G4double phiStart, G4double deltaPhi,
96 G4bool phiIsOpen, G4bool isAllBehind = false );
97 ~G4PolyconeSide() override;
98
99 G4PolyconeSide( const G4PolyconeSide& source );
100 G4PolyconeSide& operator=( const G4PolyconeSide& source );
101
102 G4bool Intersect(const G4ThreeVector& p, const G4ThreeVector& v,
103 G4bool outgoing, G4double surfTolerance,
104 G4double& distance, G4double &distFromSurface,
105 G4ThreeVector& normal, G4bool& isAllBehind) override;
106
107 G4double Distance( const G4ThreeVector& p, G4bool outgoing ) override;
108
109 EInside Inside( const G4ThreeVector& p, G4double tolerance,
110 G4double* bestDistance ) override;
111
113 G4double* bestDistance ) override;
114
115 G4double Extent( const G4ThreeVector axis ) override;
116
117 void CalculateExtent( const EAxis axis,
118 const G4VoxelLimits& voxelLimit,
119 const G4AffineTransform& tranform,
120 G4SolidExtentList& extentList ) override;
121
122 G4VCSGface* Clone() override { return new G4PolyconeSide( *this ); }
123
124 G4double SurfaceArea() override;
125 G4ThreeVector GetPointOnFace() override;
126
127 G4PolyconeSide(__void__&);
128 // Fake default constructor for usage restricted to direct object
129 // persistency for clients requiring preallocation of memory for
130 // persistifiable objects.
131
132 inline G4int GetInstanceID() const { return instanceID; }
133 // Returns the instance ID.
134
136 // Returns the private data instance manager.
137
138 protected:
139
140 G4double DistanceAway( const G4ThreeVector& p, G4bool opposite,
141 G4double& distOutside2,
142 G4double* rzNorm = nullptr );
143
144 G4double DistanceAway( const G4ThreeVector& p, G4double& distOutside2,
145 G4double* edgeRZnorm );
146
147 G4bool PointOnCone( const G4ThreeVector& hit, G4double normSign,
148 const G4ThreeVector& p,
149 const G4ThreeVector& v, G4ThreeVector& normal );
150
151 void CopyStuff( const G4PolyconeSide& source );
152
153 static void FindLineIntersect( G4double x1, G4double y1,
154 G4double tx1, G4double ty1,
155 G4double x2, G4double y2,
156 G4double tx2, G4double ty2,
157 G4double& x, G4double& y );
158
159 G4double GetPhi( const G4ThreeVector& p );
160
161 protected:
162
163 G4double r[2], z[2]; // r, z parameters, in specified order
164 G4double startPhi, // Start phi (0 to 2pi), if phiIsOpen
165 deltaPhi; // Delta phi (0 to 2pi), if phiIsOpen
166 G4bool phiIsOpen = false; // True if there is a phi slice
167 G4bool allBehind = false; // True if the entire solid is "behind" this face
168
169 G4IntersectingCone* cone = nullptr; // Our intersecting utility class
170
171 G4double rNorm, zNorm; // Normal to surface in r,z space
172 G4double rS, zS; // Unit vector along surface in r,z space
173 G4double length; // Length of face in r,z space
175 prevZS; // Unit vector along previous polyconeSide
177 nextZS; // Unit vector along next polyconeSide
178
180 zNormEdge[2]; // Normal to edges
181
183 G4ThreeVector* corners = nullptr; // The coordinates of the corners
184 // (if phiIsOpen)
185 private:
186
187 G4double kCarTolerance; // Geometrical surface thickness
188 G4double fSurfaceArea = 0.0; // Used for surface calculation
189
190 G4int instanceID;
191 // This field is used as instance ID.
192 G4GEOM_DLL static G4PlSideManager subInstanceManager;
193 // This field helps to use the class G4PlSideManager introduced above.
194};
195
196#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList) override
G4double zNormEdge[2]
G4int GetInstanceID() const
~G4PolyconeSide() override
G4ThreeVector * corners
G4double GetPhi(const G4ThreeVector &p)
static void FindLineIntersect(G4double x1, G4double y1, G4double tx1, G4double ty1, G4double x2, G4double y2, G4double tx2, G4double ty2, G4double &x, G4double &y)
G4ThreeVector GetPointOnFace() override
void CopyStuff(const G4PolyconeSide &source)
G4PolyconeSide(const G4PolyconeSideRZ *prevRZ, const G4PolyconeSideRZ *tail, const G4PolyconeSideRZ *head, const G4PolyconeSideRZ *nextRZ, G4double phiStart, G4double deltaPhi, G4bool phiIsOpen, G4bool isAllBehind=false)
G4ThreeVector Normal(const G4ThreeVector &p, G4double *bestDistance) override
G4double SurfaceArea() override
G4double DistanceAway(const G4ThreeVector &p, G4bool opposite, G4double &distOutside2, G4double *rzNorm=nullptr)
EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance) override
G4VCSGface * Clone() override
G4IntersectingCone * cone
G4bool PointOnCone(const G4ThreeVector &hit, G4double normSign, const G4ThreeVector &p, const G4ThreeVector &v, G4ThreeVector &normal)
G4double Distance(const G4ThreeVector &p, G4bool outgoing) override
G4double rNormEdge[2]
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, G4bool outgoing, G4double surfTolerance, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal, G4bool &isAllBehind) override
static const G4PlSideManager & GetSubInstanceManager()
G4double Extent(const G4ThreeVector axis) override
G4PolyconeSide & operator=(const G4PolyconeSide &source)
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
#define G4GEOM_DLL
Definition geomwdefs.hh:44