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