Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QuadrangularFacet Class Reference

#include <G4QuadrangularFacet.hh>

+ Inheritance diagram for G4QuadrangularFacet:

Public Member Functions

 G4QuadrangularFacet (const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType)
 
 G4QuadrangularFacet (const G4QuadrangularFacet &right)
 
 ~G4QuadrangularFacet ()
 
G4QuadrangularFacetoperator= (const G4QuadrangularFacet &right)
 
G4VFacetGetClone ()
 
G4ThreeVector Distance (const G4ThreeVector &p)
 
G4double Distance (const G4ThreeVector &p, G4double minDist)
 
G4double Distance (const G4ThreeVector &p, G4double minDist, const G4bool outgoing)
 
G4double Extent (const G4ThreeVector axis)
 
G4bool Intersect (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
 
G4ThreeVector GetSurfaceNormal () const
 
G4double GetArea () const
 
G4ThreeVector GetPointOnFace () const
 
G4GeometryType GetEntityType () const
 
G4bool IsDefined () const
 
G4int GetNumberOfVertices () const
 
G4ThreeVector GetVertex (G4int i) const
 
void SetVertex (G4int i, const G4ThreeVector &val)
 
void SetVertices (std::vector< G4ThreeVector > *v)
 
G4double GetRadius () const
 
G4ThreeVector GetCircumcentre () const
 
- Public Member Functions inherited from G4VFacet
 G4VFacet ()
 
virtual ~G4VFacet ()
 
G4bool operator== (const G4VFacet &right) const
 
virtual G4int GetNumberOfVertices () const =0
 
virtual G4ThreeVector GetVertex (G4int i) const =0
 
virtual void SetVertex (G4int i, const G4ThreeVector &val)=0
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetSurfaceNormal () const =0
 
virtual G4bool IsDefined () const =0
 
virtual G4ThreeVector GetCircumcentre () const =0
 
virtual G4double GetRadius () const =0
 
virtual G4VFacetGetClone ()=0
 
virtual G4double Distance (const G4ThreeVector &, G4double)=0
 
virtual G4double Distance (const G4ThreeVector &, G4double, const G4bool)=0
 
virtual G4double Extent (const G4ThreeVector)=0
 
virtual G4bool Intersect (const G4ThreeVector &, const G4ThreeVector &, const G4bool, G4double &, G4double &, G4ThreeVector &)=0
 
virtual G4double GetArea () const =0
 
virtual G4ThreeVector GetPointOnFace () const =0
 
void ApplyTranslation (const G4ThreeVector v)
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4bool IsInside (const G4ThreeVector &p) const
 
virtual G4int AllocatedMemory ()=0
 
virtual void SetVertexIndex (G4int i, G4int j)=0
 
virtual G4int GetVertexIndex (G4int i) const =0
 
virtual void SetVertices (std::vector< G4ThreeVector > *vertices)=0
 

Additional Inherited Members

- Protected Attributes inherited from G4VFacet
G4double kCarTolerance
 
- Static Protected Attributes inherited from G4VFacet
static const G4double dirTolerance = 1.0E-14
 

Detailed Description

Definition at line 61 of file G4QuadrangularFacet.hh.

Constructor & Destructor Documentation

◆ G4QuadrangularFacet() [1/2]

G4QuadrangularFacet::G4QuadrangularFacet ( const G4ThreeVector Pt0,
const G4ThreeVector vt1,
const G4ThreeVector vt2,
const G4ThreeVector vt3,
G4FacetVertexType  vertexType 
)

Definition at line 51 of file G4QuadrangularFacet.cc.

56 : G4VFacet()
57{
58 G4double delta = 1.0 * kCarTolerance; // dimension tolerance
59 G4double epsilon = 0.01 * kCarTolerance; // planarity tolerance
60
61 G4ThreeVector e1, e2, e3;
62 SetVertex(0, vt0);
63 if (vertexType == ABSOLUTE)
64 {
65 SetVertex(1, vt1);
66 SetVertex(2, vt2);
67 SetVertex(3, vt3);
68
69 e1 = vt1 - vt0;
70 e2 = vt2 - vt0;
71 e3 = vt3 - vt0;
72 }
73 else
74 {
75 SetVertex(1, vt0 + vt1);
76 SetVertex(2, vt0 + vt2);
77 SetVertex(3, vt0 + vt3);
78
79 e1 = vt1;
80 e2 = vt2;
81 e3 = vt3;
82 }
83
84 // Check length of sides and diagonals
85 //
86 G4double leng1 = e1.mag();
87 G4double leng2 = (e2-e1).mag();
88 G4double leng3 = (e3-e2).mag();
89 G4double leng4 = e3.mag();
90
91 G4double diag1 = e2.mag();
92 G4double diag2 = (e3-e1).mag();
93
94 if (leng1 <= delta || leng2 <= delta || leng3 <= delta || leng4 <= delta ||
95 diag1 <= delta || diag2 <= delta)
96 {
97 ostringstream message;
98 message << "Sides/diagonals of facet are too small." << G4endl
99 << "P0 = " << GetVertex(0) << G4endl
100 << "P1 = " << GetVertex(1) << G4endl
101 << "P2 = " << GetVertex(2) << G4endl
102 << "P3 = " << GetVertex(3) << G4endl
103 << "Side1 length (P0->P1) = " << leng1 << G4endl
104 << "Side2 length (P1->P2) = " << leng2 << G4endl
105 << "Side3 length (P2->P3) = " << leng3 << G4endl
106 << "Side4 length (P3->P0) = " << leng4 << G4endl
107 << "Diagonal1 length (P0->P2) = " << diag1 << G4endl
108 << "Diagonal2 length (P1->P3) = " << diag2;
109 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
110 "GeomSolids1001", JustWarning, message);
111 return;
112 }
113
114 // Check that vertices are not collinear
115 //
116 G4double s1 = (e1.cross(e2)).mag()*0.5;
117 G4double s2 = ((e2-e1).cross(e3-e2)).mag()*0.5;
118 G4double s3 = (e2.cross(e3)).mag()*0.5;
119 G4double s4 = (e1.cross(e3)).mag()*0.5;
120
121 G4double h1 = 2.*s1 / std::max(std::max(leng1,leng2),diag1);
122 G4double h2 = 2.*s2 / std::max(std::max(leng2,leng3),diag2);
123 G4double h3 = 2.*s3 / std::max(std::max(leng3,leng4),diag1);
124 G4double h4 = 2.*s4 / std::max(std::max(leng4,leng1),diag2);
125
126 if (h1 <= delta || h2 <= delta || h3 <= delta || h4 <= delta )
127 {
128 ostringstream message;
129 message << "Facet has three or more collinear vertices." << G4endl
130 << "P0 = " << GetVertex(0) << G4endl
131 << "P1 = " << GetVertex(1) << G4endl
132 << "P2 = " << GetVertex(2) << G4endl
133 << "P3 = " << GetVertex(3) << G4endl
134 << "Smallest heights:" << G4endl
135 << " in triangle P0-P1-P2 = " << h1 << G4endl
136 << " in triangle P1-P2-P3 = " << h2 << G4endl
137 << " in triangle P2-P3-P0 = " << h3 << G4endl
138 << " in triangle P3-P0-P1 = " << h4;
139 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
140 "GeomSolids1001", JustWarning, message);
141 return;
142 }
143
144 // Check that vertices are coplanar by computing minimal
145 // height of tetrahedron comprising of vertices
146 //
147 G4double smax = std::max( std::max(s1,s2), std::max(s3,s4) );
148 G4double hmin = 0.5 * std::fabs( e1.dot(e2.cross(e3)) ) / smax;
149 if (hmin >= epsilon)
150 {
151 ostringstream message;
152 message << "Facet is not planar." << G4endl
153 << "Disrepancy = " << hmin << G4endl
154 << "P0 = " << GetVertex(0) << G4endl
155 << "P1 = " << GetVertex(1) << G4endl
156 << "P2 = " << GetVertex(2) << G4endl
157 << "P3 = " << GetVertex(3);
158 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
159 "GeomSolids1001", JustWarning, message);
160 return;
161 }
162
163 // Check that facet is convex by computing crosspoint
164 // of diagonals
165 //
166 G4ThreeVector normal = e2.cross(e3-e1);
167 G4double s = kInfinity, t = kInfinity, magnitude2 = normal.mag2();
168 if (magnitude2 > delta*delta) // check: magnitude2 != 0.
169 {
170 s = normal.dot(e1.cross(e3-e1)) / magnitude2;
171 t = normal.dot(e1.cross(e2)) / magnitude2;
172 }
173 if (s <= 0. || s >= 1. || t <= 0. || t >= 1.)
174 {
175 ostringstream message;
176 message << "Facet is not convex." << G4endl
177 << "Parameters of crosspoint of diagonals: "
178 << s << " and " << t << G4endl
179 << "should both be within (0,1) range" << G4endl
180 << "P0 = " << GetVertex(0) << G4endl
181 << "P1 = " << GetVertex(1) << G4endl
182 << "P2 = " << GetVertex(2) << G4endl
183 << "P3 = " << GetVertex(3);
184 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
185 "GeomSolids1001", JustWarning, message);
186 return;
187 }
188
189 // Define facet
190 //
193
194 normal = normal.unit();
195 fFacet1.SetSurfaceNormal(normal);
196 fFacet2.SetSurfaceNormal(normal);
197
198 G4ThreeVector vtmp = 0.5 * (e1 + e2);
199 fCircumcentre = GetVertex(0) + vtmp;
200 G4double radiusSqr = vtmp.mag2();
201 fRadius = std::sqrt(radiusSqr);
202 // 29.02.2016 Remark by E.Tcherniaev: computation
203 // of fCircumcenter and fRadius is wrong, however
204 // it did not create any problem till now.
205 // Bizarre! Need to investigate!
206}
G4double epsilon(G4double density, G4double temperature)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
@ ABSOLUTE
Definition: G4VFacet.hh:48
#define G4endl
Definition: G4ios.hh:57
Hep3Vector unit() const
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
void SetVertex(G4int i, const G4ThreeVector &val)
G4ThreeVector GetVertex(G4int i) const
void SetSurfaceNormal(G4ThreeVector normal)
G4VFacet()
Definition: G4VFacet.cc:44
G4double kCarTolerance
Definition: G4VFacet.hh:93

