[831] | 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 | // |
---|
[850] | 27 | // $Id: G4PolyPhiFace.cc,v 1.15 2008/05/15 11:41:59 gcosmo Exp $ |
---|
| 28 | // GEANT4 tag $Name: HEAD $ |
---|
[831] | 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 | |
---|
[850] | 49 | #include "Randomize.hh" |
---|
| 50 | #include "G4TwoVector.hh" |
---|
| 51 | |
---|
[831] | 52 | // |
---|
| 53 | // Constructor |
---|
| 54 | // |
---|
| 55 | // Points r,z should be supplied in clockwise order in r,z. For example: |
---|
| 56 | // |
---|
| 57 | // [1]---------[2] ^ R |
---|
| 58 | // | | | |
---|
| 59 | // | | +--> z |
---|
| 60 | // [0]---------[3] |
---|
| 61 | // |
---|
| 62 | G4PolyPhiFace::G4PolyPhiFace( const G4ReduciblePolygon *rz, |
---|
| 63 | G4double phi, |
---|
| 64 | G4double deltaPhi, |
---|
| 65 | G4double phiOther ) |
---|
| 66 | { |
---|
| 67 | kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); |
---|
[850] | 68 | fSurfaceArea = 0.; |
---|
[831] | 69 | |
---|
| 70 | numEdges = rz->NumVertices(); |
---|
| 71 | |
---|
| 72 | rMin = rz->Amin(); |
---|
| 73 | rMax = rz->Amax(); |
---|
| 74 | zMin = rz->Bmin(); |
---|
| 75 | zMax = rz->Bmax(); |
---|
| 76 | |
---|
| 77 | // |
---|
| 78 | // Is this the "starting" phi edge of the two? |
---|
| 79 | // |
---|
| 80 | G4bool start = (phiOther > phi); |
---|
| 81 | |
---|
| 82 | // |
---|
| 83 | // Build radial vector |
---|
| 84 | // |
---|
| 85 | radial = G4ThreeVector( std::cos(phi), std::sin(phi), 0.0 ); |
---|
| 86 | |
---|
| 87 | // |
---|
| 88 | // Build normal |
---|
| 89 | // |
---|
| 90 | G4double zSign = start ? 1 : -1; |
---|
| 91 | normal = G4ThreeVector( zSign*radial.y(), -zSign*radial.x(), 0 ); |
---|
| 92 | |
---|
| 93 | // |
---|
| 94 | // Is allBehind? |
---|
| 95 | // |
---|
| 96 | allBehind = (zSign*(std::cos(phiOther)*radial.y() - std::sin(phiOther)*radial.x()) < 0); |
---|
| 97 | |
---|
| 98 | // |
---|
| 99 | // Adjacent edges |
---|
| 100 | // |
---|
| 101 | G4double midPhi = phi + (start ? +0.5 : -0.5)*deltaPhi; |
---|
| 102 | G4double cosMid = std::cos(midPhi), |
---|
| 103 | sinMid = std::sin(midPhi); |
---|
| 104 | |
---|
| 105 | // |
---|
| 106 | // Allocate corners |
---|
| 107 | // |
---|
| 108 | corners = new G4PolyPhiFaceVertex[numEdges]; |
---|
| 109 | // |
---|
| 110 | // Fill them |
---|
| 111 | // |
---|
| 112 | G4ReduciblePolygonIterator iterRZ(rz); |
---|
| 113 | |
---|
| 114 | G4PolyPhiFaceVertex *corn = corners; |
---|
[850] | 115 | G4PolyPhiFaceVertex *helper=corners; |
---|
| 116 | |
---|
[831] | 117 | iterRZ.Begin(); |
---|
| 118 | do |
---|
| 119 | { |
---|
| 120 | corn->r = iterRZ.GetA(); |
---|
| 121 | corn->z = iterRZ.GetB(); |
---|
| 122 | corn->x = corn->r*radial.x(); |
---|
| 123 | corn->y = corn->r*radial.y(); |
---|
[850] | 124 | |
---|
| 125 | // Add pointer on prev corner |
---|
| 126 | // |
---|
| 127 | if( corn == corners ) |
---|
| 128 | { corn->prev = corners+numEdges-1;} |
---|
| 129 | else |
---|
| 130 | { corn->prev = helper; } |
---|
| 131 | |
---|
| 132 | // Add pointer on next corner |
---|
| 133 | // |
---|
| 134 | if( corn < corners+numEdges-1 ) |
---|
| 135 | { corn->next = corn+1;} |
---|
| 136 | else |
---|
| 137 | { corn->next = corners; } |
---|
| 138 | |
---|
| 139 | helper = corn; |
---|
[831] | 140 | } while( ++corn, iterRZ.Next() ); |
---|
| 141 | |
---|
| 142 | // |
---|
| 143 | // Allocate edges |
---|
| 144 | // |
---|
| 145 | edges = new G4PolyPhiFaceEdge[numEdges]; |
---|
| 146 | |
---|
| 147 | // |
---|
| 148 | // Fill them |
---|
| 149 | // |
---|
| 150 | G4double rFact = std::cos(0.5*deltaPhi); |
---|
| 151 | G4double rFactNormalize = 1.0/std::sqrt(1.0+rFact*rFact); |
---|
| 152 | |
---|
| 153 | G4PolyPhiFaceVertex *prev = corners+numEdges-1, |
---|
| 154 | *here = corners; |
---|
| 155 | G4PolyPhiFaceEdge *edge = edges; |
---|
| 156 | do |
---|
| 157 | { |
---|
| 158 | G4ThreeVector sideNorm; |
---|
| 159 | |
---|
| 160 | edge->v0 = prev; |
---|
| 161 | edge->v1 = here; |
---|
| 162 | |
---|
| 163 | G4double dr = here->r - prev->r, |
---|
| 164 | dz = here->z - prev->z; |
---|
| 165 | |
---|
| 166 | edge->length = std::sqrt( dr*dr + dz*dz ); |
---|
| 167 | |
---|
| 168 | edge->tr = dr/edge->length; |
---|
| 169 | edge->tz = dz/edge->length; |
---|
| 170 | |
---|
| 171 | if ((here->r < DBL_MIN) && (prev->r < DBL_MIN)) |
---|
| 172 | { |
---|
| 173 | // |
---|
| 174 | // Sigh! Always exceptions! |
---|
| 175 | // This edge runs at r==0, so its adjoing surface is not a |
---|
| 176 | // PolyconeSide or PolyhedraSide, but the opposite PolyPhiFace. |
---|
| 177 | // |
---|
| 178 | G4double zSignOther = start ? -1 : 1; |
---|
| 179 | sideNorm = G4ThreeVector( zSignOther*std::sin(phiOther), |
---|
| 180 | -zSignOther*std::cos(phiOther), 0 ); |
---|
| 181 | } |
---|
| 182 | else |
---|
| 183 | { |
---|
| 184 | sideNorm = G4ThreeVector( edge->tz*cosMid, |
---|
| 185 | edge->tz*sinMid, |
---|
| 186 | -edge->tr*rFact ); |
---|
| 187 | sideNorm *= rFactNormalize; |
---|
| 188 | } |
---|
| 189 | sideNorm += normal; |
---|
| 190 | |
---|
| 191 | edge->norm3D = sideNorm.unit(); |
---|
| 192 | } while( edge++, prev=here, ++here < corners+numEdges ); |
---|
| 193 | |
---|
| 194 | // |
---|
| 195 | // Go back and fill in corner "normals" |
---|
| 196 | // |
---|
| 197 | G4PolyPhiFaceEdge *prevEdge = edges+numEdges-1; |
---|
| 198 | edge = edges; |
---|
| 199 | do |
---|
| 200 | { |
---|
| 201 | // |
---|
| 202 | // Calculate vertex 2D normals (on the phi surface) |
---|
| 203 | // |
---|
| 204 | G4double rPart = prevEdge->tr + edge->tr; |
---|
| 205 | G4double zPart = prevEdge->tz + edge->tz; |
---|
| 206 | G4double norm = std::sqrt( rPart*rPart + zPart*zPart ); |
---|
| 207 | G4double rNorm = +zPart/norm; |
---|
| 208 | G4double zNorm = -rPart/norm; |
---|
| 209 | |
---|
| 210 | edge->v0->rNorm = rNorm; |
---|
| 211 | edge->v0->zNorm = zNorm; |
---|
| 212 | |
---|
| 213 | // |
---|
| 214 | // Calculate the 3D normals. |
---|
| 215 | // |
---|
| 216 | // Find the vector perpendicular to the z axis |
---|
| 217 | // that defines the plane that contains the vertex normal |
---|
| 218 | // |
---|
| 219 | G4ThreeVector xyVector; |
---|
| 220 | |
---|
| 221 | if (edge->v0->r < DBL_MIN) |
---|
| 222 | { |
---|
| 223 | // |
---|
| 224 | // This is a vertex at r==0, which is a special |
---|
| 225 | // case. The normal we will construct lays in the |
---|
| 226 | // plane at the center of the phi opening. |
---|
| 227 | // |
---|
| 228 | // We also know that rNorm < 0 |
---|
| 229 | // |
---|
| 230 | G4double zSignOther = start ? -1 : 1; |
---|
| 231 | G4ThreeVector normalOther( zSignOther*std::sin(phiOther), |
---|
| 232 | -zSignOther*std::cos(phiOther), 0 ); |
---|
| 233 | |
---|
| 234 | xyVector = - normal - normalOther; |
---|
| 235 | } |
---|
| 236 | else |
---|
| 237 | { |
---|
| 238 | // |
---|
| 239 | // This is a vertex at r > 0. The plane |
---|
| 240 | // is the average of the normal and the |
---|
| 241 | // normal of the adjacent phi face |
---|
| 242 | // |
---|
| 243 | xyVector = G4ThreeVector( cosMid, sinMid, 0 ); |
---|
| 244 | if (rNorm < 0) |
---|
| 245 | xyVector -= normal; |
---|
| 246 | else |
---|
| 247 | xyVector += normal; |
---|
| 248 | } |
---|
| 249 | |
---|
| 250 | // |
---|
| 251 | // Combine it with the r/z direction from the face |
---|
| 252 | // |
---|
| 253 | edge->v0->norm3D = rNorm*xyVector.unit() + G4ThreeVector( 0, 0, zNorm ); |
---|
| 254 | } while( prevEdge=edge, ++edge < edges+numEdges ); |
---|
| 255 | |
---|
| 256 | // |
---|
| 257 | // Build point on surface |
---|
| 258 | // |
---|
| 259 | G4double rAve = 0.5*(rMax-rMin), |
---|
| 260 | zAve = 0.5*(zMax-zMin); |
---|
| 261 | surface = G4ThreeVector( rAve*radial.x(), rAve*radial.y(), zAve ); |
---|
| 262 | } |
---|
| 263 | |
---|
| 264 | |
---|
| 265 | // |
---|
| 266 | // Diagnose |
---|
| 267 | // |
---|
| 268 | // Throw an exception if something is found inconsistent with |
---|
| 269 | // the solid. |
---|
| 270 | // |
---|
| 271 | // For debugging purposes only |
---|
| 272 | // |
---|
| 273 | void G4PolyPhiFace::Diagnose( G4VSolid *owner ) |
---|
| 274 | { |
---|
| 275 | G4PolyPhiFaceVertex *corner = corners; |
---|
| 276 | do |
---|
| 277 | { |
---|
| 278 | G4ThreeVector test(corner->x, corner->y, corner->z); |
---|
| 279 | test -= 1E-6*corner->norm3D; |
---|
| 280 | |
---|
| 281 | if (owner->Inside(test) != kInside) |
---|
| 282 | G4Exception( "G4PolyPhiFace::Diagnose()", "InvalidSetup", |
---|
| 283 | FatalException, "Bad vertex normal found." ); |
---|
| 284 | } while( ++corner < corners+numEdges ); |
---|
| 285 | } |
---|
| 286 | |
---|
| 287 | |
---|
| 288 | // |
---|
| 289 | // Fake default constructor - sets only member data and allocates memory |
---|
| 290 | // for usage restricted to object persistency. |
---|
| 291 | // |
---|
| 292 | G4PolyPhiFace::G4PolyPhiFace( __void__&) |
---|
| 293 | : edges(0), corners(0) |
---|
| 294 | { |
---|
| 295 | } |
---|
| 296 | |
---|
| 297 | |
---|
| 298 | // |
---|
| 299 | // Destructor |
---|
| 300 | // |
---|
| 301 | G4PolyPhiFace::~G4PolyPhiFace() |
---|
| 302 | { |
---|
| 303 | delete [] edges; |
---|
| 304 | delete [] corners; |
---|
| 305 | } |
---|
| 306 | |
---|
| 307 | |
---|
| 308 | // |
---|
| 309 | // Copy constructor |
---|
| 310 | // |
---|
| 311 | G4PolyPhiFace::G4PolyPhiFace( const G4PolyPhiFace &source ) |
---|
| 312 | : G4VCSGface() |
---|
| 313 | { |
---|
| 314 | CopyStuff( source ); |
---|
| 315 | } |
---|
| 316 | |
---|
| 317 | |
---|
| 318 | // |
---|
| 319 | // Assignment operator |
---|
| 320 | // |
---|
| 321 | G4PolyPhiFace& G4PolyPhiFace::operator=( const G4PolyPhiFace &source ) |
---|
| 322 | { |
---|
| 323 | if (this == &source) return *this; |
---|
| 324 | |
---|
| 325 | delete [] edges; |
---|
| 326 | delete [] corners; |
---|
| 327 | |
---|
| 328 | CopyStuff( source ); |
---|
| 329 | |
---|
| 330 | return *this; |
---|
| 331 | } |
---|
| 332 | |
---|
| 333 | |
---|
| 334 | // |
---|
| 335 | // CopyStuff (protected) |
---|
| 336 | // |
---|
| 337 | void G4PolyPhiFace::CopyStuff( const G4PolyPhiFace &source ) |
---|
| 338 | { |
---|
| 339 | // |
---|
| 340 | // The simple stuff |
---|
| 341 | // |
---|
| 342 | numEdges = source.numEdges; |
---|
| 343 | normal = source.normal; |
---|
| 344 | radial = source.radial; |
---|
[850] | 345 | surface = source.surface; |
---|
[831] | 346 | rMin = source.rMin; |
---|
| 347 | rMax = source.rMax; |
---|
| 348 | zMin = source.zMin; |
---|
| 349 | zMax = source.zMax; |
---|
| 350 | allBehind = source.allBehind; |
---|
| 351 | |
---|
| 352 | kCarTolerance = source.kCarTolerance; |
---|
[850] | 353 | fSurfaceArea = source.fSurfaceArea; |
---|
[831] | 354 | |
---|
| 355 | // |
---|
| 356 | // Corner dynamic array |
---|
| 357 | // |
---|
| 358 | corners = new G4PolyPhiFaceVertex[numEdges]; |
---|
| 359 | G4PolyPhiFaceVertex *corn = corners, |
---|
| 360 | *sourceCorn = source.corners; |
---|
| 361 | do |
---|
| 362 | { |
---|
| 363 | *corn = *sourceCorn; |
---|
| 364 | } while( ++sourceCorn, ++corn < corners+numEdges ); |
---|
| 365 | |
---|
| 366 | // |
---|
| 367 | // Edge dynamic array |
---|
| 368 | // |
---|
| 369 | edges = new G4PolyPhiFaceEdge[numEdges]; |
---|
| 370 | |
---|
| 371 | G4PolyPhiFaceVertex *prev = corners+numEdges-1, |
---|
| 372 | *here = corners; |
---|
| 373 | G4PolyPhiFaceEdge *edge = edges, |
---|
| 374 | *sourceEdge = source.edges; |
---|
| 375 | do |
---|
| 376 | { |
---|
| 377 | *edge = *sourceEdge; |
---|
| 378 | edge->v0 = prev; |
---|
| 379 | edge->v1 = here; |
---|
| 380 | } while( ++sourceEdge, ++edge, prev=here, ++here < corners+numEdges ); |
---|
| 381 | } |
---|
| 382 | |
---|
| 383 | |
---|
| 384 | // |
---|
| 385 | // Intersect |
---|
| 386 | // |
---|
| 387 | G4bool G4PolyPhiFace::Intersect( const G4ThreeVector &p, |
---|
| 388 | const G4ThreeVector &v, |
---|
| 389 | G4bool outgoing, |
---|
| 390 | G4double surfTolerance, |
---|
| 391 | G4double &distance, |
---|
| 392 | G4double &distFromSurface, |
---|
| 393 | G4ThreeVector &aNormal, |
---|
| 394 | G4bool &isAllBehind ) |
---|
| 395 | { |
---|
| 396 | G4double normSign = outgoing ? +1 : -1; |
---|
| 397 | |
---|
| 398 | // |
---|
| 399 | // These don't change |
---|
| 400 | // |
---|
| 401 | isAllBehind = allBehind; |
---|
| 402 | aNormal = normal; |
---|
| 403 | |
---|
| 404 | // |
---|
| 405 | // Correct normal? Here we have straight sides, and can safely ignore |
---|
| 406 | // intersections where the dot product with the normal is zero. |
---|
| 407 | // |
---|
| 408 | G4double dotProd = normSign*normal.dot(v); |
---|
| 409 | |
---|
| 410 | if (dotProd <= 0) return false; |
---|
| 411 | |
---|
| 412 | // |
---|
| 413 | // Calculate distance to surface. If the side is too far |
---|
| 414 | // behind the point, we must reject it. |
---|
| 415 | // |
---|
| 416 | G4ThreeVector ps = p - surface; |
---|
| 417 | distFromSurface = -normSign*ps.dot(normal); |
---|
| 418 | |
---|
| 419 | if (distFromSurface < -surfTolerance) return false; |
---|
| 420 | |
---|
| 421 | // |
---|
| 422 | // Calculate precise distance to intersection with the side |
---|
| 423 | // (along the trajectory, not normal to the surface) |
---|
| 424 | // |
---|
| 425 | distance = distFromSurface/dotProd; |
---|
| 426 | |
---|
| 427 | // |
---|
| 428 | // Calculate intersection point in r,z |
---|
| 429 | // |
---|
| 430 | G4ThreeVector ip = p + distance*v; |
---|
| 431 | |
---|
| 432 | G4double r = radial.dot(ip); |
---|
| 433 | |
---|
| 434 | // |
---|
| 435 | // And is it inside the r/z extent? |
---|
| 436 | // |
---|
| 437 | return InsideEdgesExact( r, ip.z(), normSign, p, v ); |
---|
| 438 | } |
---|
| 439 | |
---|
| 440 | |
---|
| 441 | // |
---|
| 442 | // Distance |
---|
| 443 | // |
---|
| 444 | G4double G4PolyPhiFace::Distance( const G4ThreeVector &p, G4bool outgoing ) |
---|
| 445 | { |
---|
| 446 | G4double normSign = outgoing ? +1 : -1; |
---|
| 447 | // |
---|
| 448 | // Correct normal? |
---|
| 449 | // |
---|
| 450 | G4ThreeVector ps = p - surface; |
---|
| 451 | G4double distPhi = -normSign*normal.dot(ps); |
---|
| 452 | |
---|
| 453 | if (distPhi < -0.5*kCarTolerance) |
---|
| 454 | return kInfinity; |
---|
| 455 | else if (distPhi < 0) |
---|
| 456 | distPhi = 0.0; |
---|
| 457 | |
---|
| 458 | // |
---|
| 459 | // Calculate projected point in r,z |
---|
| 460 | // |
---|
| 461 | G4double r = radial.dot(p); |
---|
| 462 | |
---|
| 463 | // |
---|
| 464 | // Are we inside the face? |
---|
| 465 | // |
---|
| 466 | G4double distRZ2; |
---|
| 467 | |
---|
| 468 | if (InsideEdges( r, p.z(), &distRZ2, 0 )) |
---|
| 469 | { |
---|
| 470 | // |
---|
| 471 | // Yup, answer is just distPhi |
---|
| 472 | // |
---|
| 473 | return distPhi; |
---|
| 474 | } |
---|
| 475 | else |
---|
| 476 | { |
---|
| 477 | // |
---|
| 478 | // Nope. Penalize by distance out |
---|
| 479 | // |
---|
| 480 | return std::sqrt( distPhi*distPhi + distRZ2 ); |
---|
| 481 | } |
---|
| 482 | } |
---|
| 483 | |
---|
| 484 | |
---|
| 485 | // |
---|
| 486 | // Inside |
---|
| 487 | // |
---|
| 488 | EInside G4PolyPhiFace::Inside( const G4ThreeVector &p, |
---|
| 489 | G4double tolerance, |
---|
| 490 | G4double *bestDistance ) |
---|
| 491 | { |
---|
| 492 | // |
---|
| 493 | // Get distance along phi, which if negative means the point |
---|
| 494 | // is nominally inside the shape. |
---|
| 495 | // |
---|
| 496 | G4ThreeVector ps = p - surface; |
---|
| 497 | G4double distPhi = normal.dot(ps); |
---|
| 498 | |
---|
| 499 | // |
---|
| 500 | // Calculate projected point in r,z |
---|
| 501 | // |
---|
| 502 | G4double r = radial.dot(p); |
---|
| 503 | |
---|
| 504 | // |
---|
| 505 | // Are we inside the face? |
---|
| 506 | // |
---|
| 507 | G4double distRZ2; |
---|
| 508 | G4PolyPhiFaceVertex *base3Dnorm; |
---|
| 509 | G4ThreeVector *head3Dnorm; |
---|
| 510 | |
---|
| 511 | if (InsideEdges( r, p.z(), &distRZ2, &base3Dnorm, &head3Dnorm )) |
---|
| 512 | { |
---|
| 513 | // |
---|
| 514 | // Looks like we're inside. Distance is distance in phi. |
---|
| 515 | // |
---|
| 516 | *bestDistance = std::fabs(distPhi); |
---|
| 517 | |
---|
| 518 | // |
---|
| 519 | // Use distPhi to decide fate |
---|
| 520 | // |
---|
| 521 | if (distPhi < -tolerance) return kInside; |
---|
| 522 | if (distPhi < tolerance) return kSurface; |
---|
| 523 | return kOutside; |
---|
| 524 | } |
---|
| 525 | else |
---|
| 526 | { |
---|
| 527 | // |
---|
| 528 | // We're outside the extent of the face, |
---|
| 529 | // so the distance is penalized by distance from edges in RZ |
---|
| 530 | // |
---|
| 531 | *bestDistance = std::sqrt( distPhi*distPhi + distRZ2 ); |
---|
| 532 | |
---|
| 533 | // |
---|
| 534 | // Use edge normal to decide fate |
---|
| 535 | // |
---|
| 536 | G4ThreeVector cc( base3Dnorm->r*radial.x(), |
---|
| 537 | base3Dnorm->r*radial.y(), |
---|
| 538 | base3Dnorm->z ); |
---|
| 539 | cc = p - cc; |
---|
| 540 | G4double normDist = head3Dnorm->dot(cc); |
---|
| 541 | if ( distRZ2 > tolerance*tolerance ) |
---|
| 542 | { |
---|
| 543 | // |
---|
| 544 | // We're far enough away that kSurface is not possible |
---|
| 545 | // |
---|
| 546 | return normDist < 0 ? kInside : kOutside; |
---|
| 547 | } |
---|
| 548 | |
---|
| 549 | if (normDist < -tolerance) return kInside; |
---|
| 550 | if (normDist < tolerance) return kSurface; |
---|
| 551 | return kOutside; |
---|
| 552 | } |
---|
| 553 | } |
---|
| 554 | |
---|
| 555 | |
---|
| 556 | // |
---|
| 557 | // Normal |
---|
| 558 | // |
---|
| 559 | // This virtual member is simple for our planer shape, |
---|
| 560 | // which has only one normal |
---|
| 561 | // |
---|
| 562 | G4ThreeVector G4PolyPhiFace::Normal( const G4ThreeVector &p, |
---|
| 563 | G4double *bestDistance ) |
---|
| 564 | { |
---|
| 565 | // |
---|
| 566 | // Get distance along phi, which if negative means the point |
---|
| 567 | // is nominally inside the shape. |
---|
| 568 | // |
---|
| 569 | G4double distPhi = normal.dot(p); |
---|
| 570 | |
---|
| 571 | // |
---|
| 572 | // Calculate projected point in r,z |
---|
| 573 | // |
---|
| 574 | G4double r = radial.dot(p); |
---|
| 575 | |
---|
| 576 | // |
---|
| 577 | // Are we inside the face? |
---|
| 578 | // |
---|
| 579 | G4double distRZ2; |
---|
| 580 | |
---|
| 581 | if (InsideEdges( r, p.z(), &distRZ2, 0 )) |
---|
| 582 | { |
---|
| 583 | // |
---|
| 584 | // Yup, answer is just distPhi |
---|
| 585 | // |
---|
| 586 | *bestDistance = std::fabs(distPhi); |
---|
| 587 | } |
---|
| 588 | else |
---|
| 589 | { |
---|
| 590 | // |
---|
| 591 | // Nope. Penalize by distance out |
---|
| 592 | // |
---|
| 593 | *bestDistance = std::sqrt( distPhi*distPhi + distRZ2 ); |
---|
| 594 | } |
---|
| 595 | |
---|
| 596 | return normal; |
---|
| 597 | } |
---|
| 598 | |
---|
| 599 | |
---|
| 600 | // |
---|
| 601 | // Extent |
---|
| 602 | // |
---|
| 603 | // This actually isn't needed by polycone or polyhedra... |
---|
| 604 | // |
---|
| 605 | G4double G4PolyPhiFace::Extent( const G4ThreeVector axis ) |
---|
| 606 | { |
---|
| 607 | G4double max = -kInfinity; |
---|
| 608 | |
---|
| 609 | G4PolyPhiFaceVertex *corner = corners; |
---|
| 610 | do |
---|
| 611 | { |
---|
| 612 | G4double here = axis.x()*corner->r*radial.x() |
---|
| 613 | + axis.y()*corner->r*radial.y() |
---|
| 614 | + axis.z()*corner->z; |
---|
| 615 | if (here > max) max = here; |
---|
| 616 | } while( ++corner < corners + numEdges ); |
---|
| 617 | |
---|
| 618 | return max; |
---|
| 619 | } |
---|
| 620 | |
---|
| 621 | |
---|
| 622 | // |
---|
| 623 | // CalculateExtent |
---|
| 624 | // |
---|
| 625 | // See notes in G4VCSGface |
---|
| 626 | // |
---|
| 627 | void G4PolyPhiFace::CalculateExtent( const EAxis axis, |
---|
| 628 | const G4VoxelLimits &voxelLimit, |
---|
| 629 | const G4AffineTransform &transform, |
---|
| 630 | G4SolidExtentList &extentList ) |
---|
| 631 | { |
---|
| 632 | // |
---|
| 633 | // Construct a (sometimes big) clippable polygon, |
---|
| 634 | // |
---|
| 635 | // Perform the necessary transformations while doing so |
---|
| 636 | // |
---|
| 637 | G4ClippablePolygon polygon; |
---|
| 638 | |
---|
| 639 | G4PolyPhiFaceVertex *corner = corners; |
---|
| 640 | do |
---|
| 641 | { |
---|
| 642 | G4ThreeVector point( 0, 0, corner->z ); |
---|
| 643 | point += radial*corner->r; |
---|
| 644 | |
---|
| 645 | polygon.AddVertexInOrder( transform.TransformPoint( point ) ); |
---|
| 646 | } while( ++corner < corners + numEdges ); |
---|
| 647 | |
---|
| 648 | // |
---|
| 649 | // Clip away |
---|
| 650 | // |
---|
| 651 | if (polygon.PartialClip( voxelLimit, axis )) |
---|
| 652 | { |
---|
| 653 | // |
---|
| 654 | // Add it to the list |
---|
| 655 | // |
---|
| 656 | polygon.SetNormal( transform.TransformAxis(normal) ); |
---|
| 657 | extentList.AddSurface( polygon ); |
---|
| 658 | } |
---|
| 659 | } |
---|
| 660 | |
---|
| 661 | |
---|
| 662 | // |
---|
| 663 | //------------------------------------------------------- |
---|
| 664 | |
---|
| 665 | |
---|
| 666 | // |
---|
| 667 | // InsideEdgesExact |
---|
| 668 | // |
---|
| 669 | // Decide if the point in r,z is inside the edges of our face, |
---|
| 670 | // **but** do so consistently with other faces. |
---|
| 671 | // |
---|
| 672 | // This routine has functionality similar to InsideEdges, but uses |
---|
| 673 | // an algorithm to decide if a trajectory falls inside or outside the |
---|
| 674 | // face that uses only the trajectory p,v values and the three dimensional |
---|
| 675 | // points representing the edges of the polygon. The objective is to plug up |
---|
| 676 | // any leaks between touching G4PolyPhiFaces (at r==0) and any other face |
---|
| 677 | // that uses the same convention. |
---|
| 678 | // |
---|
| 679 | // See: "Computational Geometry in C (Second Edition)" |
---|
| 680 | // http://cs.smith.edu/~orourke/ |
---|
| 681 | // |
---|
| 682 | G4bool G4PolyPhiFace::InsideEdgesExact( G4double r, G4double z, |
---|
| 683 | G4double normSign, |
---|
| 684 | const G4ThreeVector &p, |
---|
| 685 | const G4ThreeVector &v ) |
---|
| 686 | { |
---|
| 687 | // |
---|
| 688 | // Quick check of extent |
---|
| 689 | // |
---|
| 690 | if ( (r < rMin-kCarTolerance) |
---|
| 691 | || (r > rMax+kCarTolerance) ) return false; |
---|
| 692 | |
---|
| 693 | if ( (z < zMin-kCarTolerance) |
---|
| 694 | || (z > zMax+kCarTolerance) ) return false; |
---|
| 695 | |
---|
| 696 | // |
---|
| 697 | // Exact check: loop over all vertices |
---|
| 698 | // |
---|
| 699 | G4double qx = p.x() + v.x(), |
---|
| 700 | qy = p.y() + v.y(), |
---|
| 701 | qz = p.z() + v.z(); |
---|
| 702 | |
---|
| 703 | G4int answer = 0; |
---|
| 704 | G4PolyPhiFaceVertex *corn = corners, |
---|
| 705 | *prev = corners+numEdges-1; |
---|
| 706 | |
---|
| 707 | G4double cornZ, prevZ; |
---|
| 708 | |
---|
| 709 | prevZ = ExactZOrder( z, qx, qy, qz, v, normSign, prev ); |
---|
| 710 | do |
---|
| 711 | { |
---|
| 712 | // |
---|
| 713 | // Get z order of this vertex, and compare to previous vertex |
---|
| 714 | // |
---|
| 715 | cornZ = ExactZOrder( z, qx, qy, qz, v, normSign, corn ); |
---|
| 716 | |
---|
| 717 | if (cornZ < 0) |
---|
| 718 | { |
---|
| 719 | if (prevZ < 0) continue; |
---|
| 720 | } |
---|
| 721 | else if (cornZ > 0) |
---|
| 722 | { |
---|
| 723 | if (prevZ > 0) continue; |
---|
| 724 | } |
---|
| 725 | else |
---|
| 726 | { |
---|
| 727 | // |
---|
| 728 | // By chance, we overlap exactly (within precision) with |
---|
| 729 | // the current vertex. Continue if the same happened previously |
---|
| 730 | // (e.g. the previous vertex had the same z value) |
---|
| 731 | // |
---|
| 732 | if (prevZ == 0) continue; |
---|
| 733 | |
---|
| 734 | // |
---|
| 735 | // Otherwise, to decide what to do, we need to know what is |
---|
| 736 | // coming up next. Specifically, we need to find the next vertex |
---|
| 737 | // with a non-zero z order. |
---|
| 738 | // |
---|
| 739 | // One might worry about infinite loops, but the above conditional |
---|
| 740 | // should prevent it |
---|
| 741 | // |
---|
| 742 | G4PolyPhiFaceVertex *next = corn; |
---|
| 743 | G4double nextZ; |
---|
| 744 | do |
---|
| 745 | { |
---|
| 746 | next++; |
---|
| 747 | if (next == corners+numEdges) next = corners; |
---|
| 748 | |
---|
| 749 | nextZ = ExactZOrder( z, qx, qy, qz, v, normSign, next ); |
---|
| 750 | } while( nextZ == 0 ); |
---|
| 751 | |
---|
| 752 | // |
---|
| 753 | // If we won't be changing direction, go to the next vertex |
---|
| 754 | // |
---|
| 755 | if (nextZ*prevZ < 0) continue; |
---|
| 756 | } |
---|
| 757 | |
---|
| 758 | |
---|
| 759 | // |
---|
| 760 | // We overlap in z with the side of the face that stretches from |
---|
| 761 | // vertex "prev" to "corn". On which side (left or right) do |
---|
| 762 | // we lay with respect to this segment? |
---|
| 763 | // |
---|
| 764 | G4ThreeVector qa( qx - prev->x, qy - prev->y, qz - prev->z ), |
---|
| 765 | qb( qx - corn->x, qy - corn->y, qz - corn->z ); |
---|
| 766 | |
---|
| 767 | G4double aboveOrBelow = normSign*qa.cross(qb).dot(v); |
---|
| 768 | |
---|
| 769 | if (aboveOrBelow > 0) |
---|
| 770 | answer++; |
---|
| 771 | else if (aboveOrBelow < 0) |
---|
| 772 | answer--; |
---|
| 773 | else |
---|
| 774 | { |
---|
| 775 | // |
---|
| 776 | // A precisely zero answer here means we exactly |
---|
| 777 | // intersect (within roundoff) the edge of the face. |
---|
| 778 | // Return true in this case. |
---|
| 779 | // |
---|
| 780 | return true; |
---|
| 781 | } |
---|
| 782 | } while( prevZ = cornZ, prev=corn, ++corn < corners+numEdges ); |
---|
| 783 | |
---|
| 784 | // G4int fanswer = std::abs(answer); |
---|
| 785 | // if (fanswer==1 || fanswer>2) { |
---|
| 786 | // G4cerr << "G4PolyPhiFace::InsideEdgesExact: answer is " |
---|
| 787 | // << answer << G4endl; |
---|
| 788 | // } |
---|
| 789 | |
---|
| 790 | return answer!=0; |
---|
| 791 | } |
---|
| 792 | |
---|
| 793 | |
---|
| 794 | // |
---|
| 795 | // InsideEdges (don't care aboud distance) |
---|
| 796 | // |
---|
| 797 | // Decide if the point in r,z is inside the edges of our face |
---|
| 798 | // |
---|
| 799 | // This routine can be made a zillion times quicker by implementing |
---|
| 800 | // better code, for example: |
---|
| 801 | // |
---|
| 802 | // int pnpoly(int npol, float *xp, float *yp, float x, float y) |
---|
| 803 | // { |
---|
| 804 | // int i, j, c = 0; |
---|
| 805 | // for (i = 0, j = npol-1; i < npol; j = i++) { |
---|
| 806 | // if ((((yp[i]<=y) && (y<yp[j])) || |
---|
| 807 | // ((yp[j]<=y) && (y<yp[i]))) && |
---|
| 808 | // (x < (xp[j] - xp[i]) * (y - yp[i]) / (yp[j] - yp[i]) + xp[i])) |
---|
| 809 | // |
---|
| 810 | // c = !c; |
---|
| 811 | // } |
---|
| 812 | // return c; |
---|
| 813 | // } |
---|
| 814 | // |
---|
| 815 | // See "Point in Polyon Strategies", Eric Haines [Graphic Gems IV] pp. 24-46 |
---|
| 816 | // |
---|
| 817 | // My algorithm below is rather unique, but is based on code needed to |
---|
| 818 | // calculate the distance to the shape. I left it in here because ... |
---|
| 819 | // well ... to test it better. |
---|
| 820 | // |
---|
| 821 | G4bool G4PolyPhiFace::InsideEdges( G4double r, G4double z ) |
---|
| 822 | { |
---|
| 823 | // |
---|
| 824 | // Quick check of extent |
---|
| 825 | // |
---|
| 826 | if ( r < rMin || r > rMax ) return false; |
---|
| 827 | if ( z < zMin || z > zMax ) return false; |
---|
| 828 | |
---|
| 829 | // |
---|
| 830 | // More thorough check |
---|
| 831 | // |
---|
| 832 | G4double notUsed; |
---|
| 833 | |
---|
| 834 | return InsideEdges( r, z, ¬Used, 0 ); |
---|
| 835 | } |
---|
| 836 | |
---|
| 837 | |
---|
| 838 | // |
---|
| 839 | // InsideEdges (care about distance) |
---|
| 840 | // |
---|
| 841 | // Decide if the point in r,z is inside the edges of our face |
---|
| 842 | // |
---|
| 843 | G4bool G4PolyPhiFace::InsideEdges( G4double r, G4double z, |
---|
| 844 | G4double *bestDist2, |
---|
| 845 | G4PolyPhiFaceVertex **base3Dnorm, |
---|
| 846 | G4ThreeVector **head3Dnorm ) |
---|
| 847 | { |
---|
| 848 | G4double bestDistance2 = kInfinity; |
---|
| 849 | G4bool answer = 0; |
---|
| 850 | |
---|
| 851 | G4PolyPhiFaceEdge *edge = edges; |
---|
| 852 | do |
---|
| 853 | { |
---|
| 854 | G4PolyPhiFaceVertex *testMe; |
---|
| 855 | // |
---|
| 856 | // Get distance perpendicular to the edge |
---|
| 857 | // |
---|
| 858 | G4double dr = (r-edge->v0->r), dz = (z-edge->v0->z); |
---|
| 859 | |
---|
| 860 | G4double distOut = dr*edge->tz - dz*edge->tr; |
---|
| 861 | G4double distance2 = distOut*distOut; |
---|
| 862 | if (distance2 > bestDistance2) continue; // No hope! |
---|
| 863 | |
---|
| 864 | // |
---|
| 865 | // Check to see if normal intersects edge within the edge's boundary |
---|
| 866 | // |
---|
| 867 | G4double s = dr*edge->tr + dz*edge->tz; |
---|
| 868 | |
---|
| 869 | // |
---|
| 870 | // If it doesn't, penalize distance2 appropriately |
---|
| 871 | // |
---|
| 872 | if (s < 0) |
---|
| 873 | { |
---|
| 874 | distance2 += s*s; |
---|
| 875 | testMe = edge->v0; |
---|
| 876 | } |
---|
| 877 | else if (s > edge->length) |
---|
| 878 | { |
---|
| 879 | G4double s2 = s-edge->length; |
---|
| 880 | distance2 += s2*s2; |
---|
| 881 | testMe = edge->v1; |
---|
| 882 | } |
---|
| 883 | else |
---|
| 884 | { |
---|
| 885 | testMe = 0; |
---|
| 886 | } |
---|
| 887 | |
---|
| 888 | // |
---|
| 889 | // Closest edge so far? |
---|
| 890 | // |
---|
| 891 | if (distance2 < bestDistance2) |
---|
| 892 | { |
---|
| 893 | bestDistance2 = distance2; |
---|
| 894 | if (testMe) |
---|
| 895 | { |
---|
| 896 | G4double distNorm = dr*testMe->rNorm + dz*testMe->zNorm; |
---|
| 897 | answer = (distNorm <= 0); |
---|
| 898 | if (base3Dnorm) |
---|
| 899 | { |
---|
| 900 | *base3Dnorm = testMe; |
---|
| 901 | *head3Dnorm = &testMe->norm3D; |
---|
| 902 | } |
---|
| 903 | } |
---|
| 904 | else |
---|
| 905 | { |
---|
| 906 | answer = (distOut <= 0); |
---|
| 907 | if (base3Dnorm) |
---|
| 908 | { |
---|
| 909 | *base3Dnorm = edge->v0; |
---|
| 910 | *head3Dnorm = &edge->norm3D; |
---|
| 911 | } |
---|
| 912 | } |
---|
| 913 | } |
---|
| 914 | } while( ++edge < edges + numEdges ); |
---|
| 915 | |
---|
| 916 | *bestDist2 = bestDistance2; |
---|
| 917 | return answer; |
---|
| 918 | } |
---|
[850] | 919 | |
---|
| 920 | // |
---|
| 921 | // Calculation of Surface Area of a Triangle |
---|
| 922 | // In the same time Random Point in Triangle is given |
---|
| 923 | // |
---|
| 924 | G4double G4PolyPhiFace::SurfaceTriangle( G4ThreeVector p1, |
---|
| 925 | G4ThreeVector p2, |
---|
| 926 | G4ThreeVector p3, |
---|
| 927 | G4ThreeVector *p4 ) |
---|
| 928 | { |
---|
| 929 | G4ThreeVector v, w; |
---|
| 930 | |
---|
| 931 | v = p3 - p1; |
---|
| 932 | w = p1 - p2; |
---|
| 933 | G4double lambda1 = G4UniformRand(); |
---|
| 934 | G4double lambda2 = lambda1*G4UniformRand(); |
---|
| 935 | |
---|
| 936 | *p4=p2 + lambda1*w + lambda2*v; |
---|
| 937 | return 0.5*(v.cross(w)).mag(); |
---|
| 938 | } |
---|
| 939 | |
---|
| 940 | // |
---|
| 941 | // Compute surface area |
---|
| 942 | // |
---|
| 943 | G4double G4PolyPhiFace::SurfaceArea() |
---|
| 944 | { |
---|
| 945 | if ( fSurfaceArea==0. ) { Triangulate(); } |
---|
| 946 | return fSurfaceArea; |
---|
| 947 | } |
---|
| 948 | |
---|
| 949 | // |
---|
| 950 | // Return random point on face |
---|
| 951 | // |
---|
| 952 | G4ThreeVector G4PolyPhiFace::GetPointOnFace() |
---|
| 953 | { |
---|
| 954 | Triangulate(); |
---|
| 955 | return surface_point; |
---|
| 956 | } |
---|
| 957 | |
---|
| 958 | // |
---|
| 959 | // Auxiliary Functions used for Finding the PointOnFace using Triangulation |
---|
| 960 | // |
---|
| 961 | |
---|
| 962 | // |
---|
| 963 | // Calculation of 2*Area of Triangle with Sign |
---|
| 964 | // |
---|
| 965 | G4double G4PolyPhiFace::Area2( G4TwoVector a, |
---|
| 966 | G4TwoVector b, |
---|
| 967 | G4TwoVector c ) |
---|
| 968 | { |
---|
| 969 | return ((b.x()-a.x())*(c.y()-a.y())- |
---|
| 970 | (c.x()-a.x())*(b.y()-a.y())); |
---|
| 971 | } |
---|
| 972 | |
---|
| 973 | // |
---|
| 974 | // Boolean function for sign of Surface |
---|
| 975 | // |
---|
| 976 | G4bool G4PolyPhiFace::Left( G4TwoVector a, |
---|
| 977 | G4TwoVector b, |
---|
| 978 | G4TwoVector c ) |
---|
| 979 | { |
---|
| 980 | return Area2(a,b,c)>0; |
---|
| 981 | } |
---|
| 982 | |
---|
| 983 | // |
---|
| 984 | // Boolean function for sign of Surface |
---|
| 985 | // |
---|
| 986 | G4bool G4PolyPhiFace::LeftOn( G4TwoVector a, |
---|
| 987 | G4TwoVector b, |
---|
| 988 | G4TwoVector c ) |
---|
| 989 | { |
---|
| 990 | return Area2(a,b,c)>=0; |
---|
| 991 | } |
---|
| 992 | |
---|
| 993 | // |
---|
| 994 | // Boolean function for sign of Surface |
---|
| 995 | // |
---|
| 996 | G4bool G4PolyPhiFace::Collinear( G4TwoVector a, |
---|
| 997 | G4TwoVector b, |
---|
| 998 | G4TwoVector c ) |
---|
| 999 | { |
---|
| 1000 | return Area2(a,b,c)==0; |
---|
| 1001 | } |
---|
| 1002 | |
---|
| 1003 | // |
---|
| 1004 | // Boolean function for finding "Proper" Intersection |
---|
| 1005 | // That means Intersection of two lines segments (a,b) and (c,d) |
---|
| 1006 | // |
---|
| 1007 | G4bool G4PolyPhiFace::IntersectProp( G4TwoVector a, |
---|
| 1008 | G4TwoVector b, |
---|
| 1009 | G4TwoVector c, G4TwoVector d ) |
---|
| 1010 | { |
---|
| 1011 | if( Collinear(a,b,c) || Collinear(a,b,d)|| |
---|
| 1012 | Collinear(c,d,a) || Collinear(c,d,b) ) { return false; } |
---|
| 1013 | |
---|
| 1014 | G4bool Positive; |
---|
| 1015 | Positive = !(Left(a,b,c))^!(Left(a,b,d)); |
---|
| 1016 | return Positive && (!Left(c,d,a)^!Left(c,d,b)); |
---|
| 1017 | } |
---|
| 1018 | |
---|
| 1019 | // |
---|
| 1020 | // Boolean function for determining if Point c is between a and b |
---|
| 1021 | // For the tree points(a,b,c) on the same line |
---|
| 1022 | // |
---|
| 1023 | G4bool G4PolyPhiFace::Between( G4TwoVector a, G4TwoVector b, G4TwoVector c ) |
---|
| 1024 | { |
---|
| 1025 | if( !Collinear(a,b,c) ) { return false; } |
---|
| 1026 | |
---|
| 1027 | if(a.x()!=b.x()) |
---|
| 1028 | { |
---|
| 1029 | return ((a.x()<=c.x())&&(c.x()<=b.x()))|| |
---|
| 1030 | ((a.x()>=c.x())&&(c.x()>=b.x())); |
---|
| 1031 | } |
---|
| 1032 | else |
---|
| 1033 | { |
---|
| 1034 | return ((a.y()<=c.y())&&(c.y()<=b.y()))|| |
---|
| 1035 | ((a.y()>=c.y())&&(c.y()>=b.y())); |
---|
| 1036 | } |
---|
| 1037 | } |
---|
| 1038 | |
---|
| 1039 | // |
---|
| 1040 | // Boolean function for finding Intersection "Proper" or not |
---|
| 1041 | // Between two line segments (a,b) and (c,d) |
---|
| 1042 | // |
---|
| 1043 | G4bool G4PolyPhiFace::Intersect( G4TwoVector a, |
---|
| 1044 | G4TwoVector b, |
---|
| 1045 | G4TwoVector c, G4TwoVector d ) |
---|
| 1046 | { |
---|
| 1047 | if( IntersectProp(a,b,c,d) ) |
---|
| 1048 | { return true; } |
---|
| 1049 | else if( Between(a,b,c)|| |
---|
| 1050 | Between(a,b,d)|| |
---|
| 1051 | Between(c,d,a)|| |
---|
| 1052 | Between(c,d,b) ) |
---|
| 1053 | { return true; } |
---|
| 1054 | else |
---|
| 1055 | { return false; } |
---|
| 1056 | } |
---|
| 1057 | |
---|
| 1058 | // |
---|
| 1059 | // Boolean Diagonalie help to determine |
---|
| 1060 | // if diagonal s of segment (a,b) is convex or reflex |
---|
| 1061 | // |
---|
| 1062 | G4bool G4PolyPhiFace::Diagonalie( G4PolyPhiFaceVertex *a, |
---|
| 1063 | G4PolyPhiFaceVertex *b ) |
---|
| 1064 | { |
---|
| 1065 | G4PolyPhiFaceVertex *corner = triangles; |
---|
| 1066 | G4PolyPhiFaceVertex *corner_next=triangles; |
---|
| 1067 | |
---|
| 1068 | // For each Edge (corner,corner_next) |
---|
| 1069 | do |
---|
| 1070 | { |
---|
| 1071 | corner_next=corner->next; |
---|
| 1072 | |
---|
| 1073 | // Skip edges incident to a of b |
---|
| 1074 | // |
---|
| 1075 | if( (corner!=a)&&(corner_next!=a) |
---|
| 1076 | &&(corner!=b)&&(corner_next!=b) ) |
---|
| 1077 | { |
---|
| 1078 | G4TwoVector rz1,rz2,rz3,rz4; |
---|
| 1079 | rz1 = G4TwoVector(a->r,a->z); |
---|
| 1080 | rz2 = G4TwoVector(b->r,b->z); |
---|
| 1081 | rz3 = G4TwoVector(corner->r,corner->z); |
---|
| 1082 | rz4 = G4TwoVector(corner_next->r,corner_next->z); |
---|
| 1083 | if( Intersect(rz1,rz2,rz3,rz4) ) { return false; } |
---|
| 1084 | } |
---|
| 1085 | corner=corner->next; |
---|
| 1086 | |
---|
| 1087 | } while( corner != triangles ); |
---|
| 1088 | |
---|
| 1089 | return true; |
---|
| 1090 | } |
---|
| 1091 | |
---|
| 1092 | // |
---|
| 1093 | // Boolean function that determine if b is Inside Cone (a0,a,a1) |
---|
| 1094 | // being a the center of the Cone |
---|
| 1095 | // |
---|
| 1096 | G4bool G4PolyPhiFace::InCone( G4PolyPhiFaceVertex *a, G4PolyPhiFaceVertex *b ) |
---|
| 1097 | { |
---|
| 1098 | // a0,a and a1 are consecutive vertices |
---|
| 1099 | // |
---|
| 1100 | G4PolyPhiFaceVertex *a0,*a1; |
---|
| 1101 | a1=a->next; |
---|
| 1102 | a0=a->prev; |
---|
| 1103 | |
---|
| 1104 | G4TwoVector arz,arz0,arz1,brz; |
---|
| 1105 | arz=G4TwoVector(a->r,a->z);arz0=G4TwoVector(a0->r,a0->z); |
---|
| 1106 | arz1=G4TwoVector(a1->r,a1->z);brz=G4TwoVector(b->r,b->z); |
---|
| 1107 | |
---|
| 1108 | |
---|
| 1109 | if(LeftOn(arz,arz1,arz0)) // If a is convex vertex |
---|
| 1110 | { |
---|
| 1111 | return Left(arz,brz,arz0)&&Left(brz,arz,arz1); |
---|
| 1112 | } |
---|
| 1113 | else // Else a is reflex |
---|
| 1114 | { |
---|
| 1115 | return !( LeftOn(arz,brz,arz1)&&LeftOn(brz,arz,arz0)); |
---|
| 1116 | } |
---|
| 1117 | } |
---|
| 1118 | |
---|
| 1119 | // |
---|
| 1120 | // Boolean function finding if Diagonal is possible |
---|
| 1121 | // inside Polycone or PolyHedra |
---|
| 1122 | // |
---|
| 1123 | G4bool G4PolyPhiFace::Diagonal( G4PolyPhiFaceVertex *a, G4PolyPhiFaceVertex *b ) |
---|
| 1124 | { |
---|
| 1125 | return InCone(a,b) && InCone(b,a) && Diagonalie(a,b); |
---|
| 1126 | } |
---|
| 1127 | |
---|
| 1128 | // |
---|
| 1129 | // Initialisation for Triangulisation by ear tips |
---|
| 1130 | // For details see "Computational Geometry in C" by Joseph O'Rourke |
---|
| 1131 | // |
---|
| 1132 | void G4PolyPhiFace::EarInit() |
---|
| 1133 | { |
---|
| 1134 | G4PolyPhiFaceVertex *corner = triangles; |
---|
| 1135 | G4PolyPhiFaceVertex *c_prev,*c_next; |
---|
| 1136 | |
---|
| 1137 | do |
---|
| 1138 | { |
---|
| 1139 | // We need to determine three consecutive vertices |
---|
| 1140 | // |
---|
| 1141 | c_next=corner->next; |
---|
| 1142 | c_prev=corner->prev; |
---|
| 1143 | |
---|
| 1144 | // Calculation of ears |
---|
| 1145 | // |
---|
| 1146 | corner->ear=Diagonal(c_prev,c_next); |
---|
| 1147 | corner=corner->next; |
---|
| 1148 | |
---|
| 1149 | } while( corner!=triangles ); |
---|
| 1150 | } |
---|
| 1151 | |
---|
| 1152 | // |
---|
| 1153 | // Triangulisation by ear tips for Polycone or Polyhedra |
---|
| 1154 | // For details see "Computational Geometry in C" by Joseph O'Rourke |
---|
| 1155 | // |
---|
| 1156 | void G4PolyPhiFace::Triangulate() |
---|
| 1157 | { |
---|
| 1158 | // The copy of Polycone is made and this copy is reordered in order to |
---|
| 1159 | // have a list of triangles. This list is used for GetPointOnFace(). |
---|
| 1160 | |
---|
| 1161 | G4PolyPhiFaceVertex *tri_help = new G4PolyPhiFaceVertex[numEdges]; |
---|
| 1162 | triangles = tri_help; |
---|
| 1163 | G4PolyPhiFaceVertex *triang = triangles; |
---|
| 1164 | |
---|
| 1165 | std::vector<G4double> areas; |
---|
| 1166 | std::vector<G4ThreeVector> points; |
---|
| 1167 | G4double area=0.; |
---|
| 1168 | G4PolyPhiFaceVertex *v0,*v1,*v2,*v3,*v4; |
---|
| 1169 | v2=triangles; |
---|
| 1170 | |
---|
| 1171 | // Make copy for prev/next for triang=corners |
---|
| 1172 | // |
---|
| 1173 | G4PolyPhiFaceVertex *helper = corners; |
---|
| 1174 | G4PolyPhiFaceVertex *helper2 = corners; |
---|
| 1175 | do |
---|
| 1176 | { |
---|
| 1177 | triang->r = helper->r; |
---|
| 1178 | triang->z = helper->z; |
---|
| 1179 | triang->x = helper->x; |
---|
| 1180 | triang->y= helper->y; |
---|
| 1181 | |
---|
| 1182 | // add pointer on prev corner |
---|
| 1183 | // |
---|
| 1184 | if( helper==corners ) |
---|
| 1185 | { triang->prev=triangles+numEdges-1; } |
---|
| 1186 | else |
---|
| 1187 | { triang->prev=helper2; } |
---|
| 1188 | |
---|
| 1189 | // add pointer on next corner |
---|
| 1190 | // |
---|
| 1191 | if( helper<corners+numEdges-1 ) |
---|
| 1192 | { triang->next=triang+1; } |
---|
| 1193 | else |
---|
| 1194 | { triang->next=triangles; } |
---|
| 1195 | helper2=triang; |
---|
| 1196 | helper=helper->next; |
---|
| 1197 | triang=triang->next; |
---|
| 1198 | |
---|
| 1199 | } while( helper!=corners ); |
---|
| 1200 | |
---|
| 1201 | EarInit(); |
---|
| 1202 | |
---|
| 1203 | G4int n=numEdges; |
---|
| 1204 | G4int i=0; |
---|
| 1205 | G4ThreeVector p1,p2,p3,p4; |
---|
| 1206 | const G4int max_n_loops=numEdges*10000; // protection against infinite loop |
---|
| 1207 | |
---|
| 1208 | // Each step of outer loop removes one ear |
---|
| 1209 | // |
---|
| 1210 | while(n>3) // Inner loop searches for one ear |
---|
| 1211 | { |
---|
| 1212 | v2=triangles; |
---|
| 1213 | do |
---|
| 1214 | { |
---|
| 1215 | if(v2->ear) // Ear found. Fill variables |
---|
| 1216 | { |
---|
| 1217 | // (v1,v3) is diagonal |
---|
| 1218 | // |
---|
| 1219 | v3=v2->next; v4=v3->next; |
---|
| 1220 | v1=v2->prev; v0=v1->prev; |
---|
| 1221 | |
---|
| 1222 | // Calculate areas and points |
---|
| 1223 | |
---|
| 1224 | p1=G4ThreeVector((v2)->x,(v2)->y,(v2)->z); |
---|
| 1225 | p2=G4ThreeVector((v1)->x,(v1)->y,(v1)->z); |
---|
| 1226 | p3=G4ThreeVector((v3)->x,(v3)->y,(v3)->z); |
---|
| 1227 | |
---|
| 1228 | G4double result1 = SurfaceTriangle(p1,p2,p3,&p4 ); |
---|
| 1229 | points.push_back(p4); |
---|
| 1230 | areas.push_back(result1); |
---|
| 1231 | area=area+result1; |
---|
| 1232 | |
---|
| 1233 | // Update earity of diagonal endpoints |
---|
| 1234 | // |
---|
| 1235 | v1->ear=Diagonal(v0,v3); |
---|
| 1236 | v3->ear=Diagonal(v1,v4); |
---|
| 1237 | |
---|
| 1238 | // Cut off the ear v2 |
---|
| 1239 | // Has to be done for a copy and not for real PolyPhiFace |
---|
| 1240 | // |
---|
| 1241 | v1->next=v3; |
---|
| 1242 | v3->prev=v1; |
---|
| 1243 | triangles=v3; // In case the head was v2 |
---|
| 1244 | n--; |
---|
| 1245 | |
---|
| 1246 | break; // out of inner loop |
---|
| 1247 | } // end if ear found |
---|
| 1248 | |
---|
| 1249 | v2=v2->next; |
---|
| 1250 | |
---|
| 1251 | } while( v2!=triangles ); |
---|
| 1252 | |
---|
| 1253 | i++; |
---|
| 1254 | if(i>=max_n_loops) |
---|
| 1255 | { |
---|
| 1256 | G4Exception( "G4PolyPhiFace::Triangulation()", |
---|
| 1257 | "Bad_Definition_of_Solid", FatalException, |
---|
| 1258 | "Maximum number of steps is reached for triangulation!" ); |
---|
| 1259 | } |
---|
| 1260 | } // end outer while loop |
---|
| 1261 | |
---|
| 1262 | if(v2->next) |
---|
| 1263 | { |
---|
| 1264 | // add last triangle |
---|
| 1265 | // |
---|
| 1266 | v2=v2->next; |
---|
| 1267 | p1=G4ThreeVector((v2)->x,(v2)->y,(v2)->z); |
---|
| 1268 | p2=G4ThreeVector((v2->next)->x,(v2->next)->y,(v2->next)->z); |
---|
| 1269 | p3=G4ThreeVector((v2->prev)->x,(v2->prev)->y,(v2->prev)->z); |
---|
| 1270 | G4double result1 = SurfaceTriangle(p1,p2,p3,&p4 ); |
---|
| 1271 | points.push_back(p4); |
---|
| 1272 | areas.push_back(result1); |
---|
| 1273 | area=area+result1; |
---|
| 1274 | } |
---|
| 1275 | |
---|
| 1276 | // Surface Area is stored |
---|
| 1277 | // |
---|
| 1278 | fSurfaceArea = area; |
---|
| 1279 | |
---|
| 1280 | // Second Step: choose randomly one surface |
---|
| 1281 | // |
---|
| 1282 | G4double chose = area*G4UniformRand(); |
---|
| 1283 | |
---|
| 1284 | // Third Step: Get a point on choosen surface |
---|
| 1285 | // |
---|
| 1286 | G4double Achose1, Achose2; |
---|
| 1287 | Achose1=0; Achose2=0.; |
---|
| 1288 | i=0; |
---|
| 1289 | do |
---|
| 1290 | { |
---|
| 1291 | Achose2+=areas[i]; |
---|
| 1292 | if(chose>=Achose1 && chose<Achose2) |
---|
| 1293 | { |
---|
| 1294 | G4ThreeVector point; |
---|
| 1295 | point=points[i] ; |
---|
| 1296 | surface_point=point; |
---|
| 1297 | break; |
---|
| 1298 | } |
---|
| 1299 | i++; Achose1=Achose2; |
---|
| 1300 | } while( i<numEdges-2 ); |
---|
| 1301 | |
---|
| 1302 | delete [] tri_help; |
---|
| 1303 | } |
---|