Geant4 Cross Reference |
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 class implementation. 28 // 29 // 31 October 2004 P R Truscott, QinetiQ Ltd, UK - Created. 30 // 12 October 2012 M Gayer, CERN 31 // New implementation reducing memory requirements by 50%, 32 // and considerable CPU speedup together with the new 33 // implementation of G4TessellatedSolid. 34 // 29 February 2016 E Tcherniaev, CERN 35 // Added exhaustive tests to catch various problems with a 36 // quadrangular facet: collinear vertices, non planar surface, 37 // degenerate, concave or self intersecting quadrilateral. 38 // -------------------------------------------------------------------- 39 40 #include "G4QuadrangularFacet.hh" 41 #include "geomdefs.hh" 42 #include "Randomize.hh" 43 44 using namespace std; 45 46 /////////////////////////////////////////////////////////////////////////////// 47 // 48 // Constructing two adjacent G4TriangularFacet 49 // Not efficient, but practical... 50 // 51 G4QuadrangularFacet::G4QuadrangularFacet (const G4ThreeVector& vt0, 52 const G4ThreeVector& vt1, 53 const G4ThreeVector& vt2, 54 const G4ThreeVector& vt3, 55 G4FacetVertexType vertexType) 56 { 57 G4double delta = 1.0 * kCarTolerance; // dimension tolerance 58 G4double epsilon = 0.01 * kCarTolerance; // planarity tolerance 59 60 G4ThreeVector e1, e2, e3; 61 SetVertex(0, vt0); 62 if (vertexType == ABSOLUTE) 63 { 64 SetVertex(1, vt1); 65 SetVertex(2, vt2); 66 SetVertex(3, vt3); 67 68 e1 = vt1 - vt0; 69 e2 = vt2 - vt0; 70 e3 = vt3 - vt0; 71 } 72 else 73 { 74 SetVertex(1, vt0 + vt1); 75 SetVertex(2, vt0 + vt2); 76 SetVertex(3, vt0 + vt3); 77 78 e1 = vt1; 79 e2 = vt2; 80 e3 = vt3; 81 } 82 83 // Check length of sides and diagonals 84 // 85 G4double leng1 = e1.mag(); 86 G4double leng2 = (e2-e1).mag(); 87 G4double leng3 = (e3-e2).mag(); 88 G4double leng4 = e3.mag(); 89 90 G4double diag1 = e2.mag(); 91 G4double diag2 = (e3-e1).mag(); 92 93 if (leng1 <= delta || leng2 <= delta || leng3 <= delta || leng4 <= delta || 94 diag1 <= delta || diag2 <= delta) 95 { 96 ostringstream message; 97 message << "Sides/diagonals of facet are too small." << G4endl 98 << "P0 = " << GetVertex(0) << G4endl 99 << "P1 = " << GetVertex(1) << G4endl 100 << "P2 = " << GetVertex(2) << G4endl 101 << "P3 = " << GetVertex(3) << G4endl 102 << "Side1 length (P0->P1) = " << leng1 << G4endl 103 << "Side2 length (P1->P2) = " << leng2 << G4endl 104 << "Side3 length (P2->P3) = " << leng3 << G4endl 105 << "Side4 length (P3->P0) = " << leng4 << G4endl 106 << "Diagonal1 length (P0->P2) = " << diag1 << G4endl 107 << "Diagonal2 length (P1->P3) = " << diag2; 108 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()", 109 "GeomSolids1001", JustWarning, message); 110 return; 111 } 112 113 // Check that vertices are not collinear 114 // 115 G4double s1 = (e1.cross(e2)).mag()*0.5; 116 G4double s2 = ((e2-e1).cross(e3-e2)).mag()*0.5; 117 G4double s3 = (e2.cross(e3)).mag()*0.5; 118 G4double s4 = (e1.cross(e3)).mag()*0.5; 119 120 G4double h1 = 2.*s1 / std::max(std::max(leng1,leng2),diag1); 121 G4double h2 = 2.*s2 / std::max(std::max(leng2,leng3),diag2); 122 G4double h3 = 2.*s3 / std::max(std::max(leng3,leng4),diag1); 123 G4double h4 = 2.*s4 / std::max(std::max(leng4,leng1),diag2); 124 125 if (h1 <= delta || h2 <= delta || h3 <= delta || h4 <= delta ) 126 { 127 ostringstream message; 128 message << "Facet has three or more collinear vertices." << G4endl 129 << "P0 = " << GetVertex(0) << G4endl 130 << "P1 = " << GetVertex(1) << G4endl 131 << "P2 = " << GetVertex(2) << G4endl 132 << "P3 = " << GetVertex(3) << G4endl 133 << "Smallest heights:" << G4endl 134 << " in triangle P0-P1-P2 = " << h1 << G4endl 135 << " in triangle P1-P2-P3 = " << h2 << G4endl 136 << " in triangle P2-P3-P0 = " << h3 << G4endl 137 << " in triangle P3-P0-P1 = " << h4; 138 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()", 139 "GeomSolids1001", JustWarning, message); 140 return; 141 } 142 143 // Check that vertices are coplanar by computing minimal 144 // height of tetrahedron comprising of vertices 145 // 146 G4double smax = std::max( std::max(s1,s2), std::max(s3,s4) ); 147 G4double hmin = 0.5 * std::fabs( e1.dot(e2.cross(e3)) ) / smax; 148 if (hmin >= epsilon) 149 { 150 ostringstream message; 151 message << "Facet is not planar." << G4endl 152 << "Disrepancy = " << hmin << G4endl 153 << "P0 = " << GetVertex(0) << G4endl 154 << "P1 = " << GetVertex(1) << G4endl 155 << "P2 = " << GetVertex(2) << G4endl 156 << "P3 = " << GetVertex(3); 157 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()", 158 "GeomSolids1001", JustWarning, message); 159 return; 160 } 161 162 // Check that facet is convex by computing crosspoint 163 // of diagonals 164 // 165 G4ThreeVector normal = e2.cross(e3-e1); 166 G4double s = kInfinity, t = kInfinity, magnitude2 = normal.mag2(); 167 if (magnitude2 > delta*delta) // check: magnitude2 != 0. 168 { 169 s = normal.dot(e1.cross(e3-e1)) / magnitude2; 170 t = normal.dot(e1.cross(e2)) / magnitude2; 171 } 172 if (s <= 0. || s >= 1. || t <= 0. || t >= 1.) 173 { 174 ostringstream message; 175 message << "Facet is not convex." << G4endl 176 << "Parameters of crosspoint of diagonals: " 177 << s << " and " << t << G4endl 178 << "should both be within (0,1) range" << G4endl 179 << "P0 = " << GetVertex(0) << G4endl 180 << "P1 = " << GetVertex(1) << G4endl 181 << "P2 = " << GetVertex(2) << G4endl 182 << "P3 = " << GetVertex(3); 183 G4Exception("G4QuadrangularFacet::G4QuadrangularFacet()", 184 "GeomSolids1001", JustWarning, message); 185 return; 186 } 187 188 // Define facet 189 // 190 fFacet1 = G4TriangularFacet(GetVertex(0),GetVertex(1),GetVertex(2),ABSOLUTE); 191 fFacet2 = G4TriangularFacet(GetVertex(0),GetVertex(2),GetVertex(3),ABSOLUTE); 192 193 normal = normal.unit(); 194 fFacet1.SetSurfaceNormal(normal); 195 fFacet2.SetSurfaceNormal(normal); 196 197 G4ThreeVector vtmp = 0.5 * (e1 + e2); 198 fCircumcentre = GetVertex(0) + vtmp; 199 G4double radiusSqr = vtmp.mag2(); 200 fRadius = std::sqrt(radiusSqr); 201 // 29.02.2016 Remark by E.Tcherniaev: computation 202 // of fCircumcenter and fRadius is wrong, however 203 // it did not create any problem till now. 204 // Bizarre! Need to investigate! 205 } 206 207 /////////////////////////////////////////////////////////////////////////////// 208 // 209 G4QuadrangularFacet::~G4QuadrangularFacet () = default; 210 211 /////////////////////////////////////////////////////////////////////////////// 212 // 213 G4QuadrangularFacet::G4QuadrangularFacet (const G4QuadrangularFacet& rhs) 214 : G4VFacet(rhs) 215 { 216 fFacet1 = rhs.fFacet1; 217 fFacet2 = rhs.fFacet2; 218 fRadius = 0.0; 219 } 220 221 /////////////////////////////////////////////////////////////////////////////// 222 // 223 G4QuadrangularFacet & 224 G4QuadrangularFacet::operator=(const G4QuadrangularFacet& rhs) 225 { 226 if (this == &rhs) return *this; 227 228 fFacet1 = rhs.fFacet1; 229 fFacet2 = rhs.fFacet2; 230 fRadius = 0.0; 231 232 return *this; 233 } 234 235 /////////////////////////////////////////////////////////////////////////////// 236 // 237 G4VFacet* G4QuadrangularFacet::GetClone () 238 { 239 auto c = new G4QuadrangularFacet (GetVertex(0), GetVertex(1), 240 GetVertex(2), GetVertex(3), ABSOLUTE); 241 return c; 242 } 243 244 /////////////////////////////////////////////////////////////////////////////// 245 // 246 G4ThreeVector G4QuadrangularFacet::Distance (const G4ThreeVector& p) 247 { 248 G4ThreeVector v1 = fFacet1.Distance(p); 249 G4ThreeVector v2 = fFacet2.Distance(p); 250 251 if (v1.mag2() < v2.mag2()) return v1; 252 else return v2; 253 } 254 255 /////////////////////////////////////////////////////////////////////////////// 256 // 257 G4double G4QuadrangularFacet::Distance (const G4ThreeVector& p, 258 G4double) 259 { 260 G4double dist = Distance(p).mag(); 261 return dist; 262 } 263 264 /////////////////////////////////////////////////////////////////////////////// 265 // 266 G4double G4QuadrangularFacet::Distance (const G4ThreeVector& p, G4double, 267 const G4bool outgoing) 268 { 269 G4double dist; 270 271 G4ThreeVector v = Distance(p); 272 G4double dir = v.dot(GetSurfaceNormal()); 273 if ( ((dir > dirTolerance) && (!outgoing)) 274 || ((dir < -dirTolerance) && outgoing)) 275 dist = kInfinity; 276 else 277 dist = v.mag(); 278 return dist; 279 } 280 281 /////////////////////////////////////////////////////////////////////////////// 282 // 283 G4double G4QuadrangularFacet::Extent (const G4ThreeVector axis) 284 { 285 G4double ss = 0; 286 287 for (G4int i = 0; i <= 3; ++i) 288 { 289 G4double sp = GetVertex(i).dot(axis); 290 if (sp > ss) ss = sp; 291 } 292 return ss; 293 } 294 295 /////////////////////////////////////////////////////////////////////////////// 296 // 297 G4bool G4QuadrangularFacet::Intersect (const G4ThreeVector& p, 298 const G4ThreeVector& v, 299 G4bool outgoing, 300 G4double& distance, 301 G4double& distFromSurface, 302 G4ThreeVector& normal) 303 { 304 G4bool intersect = 305 fFacet1.Intersect(p,v,outgoing,distance,distFromSurface,normal); 306 if (!intersect) intersect = 307 fFacet2.Intersect(p,v,outgoing,distance,distFromSurface,normal); 308 if (!intersect) 309 { 310 distance = distFromSurface = kInfinity; 311 normal.set(0,0,0); 312 } 313 return intersect; 314 } 315 316 /////////////////////////////////////////////////////////////////////////////// 317 // 318 // Auxiliary method to get a uniform random point on the facet 319 // 320 G4ThreeVector G4QuadrangularFacet::GetPointOnFace() const 321 { 322 G4double s1 = fFacet1.GetArea(); 323 G4double s2 = fFacet2.GetArea(); 324 return ((s1+s2)*G4UniformRand() < s1) ? 325 fFacet1.GetPointOnFace() : fFacet2.GetPointOnFace(); 326 } 327 328 /////////////////////////////////////////////////////////////////////////////// 329 // 330 // Auxiliary method for returning the surface area 331 // 332 G4double G4QuadrangularFacet::GetArea() const 333 { 334 G4double area = fFacet1.GetArea() + fFacet2.GetArea(); 335 return area; 336 } 337 338 /////////////////////////////////////////////////////////////////////////////// 339 // 340 G4String G4QuadrangularFacet::GetEntityType () const 341 { 342 return "G4QuadrangularFacet"; 343 } 344 345 /////////////////////////////////////////////////////////////////////////////// 346 // 347 G4ThreeVector G4QuadrangularFacet::GetSurfaceNormal () const 348 { 349 return fFacet1.GetSurfaceNormal(); 350 } 351