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