◆ G4QuadrangularFacet() [2/2]

G4QuadrangularFacet::G4QuadrangularFacet ( const G4QuadrangularFacet right)

Definition at line 216 of file G4QuadrangularFacet.cc.

217 : G4VFacet(rhs)
218{
219 fFacet1 = rhs.fFacet1;
220 fFacet2 = rhs.fFacet2;
221 fRadius = 0.0;
222}

◆ ~G4QuadrangularFacet()

G4QuadrangularFacet::~G4QuadrangularFacet ( )

Definition at line 210 of file G4QuadrangularFacet.cc.

211{
212}

Member Function Documentation

◆ Distance() [1/3]

G4ThreeVector G4QuadrangularFacet::Distance ( const G4ThreeVector p)

Definition at line 250 of file G4QuadrangularFacet.cc.

251{
252 G4ThreeVector v1 = fFacet1.Distance(p);
253 G4ThreeVector v2 = fFacet2.Distance(p);
254
255 if (v1.mag2() < v2.mag2()) return v1;
256 else return v2;
257}
G4ThreeVector Distance(const G4ThreeVector &p)

Referenced by Distance().

◆ Distance() [2/3]

G4double G4QuadrangularFacet::Distance ( const G4ThreeVector p,
G4double  minDist 
)
virtual

