Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QuadrangularFacet.cc
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//
28// $Id$
29//
30// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
31//
32// CHANGE HISTORY
33// --------------
34//
35// 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
36// 12 October 2012, M Gayer, CERN
37// New implementation reducing memory requirements by 50%,
38// and considerable CPU speedup together with the new
39// implementation of G4TessellatedSolid.
40//
41// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
42
44#include "geomdefs.hh"
45#include "Randomize.hh"
46
47using namespace std;
48
49///////////////////////////////////////////////////////////////////////////////
50//
51// !!!THIS IS A FUDGE!!! IT'S TWO ADJACENT G4TRIANGULARFACETS
52// --- NOT EFFICIENT BUT PRACTICAL.
53//
55 const G4ThreeVector &vt1,
56 const G4ThreeVector &vt2,
57 const G4ThreeVector &vt3,
58 G4FacetVertexType vertexType)
59{
60 G4ThreeVector e1, e2, e3;
61
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 G4double length1 = e1.mag();
84 G4double length2 = (GetVertex(2)-GetVertex(1)).mag();
85 G4double length3 = (GetVertex(3)-GetVertex(2)).mag();
86 G4double length4 = e3.mag();
87
88 G4ThreeVector normal1 = e1.cross(e2).unit();
89 G4ThreeVector normal2 = e2.cross(e3).unit();
90
91 bool isDefined = (length1 > kCarTolerance && length2 > kCarTolerance &&
92 length3 > kCarTolerance && length4 > kCarTolerance &&
93 normal1.dot(normal2) >= 0.9999999999);
94
95 if (isDefined)
96 {
97 fFacet1 = G4TriangularFacet (GetVertex(0),GetVertex(1),
99 fFacet2 = G4TriangularFacet (GetVertex(0),GetVertex(2),
101
104
105 G4ThreeVector normal12 = fFacet1.GetSurfaceNormal()
106 + fFacet2.GetSurfaceNormal();
107 G4ThreeVector normal34 = facet3.GetSurfaceNormal()
108 + facet4.GetSurfaceNormal();
109 G4ThreeVector normal = 0.25 * (normal12 + normal34);
110
111 fFacet1.SetSurfaceNormal (normal);
112 fFacet2.SetSurfaceNormal (normal);
113
114 G4ThreeVector vtmp = 0.5 * (e1 + e2);
115 fCircumcentre = GetVertex(0) + vtmp;
116 G4double radiusSqr = vtmp.mag2();
117 fRadius = std::sqrt(radiusSqr);
118 }
119 else
120 {
121 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()",
122 "GeomSolids0002", JustWarning,
123 "Length of sides of facet are too small or sides not planar.");
124 G4cerr << endl;
125 G4cerr << "P0 = " << GetVertex(0) << endl;
126 G4cerr << "P1 = " << GetVertex(1) << endl;
127 G4cerr << "P2 = " << GetVertex(2) << endl;
128 G4cerr << "P3 = " << GetVertex(3) << endl;
129 G4cerr << "Side lengths = P0->P1" << length1 << endl;
130 G4cerr << "Side lengths = P1->P2" << length2 << endl;
131 G4cerr << "Side lengths = P2->P3" << length3 << endl;
132 G4cerr << "Side lengths = P3->P0" << length4 << endl;
133 G4cerr << endl;
134 fRadius = 0.0;
135 }
136}
137
138///////////////////////////////////////////////////////////////////////////////
139//
141{
142}
143
144///////////////////////////////////////////////////////////////////////////////
145//
147 : G4VFacet(rhs)
148{
149 fFacet1 = rhs.fFacet1;
150 fFacet2 = rhs.fFacet2;
151 fRadius = 0.0;
152}
153
154///////////////////////////////////////////////////////////////////////////////
155//
158{
159 if (this == &rhs)
160 return *this;
161
162 fFacet1 = rhs.fFacet1;
163 fFacet2 = rhs.fFacet2;
164 fRadius = 0.0;
165
166 return *this;
167}
168
169///////////////////////////////////////////////////////////////////////////////
170//
172{
174 GetVertex(2), GetVertex(3),
175 ABSOLUTE);
176 return c;
177}
178
179///////////////////////////////////////////////////////////////////////////////
180//
182{
183 G4ThreeVector v1 = fFacet1.Distance(p);
184 G4ThreeVector v2 = fFacet2.Distance(p);
185
186 if (v1.mag2() < v2.mag2()) return v1;
187 else return v2;
188}
189
190///////////////////////////////////////////////////////////////////////////////
191//
193 G4double)
194{
195 G4double dist = Distance(p).mag();
196 return dist;
197}
198
199///////////////////////////////////////////////////////////////////////////////
200//
202 const G4bool outgoing)
203{
204 G4double dist;
205
206 G4ThreeVector v = Distance(p);
207 G4double dir = v.dot(GetSurfaceNormal());
208 if ( ((dir > dirTolerance) && (!outgoing))
209 || ((dir < -dirTolerance) && outgoing))
210 dist = kInfinity;
211 else
212 dist = v.mag();
213 return dist;
214}
215
216///////////////////////////////////////////////////////////////////////////////
217//
219{
220 G4double ss = 0;
221
222 for (G4int i = 0; i <= 3; ++i)
223 {
224 G4double sp = GetVertex(i).dot(axis);
225 if (sp > ss) ss = sp;
226 }
227 return ss;
228}
229
230///////////////////////////////////////////////////////////////////////////////
231//
233 const G4ThreeVector &v,
234 G4bool outgoing,
235 G4double &distance,
236 G4double &distFromSurface,
237 G4ThreeVector &normal)
238{
239 G4bool intersect =
240 fFacet1.Intersect(p,v,outgoing,distance,distFromSurface,normal);
241 if (!intersect) intersect =
242 fFacet2.Intersect(p,v,outgoing,distance,distFromSurface,normal);
243 if (!intersect)
244 {
245 distance = distFromSurface = kInfinity;
246 normal.set(0,0,0);
247 }
248 return intersect;
249}
250
251///////////////////////////////////////////////////////////////////////////////
252//
253// Auxiliary method to get a random point on surface
254//
256{
257 G4ThreeVector pr = (G4RandFlat::shoot(0.,1.) < 0.5)
258 ? fFacet1.GetPointOnFace() : fFacet2.GetPointOnFace();
259 return pr;
260}
261
262///////////////////////////////////////////////////////////////////////////////
263//
264// Auxiliary method for returning the surface area
265//
267{
268 G4double area = fFacet1.GetArea() + fFacet2.GetArea();
269 return area;
270}
271
272///////////////////////////////////////////////////////////////////////////////
273//
275{
276 return "G4QuadrangularFacet";
277}
278
279///////////////////////////////////////////////////////////////////////////////
280//
282{
283 return fFacet1.GetSurfaceNormal();
284}
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4FacetVertexType
Definition: G4VFacet.hh:56
@ ABSOLUTE
Definition: G4VFacet.hh:56
G4DLLIMPORT std::ostream G4cerr
Hep3Vector unit() const
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
G4ThreeVector Distance(const G4ThreeVector &p)
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
G4QuadrangularFacet(const G4ThreeVector &Pt0, const G4ThreeVector &vt1, const G4ThreeVector &vt2, const G4ThreeVector &vt3, G4FacetVertexType)
void SetVertex(G4int i, const G4ThreeVector &val)
G4ThreeVector GetSurfaceNormal() const
G4QuadrangularFacet & operator=(const G4QuadrangularFacet &right)
G4ThreeVector GetVertex(G4int i) const
G4GeometryType GetEntityType() const
G4double Extent(const G4ThreeVector axis)
G4ThreeVector GetPointOnFace() const
void SetSurfaceNormal(G4ThreeVector normal)
G4ThreeVector GetPointOnFace() const
G4bool Intersect(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool outgoing, G4double &distance, G4double &distFromSurface, G4ThreeVector &normal)
G4ThreeVector Distance(const G4ThreeVector &p)
G4ThreeVector GetSurfaceNormal() const
static const G4double kCarTolerance
Definition: G4VFacet.hh:102
static const G4double dirTolerance
Definition: G4VFacet.hh:101
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41