Geant4 Cross Reference |
1 // 1 2 // ******************************************* 3 // * License and Disclaimer 4 // * 5 // * The Geant4 software is copyright of th 6 // * the Geant4 Collaboration. It is provided 7 // * conditions of the Geant4 Software License 8 // * LICENSE and available at http://cern.ch/ 9 // * include a list of copyright holders. 10 // * 11 // * Neither the authors of this software syst 12 // * institutes,nor the agencies providing fin 13 // * work make any representation or warran 14 // * regarding this software system or assum 15 // * use. Please see the license in the file 16 // * for the full disclaimer and the limitatio 17 // * 18 // * This code implementation is the result 19 // * technical work of the GEANT4 collaboratio 20 // * subject to DEFCON 705 IPR conditions. 21 // * By using, copying, modifying or distri 22 // * any work based on the software) you ag 23 // * use in resulting scientific publicati 24 // * acceptance of all terms of the Geant4 Sof 25 // ******************************************* 26 // 27 // G4QuadrangularFacet class implementation. 28 // 29 // 31 October 2004 P R Truscott, QinetiQ Ltd, 30 // 12 October 2012 M Gayer, CERN 31 // New implementation reducing 32 // and considerable CPU speedu 33 // implementation of G4Tessell 34 // 29 February 2016 E Tcherniaev, CERN 35 // Added exhaustive tests to c 36 // quadrangular facet: colline 37 // degenerate, concave or self 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 (cons 52 cons 53 cons 54 cons 55 56 { 57 G4double delta = 1.0 * kCarTolerance; // 58 G4double epsilon = 0.01 * kCarTolerance; // 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 || leng 94 diag1 <= delta || diag2 <= delta) 95 { 96 ostringstream message; 97 message << "Sides/diagonals of facet are t 98 << "P0 = " << GetVertex(0) << G4en 99 << "P1 = " << GetVertex(1) << G4en 100 << "P2 = " << GetVertex(2) << G4en 101 << "P3 = " << GetVertex(3) << G4en 102 << "Side1 length (P0->P1) = " << l 103 << "Side2 length (P1->P2) = " << l 104 << "Side3 length (P2->P3) = " << l 105 << "Side4 length (P3->P0) = " << l 106 << "Diagonal1 length (P0->P2) = " 107 << "Diagonal2 length (P1->P3) = " 108 G4Exception("G4QuadrangularFacet::G4Quadra 109 "GeomSolids1001", JustWarning, 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 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(leng 121 G4double h2 = 2.*s2 / std::max(std::max(leng 122 G4double h3 = 2.*s3 / std::max(std::max(leng 123 G4double h4 = 2.*s4 / std::max(std::max(leng 124 125 if (h1 <= delta || h2 <= delta || h3 <= delt 126 { 127 ostringstream message; 128 message << "Facet has three or more collin 129 << "P0 = " << GetVertex(0) << G4en 130 << "P1 = " << GetVertex(1) << G4en 131 << "P2 = " << GetVertex(2) << G4en 132 << "P3 = " << GetVertex(3) << G4en 133 << "Smallest heights:" << G4endl 134 << " in triangle P0-P1-P2 = " << 135 << " in triangle P1-P2-P3 = " << 136 << " in triangle P2-P3-P0 = " << 137 << " in triangle P3-P0-P1 = " << 138 G4Exception("G4QuadrangularFacet::G4Quadra 139 "GeomSolids1001", JustWarning, messa 140 return; 141 } 142 143 // Check that vertices are coplanar by compu 144 // height of tetrahedron comprising of verti 145 // 146 G4double smax = std::max( std::max(s1,s2), s 147 G4double hmin = 0.5 * std::fabs( e1.dot(e2.c 148 if (hmin >= epsilon) 149 { 150 ostringstream message; 151 message << "Facet is not planar." << G4end 152 << "Disrepancy = " << hmin << G4en 153 << "P0 = " << GetVertex(0) << G4en 154 << "P1 = " << GetVertex(1) << G4en 155 << "P2 = " << GetVertex(2) << G4en 156 << "P3 = " << GetVertex(3); 157 G4Exception("G4QuadrangularFacet::G4Quadra 158 "GeomSolids1001", JustWarning, messa 159 return; 160 } 161 162 // Check that facet is convex by computing c 163 // of diagonals 164 // 165 G4ThreeVector normal = e2.cross(e3-e1); 166 G4double s = kInfinity, t = kInfinity, magni 167 if (magnitude2 > delta*delta) // check: magn 168 { 169 s = normal.dot(e1.cross(e3-e1)) / magnitud 170 t = normal.dot(e1.cross(e2)) / magnitud 171 } 172 if (s <= 0. || s >= 1. || t <= 0. || t >= 1. 173 { 174 ostringstream message; 175 message << "Facet is not convex." << G4end 176 << "Parameters of crosspoint of di 177 << s << " and " << t << G4endl 178 << "should both be within (0,1) ra 179 << "P0 = " << GetVertex(0) << G4endl 180 << "P1 = " << GetVertex(1) << G4endl 181 << "P2 = " << GetVertex(2) << G4endl 182 << "P3 = " << GetVertex(3); 183 G4Exception("G4QuadrangularFacet::G4Quadra 184 "GeomSolids1001", JustWarning, messa 185 return; 186 } 187 188 // Define facet 189 // 190 fFacet1 = G4TriangularFacet(GetVertex(0),Get 191 fFacet2 = G4TriangularFacet(GetVertex(0),Get 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: comput 202 // of fCircumcenter and fRadius is wrong, ho 203 // it did not create any problem till now. 204 // Bizarre! Need to investigate! 205 } 206 207 ////////////////////////////////////////////// 208 // 209 G4QuadrangularFacet::~G4QuadrangularFacet () = 210 211 ////////////////////////////////////////////// 212 // 213 G4QuadrangularFacet::G4QuadrangularFacet (cons 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 G4Quadran 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( 240 GetVertex( 241 return c; 242 } 243 244 ////////////////////////////////////////////// 245 // 246 G4ThreeVector G4QuadrangularFacet::Distance (c 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 258 259 { 260 G4double dist = Distance(p).mag(); 261 return dist; 262 } 263 264 ////////////////////////////////////////////// 265 // 266 G4double G4QuadrangularFacet::Distance (const 267 const 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 G4 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 G 298 const G 299 G 300 G 301 G 302 G 303 { 304 G4bool intersect = 305 fFacet1.Intersect(p,v,outgoing,distance,di 306 if (!intersect) intersect = 307 fFacet2.Intersect(p,v,outgoing,distance,di 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 po 319 // 320 G4ThreeVector G4QuadrangularFacet::GetPointOnF 321 { 322 G4double s1 = fFacet1.GetArea(); 323 G4double s2 = fFacet2.GetArea(); 324 return ((s1+s2)*G4UniformRand() < s1) ? 325 fFacet1.GetPointOnFace() : fFacet2.GetPoin 326 } 327 328 ////////////////////////////////////////////// 329 // 330 // Auxiliary method for returning the surface 331 // 332 G4double G4QuadrangularFacet::GetArea() const 333 { 334 G4double area = fFacet1.GetArea() + fFacet2. 335 return area; 336 } 337 338 ////////////////////////////////////////////// 339 // 340 G4String G4QuadrangularFacet::GetEntityType () 341 { 342 return "G4QuadrangularFacet"; 343 } 344 345 ////////////////////////////////////////////// 346 // 347 G4ThreeVector G4QuadrangularFacet::GetSurfaceN 348 { 349 return fFacet1.GetSurfaceNormal(); 350 } 351