Geant4  10.03
G4QuadrangularFacet.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: G4QuadrangularFacet.hh 95801 2016-02-25 10:59:41Z gcosmo $
28 //
29 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
30 //
31 // Class G4QuadrangularFacet
32 //
33 // Class description:
34 //
35 // The G4QuadrangularFacet class is used for the contruction of
36 // G4TessellatedSolid.
37 // It is defined by four fVertices, which shall be in the same plane and be
38 // supplied in anti-clockwise order looking from the outsider of the solid
39 // where it belongs. Its constructor
40 //
41 // G4QuadrangularFacet (const G4ThreeVector Pt0, const G4ThreeVector vt1,
42 // const G4ThreeVector vt2, const G4ThreeVector vt3,
43 // G4FacetVertexType);
44 //
45 // takes 5 parameters to define the four fVertices:
46 // 1) G4FacetvertexType = "ABSOLUTE": in this case Pt0, vt1, vt2 and vt3
47 // are the four fVertices required in anti-clockwise order when looking
48 // from the outsider.
49 // 2) G4FacetvertexType = "RELATIVE": in this case the first vertex is Pt0,
50 // the second vertex is Pt0+vt, the third vertex is Pt0+vt2 and
51 // the fourth vertex is Pt0+vt3, in anti-clockwise order when looking
52 // from the outsider.
53 
54 // CHANGE HISTORY
55 // --------------
56 //
57 // 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
58 // 12 October 2012, M Gayer, CERN, - Reviewed optimized implementation.
59 //
61 #ifndef G4QuadrangularFacet_HH
62 #define G4QuadrangularFacet_HH 1
63 
64 #include "G4VFacet.hh"
65 #include "G4Types.hh"
66 #include "G4ThreeVector.hh"
67 #include "G4TriangularFacet.hh"
68 
70 {
71  public: // with description
72 
73  G4QuadrangularFacet (const G4ThreeVector &Pt0, const G4ThreeVector &vt1,
74  const G4ThreeVector &vt2, const G4ThreeVector &vt3,
78 
80 
81  G4VFacet *GetClone ();
82 
84  G4double Distance (const G4ThreeVector &p, G4double minDist);
85  G4double Distance (const G4ThreeVector &p, G4double minDist,
86  const G4bool outgoing);
87  G4double Extent (const G4ThreeVector axis);
88  G4bool Intersect (const G4ThreeVector &p, const G4ThreeVector &v,
89  const G4bool outgoing, G4double &distance,
90  G4double &distFromSurface, G4ThreeVector &normal);
92 
93  G4double GetArea () const;
95 
97 
98  inline G4bool IsDefined () const;
99  inline G4int GetNumberOfVertices () const;
100  inline G4ThreeVector GetVertex (G4int i) const;
101  inline void SetVertex (G4int i, const G4ThreeVector &val);
102  inline void SetVertices(std::vector<G4ThreeVector> *v);
103 
104  inline G4double GetRadius () const;
105  inline G4ThreeVector GetCircumcentre () const;
106 
107  private:
108 
109  inline G4int GetVertexIndex (G4int i) const;
110  inline void SetVertexIndex (G4int i, G4int val);
111 
112  inline G4int AllocatedMemory();
113 
114  private:
115 
118 
120 };
121 
123 // Inlined Methods
125 
127 {
128  return 4;
129 }
130 
132 {
133  return i == 3 ? fFacet2.GetVertex(2) : fFacet1.GetVertex(i);
134 }
135 
136 
138 {
139  return fRadius;
140 }
141 
143 {
144  return fCircumcentre;
145 }
146 
148 {
149  switch (i)
150  {
151  case 0:
152  fFacet1.SetVertex(0, val);
153  fFacet2.SetVertex(0, val);
154  break;
155  case 1:
156  fFacet1.SetVertex(1, val);
157  break;
158  case 2:
159  fFacet1.SetVertex(2, val);
160  fFacet2.SetVertex(1, val);
161  break;
162  case 3:
163  fFacet2.SetVertex(2, val);
164  break;
165  }
166 }
167 
168 inline void G4QuadrangularFacet::SetVertices(std::vector<G4ThreeVector> *v)
169 {
170  fFacet1.SetVertices(v);
171  fFacet2.SetVertices(v);
172 }
173 
175 {
176  return fFacet1.IsDefined();
177 }
178 
180 {
181  return i == 3 ? fFacet2.GetVertexIndex(2) : fFacet1.GetVertexIndex(i);
182 }
183 
184 
186 {
187  switch (i)
188  {
189  case 0:
190  fFacet1.SetVertexIndex(0, val);
191  fFacet2.SetVertexIndex(0, val);
192  break;
193  case 1:
194  fFacet1.SetVertexIndex(1, val);
195  break;
196  case 2:
197  fFacet1.SetVertexIndex(2, val);
198  fFacet2.SetVertexIndex(1, val);
199  break;
200  case 3:
201  fFacet2.SetVertexIndex(2, val);
202  break;
203  }
204 }
205 
207 {
208  return sizeof(*this) + fFacet1.AllocatedMemory() + fFacet2.AllocatedMemory();
209 }
210 
211 #endif
G4bool IsDefined() const
G4TriangularFacet fFacet1
G4int GetVertexIndex(G4int i) const
void SetVertexIndex(G4int i, G4int j)
CLHEP::Hep3Vector G4ThreeVector
void SetVertices(std::vector< G4ThreeVector > *v)
G4double Extent(const G4ThreeVector axis)
void SetVertex(G4int i, const G4ThreeVector &val)
G4ThreeVector GetSurfaceNormal() const
int G4int
Definition: G4Types.hh:78
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType)
void SetVertices(std::vector< G4ThreeVector > *v)
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4FacetVertexType
Definition: G4VFacet.hh:56
G4int GetVertexIndex(G4int i) const
G4ThreeVector GetVertex(G4int i) const
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector GetPointOnFace() const
G4QuadrangularFacet & operator=(const G4QuadrangularFacet &right)
void SetVertexIndex(G4int i, G4int val)
void SetVertex(G4int i, const G4ThreeVector &val)
G4TriangularFacet fFacet2
G4ThreeVector GetCircumcentre() const
G4double GetArea() const
G4ThreeVector GetVertex(G4int i) const
G4double GetRadius() const
double G4double
Definition: G4Types.hh:76
G4ThreeVector Distance(const G4ThreeVector &p)
G4GeometryType GetEntityType() const
G4int GetNumberOfVertices() const