Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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+vt1, 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// Author: P.R.Truscott (QinetiQ Ltd, UK), 31.10.2004 - Created
51// M.Gayer (CERN), 12.10.2012 - Reviewed optimised 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
61/**
62 * @brief G4QuadrangularFacet defines a facet with 4 vertices, used for the
63 * contruction of G4TessellatedSolid. Vertices shall be in the same plane and
64 * be supplied in anti-clockwise order looking from the outsider of the solid
65 * where it belongs.
66 */
67
69{
70 public:
71
72 /**
73 * Constructs a facet with 4 vertices, given its parameters.
74 * @param[in] Pt0 The anchor point, first vertex.
75 * @param[in] vt1 Second vertex.
76 * @param[in] vt2 Third vertex.
77 * @param[in] vt3 Fourth vertex.
78 * @param[in] vType The positioning type for the vertices, either:
79 * "ABSOLUTE" - vertices set in anti-clockwise order
80 * when looking from the outsider.
81 * "RELATIVE" - first vertex is Pt0, second is Pt0+vt1,
82 * third vertex is Pt0+vt2 and fourth is Pt0+vt3,
83 * still in anti-clockwise order.
84 */
85 G4QuadrangularFacet (const G4ThreeVector& Pt0, const G4ThreeVector& vt1,
86 const G4ThreeVector& vt2, const G4ThreeVector& vt3,
87 G4FacetVertexType vType);
88
89 /**
90 * Default Destructor.
91 */
92 ~G4QuadrangularFacet () override = default;
93
94 /**
95 * Copy constructor and assignment operator.
96 */
99
100 /**
101 * Returns a pointer to a newly allocated duplicate copy of the facet.
102 */
103 G4VFacet* GetClone () override;
104
105 /**
106 * Determines the vector between p and the closest point on the facet to p.
107 */
109
110 /**
111 * Determines the closest distance between point p and the facet.
112 */
113 G4double Distance (const G4ThreeVector& p, G4double minDist) override;
114
115 /**
116 * Determines the distance to point 'p'. kInfinity is returned if either:
117 * (1) outgoing is TRUE and the dot product of the normal vector to the
118 * facet and the displacement vector from p to the triangle is negative.
119 * (2) outgoing is FALSE and the dot product of the normal vector to the
120 * facet and the displacement vector from p to the triangle is positive.
121 */
122 G4double Distance (const G4ThreeVector& p, G4double minDist,
123 const G4bool outgoing) override;
124
125 /**
126 * Calculates the furthest the quadrangle extends in fA particular
127 * direction defined by the vector axis.
128 */
129 G4double Extent (const G4ThreeVector axis) override;
130
131 /**
132 * Finds the next intersection when going from 'p' in the direction of 'v'.
133 * If 'outgoing' is true, only consider the face if we are going out
134 * through the face; otherwise, if false, only consider the face if we are
135 * going in through the face.
136 * @returns true if there is an intersection, false otherwise.
137 */
138 G4bool Intersect (const G4ThreeVector& p, const G4ThreeVector& v,
139 const G4bool outgoing, G4double& distance,
140 G4double& distFromSurface,
141 G4ThreeVector& normal) override;
142 /**
143 * Returns the normal vector to the face.
144 */
145 G4ThreeVector GetSurfaceNormal () const override;
146
147 /**
148 * Auxiliary method for returning the surface area.
149 */
150 G4double GetArea () const override;
151
152 /**
153 * Auxiliary method to get a uniform random point on the facet.
154 */
155 G4ThreeVector GetPointOnFace () const override;
156
157 /**
158 * Returns the type ID, "G4QuadrangularFacet" of the facet.
159 */
160 G4GeometryType GetEntityType () const override;
161
162 /**
163 * Returns true if the facet is defined.
164 */
165 inline G4bool IsDefined () const override;
166
167 /**
168 * Returns the number of vertices, i.e. 4.
169 */
170 inline G4int GetNumberOfVertices () const override;
171
172 /**
173 * Returns the vertex based on the index 'i'.
174 */
175 inline G4ThreeVector GetVertex (G4int i) const override;
176
177 /**
178 * Methods to set the vertices.
179 */
180 inline void SetVertex (G4int i, const G4ThreeVector& val) override;
181 inline void SetVertices(std::vector<G4ThreeVector>* v) override;
182
183 /**
184 * Returns the radius to the anchor point and centered on the circumcentre.
185 */
186 inline G4double GetRadius () const override;
187
188 /**
189 * Returns the circumcentre point of the facet.
190 */
191 inline G4ThreeVector GetCircumcentre () const override;
192
193 private:
194
195 /**
196 * Private accessor/setter for the vertex index.
197 */
198 inline G4int GetVertexIndex (G4int i) const override;
199 inline void SetVertexIndex (G4int i, G4int val) override;
200
201 /**
202 * Returns the allocated memory (sizeof) by the facet.
203 */
204 inline G4int AllocatedMemory() override;
205
206 private:
207
208 G4double fRadius = 0.0;
209
210 G4ThreeVector fCircumcentre;
211
212 G4TriangularFacet fFacet1, fFacet2;
213};
214
215// --------------------------------------------------------------------
216// Inline Methods
217// --------------------------------------------------------------------
218
219#include "G4QuadrangularFacet.icc"
220
221#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
G4int GetNumberOfVertices() const override
G4double GetArea() const override
~G4QuadrangularFacet() override=default
G4ThreeVector Distance(const G4ThreeVector &p)
G4ThreeVector GetSurfaceNormal() const override
G4ThreeVector GetCircumcentre() const override
G4QuadrangularFacet & operator=(const G4QuadrangularFacet &right)
G4double GetRadius() const override
G4ThreeVector GetVertex(G4int i) const override
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType vType)
G4bool IsDefined() const override
void SetVertices(std::vector< G4ThreeVector > *v) override
G4VFacet * GetClone() override
G4GeometryType GetEntityType() const override
void SetVertex(G4int i, const G4ThreeVector &val) override
G4double Extent(const G4ThreeVector axis) override
G4ThreeVector GetPointOnFace() const override
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal) override
G4TriangularFacet defines a facet with 3 vertices, used for the contruction of G4TessellatedSolid....
virtual G4int AllocatedMemory()=0
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668