Geant4 11.4.0
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 G4TessellatedSolid.
32// It is defined by three fVertices, which shall be supplied in anti-clockwise
33// order looking from the outsider of the solid where it belongs.
34// Its constructor:
35//
36// G4TriangularFacet (const G4ThreeVector Pt0, const G4ThreeVector vt1,
37// const G4ThreeVector vt2, G4FacetVertexType);
38//
39// takes 4 parameters to define the three fVertices:
40// 1) G4FacetvertexType = "ABSOLUTE": in this case Pt0, vt1 and vt2 are
41// the 3 fVertices in anti-clockwise order looking from the outsider.
42// 2) G4FacetvertexType = "RELATIVE": in this case the first vertex is Pt0,
43// the second vertex is Pt0+vt1 and the third vertex is Pt0+vt2, all
44// in anti-clockwise order when looking from the outsider.
45
46// Author: P.R.Truscott (QinetiQ Ltd, UK), 31.10.2004 - Created
47// M.Gayer (CERN), 12.10.2012 - Reviewed optimised implementation
48// --------------------------------------------------------------------
49#ifndef G4TRIANGULARFACET_HH
50#define G4TRIANGULARFACET_HH
51
52#include "G4VFacet.hh"
53#include "G4Types.hh"
54#include "G4ThreeVector.hh"
55
56#include <vector>
57#include <array>
58
59/**
60 * @brief G4TriangularFacet defines a facet with 3 vertices, used for the
61 * contruction of G4TessellatedSolid. Vertices shall be supplied in
62 * anti-clockwise order looking from the outsider of the solid where it belongs.
63 */
64
66{
67 public:
68
69 /**
70 * Default Constructor.
71 */
73
74 /**
75 * Constructs a facet with 3 vertices, given its parameters.
76 * @param[in] Pt0 The anchor point, first vertex.
77 * @param[in] vt1 Second vertex.
78 * @param[in] vt2 Third vertex.
79 * @param[in] vType The positioning type for the vertices, either:
80 * "ABSOLUTE" - vertices set in anti-clockwise order
81 * when looking from the outsider.
82 * "RELATIVE" - first vertex is Pt0, second is Pt0+vt1,
83 * and the third vertex is Pt0+vt2,
84 * still in anti-clockwise order.
85 */
86 G4TriangularFacet (const G4ThreeVector& Pt0, const G4ThreeVector& vt1,
87 const G4ThreeVector& vt2, G4FacetVertexType vType);
88
89 /**
90 * Destructor.
91 */
92 ~G4TriangularFacet () override;
93
94 /**
95 * Copy and move constructors.
96 */
98 G4TriangularFacet ( G4TriangularFacet&& right) noexcept ;
99
100 /**
101 * Assignment and move assignment operators.
102 */
104 G4TriangularFacet& operator=( G4TriangularFacet&& right) noexcept ;
105
106 /**
107 * Returns a pointer to a newly allocated duplicate copy of the facet.
108 */
109 G4VFacet* GetClone () override;
110
111 /**
112 * Generates and returns an identical facet, but with the normal vector
113 * pointing at 180 degrees.
114 */
116
117 /**
118 * Determines the vector between p and the closest point on the facet to p.
119 */
121
122 /**
123 * Determines the closest distance between point p and the facet.
124 */
125 G4double Distance (const G4ThreeVector& p, G4double minDist) override;
126
127 /**
128 * Determines the distance to point 'p'. kInfinity is returned if either:
129 * (1) outgoing is TRUE and the dot product of the normal vector to the
130 * facet and the displacement vector from p to the triangle is negative.
131 * (2) outgoing is FALSE and the dot product of the normal vector to the
132 * facet and the displacement vector from p to the triangle is positive.
133 */
134 G4double Distance (const G4ThreeVector& p, G4double minDist,
135 const G4bool outgoing) override;
136
137 /**
138 * Calculates the furthest the triangle extends in fA particular
139 * direction defined by the vector axis.
140 */
141 G4double Extent (const G4ThreeVector axis) override;
142
143 /**
144 * Finds the next intersection when going from 'p' in the direction of 'v'.
145 * If 'outgoing' is true, only consider the face if we are going out
146 * through the face; otherwise, if false, only consider the face if we are
147 * going in through the face.
148 * @returns true if there is an intersection, false otherwise.
149 */
150 G4bool Intersect (const G4ThreeVector& p, const G4ThreeVector& v,
151 const G4bool outgoing, G4double& distance,
152 G4double& distFromSurface,
153 G4ThreeVector& normal) override;
154
155 /**
156 * Auxiliary method for returning the surface area.
157 */
158 G4double GetArea () const override;
159
160 /**
161 * Auxiliary method to get a uniform random point on the facet.
162 */
163 G4ThreeVector GetPointOnFace () const override;
164
165 /**
166 * Returns/sets the normal vector to the facet.
167 */
168 G4ThreeVector GetSurfaceNormal () const override;
169 void SetSurfaceNormal (const G4ThreeVector& normal);
170
171 /**
172 * Returns the type ID, "G4TriangularFacet" of the facet.
173 */
174 G4GeometryType GetEntityType () const override;
175
176 /**
177 * Returns true if the facet is defined.
178 */
179 inline G4bool IsDefined () const override;
180
181 /**
182 * Returns the number of vertices, i.e. 3.
183 */
184 inline G4int GetNumberOfVertices () const override;
185
186 /**
187 * Returns the vertex based on the index 'i'.
188 */
189 inline G4ThreeVector GetVertex (G4int i) const override;
190
191 /**
192 * Methods to set the vertices.
193 */
194 inline void SetVertex (G4int i, const G4ThreeVector& val) override;
195 inline void SetVertices(std::vector<G4ThreeVector>* v) override;
196
197 /**
198 * Returns the radius to the anchor point and centered on the circumcentre.
199 */
200 inline G4double GetRadius () const override;
201
202 /**
203 * Returns the circumcentre point of the facet.
204 */
205 inline G4ThreeVector GetCircumcentre () const override;
206
207 /**
208 * Returns the allocated memory (sizeof) by the facet.
209 */
210 inline G4int AllocatedMemory() override;
211
212 /**
213 * Accessor/setter for the vertex index.
214 */
215 inline G4int GetVertexIndex (G4int i) const override;
216 inline void SetVertexIndex (G4int i, G4int j) override;
217
218 private:
219
220 /**
221 * Utilities for copying/moving data.
222 */
223 void CopyFrom(const G4TriangularFacet& rhs);
224 void MoveFrom(G4TriangularFacet& rhs);
225
226 G4ThreeVector fSurfaceNormal;
227 G4double fArea = 0.0;
228 G4ThreeVector fCircumcentre;
229 G4double fRadius = 0.0;
230 std::array<G4int, 3> fIndices;
231 std::vector<G4ThreeVector>* fVertices = nullptr;
232
233 G4double fA, fB, fC;
234 G4double fDet;
235 G4double fSqrDist = 0.0;
236 G4ThreeVector fE1, fE2;
237 G4bool fIsDefined = false;
238};
239
240// --------------------------------------------------------------------
241// Inline Methods
242// --------------------------------------------------------------------
243
244#include "G4TriangularFacet.icc"
245
246#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
G4String G4GeometryType
Definition G4VSolid.hh:70
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
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668