Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4TriangularFacet.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 and of QinetiQ Ltd, *
20 // * subject to DEFCON 705 IPR conditions. *
21 // * By using, copying, modifying or distributing the software (or *
22 // * any work based on the software) you agree to acknowledge its *
23 // * use in resulting scientific publications, and indicate your *
24 // * acceptance of all terms of the Geant4 Software license. *
25 // ********************************************************************
26 //
27 // $Id: G4TriangularFacet.hh 95801 2016-02-25 10:59:41Z gcosmo $
28 //
29 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
30 //
31 // Class G4TriangularFacet
32 //
33 // Class description:
34 //
35 // The G4TriangularFacet class is used for the contruction of
36 // G4TessellatedSolid.
37 // It is defined by three fVertices, which shall be supplied in anti-clockwise
38 // order looking from the outsider of the solid where it belongs.
39 // Its constructor:
40 //
41 // G4TriangularFacet (const G4ThreeVector Pt0, const G4ThreeVector vt1,
42 // const G4ThreeVector vt2, G4FacetVertexType);
43 //
44 // takes 4 parameters to define the three fVertices:
45 // 1) G4FacetvertexType = "ABSOLUTE": in this case Pt0, vt1 and vt2 are
46 // the 3 fVertices in anti-clockwise order looking from the outsider.
47 // 2) G4FacetvertexType = "RELATIVE": in this case the first vertex is Pt0,
48 // the second vertex is Pt0+vt1 and the third vertex is Pt0+vt2, all
49 // in anti-clockwise order when looking from the outsider.
50 
51 // CHANGE HISTORY
52 // --------------
53 //
54 // 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
55 // 12 October 2012, M Gayer, CERN, - Reviewed optimized implementation.
56 //
58 #ifndef G4TriangularFacet_hh
59 #define G4TriangularFacet_hh 1
60 
61 #include "G4VFacet.hh"
62 #include "G4Types.hh"
63 #include "G4ThreeVector.hh"
64 
66 {
67  public: // with desctiption
68 
71 
72  G4TriangularFacet (const G4ThreeVector &vt0, const G4ThreeVector &vt1,
73  const G4ThreeVector &vt2, G4FacetVertexType);
75 
77 
78  G4VFacet *GetClone ();
80 
82  G4double Distance (const G4ThreeVector &p, G4double minDist);
83  G4double Distance (const G4ThreeVector &p, G4double minDist,
84  const G4bool outgoing);
85  G4double Extent (const G4ThreeVector axis);
86  G4bool Intersect (const G4ThreeVector &p, const G4ThreeVector &v,
87  const G4bool outgoing, G4double &distance,
88  G4double &distFromSurface, G4ThreeVector &normal);
89  G4double GetArea () const;
91 
93  void SetSurfaceNormal (G4ThreeVector normal);
94 
96 
97  inline G4bool IsDefined () const;
98  inline G4int GetNumberOfVertices () const;
99  inline G4ThreeVector GetVertex (G4int i) const;
100  inline void SetVertex (G4int i, const G4ThreeVector &val);
101 
102  inline G4ThreeVector GetCircumcentre () const;
103  inline G4double GetRadius () const;
104 
105  inline G4int AllocatedMemory();
106 
107  inline G4int GetVertexIndex (G4int i) const;
108  inline void SetVertexIndex (G4int i, G4int j);
109  inline void SetVertices(std::vector<G4ThreeVector> *v);
110 
111  private:
112 
113  void CopyFrom(const G4TriangularFacet &rhs);
114 
115  G4ThreeVector fSurfaceNormal;
116  G4double fArea;
117  G4ThreeVector fCircumcentre;
118  G4double fRadius;
119  G4int fIndices[3];
120 
121  std::vector<G4ThreeVector> *fVertices;
122 
123  G4double fA, fB, fC;
124  G4double fDet;
125  G4double fSqrDist;
126  G4ThreeVector fE1, fE2;
127  G4bool fIsDefined;
128 };
129 
131 // Inlined Methods
133 
135 {
136  return fIsDefined;
137 }
138 
140 {
141  return 3;
142 }
143 
145 {
146  G4int indice = fIndices[i];
147  return indice < 0 ? (*fVertices)[i] : (*fVertices)[indice];
148 }
149 
151 {
152  (*fVertices)[i] = val;
153 }
154 
156 {
157  return fCircumcentre;
158 }
159 
161 {
162  return fRadius;
163 }
164 
166 {
167  G4int size = sizeof(*this);
168  size += GetNumberOfVertices() * sizeof(G4ThreeVector);
169  return size;
170 }
171 
173 {
174  if (i < 3) return fIndices[i];
175  else return 999999999;
176 }
177 
179 {
180  fIndices[i] = j;
181 }
182 
183 inline void G4TriangularFacet::SetVertices(std::vector<G4ThreeVector> *v)
184 {
185  if (fIndices[0] < 0 && fVertices)
186  {
187  delete fVertices;
188  fVertices = 0;
189  }
190  fVertices = v;
191 }
192 
193 #endif
G4bool IsDefined() const
G4double GetArea() const
G4int GetVertexIndex(G4int i) const
G4TriangularFacet & operator=(const G4TriangularFacet &right)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
void SetVertexIndex(G4int i, G4int j)
CLHEP::Hep3Vector G4ThreeVector
void SetVertices(std::vector< G4ThreeVector > *v)
const char * p
Definition: xmltok.h:285
void SetVertex(G4int i, const G4ThreeVector &val)
int G4int
Definition: G4Types.hh:78
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4FacetVertexType
Definition: G4VFacet.hh:56
G4GeometryType GetEntityType() const
bool G4bool
Definition: G4Types.hh:79
void SetSurfaceNormal(G4ThreeVector normal)
G4double Extent(const G4ThreeVector axis)
tuple v
Definition: test.py:18
G4ThreeVector GetCircumcentre() const
G4TriangularFacet * GetFlippedFacet()
G4ThreeVector GetVertex(G4int i) const
G4int GetNumberOfVertices() const
double G4double
Definition: G4Types.hh:76
G4ThreeVector GetPointOnFace() const
G4ThreeVector GetSurfaceNormal() const
G4double GetRadius() const
G4ThreeVector Distance(const G4ThreeVector &p)