| 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 | // $Id: G4TriangularFacet.cc,v 1.13 2010/04/28 16:21:21 flei Exp $
|
|---|
| 28 | // GEANT4 tag $Name: geant4-09-04-beta-01 $
|
|---|
| 29 | //
|
|---|
| 30 | // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
|---|
| 31 | //
|
|---|
| 32 | // MODULE: G4TriangularFacet.cc
|
|---|
| 33 | //
|
|---|
| 34 | // Date: 15/06/2005
|
|---|
| 35 | // Author: P R Truscott
|
|---|
| 36 | // Organisation: QinetiQ Ltd, UK
|
|---|
| 37 | // Customer: UK Ministry of Defence : RAO CRP TD Electronic Systems
|
|---|
| 38 | // Contract: C/MAT/N03517
|
|---|
| 39 | //
|
|---|
| 40 | // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
|---|
| 41 | //
|
|---|
| 42 | // CHANGE HISTORY
|
|---|
| 43 | // --------------
|
|---|
| 44 | //
|
|---|
| 45 | // 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
|
|---|
| 46 | //
|
|---|
| 47 | // 01 August 2007 P R Truscott, QinetiQ Ltd, UK
|
|---|
| 48 | // Significant modification to correct for errors and enhance
|
|---|
| 49 | // based on patches/observations kindly provided by Rickard
|
|---|
| 50 | // Holmberg
|
|---|
| 51 | //
|
|---|
| 52 | // 26 September 2007
|
|---|
| 53 | // P R Truscott, QinetiQ Ltd, UK
|
|---|
| 54 | // Further chamges implemented to the Intersect member
|
|---|
| 55 | // function to correctly treat rays nearly parallel to the
|
|---|
| 56 | // plane of the triangle.
|
|---|
| 57 | //
|
|---|
| 58 | // 12 April 2010 P R Truscott, QinetiQ, bug fixes to treat optical
|
|---|
| 59 | // photon transport, in particular internal reflection
|
|---|
| 60 | // at surface.
|
|---|
| 61 | // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
|
|---|
| 62 |
|
|---|
| 63 | #include "G4TriangularFacet.hh"
|
|---|
| 64 | #include "G4TwoVector.hh"
|
|---|
| 65 | #include "globals.hh"
|
|---|
| 66 | #include "Randomize.hh"
|
|---|
| 67 |
|
|---|
| 68 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 69 | //
|
|---|
| 70 | // Definition of triangular facet using absolute vectors to vertices.
|
|---|
| 71 | // From this for first vector is retained to define the facet location and
|
|---|
| 72 | // two relative vectors (E0 and E1) define the sides and orientation of
|
|---|
| 73 | // the outward surface normal.
|
|---|
| 74 | //
|
|---|
| 75 | G4TriangularFacet::G4TriangularFacet (const G4ThreeVector Pt0,
|
|---|
| 76 | const G4ThreeVector vt1, const G4ThreeVector vt2,
|
|---|
| 77 | G4FacetVertexType vertexType)
|
|---|
| 78 | : G4VFacet()
|
|---|
| 79 | {
|
|---|
| 80 | tGeomAlg = G4TessellatedGeometryAlgorithms::GetInstance();
|
|---|
| 81 | P0 = Pt0;
|
|---|
| 82 | nVertices = 3;
|
|---|
| 83 | if (vertexType == ABSOLUTE)
|
|---|
| 84 | {
|
|---|
| 85 | P.push_back(vt1);
|
|---|
| 86 | P.push_back(vt2);
|
|---|
| 87 |
|
|---|
| 88 | E.push_back(vt1 - P0);
|
|---|
| 89 | E.push_back(vt2 - P0);
|
|---|
| 90 | }
|
|---|
| 91 | else
|
|---|
| 92 | {
|
|---|
| 93 | P.push_back(P0 + vt1);
|
|---|
| 94 | P.push_back(P0 + vt2);
|
|---|
| 95 |
|
|---|
| 96 | E.push_back(vt1);
|
|---|
| 97 | E.push_back(vt2);
|
|---|
| 98 | }
|
|---|
| 99 |
|
|---|
| 100 | G4double Emag1 = E[0].mag();
|
|---|
| 101 | G4double Emag2 = E[1].mag();
|
|---|
| 102 | G4double Emag3 = (E[1]-E[0]).mag();
|
|---|
| 103 |
|
|---|
| 104 | if (Emag1 <= kCarTolerance || Emag2 <= kCarTolerance ||
|
|---|
| 105 | Emag3 <= kCarTolerance)
|
|---|
| 106 | {
|
|---|
| 107 | G4Exception("G4TriangularFacet::G4TriangularFacet()", "InvalidSetup",
|
|---|
| 108 | JustWarning, "Length of sides of facet are too small.");
|
|---|
| 109 | G4cerr << G4endl;
|
|---|
| 110 | G4cerr << "P0 = " << P0 << G4endl;
|
|---|
| 111 | G4cerr << "P1 = " << P[0] << G4endl;
|
|---|
| 112 | G4cerr << "P2 = " << P[1] << G4endl;
|
|---|
| 113 | G4cerr << "Side lengths = P0->P1" << Emag1 << G4endl;
|
|---|
| 114 | G4cerr << "Side lengths = P0->P2" << Emag2 << G4endl;
|
|---|
| 115 | G4cerr << "Side lengths = P1->P2" << Emag3 << G4endl;
|
|---|
| 116 | G4cerr << G4endl;
|
|---|
| 117 |
|
|---|
| 118 | isDefined = false;
|
|---|
| 119 | geometryType = "G4TriangularFacet";
|
|---|
| 120 | surfaceNormal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 121 | a = 0.0;
|
|---|
| 122 | b = 0.0;
|
|---|
| 123 | c = 0.0;
|
|---|
| 124 | det = 0.0;
|
|---|
| 125 | }
|
|---|
| 126 | else
|
|---|
| 127 | {
|
|---|
| 128 | isDefined = true;
|
|---|
| 129 | geometryType = "G4TriangularFacet";
|
|---|
| 130 | surfaceNormal = E[0].cross(E[1]).unit();
|
|---|
| 131 | a = E[0].mag2();
|
|---|
| 132 | b = E[0].dot(E[1]);
|
|---|
| 133 | c = E[1].mag2();
|
|---|
| 134 | det = std::abs(a*c - b*b);
|
|---|
| 135 |
|
|---|
| 136 | sMin = -0.5*kCarTolerance/std::sqrt(a);
|
|---|
| 137 | sMax = 1.0 - sMin;
|
|---|
| 138 | tMin = -0.5*kCarTolerance/std::sqrt(c);
|
|---|
| 139 |
|
|---|
| 140 | area = 0.5 * (E[0].cross(E[1])).mag();
|
|---|
| 141 |
|
|---|
| 142 | // G4ThreeVector vtmp = 0.25 * (E[0] + E[1]);
|
|---|
| 143 | G4double lambda0 = (a-b) * c / (8.0*area*area);
|
|---|
| 144 | G4double lambda1 = (c-b) * a / (8.0*area*area);
|
|---|
| 145 | circumcentre = P0 + lambda0*E[0] + lambda1*E[1];
|
|---|
| 146 | radiusSqr = (circumcentre-P0).mag2();
|
|---|
| 147 | radius = std::sqrt(radiusSqr);
|
|---|
| 148 |
|
|---|
| 149 | for (size_t i=0; i<3; i++) { I.push_back(0); }
|
|---|
| 150 | }
|
|---|
| 151 | }
|
|---|
| 152 |
|
|---|
| 153 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 154 | //
|
|---|
| 155 | // ~G4TriangularFacet
|
|---|
| 156 | //
|
|---|
| 157 | // A pretty boring destructor indeed!
|
|---|
| 158 | //
|
|---|
| 159 | G4TriangularFacet::~G4TriangularFacet ()
|
|---|
| 160 | {
|
|---|
| 161 | P.clear();
|
|---|
| 162 | E.clear();
|
|---|
| 163 | I.clear();
|
|---|
| 164 | }
|
|---|
| 165 |
|
|---|
| 166 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 167 | //
|
|---|
| 168 | // GetClone
|
|---|
| 169 | //
|
|---|
| 170 | // Simple member function to generate a diplicate of the triangular facet.
|
|---|
| 171 | //
|
|---|
| 172 | G4VFacet *G4TriangularFacet::GetClone ()
|
|---|
| 173 | {
|
|---|
| 174 | G4TriangularFacet *fc = new G4TriangularFacet (P0, P[0], P[1], ABSOLUTE);
|
|---|
| 175 | G4VFacet *cc = 0;
|
|---|
| 176 | cc = fc;
|
|---|
| 177 | return cc;
|
|---|
| 178 | }
|
|---|
| 179 |
|
|---|
| 180 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 181 | //
|
|---|
| 182 | // GetFlippedFacet
|
|---|
| 183 | //
|
|---|
| 184 | // Member function to generate an identical facet, but with the normal vector
|
|---|
| 185 | // pointing at 180 degrees.
|
|---|
| 186 | //
|
|---|
| 187 | G4TriangularFacet *G4TriangularFacet::GetFlippedFacet ()
|
|---|
| 188 | {
|
|---|
| 189 | G4TriangularFacet *flipped = new G4TriangularFacet (P0, P[1], P[0], ABSOLUTE);
|
|---|
| 190 | return flipped;
|
|---|
| 191 | }
|
|---|
| 192 |
|
|---|
| 193 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 194 | //
|
|---|
| 195 | // Distance (G4ThreeVector)
|
|---|
| 196 | //
|
|---|
| 197 | // Determines the vector between p and the closest point on the facet to p.
|
|---|
| 198 | // This is based on the algorithm published in "Geometric Tools for Computer
|
|---|
| 199 | // Graphics," Philip J Scheider and David H Eberly, Elsevier Science (USA),
|
|---|
| 200 | // 2003. at the time of writing, the algorithm is also available in a
|
|---|
| 201 | // technical note "Distance between point and triangle in 3D," by David Eberly
|
|---|
| 202 | // at http://www.geometrictools.com/Documentation/DistancePoint3Triangle3.pdf
|
|---|
| 203 | //
|
|---|
| 204 | // The by-product is the square-distance sqrDist, which is retained
|
|---|
| 205 | // in case needed by the other "Distance" member functions.
|
|---|
| 206 | //
|
|---|
| 207 | G4ThreeVector G4TriangularFacet::Distance (const G4ThreeVector &p)
|
|---|
| 208 | {
|
|---|
| 209 | G4ThreeVector D = P0 - p;
|
|---|
| 210 | G4double d = E[0].dot(D);
|
|---|
| 211 | G4double e = E[1].dot(D);
|
|---|
| 212 | G4double f = D.mag2();
|
|---|
| 213 | G4double s = b*e - c*d;
|
|---|
| 214 | G4double t = b*d - a*e;
|
|---|
| 215 |
|
|---|
| 216 | sqrDist = 0.0;
|
|---|
| 217 |
|
|---|
| 218 | if (s+t <= det)
|
|---|
| 219 | {
|
|---|
| 220 | if (s < 0.0)
|
|---|
| 221 | {
|
|---|
| 222 | if (t < 0.0)
|
|---|
| 223 | {
|
|---|
| 224 | //
|
|---|
| 225 | // We are in region 4.
|
|---|
| 226 | //
|
|---|
| 227 | if (d < 0.0)
|
|---|
| 228 | {
|
|---|
| 229 | t = 0.0;
|
|---|
| 230 | if (-d >= a) {s = 1.0; sqrDist = a + 2.0*d + f;}
|
|---|
| 231 | else {s = -d/a; sqrDist = d*s + f;}
|
|---|
| 232 | }
|
|---|
| 233 | else
|
|---|
| 234 | {
|
|---|
| 235 | s = 0.0;
|
|---|
| 236 | if (e >= 0.0) {t = 0.0; sqrDist = f;}
|
|---|
| 237 | else if (-e >= c) {t = 1.0; sqrDist = c + 2.0*e + f;}
|
|---|
| 238 | else {t = -e/c; sqrDist = e*t + f;}
|
|---|
| 239 | }
|
|---|
| 240 | }
|
|---|
| 241 | else
|
|---|
| 242 | {
|
|---|
| 243 | //
|
|---|
| 244 | // We are in region 3.
|
|---|
| 245 | //
|
|---|
| 246 | s = 0.0;
|
|---|
| 247 | if (e >= 0.0) {t = 0.0; sqrDist = f;}
|
|---|
| 248 | else if (-e >= c) {t = 1.0; sqrDist = c + 2.0*e + f;}
|
|---|
| 249 | else {t = -e/c; sqrDist = e*t + f;}
|
|---|
| 250 | }
|
|---|
| 251 | }
|
|---|
| 252 | else if (t < 0.0)
|
|---|
| 253 | {
|
|---|
| 254 | //
|
|---|
| 255 | // We are in region 5.
|
|---|
| 256 | //
|
|---|
| 257 | t = 0.0;
|
|---|
| 258 | if (d >= 0.0) {s = 0.0; sqrDist = f;}
|
|---|
| 259 | else if (-d >= a) {s = 1.0; sqrDist = a + 2.0*d + f;}
|
|---|
| 260 | else {s = -d/a; sqrDist = d*s + f;}
|
|---|
| 261 | }
|
|---|
| 262 | else
|
|---|
| 263 | {
|
|---|
| 264 | //
|
|---|
| 265 | // We are in region 0.
|
|---|
| 266 | //
|
|---|
| 267 | s = s / det;
|
|---|
| 268 | t = t / det;
|
|---|
| 269 | sqrDist = s*(a*s + b*t + 2.0*d) + t*(b*s + c*t + 2.0*e) + f;
|
|---|
| 270 | }
|
|---|
| 271 | }
|
|---|
| 272 | else
|
|---|
| 273 | {
|
|---|
| 274 | if (s < 0.0)
|
|---|
| 275 | {
|
|---|
| 276 | //
|
|---|
| 277 | // We are in region 2.
|
|---|
| 278 | //
|
|---|
| 279 | G4double tmp0 = b + d;
|
|---|
| 280 | G4double tmp1 = c + e;
|
|---|
| 281 | if (tmp1 > tmp0)
|
|---|
| 282 | {
|
|---|
| 283 | G4double numer = tmp1 - tmp0;
|
|---|
| 284 | G4double denom = a - 2.0*b + c;
|
|---|
| 285 | if (numer >= denom) {s = 1.0; t = 0.0; sqrDist = a + 2.0*d + f;}
|
|---|
| 286 | else
|
|---|
| 287 | {
|
|---|
| 288 | s = numer/denom;
|
|---|
| 289 | t = 1.0 - s;
|
|---|
| 290 | sqrDist = s*(a*s + b*t +2.0*d) + t*(b*s + c*t + 2.0*e) + f;
|
|---|
| 291 | }
|
|---|
| 292 | }
|
|---|
| 293 | else
|
|---|
| 294 | {
|
|---|
| 295 | s = 0.0;
|
|---|
| 296 | if (tmp1 <= 0.0) {t = 1.0; sqrDist = c + 2.0*e + f;}
|
|---|
| 297 | else if (e >= 0.0) {t = 0.0; sqrDist = f;}
|
|---|
| 298 | else {t = -e/c; sqrDist = e*t + f;}
|
|---|
| 299 | }
|
|---|
| 300 | }
|
|---|
| 301 | else if (t < 0.0)
|
|---|
| 302 | {
|
|---|
| 303 | //
|
|---|
| 304 | // We are in region 6.
|
|---|
| 305 | //
|
|---|
| 306 | G4double tmp0 = b + e;
|
|---|
| 307 | G4double tmp1 = a + d;
|
|---|
| 308 | if (tmp1 > tmp0)
|
|---|
| 309 | {
|
|---|
| 310 | G4double numer = tmp1 - tmp0;
|
|---|
| 311 | G4double denom = a - 2.0*b + c;
|
|---|
| 312 | if (numer >= denom) {t = 1.0; s = 0.0; sqrDist = c + 2.0*e + f;}
|
|---|
| 313 | else
|
|---|
| 314 | {
|
|---|
| 315 | t = numer/denom;
|
|---|
| 316 | s = 1.0 - t;
|
|---|
| 317 | sqrDist = s*(a*s + b*t +2.0*d) + t*(b*s + c*t + 2.0*e) + f;
|
|---|
| 318 | }
|
|---|
| 319 | }
|
|---|
| 320 | else
|
|---|
| 321 | {
|
|---|
| 322 | t = 0.0;
|
|---|
| 323 | if (tmp1 <= 0.0) {s = 1.0; sqrDist = a + 2.0*d + f;}
|
|---|
| 324 | else if (d >= 0.0) {s = 0.0; sqrDist = f;}
|
|---|
| 325 | else {s = -d/a; sqrDist = d*s + f;}
|
|---|
| 326 | }
|
|---|
| 327 | }
|
|---|
| 328 | else
|
|---|
| 329 | //
|
|---|
| 330 | // We are in region 1.
|
|---|
| 331 | //
|
|---|
| 332 | {
|
|---|
| 333 | G4double numer = c + e - b - d;
|
|---|
| 334 | if (numer <= 0.0)
|
|---|
| 335 | {
|
|---|
| 336 | s = 0.0;
|
|---|
| 337 | t = 1.0;
|
|---|
| 338 | sqrDist = c + 2.0*e + f;
|
|---|
| 339 | }
|
|---|
| 340 | else
|
|---|
| 341 | {
|
|---|
| 342 | G4double denom = a - 2.0*b + c;
|
|---|
| 343 | if (numer >= denom) {s = 1.0; t = 0.0; sqrDist = a + 2.0*d + f;}
|
|---|
| 344 | else
|
|---|
| 345 | {
|
|---|
| 346 | s = numer/denom;
|
|---|
| 347 | t = 1.0 - s;
|
|---|
| 348 | sqrDist = s*(a*s + b*t + 2.0*d) + t*(b*s + c*t + 2.0*e) + f;
|
|---|
| 349 | }
|
|---|
| 350 | }
|
|---|
| 351 | }
|
|---|
| 352 | }
|
|---|
| 353 | //
|
|---|
| 354 | //
|
|---|
| 355 | // Do a check for rounding errors in the distance-squared. It appears that
|
|---|
| 356 | // the conventional methods for calculating sqrDist breaks down when very near
|
|---|
| 357 | // to or at the surface (as required by transport). We'll therefore also use
|
|---|
| 358 | // the magnitude-squared of the vector displacement. (Note that I've also
|
|---|
| 359 | // tried to get around this problem by using the existing equations for
|
|---|
| 360 | //
|
|---|
| 361 | // sqrDist = function(a,b,c,d,s,t)
|
|---|
| 362 | //
|
|---|
| 363 | // and use a more accurate addition process which minimises errors and
|
|---|
| 364 | // breakdown of cummutitivity [where (A+B)+C != A+(B+C)] but this still
|
|---|
| 365 | // doesn't work.
|
|---|
| 366 | // Calculation from u = D + s*E[0] + t*E[1] is less efficient, but appears
|
|---|
| 367 | // more robust.
|
|---|
| 368 | //
|
|---|
| 369 | if (sqrDist < 0.0) { sqrDist = 0.0; }
|
|---|
| 370 | G4ThreeVector u = D + s*E[0] + t*E[1];
|
|---|
| 371 | G4double u2 = u.mag2();
|
|---|
| 372 | //
|
|---|
| 373 | //
|
|---|
| 374 | // The following (part of the roundoff error check) is from Oliver Merle's
|
|---|
| 375 | // updates.
|
|---|
| 376 | //
|
|---|
| 377 | if ( sqrDist > u2 ) sqrDist = u2;
|
|---|
| 378 |
|
|---|
| 379 | return u;
|
|---|
| 380 | }
|
|---|
| 381 |
|
|---|
| 382 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 383 | //
|
|---|
| 384 | // Distance (G4ThreeVector, G4double)
|
|---|
| 385 | //
|
|---|
| 386 | // Determines the closest distance between point p and the facet. This makes
|
|---|
| 387 | // use of G4ThreeVector G4TriangularFacet::Distance, which stores the
|
|---|
| 388 | // square of the distance in variable sqrDist. If approximate methods show
|
|---|
| 389 | // the distance is to be greater than minDist, then forget about further
|
|---|
| 390 | // computation and return a very large number.
|
|---|
| 391 | //
|
|---|
| 392 | G4double G4TriangularFacet::Distance (const G4ThreeVector &p,
|
|---|
| 393 | const G4double minDist)
|
|---|
| 394 | {
|
|---|
| 395 | //
|
|---|
| 396 | //
|
|---|
| 397 | // Start with quicky test to determine if the surface of the sphere enclosing
|
|---|
| 398 | // the triangle is any closer to p than minDist. If not, then don't bother
|
|---|
| 399 | // about more accurate test.
|
|---|
| 400 | //
|
|---|
| 401 | G4double dist = kInfinity;
|
|---|
| 402 | if ((p-circumcentre).mag()-radius < minDist)
|
|---|
| 403 | {
|
|---|
| 404 | //
|
|---|
| 405 | //
|
|---|
| 406 | // It's possible that the triangle is closer than minDist, so do more accurate
|
|---|
| 407 | // assessment.
|
|---|
| 408 | //
|
|---|
| 409 | dist = Distance(p).mag();
|
|---|
| 410 | // dist = std::sqrt(sqrDist);
|
|---|
| 411 | }
|
|---|
| 412 |
|
|---|
| 413 | return dist;
|
|---|
| 414 | }
|
|---|
| 415 |
|
|---|
| 416 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 417 | //
|
|---|
| 418 | // Distance (G4ThreeVector, G4double, G4double)
|
|---|
| 419 | //
|
|---|
| 420 | // Determine the distance to point p. kInfinity is returned if either:
|
|---|
| 421 | // (1) outgoing is TRUE and the dot product of the normal vector to the facet
|
|---|
| 422 | // and the displacement vector from p to the triangle is negative.
|
|---|
| 423 | // (2) outgoing is FALSE and the dot product of the normal vector to the facet
|
|---|
| 424 | // and the displacement vector from p to the triangle is positive.
|
|---|
| 425 | // If approximate methods show the distance is to be greater than minDist, then
|
|---|
| 426 | // forget about further computation and return a very large number.
|
|---|
| 427 | //
|
|---|
| 428 | // This method has been heavily modified thanks to the valuable comments and
|
|---|
| 429 | // corrections of Rickard Holmberg.
|
|---|
| 430 | //
|
|---|
| 431 | G4double G4TriangularFacet::Distance (const G4ThreeVector &p,
|
|---|
| 432 | const G4double minDist, const G4bool outgoing)
|
|---|
| 433 | {
|
|---|
| 434 | //
|
|---|
| 435 | //
|
|---|
| 436 | // Start with quicky test to determine if the surface of the sphere enclosing
|
|---|
| 437 | // the triangle is any closer to p than minDist. If not, then don't bother
|
|---|
| 438 | // about more accurate test.
|
|---|
| 439 | //
|
|---|
| 440 | G4double dist = kInfinity;
|
|---|
| 441 | if ((p-circumcentre).mag()-radius < minDist)
|
|---|
| 442 | {
|
|---|
| 443 | //
|
|---|
| 444 | //
|
|---|
| 445 | // It's possible that the triangle is closer than minDist, so do more accurate
|
|---|
| 446 | // assessment.
|
|---|
| 447 | //
|
|---|
| 448 | G4ThreeVector v = Distance(p);
|
|---|
| 449 | G4double dist1 = std::sqrt(sqrDist);
|
|---|
| 450 | G4double dir = v.dot(surfaceNormal);
|
|---|
| 451 | G4bool wrongSide = (dir > 0.0 && !outgoing) || (dir < 0.0 && outgoing);
|
|---|
| 452 | if (dist1 <= kCarTolerance*0.5)
|
|---|
| 453 | {
|
|---|
| 454 | //
|
|---|
| 455 | //
|
|---|
| 456 | // Point p is very close to triangle. Check if it's on the wrong side, in
|
|---|
| 457 | // which case return distance of 0.0 otherwise .
|
|---|
| 458 | //
|
|---|
| 459 | if (wrongSide) dist = 0.0;
|
|---|
| 460 | else dist = dist1;
|
|---|
| 461 | }
|
|---|
| 462 | else if (!wrongSide) dist = dist1;
|
|---|
| 463 | }
|
|---|
| 464 |
|
|---|
| 465 | return dist;
|
|---|
| 466 | }
|
|---|
| 467 |
|
|---|
| 468 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 469 | //
|
|---|
| 470 | // Extent
|
|---|
| 471 | //
|
|---|
| 472 | // Calculates the furthest the triangle extends in a particular direction
|
|---|
| 473 | // defined by the vector axis.
|
|---|
| 474 | //
|
|---|
| 475 | G4double G4TriangularFacet::Extent (const G4ThreeVector axis)
|
|---|
| 476 | {
|
|---|
| 477 | G4double s = P0.dot(axis);
|
|---|
| 478 | G4double sp = P[0].dot(axis);
|
|---|
| 479 | if (sp > s) s = sp;
|
|---|
| 480 | sp = P[1].dot(axis);
|
|---|
| 481 | if (sp > s) s = sp;
|
|---|
| 482 |
|
|---|
| 483 | return s;
|
|---|
| 484 | }
|
|---|
| 485 |
|
|---|
| 486 | ///////////////////////////////////////////////////////////////////////////////
|
|---|
| 487 | //
|
|---|
| 488 | // Intersect
|
|---|
| 489 | //
|
|---|
| 490 | // Member function to find the next intersection when going from p in the
|
|---|
| 491 | // direction of v. If:
|
|---|
| 492 | // (1) "outgoing" is TRUE, only consider the face if we are going out through
|
|---|
| 493 | // the face.
|
|---|
| 494 | // (2) "outgoing" is FALSE, only consider the face if we are going in through
|
|---|
| 495 | // the face.
|
|---|
| 496 | // Member functions returns TRUE if there is an intersection, FALSE otherwise.
|
|---|
| 497 | // Sets the distance (distance along w), distFromSurface (orthogonal distance)
|
|---|
| 498 | // and normal.
|
|---|
| 499 | //
|
|---|
| 500 | // Also considers intersections that happen with negative distance for small
|
|---|
| 501 | // distances of distFromSurface = 0.5*kCarTolerance in the wrong direction.
|
|---|
| 502 | // This is to detect kSurface without doing a full Inside(p) in
|
|---|
| 503 | // G4TessellatedSolid::Distance(p,v) calculation.
|
|---|
| 504 | //
|
|---|
| 505 | // This member function is thanks the valuable work of Rickard Holmberg. PT.
|
|---|
| 506 | // However, "gotos" are the Work of the Devil have been exorcised with
|
|---|
| 507 | // extreme prejudice!!
|
|---|
| 508 | //
|
|---|
| 509 | // IMPORTANT NOTE: These calculations are predicated on v being a unit
|
|---|
| 510 | // vector. If G4TessellatedSolid or other classes call this member function
|
|---|
| 511 | // with |v| != 1 then there will be errors.
|
|---|
| 512 | //
|
|---|
| 513 | G4bool G4TriangularFacet::Intersect (const G4ThreeVector &p,
|
|---|
| 514 | const G4ThreeVector &v, G4bool outgoing, G4double &distance,
|
|---|
| 515 | G4double &distFromSurface, G4ThreeVector &normal)
|
|---|
| 516 | {
|
|---|
| 517 | //
|
|---|
| 518 | //
|
|---|
| 519 | // Check whether the direction of the facet is consistent with the vector v
|
|---|
| 520 | // and the need to be outgoing or ingoing. If inconsistent, disregard and
|
|---|
| 521 | // return false.
|
|---|
| 522 | //
|
|---|
| 523 | G4double w = v.dot(surfaceNormal);
|
|---|
| 524 | if ((outgoing && (w <-dirTolerance)) || (!outgoing && (w > dirTolerance)))
|
|---|
| 525 | {
|
|---|
| 526 | distance = kInfinity;
|
|---|
| 527 | distFromSurface = kInfinity;
|
|---|
| 528 | normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 529 | return false;
|
|---|
| 530 | }
|
|---|
| 531 | //
|
|---|
| 532 | //
|
|---|
| 533 | // Calculate the orthogonal distance from p to the surface containing the
|
|---|
| 534 | // triangle. Then determine if we're on the right or wrong side of the
|
|---|
| 535 | // surface (at a distance greater than kCarTolerance) to be consistent with
|
|---|
| 536 | // "outgoing".
|
|---|
| 537 | //
|
|---|
| 538 | G4ThreeVector D = P0 - p;
|
|---|
| 539 | distFromSurface = D.dot(surfaceNormal);
|
|---|
| 540 | G4bool wrongSide = (outgoing && (distFromSurface < -0.5*kCarTolerance)) ||
|
|---|
| 541 | (!outgoing && (distFromSurface > 0.5*kCarTolerance));
|
|---|
| 542 | if (wrongSide)
|
|---|
| 543 | {
|
|---|
| 544 | distance = kInfinity;
|
|---|
| 545 | distFromSurface = kInfinity;
|
|---|
| 546 | normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 547 | return false;
|
|---|
| 548 | }
|
|---|
| 549 |
|
|---|
| 550 | wrongSide = (outgoing && (distFromSurface < 0.0)) ||
|
|---|
| 551 | (!outgoing && (distFromSurface > 0.0));
|
|---|
| 552 | if (wrongSide)
|
|---|
| 553 | {
|
|---|
| 554 | //
|
|---|
| 555 | //
|
|---|
| 556 | // We're slightly on the wrong side of the surface. Check if we're close
|
|---|
| 557 | // enough using a precise distance calculation.
|
|---|
| 558 | //
|
|---|
| 559 | G4ThreeVector u = Distance(p);
|
|---|
| 560 | if (std::sqrt(sqrDist) <= 0.5*kCarTolerance)
|
|---|
| 561 | {
|
|---|
| 562 | //
|
|---|
| 563 | //
|
|---|
| 564 | // We're very close. Therefore return a small negative number to pretend
|
|---|
| 565 | // we intersect.
|
|---|
| 566 | //
|
|---|
| 567 | // distance = -0.5*kCarTolerance;
|
|---|
| 568 | distance = 0.0;
|
|---|
| 569 | normal = surfaceNormal;
|
|---|
| 570 | return true;
|
|---|
| 571 | }
|
|---|
| 572 | else
|
|---|
| 573 | {
|
|---|
| 574 | //
|
|---|
| 575 | //
|
|---|
| 576 | // We're close to the surface containing the triangle, but sufficiently
|
|---|
| 577 | // far from the triangle, and on the wrong side compared to the directions
|
|---|
| 578 | // of the surface normal and v. There is no intersection.
|
|---|
| 579 | //
|
|---|
| 580 | distance = kInfinity;
|
|---|
| 581 | distFromSurface = kInfinity;
|
|---|
| 582 | normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 583 | return false;
|
|---|
| 584 | }
|
|---|
| 585 | }
|
|---|
| 586 | if (w < dirTolerance && w > -dirTolerance)
|
|---|
| 587 | {
|
|---|
| 588 | //
|
|---|
| 589 | //
|
|---|
| 590 | // The ray is within the plane of the triangle. Project the problem into 2D
|
|---|
| 591 | // in the plane of the triangle. First try to create orthogonal unit vectors
|
|---|
| 592 | // mu and nu, where mu is E[0]/|E[0]|. This is kinda like
|
|---|
| 593 | // the original algorithm due to Rickard Holmberg, but with better mathematical
|
|---|
| 594 | // justification than the original method ... however, beware Rickard's was less
|
|---|
| 595 | // time-consuming.
|
|---|
| 596 | //
|
|---|
| 597 | // Note that vprime is not a unit vector. We need to keep it unnormalised
|
|---|
| 598 | // since the values of distance along vprime (s0 and s1) for intersection with
|
|---|
| 599 | // the triangle will be used to determine if we cut the plane at the same
|
|---|
| 600 | // time.
|
|---|
| 601 | //
|
|---|
| 602 | G4ThreeVector mu = E[0].unit();
|
|---|
| 603 | G4ThreeVector nu = surfaceNormal.cross(mu);
|
|---|
| 604 | G4TwoVector pprime(p.dot(mu),p.dot(nu));
|
|---|
| 605 | G4TwoVector vprime(v.dot(mu),v.dot(nu));
|
|---|
| 606 | G4TwoVector P0prime(P0.dot(mu),P0.dot(nu));
|
|---|
| 607 | G4TwoVector E0prime(E[0].mag(),0.0);
|
|---|
| 608 | G4TwoVector E1prime(E[1].dot(mu),E[1].dot(nu));
|
|---|
| 609 |
|
|---|
| 610 | G4TwoVector loc[2];
|
|---|
| 611 | if ( tGeomAlg->IntersectLineAndTriangle2D(pprime,vprime,P0prime,
|
|---|
| 612 | E0prime,E1prime,loc) )
|
|---|
| 613 | {
|
|---|
| 614 | //
|
|---|
| 615 | //
|
|---|
| 616 | // There is an intersection between the line and triangle in 2D. Now check
|
|---|
| 617 | // which part of the line intersects with the plane containing the triangle
|
|---|
| 618 | // in 3D.
|
|---|
| 619 | //
|
|---|
| 620 | G4double vprimemag = vprime.mag();
|
|---|
| 621 | G4double s0 = (loc[0] - pprime).mag()/vprimemag;
|
|---|
| 622 | G4double s1 = (loc[1] - pprime).mag()/vprimemag;
|
|---|
| 623 | G4double normDist0 = surfaceNormal.dot(s0*v) - distFromSurface;
|
|---|
| 624 | G4double normDist1 = surfaceNormal.dot(s1*v) - distFromSurface;
|
|---|
| 625 |
|
|---|
| 626 | if ((normDist0 < 0.0 && normDist1 < 0.0) ||
|
|---|
| 627 | (normDist0 > 0.0 && normDist1 > 0.0))
|
|---|
| 628 | {
|
|---|
| 629 | distance = kInfinity;
|
|---|
| 630 | distFromSurface = kInfinity;
|
|---|
| 631 | normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 632 | return false;
|
|---|
| 633 | }
|
|---|
| 634 | else
|
|---|
| 635 | {
|
|---|
| 636 | G4double dnormDist = normDist1-normDist0;
|
|---|
| 637 | if (std::abs(dnormDist) < DBL_EPSILON)
|
|---|
| 638 | {
|
|---|
| 639 | distance = s0;
|
|---|
| 640 | normal = surfaceNormal;
|
|---|
| 641 | if (!outgoing) distFromSurface = -distFromSurface;
|
|---|
| 642 | return true;
|
|---|
| 643 | }
|
|---|
| 644 | else
|
|---|
| 645 | {
|
|---|
| 646 | distance = s0 - normDist0*(s1-s0)/dnormDist;
|
|---|
| 647 | normal = surfaceNormal;
|
|---|
| 648 | if (!outgoing) distFromSurface = -distFromSurface;
|
|---|
| 649 | return true;
|
|---|
| 650 | }
|
|---|
| 651 | }
|
|---|
| 652 |
|
|---|
| 653 | // G4ThreeVector dloc = loc1 - loc0;
|
|---|
| 654 | // G4ThreeVector dlocXv = dloc.cross(v);
|
|---|
| 655 | // G4double dlocXvmag = dlocXv.mag();
|
|---|
| 656 | // if (dloc.mag() <= 0.5*kCarTolerance || dlocXvmag <= DBL_EPSILON)
|
|---|
| 657 | // {
|
|---|
| 658 | // distance = loc0.mag();
|
|---|
| 659 | // normal = surfaceNormal;
|
|---|
| 660 | // if (!outgoing) distFromSurface = -distFromSurface;
|
|---|
| 661 | // return true;
|
|---|
| 662 | // }
|
|---|
| 663 |
|
|---|
| 664 | // G4ThreeVector loc0Xv = loc0.cross(v);
|
|---|
| 665 | // G4ThreeVector loc1Xv = loc1.cross(v);
|
|---|
| 666 | // G4double sameDir = -loc0Xv.dot(loc1Xv);
|
|---|
| 667 | // if (sameDir < 0.0)
|
|---|
| 668 | // {
|
|---|
| 669 | // distance = kInfinity;
|
|---|
| 670 | // distFromSurface = kInfinity;
|
|---|
| 671 | // normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 672 | // return false;
|
|---|
| 673 | // }
|
|---|
| 674 | // else
|
|---|
| 675 | // {
|
|---|
| 676 | // distance = loc0.mag() + loc0Xv.mag() * dloc.mag()/dlocXvmag;
|
|---|
| 677 | // normal = surfaceNormal;
|
|---|
| 678 | // if (!outgoing) distFromSurface = -distFromSurface;
|
|---|
| 679 | // return true;
|
|---|
| 680 | // }
|
|---|
| 681 | }
|
|---|
| 682 | else
|
|---|
| 683 | {
|
|---|
| 684 | distance = kInfinity;
|
|---|
| 685 | distFromSurface = kInfinity;
|
|---|
| 686 | normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 687 | return false;
|
|---|
| 688 | }
|
|---|
| 689 | }
|
|---|
| 690 | //
|
|---|
| 691 | //
|
|---|
| 692 | // Use conventional algorithm to determine the whether there is an
|
|---|
| 693 | // intersection. This involves determining the point of intersection of the
|
|---|
| 694 | // line with the plane containing the triangle, and then calculating if the
|
|---|
| 695 | // point is within the triangle.
|
|---|
| 696 | //
|
|---|
| 697 | distance = distFromSurface / w;
|
|---|
| 698 | G4ThreeVector pp = p + v*distance;
|
|---|
| 699 | G4ThreeVector DD = P0 - pp;
|
|---|
| 700 | G4double d = E[0].dot(DD);
|
|---|
| 701 | G4double e = E[1].dot(DD);
|
|---|
| 702 | G4double s = b*e - c*d;
|
|---|
| 703 | G4double t = b*d - a*e;
|
|---|
| 704 |
|
|---|
| 705 | if (s < 0.0 || t < 0.0 || s+t > det)
|
|---|
| 706 | {
|
|---|
| 707 | //
|
|---|
| 708 | //
|
|---|
| 709 | // The intersection is outside of the triangle.
|
|---|
| 710 | //
|
|---|
| 711 | distance = kInfinity;
|
|---|
| 712 | distFromSurface = kInfinity;
|
|---|
| 713 | normal = G4ThreeVector(0.0,0.0,0.0);
|
|---|
| 714 | return false;
|
|---|
| 715 | }
|
|---|
| 716 | else
|
|---|
| 717 | {
|
|---|
| 718 | //
|
|---|
| 719 | //
|
|---|
| 720 | // There is an intersection. Now we only need to set the surface normal.
|
|---|
| 721 | //
|
|---|
| 722 | normal = surfaceNormal;
|
|---|
| 723 | if (!outgoing) distFromSurface = -distFromSurface;
|
|---|
| 724 | return true;
|
|---|
| 725 | }
|
|---|
| 726 | }
|
|---|
| 727 |
|
|---|
| 728 | ////////////////////////////////////////////////////////////////////////
|
|---|
| 729 | //
|
|---|
| 730 | // GetPointOnFace
|
|---|
| 731 | //
|
|---|
| 732 | // Auxiliary method for get a random point on surface
|
|---|
| 733 |
|
|---|
| 734 | G4ThreeVector G4TriangularFacet::GetPointOnFace() const
|
|---|
| 735 | {
|
|---|
| 736 | G4double alpha = CLHEP::RandFlat::shoot(0.,1.);
|
|---|
| 737 | G4double beta = CLHEP::RandFlat::shoot(0.,1);
|
|---|
| 738 | G4double lambda1=alpha*beta;
|
|---|
| 739 | G4double lambda0=alpha-lambda1;
|
|---|
| 740 |
|
|---|
| 741 | return (P0 + lambda0*E[0] + lambda1*E[1]);
|
|---|
| 742 | }
|
|---|
| 743 |
|
|---|
| 744 | ////////////////////////////////////////////////////////////////////////
|
|---|
| 745 | //
|
|---|
| 746 | // GetArea
|
|---|
| 747 | //
|
|---|
| 748 | // Auxiliary method for returning the surface area
|
|---|
| 749 |
|
|---|
| 750 | G4double G4TriangularFacet::GetArea()
|
|---|
| 751 | {
|
|---|
| 752 | return area;
|
|---|
| 753 | }
|
|---|
| 754 | ////////////////////////////////////////////////////////////////////////
|
|---|
| 755 | //
|
|---|