Implements G4VFacet.

Definition at line 261 of file G4QuadrangularFacet.cc.

263{
264 G4double dist = Distance(p).mag();
265 return dist;
266}
G4ThreeVector Distance(const G4ThreeVector &p)

◆ Distance() [3/3]

G4double G4QuadrangularFacet::Distance ( const G4ThreeVector p,
G4double  minDist,
const G4bool  outgoing 
)
virtual

Implements G4VFacet.

Definition at line 270 of file G4QuadrangularFacet.cc.

272{
273 G4double dist;
274
275 G4ThreeVector v = Distance(p);
276 G4double dir = v.dot(GetSurfaceNormal());
277 if ( ((dir > dirTolerance) && (!outgoing))
278 || ((dir < -dirTolerance) && outgoing))
279 dist = kInfinity;
280 else
281 dist = v.mag();
282 return dist;
283}
G4ThreeVector GetSurfaceNormal() const
static const G4double dirTolerance
Definition: G4VFacet.hh:92

◆ Extent()

G4double G4QuadrangularFacet::Extent ( const G4ThreeVector  axis)
virtual

Implements G4VFacet.

Definition at line 287 of file G4QuadrangularFacet.cc.

288{
289 G4double ss = 0;
290
291 for (G4int i = 0; i <= 3; ++i)
292 {
293 G4double sp = GetVertex(i).dot(axis);
294 if (sp > ss) ss = sp;
295 }
296 return ss;
297}
int G4int
Definition: G4Types.hh:85

◆ GetArea()

G4double G4QuadrangularFacet::GetArea ( ) const
virtual

Implements G4VFacet.

Definition at line 336 of file G4QuadrangularFacet.cc.

337{
338 G4double area = fFacet1.GetArea() + fFacet2.GetArea();
339 return area;
340}
G4double GetArea() const

◆ GetCircumcentre()

G4ThreeVector G4QuadrangularFacet::GetCircumcentre ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 135 of file G4QuadrangularFacet.hh.

136{
137 return fCircumcentre;
138}

◆ GetClone()

