| 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. *
|
|---|
| 20 | // * By using, copying, modifying or distributing the software (or *
|
|---|
| 21 | // * any work based on the software) you agree to acknowledge its *
|
|---|
| 22 | // * use in resulting scientific publications, and indicate your *
|
|---|
| 23 | // * acceptance of all terms of the Geant4 Software license. *
|
|---|
| 24 | // ********************************************************************
|
|---|
| 25 | //
|
|---|
| 26 | //
|
|---|
| 27 | // $Id: G4PolyPhiFace.cc,v 1.13 2007/07/19 12:57:14 gcosmo Exp $
|
|---|
| 28 | // GEANT4 tag $Name: $
|
|---|
| 29 | //
|
|---|
| 30 | //
|
|---|
| 31 | // --------------------------------------------------------------------
|
|---|
| 32 | // GEANT 4 class source file
|
|---|
| 33 | //
|
|---|
| 34 | //
|
|---|
| 35 | // G4PolyPhiFace.cc
|
|---|
| 36 | //
|
|---|
| 37 | // Implementation of the face that bounds a polycone or polyhedra at
|
|---|
| 38 | // its phi opening.
|
|---|
| 39 | //
|
|---|
| 40 | // --------------------------------------------------------------------
|
|---|
| 41 |
|
|---|
| 42 | #include "G4PolyPhiFace.hh"
|
|---|
| 43 | #include "G4ClippablePolygon.hh"
|
|---|
| 44 | #include "G4ReduciblePolygon.hh"
|
|---|
| 45 | #include "G4AffineTransform.hh"
|
|---|
| 46 | #include "G4SolidExtentList.hh"
|
|---|
| 47 | #include "G4GeometryTolerance.hh"
|
|---|
| 48 |
|
|---|
| 49 | //
|
|---|
| 50 | // Constructor
|
|---|
| 51 | //
|
|---|
| 52 | // Points r,z should be supplied in clockwise order in r,z. For example:
|
|---|
| 53 | //
|
|---|
| 54 | // [1]---------[2] ^ R
|
|---|
| 55 | // | | |
|
|---|
| 56 | // | | +--> z
|
|---|
| 57 | // [0]---------[3]
|
|---|
| 58 | //
|
|---|
| 59 | G4PolyPhiFace::G4PolyPhiFace( const G4ReduciblePolygon *rz,
|
|---|
| 60 | G4double phi,
|
|---|
| 61 | G4double deltaPhi,
|
|---|
| 62 | G4double phiOther )
|
|---|
| 63 | {
|
|---|
| 64 | kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
|
|---|
| 65 |
|
|---|
| 66 | numEdges = rz->NumVertices();
|
|---|
| 67 |
|
|---|
| 68 | rMin = rz->Amin();
|
|---|
| 69 | rMax = rz->Amax();
|
|---|
| 70 | zMin = rz->Bmin();
|
|---|
| 71 | zMax = rz->Bmax();
|
|---|
| 72 |
|
|---|
| 73 | //
|
|---|
| 74 | // Is this the "starting" phi edge of the two?
|
|---|
| 75 | //
|
|---|
| 76 | G4bool start = (phiOther > phi);
|
|---|
| 77 |
|
|---|
| 78 | //
|
|---|
| 79 | // Build radial vector
|
|---|
| 80 | //
|
|---|
| 81 | radial = G4ThreeVector( std::cos(phi), std::sin(phi), 0.0 );
|
|---|
| 82 |
|
|---|
| 83 | //
|
|---|
| 84 | // Build normal
|
|---|
| 85 | //
|
|---|
| 86 | G4double zSign = start ? 1 : -1;
|
|---|
| 87 | normal = G4ThreeVector( zSign*radial.y(), -zSign*radial.x(), 0 );
|
|---|
| 88 |
|
|---|
| 89 | //
|
|---|
| 90 | // Is allBehind?
|
|---|
| 91 | //
|
|---|
| 92 | allBehind = (zSign*(std::cos(phiOther)*radial.y() - std::sin(phiOther)*radial.x()) < 0);
|
|---|
| 93 |
|
|---|
| 94 | //
|
|---|
| 95 | // Adjacent edges
|
|---|
| 96 | //
|
|---|
| 97 | G4double midPhi = phi + (start ? +0.5 : -0.5)*deltaPhi;
|
|---|
| 98 | G4double cosMid = std::cos(midPhi),
|
|---|
| 99 | sinMid = std::sin(midPhi);
|
|---|
| 100 |
|
|---|
| 101 | //
|
|---|
| 102 | // Allocate corners
|
|---|
| 103 | //
|
|---|
| 104 | corners = new G4PolyPhiFaceVertex[numEdges];
|
|---|
| 105 |
|
|---|
| 106 | //
|
|---|
| 107 | // Fill them
|
|---|
| 108 | //
|
|---|
| 109 | G4ReduciblePolygonIterator iterRZ(rz);
|
|---|
| 110 |
|
|---|
| 111 | G4PolyPhiFaceVertex *corn = corners;
|
|---|
| 112 | iterRZ.Begin();
|
|---|
| 113 | do
|
|---|
| 114 | {
|
|---|
| 115 | corn->r = iterRZ.GetA();
|
|---|
| 116 | corn->z = iterRZ.GetB();
|
|---|
| 117 | corn->x = corn->r*radial.x();
|
|---|
| 118 | corn->y = corn->r*radial.y();
|
|---|
| 119 | } while( ++corn, iterRZ.Next() );
|
|---|
| 120 |
|
|---|
| 121 | //
|
|---|
| 122 | // Allocate edges
|
|---|
| 123 | //
|
|---|
| 124 | edges = new G4PolyPhiFaceEdge[numEdges];
|
|---|
| 125 |
|
|---|
| 126 | //
|
|---|
| 127 | // Fill them
|
|---|
| 128 | //
|
|---|
| 129 | G4double rFact = std::cos(0.5*deltaPhi);
|
|---|
| 130 | G4double rFactNormalize = 1.0/std::sqrt(1.0+rFact*rFact);
|
|---|
| 131 |
|
|---|
| 132 | G4PolyPhiFaceVertex *prev = corners+numEdges-1,
|
|---|
| 133 | *here = corners;
|
|---|
| 134 | G4PolyPhiFaceEdge *edge = edges;
|
|---|
| 135 | do
|
|---|
| 136 | {
|
|---|
| 137 | G4ThreeVector sideNorm;
|
|---|
| 138 |
|
|---|
| 139 | edge->v0 = prev;
|
|---|
| 140 | edge->v1 = here;
|
|---|
| 141 |
|
|---|
| 142 | G4double dr = here->r - prev->r,
|
|---|
| 143 | dz = here->z - prev->z;
|
|---|
| 144 |
|
|---|
| 145 | edge->length = std::sqrt( dr*dr + dz*dz );
|
|---|
| 146 |
|
|---|
| 147 | edge->tr = dr/edge->length;
|
|---|
| 148 | edge->tz = dz/edge->length;
|
|---|
| 149 |
|
|---|
| 150 | if ((here->r < DBL_MIN) && (prev->r < DBL_MIN))
|
|---|
| 151 | {
|
|---|
| 152 | //
|
|---|
| 153 | // Sigh! Always exceptions!
|
|---|
| 154 | // This edge runs at r==0, so its adjoing surface is not a
|
|---|
| 155 | // PolyconeSide or PolyhedraSide, but the opposite PolyPhiFace.
|
|---|
| 156 | //
|
|---|
| 157 | G4double zSignOther = start ? -1 : 1;
|
|---|
| 158 | sideNorm = G4ThreeVector( zSignOther*std::sin(phiOther),
|
|---|
| 159 | -zSignOther*std::cos(phiOther), 0 );
|
|---|
| 160 | }
|
|---|
| 161 | else
|
|---|
| 162 | {
|
|---|
| 163 | sideNorm = G4ThreeVector( edge->tz*cosMid,
|
|---|
| 164 | edge->tz*sinMid,
|
|---|
| 165 | -edge->tr*rFact );
|
|---|
| 166 | sideNorm *= rFactNormalize;
|
|---|
| 167 | }
|
|---|
| 168 | sideNorm += normal;
|
|---|
| 169 |
|
|---|
| 170 | edge->norm3D = sideNorm.unit();
|
|---|
| 171 | } while( edge++, prev=here, ++here < corners+numEdges );
|
|---|
| 172 |
|
|---|
| 173 | //
|
|---|
| 174 | // Go back and fill in corner "normals"
|
|---|
| 175 | //
|
|---|
| 176 | G4PolyPhiFaceEdge *prevEdge = edges+numEdges-1;
|
|---|
| 177 | edge = edges;
|
|---|
| 178 | do
|
|---|
| 179 | {
|
|---|
| 180 | //
|
|---|
| 181 | // Calculate vertex 2D normals (on the phi surface)
|
|---|
| 182 | //
|
|---|
| 183 | G4double rPart = prevEdge->tr + edge->tr;
|
|---|
| 184 | G4double zPart = prevEdge->tz + edge->tz;
|
|---|
| 185 | G4double norm = std::sqrt( rPart*rPart + zPart*zPart );
|
|---|
| 186 | G4double rNorm = +zPart/norm;
|
|---|
| 187 | G4double zNorm = -rPart/norm;
|
|---|
| 188 |
|
|---|
| 189 | edge->v0->rNorm = rNorm;
|
|---|
| 190 | edge->v0->zNorm = zNorm;
|
|---|
| 191 |
|
|---|
| 192 | //
|
|---|
| 193 | // Calculate the 3D normals.
|
|---|
| 194 | //
|
|---|
| 195 | // Find the vector perpendicular to the z axis
|
|---|
| 196 | // that defines the plane that contains the vertex normal
|
|---|
| 197 | //
|
|---|
| 198 | G4ThreeVector xyVector;
|
|---|
| 199 |
|
|---|
| 200 | if (edge->v0->r < DBL_MIN)
|
|---|
| 201 | {
|
|---|
| 202 | //
|
|---|
| 203 | // This is a vertex at r==0, which is a special
|
|---|
| 204 | // case. The normal we will construct lays in the
|
|---|
| 205 | // plane at the center of the phi opening.
|
|---|
| 206 | //
|
|---|
| 207 | // We also know that rNorm < 0
|
|---|
| 208 | //
|
|---|
| 209 | G4double zSignOther = start ? -1 : 1;
|
|---|
| 210 | G4ThreeVector normalOther( zSignOther*std::sin(phiOther),
|
|---|
| 211 | -zSignOther*std::cos(phiOther), 0 );
|
|---|
| 212 |
|
|---|
| 213 | xyVector = - normal - normalOther;
|
|---|
| 214 | }
|
|---|
| 215 | else
|
|---|
| 216 | {
|
|---|
| 217 | //
|
|---|
| 218 | // This is a vertex at r > 0. The plane
|
|---|
| 219 | // is the average of the normal and the
|
|---|
| 220 | // normal of the adjacent phi face
|
|---|
| 221 | //
|
|---|
| 222 | xyVector = G4ThreeVector( cosMid, sinMid, 0 );
|
|---|
| 223 | if (rNorm < 0)
|
|---|
| 224 | xyVector -= normal;
|
|---|
| 225 | else
|
|---|
| 226 | xyVector += normal;
|
|---|
| 227 | }
|
|---|
| 228 |
|
|---|
| 229 | //
|
|---|
| 230 | // Combine it with the r/z direction from the face
|
|---|
| 231 | //
|
|---|
| 232 | edge->v0->norm3D = rNorm*xyVector.unit() + G4ThreeVector( 0, 0, zNorm );
|
|---|
| 233 | } while( prevEdge=edge, ++edge < edges+numEdges );
|
|---|
| 234 |
|
|---|
| 235 | //
|
|---|
| 236 | // Build point on surface
|
|---|
| 237 | //
|
|---|
| 238 | G4double rAve = 0.5*(rMax-rMin),
|
|---|
| 239 | zAve = 0.5*(zMax-zMin);
|
|---|
| 240 | surface = G4ThreeVector( rAve*radial.x(), rAve*radial.y(), zAve );
|
|---|
| 241 | }
|
|---|
| 242 |
|
|---|
| 243 |
|
|---|
| 244 | //
|
|---|
| 245 | // Diagnose
|
|---|
| 246 | //
|
|---|
| 247 | // Throw an exception if something is found inconsistent with
|
|---|
| 248 | // the solid.
|
|---|
| 249 | //
|
|---|
| 250 | // For debugging purposes only
|
|---|
| 251 | //
|
|---|
| 252 | void G4PolyPhiFace::Diagnose( G4VSolid *owner )
|
|---|
| 253 | {
|
|---|
| 254 | G4PolyPhiFaceVertex *corner = corners;
|
|---|
| 255 | do
|
|---|
| 256 | {
|
|---|
| 257 | G4ThreeVector test(corner->x, corner->y, corner->z);
|
|---|
| 258 | test -= 1E-6*corner->norm3D;
|
|---|
| 259 |
|
|---|
| 260 | if (owner->Inside(test) != kInside)
|
|---|
| 261 | G4Exception( "G4PolyPhiFace::Diagnose()", "InvalidSetup",
|
|---|
| 262 | FatalException, "Bad vertex normal found." );
|
|---|
| 263 | } while( ++corner < corners+numEdges );
|
|---|
| 264 | }
|
|---|
| 265 |
|
|---|
| 266 |
|
|---|
| 267 | //
|
|---|
| 268 | // Fake default constructor - sets only member data and allocates memory
|
|---|
| 269 | // for usage restricted to object persistency.
|
|---|
| 270 | //
|
|---|
| 271 | G4PolyPhiFace::G4PolyPhiFace( __void__&)
|
|---|
| 272 | : edges(0), corners(0)
|
|---|
| 273 | {
|
|---|
| 274 | }
|
|---|
| 275 |
|
|---|
| 276 |
|
|---|
| 277 | //
|
|---|
| 278 | // Destructor
|
|---|
| 279 | //
|
|---|
| 280 | G4PolyPhiFace::~G4PolyPhiFace()
|
|---|
| 281 | {
|
|---|
| 282 | delete [] edges;
|
|---|
| 283 | delete [] corners;
|
|---|
| 284 | }
|
|---|
| 285 |
|
|---|
| 286 |
|
|---|
| 287 | //
|
|---|
| 288 | // Copy constructor
|
|---|
| 289 | //
|
|---|
| 290 | G4PolyPhiFace::G4PolyPhiFace( const G4PolyPhiFace &source )
|
|---|
| 291 | : G4VCSGface()
|
|---|
| 292 | {
|
|---|
| 293 | CopyStuff( source );
|
|---|
| 294 | }
|
|---|
| 295 |
|
|---|
| 296 |
|
|---|
| 297 | //
|
|---|
| 298 | // Assignment operator
|
|---|
| 299 | //
|
|---|
| 300 | G4PolyPhiFace& G4PolyPhiFace::operator=( const G4PolyPhiFace &source )
|
|---|
| 301 | {
|
|---|
| 302 | if (this == &source) return *this;
|
|---|
| 303 |
|
|---|
| 304 | delete [] edges;
|
|---|
| 305 | delete [] corners;
|
|---|
| 306 |
|
|---|
| 307 | CopyStuff( source );
|
|---|
| 308 |
|
|---|
| 309 | return *this;
|
|---|
| 310 | }
|
|---|
| 311 |
|
|---|
| 312 |
|
|---|
| 313 | //
|
|---|
| 314 | // CopyStuff (protected)
|
|---|
| 315 | //
|
|---|
| 316 | void G4PolyPhiFace::CopyStuff( const G4PolyPhiFace &source )
|
|---|
| 317 | {
|
|---|
| 318 | //
|
|---|
| 319 | // The simple stuff
|
|---|
| 320 | //
|
|---|
| 321 | numEdges = source.numEdges;
|
|---|
| 322 | normal = source.normal;
|
|---|
| 323 | radial = source.radial;
|
|---|
| 324 | surface = source.surface;
|
|---|
| 325 | rMin = source.rMin;
|
|---|
| 326 | rMax = source.rMax;
|
|---|
| 327 | zMin = source.zMin;
|
|---|
| 328 | zMax = source.zMax;
|
|---|
| 329 | allBehind = source.allBehind;
|
|---|
| 330 |
|
|---|
| 331 | kCarTolerance = source.kCarTolerance;
|
|---|
| 332 |
|
|---|
| 333 | //
|
|---|
| 334 | // Corner dynamic array
|
|---|
| 335 | //
|
|---|
| 336 | corners = new G4PolyPhiFaceVertex[numEdges];
|
|---|
| 337 | G4PolyPhiFaceVertex *corn = corners,
|
|---|
| 338 | *sourceCorn = source.corners;
|
|---|
| 339 | do
|
|---|
| 340 | {
|
|---|
| 341 | *corn = *sourceCorn;
|
|---|
| 342 | } while( ++sourceCorn, ++corn < corners+numEdges );
|
|---|
| 343 |
|
|---|
| 344 | //
|
|---|
| 345 | // Edge dynamic array
|
|---|
| 346 | //
|
|---|
| 347 | edges = new G4PolyPhiFaceEdge[numEdges];
|
|---|
| 348 |
|
|---|
| 349 | G4PolyPhiFaceVertex *prev = corners+numEdges-1,
|
|---|
| 350 | *here = corners;
|
|---|
| 351 | G4PolyPhiFaceEdge *edge = edges,
|
|---|
| 352 | *sourceEdge = source.edges;
|
|---|
| 353 | do
|
|---|
| 354 | {
|
|---|
| 355 | *edge = *sourceEdge;
|
|---|
| 356 | edge->v0 = prev;
|
|---|
| 357 | edge->v1 = here;
|
|---|
| 358 | } while( ++sourceEdge, ++edge, prev=here, ++here < corners+numEdges );
|
|---|
| 359 | }
|
|---|
| 360 |
|
|---|
| 361 |
|
|---|
| 362 | //
|
|---|
| 363 | // Intersect
|
|---|
| 364 | //
|
|---|
| 365 | G4bool G4PolyPhiFace::Intersect( const G4ThreeVector &p,
|
|---|
| 366 | const G4ThreeVector &v,
|
|---|
| 367 | G4bool outgoing,
|
|---|
| 368 | G4double surfTolerance,
|
|---|
| 369 | G4double &distance,
|
|---|
| 370 | G4double &distFromSurface,
|
|---|
| 371 | G4ThreeVector &aNormal,
|
|---|
| 372 | G4bool &isAllBehind )
|
|---|
| 373 | {
|
|---|
| 374 | G4double normSign = outgoing ? +1 : -1;
|
|---|
| 375 |
|
|---|
| 376 | //
|
|---|
| 377 | // These don't change
|
|---|
| 378 | //
|
|---|
| 379 | isAllBehind = allBehind;
|
|---|
| 380 | aNormal = normal;
|
|---|
| 381 |
|
|---|
| 382 | //
|
|---|
| 383 | // Correct normal? Here we have straight sides, and can safely ignore
|
|---|
| 384 | // intersections where the dot product with the normal is zero.
|
|---|
| 385 | //
|
|---|
| 386 | G4double dotProd = normSign*normal.dot(v);
|
|---|
| 387 |
|
|---|
| 388 | if (dotProd <= 0) return false;
|
|---|
| 389 |
|
|---|
| 390 | //
|
|---|
| 391 | // Calculate distance to surface. If the side is too far
|
|---|
| 392 | // behind the point, we must reject it.
|
|---|
| 393 | //
|
|---|
| 394 | G4ThreeVector ps = p - surface;
|
|---|
| 395 | distFromSurface = -normSign*ps.dot(normal);
|
|---|
| 396 |
|
|---|
| 397 | if (distFromSurface < -surfTolerance) return false;
|
|---|
| 398 |
|
|---|
| 399 | //
|
|---|
| 400 | // Calculate precise distance to intersection with the side
|
|---|
| 401 | // (along the trajectory, not normal to the surface)
|
|---|
| 402 | //
|
|---|
| 403 | distance = distFromSurface/dotProd;
|
|---|
| 404 |
|
|---|
| 405 | //
|
|---|
| 406 | // Calculate intersection point in r,z
|
|---|
| 407 | //
|
|---|
| 408 | G4ThreeVector ip = p + distance*v;
|
|---|
| 409 |
|
|---|
| 410 | G4double r = radial.dot(ip);
|
|---|
| 411 |
|
|---|
| 412 | //
|
|---|
| 413 | // And is it inside the r/z extent?
|
|---|
| 414 | //
|
|---|
| 415 | return InsideEdgesExact( r, ip.z(), normSign, p, v );
|
|---|
| 416 | }
|
|---|
| 417 |
|
|---|
| 418 |
|
|---|
| 419 | //
|
|---|
| 420 | // Distance
|
|---|
| 421 | //
|
|---|
| 422 | G4double G4PolyPhiFace::Distance( const G4ThreeVector &p, G4bool outgoing )
|
|---|
| 423 | {
|
|---|
| 424 | G4double normSign = outgoing ? +1 : -1;
|
|---|
| 425 | //
|
|---|
| 426 | // Correct normal?
|
|---|
| 427 | //
|
|---|
| 428 | G4ThreeVector ps = p - surface;
|
|---|
| 429 | G4double distPhi = -normSign*normal.dot(ps);
|
|---|
| 430 |
|
|---|
| 431 | if (distPhi < -0.5*kCarTolerance)
|
|---|
| 432 | return kInfinity;
|
|---|
| 433 | else if (distPhi < 0)
|
|---|
| 434 | distPhi = 0.0;
|
|---|
| 435 |
|
|---|
| 436 | //
|
|---|
| 437 | // Calculate projected point in r,z
|
|---|
| 438 | //
|
|---|
| 439 | G4double r = radial.dot(p);
|
|---|
| 440 |
|
|---|
| 441 | //
|
|---|
| 442 | // Are we inside the face?
|
|---|
| 443 | //
|
|---|
| 444 | G4double distRZ2;
|
|---|
| 445 |
|
|---|
| 446 | if (InsideEdges( r, p.z(), &distRZ2, 0 ))
|
|---|
| 447 | {
|
|---|
| 448 | //
|
|---|
| 449 | // Yup, answer is just distPhi
|
|---|
| 450 | //
|
|---|
| 451 | return distPhi;
|
|---|
| 452 | }
|
|---|
| 453 | else
|
|---|
| 454 | {
|
|---|
| 455 | //
|
|---|
| 456 | // Nope. Penalize by distance out
|
|---|
| 457 | //
|
|---|
| 458 | return std::sqrt( distPhi*distPhi + distRZ2 );
|
|---|
| 459 | }
|
|---|
| 460 | }
|
|---|
| 461 |
|
|---|
| 462 |
|
|---|
| 463 | //
|
|---|
| 464 | // Inside
|
|---|
| 465 | //
|
|---|
| 466 | EInside G4PolyPhiFace::Inside( const G4ThreeVector &p,
|
|---|
| 467 | G4double tolerance,
|
|---|
| 468 | G4double *bestDistance )
|
|---|
| 469 | {
|
|---|
| 470 | //
|
|---|
| 471 | // Get distance along phi, which if negative means the point
|
|---|
| 472 | // is nominally inside the shape.
|
|---|
| 473 | //
|
|---|
| 474 | G4ThreeVector ps = p - surface;
|
|---|
| 475 | G4double distPhi = normal.dot(ps);
|
|---|
| 476 |
|
|---|
| 477 | //
|
|---|
| 478 | // Calculate projected point in r,z
|
|---|
| 479 | //
|
|---|
| 480 | G4double r = radial.dot(p);
|
|---|
| 481 |
|
|---|
| 482 | //
|
|---|
| 483 | // Are we inside the face?
|
|---|
| 484 | //
|
|---|
| 485 | G4double distRZ2;
|
|---|
| 486 | G4PolyPhiFaceVertex *base3Dnorm;
|
|---|
| 487 | G4ThreeVector *head3Dnorm;
|
|---|
| 488 |
|
|---|
| 489 | if (InsideEdges( r, p.z(), &distRZ2, &base3Dnorm, &head3Dnorm ))
|
|---|
| 490 | {
|
|---|
| 491 | //
|
|---|
| 492 | // Looks like we're inside. Distance is distance in phi.
|
|---|
| 493 | //
|
|---|
| 494 | *bestDistance = std::fabs(distPhi);
|
|---|
| 495 |
|
|---|
| 496 | //
|
|---|
| 497 | // Use distPhi to decide fate
|
|---|
| 498 | //
|
|---|
| 499 | if (distPhi < -tolerance) return kInside;
|
|---|
| 500 | if (distPhi < tolerance) return kSurface;
|
|---|
| 501 | return kOutside;
|
|---|
| 502 | }
|
|---|
| 503 | else
|
|---|
| 504 | {
|
|---|
| 505 | //
|
|---|
| 506 | // We're outside the extent of the face,
|
|---|
| 507 | // so the distance is penalized by distance from edges in RZ
|
|---|
| 508 | //
|
|---|
| 509 | *bestDistance = std::sqrt( distPhi*distPhi + distRZ2 );
|
|---|
| 510 |
|
|---|
| 511 | //
|
|---|
| 512 | // Use edge normal to decide fate
|
|---|
| 513 | //
|
|---|
| 514 | G4ThreeVector cc( base3Dnorm->r*radial.x(),
|
|---|
| 515 | base3Dnorm->r*radial.y(),
|
|---|
| 516 | base3Dnorm->z );
|
|---|
| 517 | cc = p - cc;
|
|---|
| 518 | G4double normDist = head3Dnorm->dot(cc);
|
|---|
| 519 | if ( distRZ2 > tolerance*tolerance )
|
|---|
| 520 | {
|
|---|
| 521 | //
|
|---|
| 522 | // We're far enough away that kSurface is not possible
|
|---|
| 523 | //
|
|---|
| 524 | return normDist < 0 ? kInside : kOutside;
|
|---|
| 525 | }
|
|---|
| 526 |
|
|---|
| 527 | if (normDist < -tolerance) return kInside;
|
|---|
| 528 | if (normDist < tolerance) return kSurface;
|
|---|
| 529 | return kOutside;
|
|---|
| 530 | }
|
|---|
| 531 | }
|
|---|
| 532 |
|
|---|
| 533 |
|
|---|
| 534 | //
|
|---|
| 535 | // Normal
|
|---|
| 536 | //
|
|---|
| 537 | // This virtual member is simple for our planer shape,
|
|---|
| 538 | // which has only one normal
|
|---|
| 539 | //
|
|---|
| 540 | G4ThreeVector G4PolyPhiFace::Normal( const G4ThreeVector &p,
|
|---|
| 541 | G4double *bestDistance )
|
|---|
| 542 | {
|
|---|
| 543 | //
|
|---|
| 544 | // Get distance along phi, which if negative means the point
|
|---|
| 545 | // is nominally inside the shape.
|
|---|
| 546 | //
|
|---|
| 547 | G4double distPhi = normal.dot(p);
|
|---|
| 548 |
|
|---|
| 549 | //
|
|---|
| 550 | // Calculate projected point in r,z
|
|---|
| 551 | //
|
|---|
| 552 | G4double r = radial.dot(p);
|
|---|
| 553 |
|
|---|
| 554 | //
|
|---|
| 555 | // Are we inside the face?
|
|---|
| 556 | //
|
|---|
| 557 | G4double distRZ2;
|
|---|
| 558 |
|
|---|
| 559 | if (InsideEdges( r, p.z(), &distRZ2, 0 ))
|
|---|
| 560 | {
|
|---|
| 561 | //
|
|---|
| 562 | // Yup, answer is just distPhi
|
|---|
| 563 | //
|
|---|
| 564 | *bestDistance = std::fabs(distPhi);
|
|---|
| 565 | }
|
|---|
| 566 | else
|
|---|
| 567 | {
|
|---|
| 568 | //
|
|---|
| 569 | // Nope. Penalize by distance out
|
|---|
| 570 | //
|
|---|
| 571 | *bestDistance = std::sqrt( distPhi*distPhi + distRZ2 );
|
|---|
| 572 | }
|
|---|
| 573 |
|
|---|
| 574 | return normal;
|
|---|
| 575 | }
|
|---|
| 576 |
|
|---|
| 577 |
|
|---|
| 578 | //
|
|---|
| 579 | // Extent
|
|---|
| 580 | //
|
|---|
| 581 | // This actually isn't needed by polycone or polyhedra...
|
|---|
| 582 | //
|
|---|
| 583 | G4double G4PolyPhiFace::Extent( const G4ThreeVector axis )
|
|---|
| 584 | {
|
|---|
| 585 | G4double max = -kInfinity;
|
|---|
| 586 |
|
|---|
| 587 | G4PolyPhiFaceVertex *corner = corners;
|
|---|
| 588 | do
|
|---|
| 589 | {
|
|---|
| 590 | G4double here = axis.x()*corner->r*radial.x()
|
|---|
| 591 | + axis.y()*corner->r*radial.y()
|
|---|
| 592 | + axis.z()*corner->z;
|
|---|
| 593 | if (here > max) max = here;
|
|---|
| 594 | } while( ++corner < corners + numEdges );
|
|---|
| 595 |
|
|---|
| 596 | return max;
|
|---|
| 597 | }
|
|---|
| 598 |
|
|---|
| 599 |
|
|---|
| 600 | //
|
|---|
| 601 | // CalculateExtent
|
|---|
| 602 | //
|
|---|
| 603 | // See notes in G4VCSGface
|
|---|
| 604 | //
|
|---|
| 605 | void G4PolyPhiFace::CalculateExtent( const EAxis axis,
|
|---|
| 606 | const G4VoxelLimits &voxelLimit,
|
|---|
| 607 | const G4AffineTransform &transform,
|
|---|
| 608 | G4SolidExtentList &extentList )
|
|---|
| 609 | {
|
|---|
| 610 | //
|
|---|
| 611 | // Construct a (sometimes big) clippable polygon,
|
|---|
| 612 | //
|
|---|
| 613 | // Perform the necessary transformations while doing so
|
|---|
| 614 | //
|
|---|
| 615 | G4ClippablePolygon polygon;
|
|---|
| 616 |
|
|---|
| 617 | G4PolyPhiFaceVertex *corner = corners;
|
|---|
| 618 | do
|
|---|
| 619 | {
|
|---|
| 620 | G4ThreeVector point( 0, 0, corner->z );
|
|---|
| 621 | point += radial*corner->r;
|
|---|
| 622 |
|
|---|
| 623 | polygon.AddVertexInOrder( transform.TransformPoint( point ) );
|
|---|
| 624 | } while( ++corner < corners + numEdges );
|
|---|
| 625 |
|
|---|
| 626 | //
|
|---|
| 627 | // Clip away
|
|---|
| 628 | //
|
|---|
| 629 | if (polygon.PartialClip( voxelLimit, axis ))
|
|---|
| 630 | {
|
|---|
| 631 | //
|
|---|
| 632 | // Add it to the list
|
|---|
| 633 | //
|
|---|
| 634 | polygon.SetNormal( transform.TransformAxis(normal) );
|
|---|
| 635 | extentList.AddSurface( polygon );
|
|---|
| 636 | }
|
|---|
| 637 | }
|
|---|
| 638 |
|
|---|
| 639 |
|
|---|
| 640 | //
|
|---|
| 641 | //-------------------------------------------------------
|
|---|
| 642 |
|
|---|
| 643 |
|
|---|
| 644 | //
|
|---|
| 645 | // InsideEdgesExact
|
|---|
| 646 | //
|
|---|
| 647 | // Decide if the point in r,z is inside the edges of our face,
|
|---|
| 648 | // **but** do so consistently with other faces.
|
|---|
| 649 | //
|
|---|
| 650 | // This routine has functionality similar to InsideEdges, but uses
|
|---|
| 651 | // an algorithm to decide if a trajectory falls inside or outside the
|
|---|
| 652 | // face that uses only the trajectory p,v values and the three dimensional
|
|---|
| 653 | // points representing the edges of the polygon. The objective is to plug up
|
|---|
| 654 | // any leaks between touching G4PolyPhiFaces (at r==0) and any other face
|
|---|
| 655 | // that uses the same convention.
|
|---|
| 656 | //
|
|---|
| 657 | // See: "Computational Geometry in C (Second Edition)"
|
|---|
| 658 | // http://cs.smith.edu/~orourke/
|
|---|
| 659 | //
|
|---|
| 660 | G4bool G4PolyPhiFace::InsideEdgesExact( G4double r, G4double z,
|
|---|
| 661 | G4double normSign,
|
|---|
| 662 | const G4ThreeVector &p,
|
|---|
| 663 | const G4ThreeVector &v )
|
|---|
| 664 | {
|
|---|
| 665 | //
|
|---|
| 666 | // Quick check of extent
|
|---|
| 667 | //
|
|---|
| 668 | if ( (r < rMin-kCarTolerance)
|
|---|
| 669 | || (r > rMax+kCarTolerance) ) return false;
|
|---|
| 670 |
|
|---|
| 671 | if ( (z < zMin-kCarTolerance)
|
|---|
| 672 | || (z > zMax+kCarTolerance) ) return false;
|
|---|
| 673 |
|
|---|
| 674 | //
|
|---|
| 675 | // Exact check: loop over all vertices
|
|---|
| 676 | //
|
|---|
| 677 | G4double qx = p.x() + v.x(),
|
|---|
| 678 | qy = p.y() + v.y(),
|
|---|
| 679 | qz = p.z() + v.z();
|
|---|
| 680 |
|
|---|
| 681 | G4int answer = 0;
|
|---|
| 682 | G4PolyPhiFaceVertex *corn = corners,
|
|---|
| 683 | *prev = corners+numEdges-1;
|
|---|
| 684 |
|
|---|
| 685 | G4double cornZ, prevZ;
|
|---|
| 686 |
|
|---|
| 687 | prevZ = ExactZOrder( z, qx, qy, qz, v, normSign, prev );
|
|---|
| 688 | do
|
|---|
| 689 | {
|
|---|
| 690 | //
|
|---|
| 691 | // Get z order of this vertex, and compare to previous vertex
|
|---|
| 692 | //
|
|---|
| 693 | cornZ = ExactZOrder( z, qx, qy, qz, v, normSign, corn );
|
|---|
| 694 |
|
|---|
| 695 | if (cornZ < 0)
|
|---|
| 696 | {
|
|---|
| 697 | if (prevZ < 0) continue;
|
|---|
| 698 | }
|
|---|
| 699 | else if (cornZ > 0)
|
|---|
| 700 | {
|
|---|
| 701 | if (prevZ > 0) continue;
|
|---|
| 702 | }
|
|---|
| 703 | else
|
|---|
| 704 | {
|
|---|
| 705 | //
|
|---|
| 706 | // By chance, we overlap exactly (within precision) with
|
|---|
| 707 | // the current vertex. Continue if the same happened previously
|
|---|
| 708 | // (e.g. the previous vertex had the same z value)
|
|---|
| 709 | //
|
|---|
| 710 | if (prevZ == 0) continue;
|
|---|
| 711 |
|
|---|
| 712 | //
|
|---|
| 713 | // Otherwise, to decide what to do, we need to know what is
|
|---|
| 714 | // coming up next. Specifically, we need to find the next vertex
|
|---|
| 715 | // with a non-zero z order.
|
|---|
| 716 | //
|
|---|
| 717 | // One might worry about infinite loops, but the above conditional
|
|---|
| 718 | // should prevent it
|
|---|
| 719 | //
|
|---|
| 720 | G4PolyPhiFaceVertex *next = corn;
|
|---|
| 721 | G4double nextZ;
|
|---|
| 722 | do
|
|---|
| 723 | {
|
|---|
| 724 | next++;
|
|---|
| 725 | if (next == corners+numEdges) next = corners;
|
|---|
| 726 |
|
|---|
| 727 | nextZ = ExactZOrder( z, qx, qy, qz, v, normSign, next );
|
|---|
| 728 | } while( nextZ == 0 );
|
|---|
| 729 |
|
|---|
| 730 | //
|
|---|
| 731 | // If we won't be changing direction, go to the next vertex
|
|---|
| 732 | //
|
|---|
| 733 | if (nextZ*prevZ < 0) continue;
|
|---|
| 734 | }
|
|---|
| 735 |
|
|---|
| 736 |
|
|---|
| 737 | //
|
|---|
| 738 | // We overlap in z with the side of the face that stretches from
|
|---|
| 739 | // vertex "prev" to "corn". On which side (left or right) do
|
|---|
| 740 | // we lay with respect to this segment?
|
|---|
| 741 | //
|
|---|
| 742 | G4ThreeVector qa( qx - prev->x, qy - prev->y, qz - prev->z ),
|
|---|
| 743 | qb( qx - corn->x, qy - corn->y, qz - corn->z );
|
|---|
| 744 |
|
|---|
| 745 | G4double aboveOrBelow = normSign*qa.cross(qb).dot(v);
|
|---|
| 746 |
|
|---|
| 747 | if (aboveOrBelow > 0)
|
|---|
| 748 | answer++;
|
|---|
| 749 | else if (aboveOrBelow < 0)
|
|---|
| 750 | answer--;
|
|---|
| 751 | else
|
|---|
| 752 | {
|
|---|
| 753 | //
|
|---|
| 754 | // A precisely zero answer here means we exactly
|
|---|
| 755 | // intersect (within roundoff) the edge of the face.
|
|---|
| 756 | // Return true in this case.
|
|---|
| 757 | //
|
|---|
| 758 | return true;
|
|---|
| 759 | }
|
|---|
| 760 | } while( prevZ = cornZ, prev=corn, ++corn < corners+numEdges );
|
|---|
| 761 |
|
|---|
| 762 | // G4int fanswer = std::abs(answer);
|
|---|
| 763 | // if (fanswer==1 || fanswer>2) {
|
|---|
| 764 | // G4cerr << "G4PolyPhiFace::InsideEdgesExact: answer is "
|
|---|
| 765 | // << answer << G4endl;
|
|---|
| 766 | // }
|
|---|
| 767 |
|
|---|
| 768 | return answer!=0;
|
|---|
| 769 | }
|
|---|
| 770 |
|
|---|
| 771 |
|
|---|
| 772 | //
|
|---|
| 773 | // InsideEdges (don't care aboud distance)
|
|---|
| 774 | //
|
|---|
| 775 | // Decide if the point in r,z is inside the edges of our face
|
|---|
| 776 | //
|
|---|
| 777 | // This routine can be made a zillion times quicker by implementing
|
|---|
| 778 | // better code, for example:
|
|---|
| 779 | //
|
|---|
| 780 | // int pnpoly(int npol, float *xp, float *yp, float x, float y)
|
|---|
| 781 | // {
|
|---|
| 782 | // int i, j, c = 0;
|
|---|
| 783 | // for (i = 0, j = npol-1; i < npol; j = i++) {
|
|---|
| 784 | // if ((((yp[i]<=y) && (y<yp[j])) ||
|
|---|
| 785 | // ((yp[j]<=y) && (y<yp[i]))) &&
|
|---|
| 786 | // (x < (xp[j] - xp[i]) * (y - yp[i]) / (yp[j] - yp[i]) + xp[i]))
|
|---|
| 787 | //
|
|---|
| 788 | // c = !c;
|
|---|
| 789 | // }
|
|---|
| 790 | // return c;
|
|---|
| 791 | // }
|
|---|
| 792 | //
|
|---|
| 793 | // See "Point in Polyon Strategies", Eric Haines [Graphic Gems IV] pp. 24-46
|
|---|
| 794 | //
|
|---|
| 795 | // My algorithm below is rather unique, but is based on code needed to
|
|---|
| 796 | // calculate the distance to the shape. I left it in here because ...
|
|---|
| 797 | // well ... to test it better.
|
|---|
| 798 | //
|
|---|
| 799 | G4bool G4PolyPhiFace::InsideEdges( G4double r, G4double z )
|
|---|
| 800 | {
|
|---|
| 801 | //
|
|---|
| 802 | // Quick check of extent
|
|---|
| 803 | //
|
|---|
| 804 | if ( r < rMin || r > rMax ) return false;
|
|---|
| 805 | if ( z < zMin || z > zMax ) return false;
|
|---|
| 806 |
|
|---|
| 807 | //
|
|---|
| 808 | // More thorough check
|
|---|
| 809 | //
|
|---|
| 810 | G4double notUsed;
|
|---|
| 811 |
|
|---|
| 812 | return InsideEdges( r, z, ¬Used, 0 );
|
|---|
| 813 | }
|
|---|
| 814 |
|
|---|
| 815 |
|
|---|
| 816 | //
|
|---|
| 817 | // InsideEdges (care about distance)
|
|---|
| 818 | //
|
|---|
| 819 | // Decide if the point in r,z is inside the edges of our face
|
|---|
| 820 | //
|
|---|
| 821 | G4bool G4PolyPhiFace::InsideEdges( G4double r, G4double z,
|
|---|
| 822 | G4double *bestDist2,
|
|---|
| 823 | G4PolyPhiFaceVertex **base3Dnorm,
|
|---|
| 824 | G4ThreeVector **head3Dnorm )
|
|---|
| 825 | {
|
|---|
| 826 | G4double bestDistance2 = kInfinity;
|
|---|
| 827 | G4bool answer = 0;
|
|---|
| 828 |
|
|---|
| 829 | G4PolyPhiFaceEdge *edge = edges;
|
|---|
| 830 | do
|
|---|
| 831 | {
|
|---|
| 832 | G4PolyPhiFaceVertex *testMe;
|
|---|
| 833 | //
|
|---|
| 834 | // Get distance perpendicular to the edge
|
|---|
| 835 | //
|
|---|
| 836 | G4double dr = (r-edge->v0->r), dz = (z-edge->v0->z);
|
|---|
| 837 |
|
|---|
| 838 | G4double distOut = dr*edge->tz - dz*edge->tr;
|
|---|
| 839 | G4double distance2 = distOut*distOut;
|
|---|
| 840 | if (distance2 > bestDistance2) continue; // No hope!
|
|---|
| 841 |
|
|---|
| 842 | //
|
|---|
| 843 | // Check to see if normal intersects edge within the edge's boundary
|
|---|
| 844 | //
|
|---|
| 845 | G4double s = dr*edge->tr + dz*edge->tz;
|
|---|
| 846 |
|
|---|
| 847 | //
|
|---|
| 848 | // If it doesn't, penalize distance2 appropriately
|
|---|
| 849 | //
|
|---|
| 850 | if (s < 0)
|
|---|
| 851 | {
|
|---|
| 852 | distance2 += s*s;
|
|---|
| 853 | testMe = edge->v0;
|
|---|
| 854 | }
|
|---|
| 855 | else if (s > edge->length)
|
|---|
| 856 | {
|
|---|
| 857 | G4double s2 = s-edge->length;
|
|---|
| 858 | distance2 += s2*s2;
|
|---|
| 859 | testMe = edge->v1;
|
|---|
| 860 | }
|
|---|
| 861 | else
|
|---|
| 862 | {
|
|---|
| 863 | testMe = 0;
|
|---|
| 864 | }
|
|---|
| 865 |
|
|---|
| 866 | //
|
|---|
| 867 | // Closest edge so far?
|
|---|
| 868 | //
|
|---|
| 869 | if (distance2 < bestDistance2)
|
|---|
| 870 | {
|
|---|
| 871 | bestDistance2 = distance2;
|
|---|
| 872 | if (testMe)
|
|---|
| 873 | {
|
|---|
| 874 | G4double distNorm = dr*testMe->rNorm + dz*testMe->zNorm;
|
|---|
| 875 | answer = (distNorm <= 0);
|
|---|
| 876 | if (base3Dnorm)
|
|---|
| 877 | {
|
|---|
| 878 | *base3Dnorm = testMe;
|
|---|
| 879 | *head3Dnorm = &testMe->norm3D;
|
|---|
| 880 | }
|
|---|
| 881 | }
|
|---|
| 882 | else
|
|---|
| 883 | {
|
|---|
| 884 | answer = (distOut <= 0);
|
|---|
| 885 | if (base3Dnorm)
|
|---|
| 886 | {
|
|---|
| 887 | *base3Dnorm = edge->v0;
|
|---|
| 888 | *head3Dnorm = &edge->norm3D;
|
|---|
| 889 | }
|
|---|
| 890 | }
|
|---|
| 891 | }
|
|---|
| 892 | } while( ++edge < edges + numEdges );
|
|---|
| 893 |
|
|---|
| 894 | *bestDist2 = bestDistance2;
|
|---|
| 895 | return answer;
|
|---|
| 896 | }
|
|---|