Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VTwistSurface.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// G4VTwistSurface
27//
28// Class description:
29//
30// Abstract base class for boundary surface of G4VSolid.
31
32// 01-Aug-2002 - Kotoyo Hoshina (hoshina@hepburn.s.chiba-u.ac.jp), created.
33// 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), Integration in Geant4
34// from original version in Jupiter-2.5.02 application.
35// --------------------------------------------------------------------
36#ifndef G4VTWISTSURFACE_HH
37#define G4VTWISTSURFACE_HH
38
40
41#include "G4VSolid.hh"
42#include "geomdefs.hh"
43
44#include "G4RotationMatrix.hh"
45
46#define G4VSURFACENXX 10
47
49{
50 public: // without description
51
54
55 public: // with description
56
57 G4VTwistSurface (const G4String& name);
58 G4VTwistSurface (const G4String& name,
59 const G4RotationMatrix& rot,
60 const G4ThreeVector& tlate,
61 G4int handedness,
62 const EAxis axis1,
63 const EAxis axis2,
64 G4double axis0min = -kInfinity,
65 G4double axis1min = -kInfinity,
66 G4double axis0max = kInfinity,
67 G4double axis1max = kInfinity);
68
69 virtual ~G4VTwistSurface();
70
71 virtual G4int AmIOnLeftSide(const G4ThreeVector& me,
72 const G4ThreeVector& vec,
73 G4bool withTol = true);
74
75 virtual G4double DistanceToBoundary( G4int areacode,
76 G4ThreeVector& xx,
77 const G4ThreeVector& p) ;
78
79 virtual G4double DistanceToIn(const G4ThreeVector& gp,
80 const G4ThreeVector& gv,
81 G4ThreeVector& gxxbest);
82 virtual G4double DistanceToOut(const G4ThreeVector& gp,
83 const G4ThreeVector& gv,
84 G4ThreeVector& gxxbest);
85 virtual G4double DistanceTo(const G4ThreeVector& gp,
86 G4ThreeVector& gxx);
87
89 const G4ThreeVector& gv,
90 G4ThreeVector gxx[],
91 G4double distance[],
92 G4int areacode[],
93 G4bool isvalid[],
94 EValidate validate=kValidateWithTol) = 0;
95
97 G4ThreeVector gxx[],
98 G4double distance[],
99 G4int areacode[]) = 0;
100
101 void DebugPrint() const;
102
103 // get methods
104
105 virtual G4ThreeVector GetNormal(const G4ThreeVector& xx,G4bool isGlobal) = 0;
106
107 virtual G4String GetName() const { return fName; }
108 virtual void GetBoundaryParameters(const G4int& areacode,
109 G4ThreeVector& d,
110 G4ThreeVector& x0,
111 G4int& boundarytype) const;
112 virtual G4ThreeVector GetBoundaryAtPZ(G4int areacode,
113 const G4ThreeVector& p) const;
114
116 const G4ThreeVector& v,
117 const G4ThreeVector& x0,
118 const G4ThreeVector& n0,
119 G4ThreeVector& xx);
120
122 const G4ThreeVector& x0,
123 const G4ThreeVector& n0,
124 G4ThreeVector& xx);
125
127 const G4ThreeVector& x0,
128 const G4ThreeVector& t1,
129 const G4ThreeVector& t2,
130 G4ThreeVector& xx,
131 G4ThreeVector& n);
132
134 const G4ThreeVector& x0,
135 const G4ThreeVector& d,
136 G4ThreeVector& xx);
137
138 inline G4bool IsAxis0 (G4int areacode) const;
139 inline G4bool IsAxis1 (G4int areacode) const;
140 inline G4bool IsOutside (G4int areacode) const;
141 inline G4bool IsInside (G4int areacode, G4bool testbitmode = false) const;
142 inline G4bool IsBoundary (G4int areacode, G4bool testbitmode = false) const;
143 inline G4bool IsCorner (G4int areacode, G4bool testbitmode = false) const;
144 inline G4bool IsValidNorm() const { return fIsValidNorm; }
145 G4bool IsSameBoundary (G4VTwistSurface* surface1, G4int areacode1,
146 G4VTwistSurface* surface2, G4int areacode2 ) const;
147 inline G4int GetAxisType(G4int areacode, G4int whichaxis) const;
148
153
154 // set methods
155
156 inline void SetAxis(G4int i, const EAxis axis) { fAxis[i] = axis; }
157 inline void SetNeighbours(G4VTwistSurface* ax0min, G4VTwistSurface* ax1min,
158 G4VTwistSurface* ax0max, G4VTwistSurface* ax1max);
159
161 G4bool isGlobal = false ) = 0 ;
164 virtual G4double GetSurfaceArea() = 0 ;
165 virtual void GetFacets(G4int m, G4int n, G4double xyz[][3],
166 G4int faces[][4], G4int iside) = 0 ;
167 G4int GetNode( G4int i, G4int j, G4int m, G4int n, G4int iside ) ;
168 G4int GetFace( G4int i, G4int j, G4int m, G4int n, G4int iside ) ;
170 G4int number, G4int orientation) ;
171
172 public: // without description
173
174 G4VTwistSurface(__void__&);
175 // Fake default constructor for usage restricted to direct object
176 // persistency for clients requiring preallocation of memory for
177 // persistifiable objects.
178
179 protected: // with description
180
181 // get methods
182
183 inline G4VTwistSurface** GetNeighbours() { return fNeighbours; }
184 inline G4int GetNeighbours(G4int areacode, G4VTwistSurface* surfaces[]);
185 inline G4ThreeVector GetCorner(G4int areacode) const;
186 void GetBoundaryAxis(G4int areacode, EAxis axis[]) const;
187 void GetBoundaryLimit(G4int areacode, G4double limit[]) const;
188 virtual G4int GetAreaCode(const G4ThreeVector& xx, G4bool withtol=true) = 0;
189
190 // set methods
191
192 virtual void SetBoundary(const G4int& axiscode,
193 const G4ThreeVector& direction,
194 const G4ThreeVector& x0,
195 const G4int& boundarytype);
196 // areacode must be one of them:
197 // sAxis0 & sAxisMin, sAxis0 & sAxisMax,
198 // sAxis1 & sAxisMin, sAxis1 & sAxisMax.
199 // boundarytype represents the shape of locus
200 // from the start point to end point of boundary.
201 // ex.
202 // sAxisRho = linear line which start point is fixed at origin.
203 // sAxisPhi = part of circle which center placed at the origin.
204
205 void SetCorner(G4int areacode, G4double x, G4double y, G4double z);
206
207 private:
208
209 virtual void SetBoundaries() = 0;
210 virtual void SetCorners() = 0;
211
212 // data members ---------------------------------------------------------
213
214 public:
215
216 static const G4int sOutside ;
217 static const G4int sInside ;
218 static const G4int sBoundary;
219 static const G4int sCorner;
220 static const G4int sC0Min1Min;
221 static const G4int sC0Max1Min;
222 static const G4int sC0Max1Max;
223 static const G4int sC0Min1Max;
224 static const G4int sAxisMin;
225 static const G4int sAxisMax;
226 static const G4int sAxisX;
227 static const G4int sAxisY;
228 static const G4int sAxisZ;
229 static const G4int sAxisRho;
230 static const G4int sAxisPhi;
231 static const G4int sAxis0;
232 static const G4int sAxis1;
233 static const G4int sSizeMask;
234 static const G4int sAxisMask;
235 static const G4int sAreaMask;
236
237 protected:
238
240 {
241 public:
242
244 virtual ~CurrentStatus();
245
246 inline G4ThreeVector GetXX(G4int i) const { return fXX[i]; }
247 inline G4double GetDistance(G4int i) const { return fDistance[i]; }
248 inline G4int GetAreacode(G4int i) const { return fAreacode[i]; }
249 inline G4int GetNXX() const { return fNXX; }
250 inline G4bool IsDone() const { return fDone; }
251 inline G4bool IsValid(G4int i) const { return fIsValid[i]; }
252
253 void SetCurrentStatus(G4int i,
254 G4ThreeVector& xx,
255 G4double& dist,
256 G4int& areacode,
257 G4bool& isvalid,
258 G4int nxx,
259 EValidate validate,
260 const G4ThreeVector* p,
261 const G4ThreeVector* v = nullptr);
262
263 void ResetfDone(EValidate validate,
264 const G4ThreeVector* p,
265 const G4ThreeVector* v = nullptr);
266
267
268 void DebugPrint() const;
269
270 private:
271
272 G4double fDistance[G4VSURFACENXX];
274 G4int fAreacode[G4VSURFACENXX];
275 G4bool fIsValid[G4VSURFACENXX];
276 G4int fNXX;
277 G4ThreeVector fLastp;
278 G4ThreeVector fLastv;
279 EValidate fLastValidate;
280 G4bool fDone;
281 };
282
283 class Boundary
284 {
285 public:
286 Boundary();
287 virtual ~Boundary();
288
289 void SetFields(const G4int& areacode,
290 const G4ThreeVector& d,
291 const G4ThreeVector& x0,
292 const G4int& boundarytype);
293
294 G4bool IsEmpty() const;
295
296 G4bool GetBoundaryParameters(const G4int& areacode,
297 G4ThreeVector& d,
298 G4ThreeVector& x0,
299 G4int& boundarytype) const;
300
301 private:
302 G4int fBoundaryAcode;
303 G4ThreeVector fBoundaryDirection;
304 G4ThreeVector fBoundaryX0;
305 G4int fBoundaryType;
306 };
307
317 {
318 public:
321 };
325
326 private:
327
328 G4VTwistSurface* fNeighbours[4]; // {0,1,2,3} = sAxis0min, sAxis1min,
329 // sAxis0max, sAxis1max
330 G4ThreeVector fCorners[4]; // corners of the surface in local coordinate
331 Boundary fBoundaries[4]; // boundaries of the surface.
332 G4String fName;
333
334 class G4SurfSideQuery
335 {
336 public:
337 G4ThreeVector me;
338 G4ThreeVector vec;
339 G4bool withTol;
340 G4int amIOnLeftSide;
341 };
342 G4SurfSideQuery fAmIOnLeftSide;
343};
344
345//========================================================
346// inline functions
347//========================================================
348
350{
351 G4double phi ; // parameter phi
352 G4double u ; // parameter u
353 G4ThreeVector xx ; // intersection point in cartesian
354 G4double distance ; // distance to intersection
355 G4int areacode ; // the areacode of the intersection
356 G4bool isvalid ; // valid intersection ??
357
358};
359
360inline
362{
363 return a.distance < b.distance ;
364}
365
366inline
368{
369 return ( ( a.xx - b.xx ).mag() < 1E-9*CLHEP::mm ) ;
370}
371
372#include "G4VTwistSurface.icc"
373
374#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4bool DistanceSort(const Intersection &a, const Intersection &b)
G4bool EqualIntersection(const Intersection &a, const Intersection &b)
#define G4VSURFACENXX
void SetFields(const G4int &areacode, const G4ThreeVector &d, const G4ThreeVector &x0, const G4int &boundarytype)
G4bool GetBoundaryParameters(const G4int &areacode, G4ThreeVector &d, G4ThreeVector &x0, G4int &boundarytype) const
G4int GetAreacode(G4int i) const
G4double GetDistance(G4int i) const
G4bool IsValid(G4int i) const
void SetCurrentStatus(G4int i, G4ThreeVector &xx, G4double &dist, G4int &areacode, G4bool &isvalid, G4int nxx, EValidate validate, const G4ThreeVector *p, const G4ThreeVector *v=nullptr)
G4ThreeVector GetXX(G4int i) const
void ResetfDone(EValidate validate, const G4ThreeVector *p, const G4ThreeVector *v=nullptr)
virtual G4int AmIOnLeftSide(const G4ThreeVector &me, const G4ThreeVector &vec, G4bool withTol=true)
static const G4int sAxisMask
static const G4int sC0Min1Min
static const G4int sC0Min1Max
virtual void GetFacets(G4int m, G4int n, G4double xyz[][3], G4int faces[][4], G4int iside)=0
G4bool IsAxis1(G4int areacode) const
G4double DistanceToPlane(const G4ThreeVector &p, const G4ThreeVector &x0, const G4ThreeVector &n0, G4ThreeVector &xx)
G4int GetNode(G4int i, G4int j, G4int m, G4int n, G4int iside)
static const G4int sOutside
void SetAxis(G4int i, const EAxis axis)
void SetNeighbours(G4VTwistSurface *ax0min, G4VTwistSurface *ax1min, G4VTwistSurface *ax0max, G4VTwistSurface *ax1max)
virtual ~G4VTwistSurface()
G4ThreeVector ComputeGlobalDirection(const G4ThreeVector &lp) const
G4bool IsCorner(G4int areacode, G4bool testbitmode=false) const
static const G4int sAxisMax
static const G4int sAxis0
virtual G4double GetBoundaryMin(G4double)=0
G4int GetFace(G4int i, G4int j, G4int m, G4int n, G4int iside)
G4double fAxisMax[2]
G4bool IsValidNorm() const
void GetBoundaryAxis(G4int areacode, EAxis axis[]) const
virtual G4double DistanceTo(const G4ThreeVector &gp, G4ThreeVector &gxx)
G4bool IsAxis0(G4int areacode) const
G4RotationMatrix fRot
G4int GetEdgeVisibility(G4int i, G4int j, G4int m, G4int n, G4int number, G4int orientation)
virtual G4int DistanceToSurface(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector gxx[], G4double distance[], G4int areacode[], G4bool isvalid[], EValidate validate=kValidateWithTol)=0
G4VTwistSurface ** GetNeighbours()
G4double DistanceToLine(const G4ThreeVector &p, const G4ThreeVector &x0, const G4ThreeVector &d, G4ThreeVector &xx)
G4ThreeVector ComputeLocalDirection(const G4ThreeVector &gp) const
void GetBoundaryLimit(G4int areacode, G4double limit[]) const
void DebugPrint() const
G4int GetAxisType(G4int areacode, G4int whichaxis) const
static const G4int sAxisPhi
static const G4int sAxisMin
static const G4int sC0Max1Max
virtual G4double DistanceToOut(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
virtual G4double DistanceToIn(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
static const G4int sAxis1
virtual G4ThreeVector GetBoundaryAtPZ(G4int areacode, const G4ThreeVector &p) const
G4bool IsInside(G4int areacode, G4bool testbitmode=false) const
virtual G4ThreeVector SurfacePoint(G4double, G4double, G4bool isGlobal=false)=0
G4ThreeVector fTrans
virtual void SetBoundary(const G4int &axiscode, const G4ThreeVector &direction, const G4ThreeVector &x0, const G4int &boundarytype)
G4ThreeVector ComputeLocalPoint(const G4ThreeVector &gp) const
virtual G4double GetBoundaryMax(G4double)=0
static const G4int sAxisRho
void SetCorner(G4int areacode, G4double x, G4double y, G4double z)
G4ThreeVector GetCorner(G4int areacode) const
virtual G4int GetAreaCode(const G4ThreeVector &xx, G4bool withtol=true)=0
static const G4int sBoundary
static const G4int sAxisZ
G4bool IsOutside(G4int areacode) const
virtual G4double DistanceToBoundary(G4int areacode, G4ThreeVector &xx, const G4ThreeVector &p)
G4double fAxisMin[2]
virtual G4int DistanceToSurface(const G4ThreeVector &gp, G4ThreeVector gxx[], G4double distance[], G4int areacode[])=0
static const G4int sCorner
static const G4int sC0Max1Min
static const G4int sInside
G4bool IsSameBoundary(G4VTwistSurface *surface1, G4int areacode1, G4VTwistSurface *surface2, G4int areacode2) const
virtual G4String GetName() const
CurrentStatus fCurStatWithV
virtual G4ThreeVector GetNormal(const G4ThreeVector &xx, G4bool isGlobal)=0
G4int GetNeighbours(G4int areacode, G4VTwistSurface *surfaces[])
static const G4int sAxisY
static const G4int sSizeMask
static const G4int sAxisX
G4double DistanceToPlaneWithV(const G4ThreeVector &p, const G4ThreeVector &v, const G4ThreeVector &x0, const G4ThreeVector &n0, G4ThreeVector &xx)
static const G4int sAreaMask
G4bool IsBoundary(G4int areacode, G4bool testbitmode=false) const
G4ThreeVector ComputeGlobalPoint(const G4ThreeVector &lp) const
G4double DistanceToPlane(const G4ThreeVector &p, const G4ThreeVector &x0, const G4ThreeVector &t1, const G4ThreeVector &t2, G4ThreeVector &xx, G4ThreeVector &n)
G4SurfCurNormal fCurrentNormal
virtual G4double GetSurfaceArea()=0
virtual void GetBoundaryParameters(const G4int &areacode, G4ThreeVector &d, G4ThreeVector &x0, G4int &boundarytype) const
CurrentStatus fCurStat
EAxis
Definition: geomdefs.hh:54
G4ThreeVector xx