Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// G4TriangularFacet
28//
29// Class description:
30//
31// The G4TriangularFacet class is used for the contruction of
32// G4TessellatedSolid.
33// It is defined by three fVertices, which shall be supplied in anti-clockwise
34// order looking from the outsider of the solid where it belongs.
35// Its constructor:
36//
37// G4TriangularFacet (const G4ThreeVector Pt0, const G4ThreeVector vt1,
38// const G4ThreeVector vt2, G4FacetVertexType);
39//
40// takes 4 parameters to define the three fVertices:
41// 1) G4FacetvertexType = "ABSOLUTE": in this case Pt0, vt1 and vt2 are
42// the 3 fVertices in anti-clockwise order looking from the outsider.
43// 2) G4FacetvertexType = "RELATIVE": in this case the first vertex is Pt0,
44// the second vertex is Pt0+vt1 and the third vertex is Pt0+vt2, all
45// in anti-clockwise order when looking from the outsider.
46
47// 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
48// 12 October 2012, M Gayer, CERN, - Reviewed optimized implementation.
49// --------------------------------------------------------------------
50#ifndef G4TRIANGULARFACET_HH
51#define G4TRIANGULARFACET_HH 1
52
53#include "G4VFacet.hh"
54#include "G4Types.hh"
55#include "G4ThreeVector.hh"
56
57#include <vector>
58#include <array>
59
61{
62 public:
63
65 ~G4TriangularFacet () override;
66
67 G4TriangularFacet (const G4ThreeVector& vt0, const G4ThreeVector& vt1,
70 G4TriangularFacet ( G4TriangularFacet&& right) noexcept ;
71
73 G4TriangularFacet& operator=( G4TriangularFacet&& right) noexcept ;
74
75 G4VFacet* GetClone () override;
77
79 G4double Distance (const G4ThreeVector& p, G4double minDist) override;
80 G4double Distance (const G4ThreeVector& p, G4double minDist,
81 const G4bool outgoing) override;
82 G4double Extent (const G4ThreeVector axis) override;
83 G4bool Intersect (const G4ThreeVector& p, const G4ThreeVector& v,
84 const G4bool outgoing, G4double& distance,
85 G4double& distFromSurface,
86 G4ThreeVector& normal) override;
87 G4double GetArea () const override;
88 G4ThreeVector GetPointOnFace () const override;
89
90 G4ThreeVector GetSurfaceNormal () const override;
91 void SetSurfaceNormal (const G4ThreeVector& normal);
92
93 G4GeometryType GetEntityType () const override;
94
95 inline G4bool IsDefined () const override;
96 inline G4int GetNumberOfVertices () const override;
97 inline G4ThreeVector GetVertex (G4int i) const override;
98 inline void SetVertex (G4int i, const G4ThreeVector& val) override;
99
100 inline G4ThreeVector GetCircumcentre () const override;
101 inline G4double GetRadius () const override;
102
103 inline G4int AllocatedMemory() override;
104
105 inline G4int GetVertexIndex (G4int i) const override;
106 inline void SetVertexIndex (G4int i, G4int j) override;
107 inline void SetVertices(std::vector<G4ThreeVector>* v) override;
108
109 private:
110
111 void CopyFrom(const G4TriangularFacet& rhs);
112 void MoveFrom(G4TriangularFacet& rhs);
113
114 G4ThreeVector fSurfaceNormal;
115 G4double fArea = 0.0;
116 G4ThreeVector fCircumcentre;
117 G4double fRadius = 0.0;
118 std::array<G4int, 3> fIndices;
119 std::vector<G4ThreeVector>* fVertices = nullptr;
120
121 G4double fA, fB, fC;
122 G4double fDet;
123 G4double fSqrDist = 0.0;
124 G4ThreeVector fE1, fE2;
125 G4bool fIsDefined = false;
126};
127
128// --------------------------------------------------------------------
129// Inlined Methods
130// --------------------------------------------------------------------
131
133{
134 return fIsDefined;
135}
136
138{
139 return 3;
140}
141
143{
144 G4int indice = fIndices[i];
145 return indice < 0 ? (*fVertices)[i] : (*fVertices)[indice];
146}
147
149{
150 (*fVertices)[i] = val;
151}
152
154{
155 return fCircumcentre;
156}
157
159{
160 return fRadius;
161}
162
164{
165 G4int size = sizeof(*this);
166 size += GetNumberOfVertices() * sizeof(G4ThreeVector);
167 return size;
168}
169
171{
172 if (i < 3) return fIndices[i];
173 else return 999999999;
174}
175
177{
178 fIndices[i] = j;
179}
180
181inline void G4TriangularFacet::SetVertices(std::vector<G4ThreeVector>* v)
182{
183 if (fIndices[0] < 0 && (fVertices != nullptr))
184 {
185 delete fVertices;
186 fVertices = nullptr;
187 }
188 fVertices = v;
189}
190
191#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4FacetVertexType
Definition G4VFacet.hh:48
G4GeometryType GetEntityType() const override
G4int GetVertexIndex(G4int i) const override
G4double Extent(const G4ThreeVector axis) override
G4double GetRadius() const override
void SetVertex(G4int i, const G4ThreeVector &val) override
G4ThreeVector GetSurfaceNormal() const override
G4ThreeVector GetVertex(G4int i) const override
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal) override
G4ThreeVector GetCircumcentre() const override
G4TriangularFacet & operator=(const G4TriangularFacet &right)
void SetVertices(std::vector< G4ThreeVector > *v) override
G4ThreeVector GetPointOnFace() const override
G4bool IsDefined() const override
G4TriangularFacet * GetFlippedFacet()
G4ThreeVector Distance(const G4ThreeVector &p)
G4VFacet * GetClone() override
G4int GetNumberOfVertices() const override
void SetSurfaceNormal(const G4ThreeVector &normal)
void SetVertexIndex(G4int i, G4int j) override
G4double GetArea() const override
G4int AllocatedMemory() override