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