G4VFacet * G4QuadrangularFacet::GetClone ( )
virtual

Implements G4VFacet.

Definition at line 240 of file G4QuadrangularFacet.cc.

241{
243 GetVertex(2), GetVertex(3),
244 ABSOLUTE);
245 return c;
246}

◆ GetEntityType()

G4String G4QuadrangularFacet::GetEntityType ( ) const
virtual

Implements G4VFacet.

Definition at line 344 of file G4QuadrangularFacet.cc.

345{
346 return "G4QuadrangularFacet";
347}

◆ GetNumberOfVertices()

G4int G4QuadrangularFacet::GetNumberOfVertices ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 119 of file G4QuadrangularFacet.hh.

120{
121 return 4;
122}

◆ GetPointOnFace()

G4ThreeVector G4QuadrangularFacet::GetPointOnFace ( ) const
virtual

Implements G4VFacet.

Definition at line 324 of file G4QuadrangularFacet.cc.

325{
326 G4double s1 = fFacet1.GetArea();
327 G4double s2 = fFacet2.GetArea();
328 return ((s1+s2)*G4UniformRand() < s1) ?
329 fFacet1.GetPointOnFace() : fFacet2.GetPointOnFace();
330}
#define G4UniformRand()
Definition: Randomize.hh:52
G4ThreeVector GetPointOnFace() const

◆ GetRadius()

G4double G4QuadrangularFacet::GetRadius ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 130 of file G4QuadrangularFacet.hh.

131{
132 return fRadius;
133}

◆ GetSurfaceNormal()

G4ThreeVector G4QuadrangularFacet::GetSurfaceNormal ( ) const
virtual

Implements G4VFacet.

Definition at line 351 of file G4QuadrangularFacet.cc.

352{
353 return fFacet1.GetSurfaceNormal();
354}
G4ThreeVector GetSurfaceNormal() const

Referenced by Distance().

◆ GetVertex()

G4ThreeVector G4QuadrangularFacet::GetVertex ( G4int  i) const
inlinevirtual

Implements G4VFacet.

Definition at line 124 of file G4QuadrangularFacet.hh.

125{
126 return i == 3 ? fFacet2.GetVertex(2) : fFacet1.GetVertex(i);
127}
G4ThreeVector GetVertex(G4int i) const

Referenced by Extent(), G4QuadrangularFacet(), and GetClone().

◆ Intersect()

G4bool G4QuadrangularFacet::Intersect ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  outgoing,
G4double distance,
G4double distFromSurface,
G4ThreeVector normal 
)
virtual

Implements G4VFacet.

Definition at line 301 of file G4QuadrangularFacet.cc.

307{
308 G4bool intersect =
309 fFacet1.Intersect(p,v,outgoing,distance,distFromSurface,normal);
310 if (!intersect) intersect =
311 fFacet2.Intersect(p,v,outgoing,distance,distFromSurface,normal);
312 if (!intersect)
313 {
314 distance = distFromSurface = kInfinity;
315 normal.set(0,0,0);
316 }
317 return intersect;
318}
bool G4bool
Definition: G4Types.hh:86
void set(double x, double y, double z)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)

◆ IsDefined()

G4bool G4QuadrangularFacet::IsDefined ( ) const
inlinevirtual

Implements G4VFacet.

Definition at line 167 of file G4QuadrangularFacet.hh.

168{
169 return fFacet1.IsDefined();
170}
G4bool IsDefined() const

◆ operator=()

G4QuadrangularFacet & G4QuadrangularFacet::operator= ( const G4QuadrangularFacet right)

Definition at line 227 of file G4QuadrangularFacet.cc.

228{
229 if (this == &rhs) return *this;
230
231 fFacet1 = rhs.fFacet1;
232 fFacet2 = rhs.fFacet2;
233 fRadius = 0.0;
234
235 return *this;
236}

◆ SetVertex()

void G4QuadrangularFacet::SetVertex ( G4int  i,
const G4ThreeVector val 
)
inlinevirtual

Implements G4VFacet.

Definition at line 140 of file G4QuadrangularFacet.hh.

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}
void SetVertex(G4int i, const G4ThreeVector &val)

Referenced by G4QuadrangularFacet().

◆ SetVertices()

void G4QuadrangularFacet::SetVertices ( std::vector< G4ThreeVector > *  v)
inlinevirtual

Implements G4VFacet.

Definition at line 161 of file G4QuadrangularFacet.hh.

162{
163 fFacet1.SetVertices(v);
164 fFacet2.SetVertices(v);
165}
void SetVertices(std::vector< G4ThreeVector > *v)

The documentation for this class was generated from the following files: