Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 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 // G4UAdapter 27 // 28 // Class description: 29 // 30 // Utility class for adapting VecGeom solids API to Geant4 solids. 31 // NOTE: Using protected inheritance since the Adapter is supposed to 32 // be a G4VSolid "implemented-in-terms-of" the VecGeom UnplacedVolume_t. 33 // The choice of protected vs private is due to the fact that we want 34 // to propagate functions further down in the inheritance hierarchy. 35 36 // Author: 37 // 17.05.17 G.Cosmo: Adapted for G4VSolid from original G4USolids bridge 38 // class and the USolidsAdapter class in VecGeom. 39 // ------------------------------------------------------------------------ 40 #ifndef G4UADAPTER_HH 41 #define G4UADAPTER_HH 42 43 #include "G4ThreeVector.hh" 44 #include "G4VSolid.hh" 45 46 // Required for inline visualization adapter functions 47 // 48 #include "G4AffineTransform.hh" 49 #include "G4VoxelLimits.hh" 50 #include "G4VGraphicsScene.hh" 51 #include "G4Polyhedron.hh" 52 #include "G4VisExtent.hh" 53 #include "G4BoundingEnvelope.hh" 54 #include "G4AutoLock.hh" 55 56 #include "G4GeomTypes.hh" 57 58 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) ) 59 60 #include <VecGeom/base/Global.h> 61 #include <VecGeom/base/Vector3D.h> 62 63 class G4VPVParameterisation; 64 65 template <class UnplacedVolume_t> 66 class G4UAdapter : public G4VSolid, protected UnplacedVolume_t 67 { 68 public: 69 70 using U3Vector = vecgeom::Vector3D<G4double>; 71 72 using UnplacedVolume_t::operator delete; 73 using UnplacedVolume_t::operator new; 74 // VecGeom volumes have special delete/new ("AlignedBase") 75 // and we need to make these functions public again 76 77 G4UAdapter(const G4String& name) 78 : G4VSolid(name) 79 { kHalfTolerance = 0.5*kCarTolerance; } 80 81 template <typename... T> 82 G4UAdapter(const G4String& name, const T &... params) 83 : G4VSolid(name), UnplacedVolume_t(params...) 84 { kHalfTolerance = 0.5*kCarTolerance; } 85 86 virtual ~G4UAdapter(); 87 88 G4bool operator==(const G4UAdapter& s) const; 89 // Return true only if addresses are the same. 90 91 virtual G4bool CalculateExtent(const EAxis pAxis, 92 const G4VoxelLimits& pVoxelLimit, 93 const G4AffineTransform& pTransform, 94 G4double& pMin, G4double& pMax) const override; 95 // Calculate the minimum and maximum extent of the solid, when under the 96 // specified transform, and within the specified limits. If the solid 97 // is not intersected by the region, return false, else return true. 98 99 virtual EInside Inside(const G4ThreeVector& p) const override; 100 // Returns kOutside if the point at offset p is outside the shapes 101 // boundaries plus Tolerance/2, kSurface if the point is <= Tolerance/2 102 // from a surface, otherwise kInside. 103 104 virtual G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const override; 105 // Returns the outwards pointing unit normal of the shape for the 106 // surface closest to the point at offset p. 107 108 virtual G4double DistanceToIn(const G4ThreeVector& p, 109 const G4ThreeVector& v) const override; 110 // Return the distance along the normalised vector v to the shape, 111 // from the point at offset p. If there is no intersection, return 112 // kInfinity. The first intersection resulting from `leaving' a 113 // surface/volume is discarded. Hence, it is tolerant of points on 114 // the surface of the shape. 115 116 virtual G4double DistanceToIn(const G4ThreeVector& p) const override; 117 // Calculate the distance to the nearest surface of a shape from an 118 // outside point. The distance can be an underestimate. 119 120 virtual G4double DistanceToOut(const G4ThreeVector& p, 121 const G4ThreeVector& v, 122 const G4bool calcNorm = false, 123 G4bool* validNorm = 0, 124 G4ThreeVector* n = 0) const override; 125 // Return the distance along the normalised vector v to the shape, 126 // from a point at an offset p inside or on the surface of the shape. 127 // Intersections with surfaces, when the point is < Tolerance/2 from a 128 // surface must be ignored. 129 // If calcNorm==true: 130 // validNorm set true if the solid lies entirely behind or on the 131 // exiting surface. 132 // n set to exiting outwards normal vector (undefined Magnitude). 133 // validNorm set to false if the solid does not lie entirely behind 134 // or on the exiting surface 135 // If calcNorm==false: 136 // validNorm and n are unused. 137 // 138 // Must be called as solid.DistanceToOut(p,v) or by specifying all 139 // the parameters. 140 141 virtual G4double DistanceToOut(const G4ThreeVector& p) const override; 142 // Calculate the distance to the nearest surface of a shape from an 143 // inside point. The distance can be an underestimate. 144 145 virtual void ComputeDimensions(G4VPVParameterisation* p, 146 const G4int n, 147 const G4VPhysicalVolume* pRep) override; 148 // Throw exception if ComputeDimensions called from an illegal 149 // derived class. 150 151 virtual G4double GetCubicVolume() override; 152 // Returns an estimation of the solid volume in internal units. 153 // This method may be overloaded by derived classes to compute the 154 // exact geometrical quantity for solids where this is possible, 155 // or anyway to cache the computed value. 156 // Note: the computed value is NOT cached. 157 158 virtual G4double GetSurfaceArea() override; 159 // Return an estimation of the solid surface area in internal units. 160 // This method may be overloaded by derived classes to compute the 161 // exact geometrical quantity for solids where this is possible, 162 // or anyway to cache the computed value. 163 // Note: the computed value is NOT cached. 164 165 virtual G4ThreeVector GetPointOnSurface() const override; 166 // Returns a random point located on the surface of the solid. 167 168 virtual G4int GetNumOfConstituents() const override; 169 // Returns the number of constituents used for construction of the solid. 170 // For non-Boolean solids the return value is one. 171 172 virtual G4bool IsFaceted() const override; 173 // Returns true if the solid has only planar faces, false otherwise. 174 175 virtual G4GeometryType GetEntityType() const override; 176 // Provide identification of the class of an object. 177 // (required for persistency) 178 179 virtual G4VSolid* Clone() const override; 180 // Returns a pointer of a dynamically allocated copy of the solid. 181 // Returns NULL pointer with warning in case the concrete solid does not 182 // implement this method. The caller has responsibility for ownership. 183 184 virtual std::ostream& StreamInfo(std::ostream& os) const override; 185 // Dumps contents of the solid to a stream. 186 187 virtual void DescribeYourselfTo(G4VGraphicsScene& scene) const override; 188 // A "double dispatch" function which identifies the solid 189 // to the graphics scene for visualization. 190 191 virtual G4VisExtent GetExtent() const override; 192 // Provide extent (bounding box) as possible hint to the graphics view. 193 virtual G4Polyhedron* CreatePolyhedron() const override; 194 // Create Polyhedron used for Visualisation 195 virtual G4Polyhedron* GetPolyhedron() const override; 196 // Smart access function - creates on request and stores for future 197 // access. A null pointer means "not available". 198 199 public: // without description 200 201 G4UAdapter(__void__&); 202 // Fake default constructor for usage restricted to direct object 203 // persistency for clients requiring preallocation of memory for 204 // persistifiable objects. 205 206 G4UAdapter(const G4UAdapter& rhs); 207 G4UAdapter& operator=(const G4UAdapter& rhs); 208 // Copy constructor and assignment operator. 209 210 public: // VecGeom overridden methods 211 212 vecgeom::Precision 213 DistanceToOut(U3Vector const& position, U3Vector const& direction, 214 vecgeom::Precision stepMax = kInfinity) const override 215 { 216 return UnplacedVolume_t::DistanceToOut(position, direction, stepMax); 217 } 218 219 vecgeom::EnumInside 220 Inside(U3Vector const& aPoint) const override 221 { 222 return UnplacedVolume_t::Inside(aPoint); 223 } 224 225 vecgeom::Precision 226 DistanceToIn(U3Vector const& position, U3Vector const& direction, 227 const vecgeom::Precision step_max = kInfinity) const override 228 { 229 return UnplacedVolume_t::DistanceToIn(position, direction, step_max); 230 } 231 232 G4bool Normal(U3Vector const& aPoint, U3Vector& aNormal) const override 233 { 234 return UnplacedVolume_t::Normal(aPoint, aNormal); 235 } 236 237 void Extent(U3Vector& aMin, U3Vector& aMax) const override 238 { 239 return UnplacedVolume_t::Extent(aMin, aMax); 240 } 241 242 U3Vector SamplePointOnSurface() const override 243 { 244 return UnplacedVolume_t::SamplePointOnSurface(); 245 } 246 247 protected: // data 248 249 mutable G4bool fRebuildPolyhedron = false; 250 mutable G4Polyhedron* fPolyhedron = nullptr; 251 252 G4double kHalfTolerance; // Cached geometrical tolerance 253 254 using UnplacedVolume_t::DistanceToOut; 255 using UnplacedVolume_t::DistanceToIn; 256 }; 257 258 // Inline implementations 259 260 template <class UnplacedVolume_t> 261 G4UAdapter<UnplacedVolume_t>::G4UAdapter(__void__& a) 262 : G4VSolid(a), UnplacedVolume_t(*this), 263 kHalfTolerance(0.5*kCarTolerance) 264 { 265 } 266 267 template <class UnplacedVolume_t> 268 G4UAdapter<UnplacedVolume_t>::~G4UAdapter() 269 { 270 delete fPolyhedron; fPolyhedron = nullptr; 271 } 272 273 template <class UnplacedVolume_t> 274 G4bool G4UAdapter<UnplacedVolume_t>:: 275 operator==(const G4UAdapter& rhs) const 276 { 277 return (this == &rhs) ? true : false; 278 } 279 280 template <class UnplacedVolume_t> 281 G4UAdapter<UnplacedVolume_t>:: 282 G4UAdapter(const G4UAdapter& rhs) 283 : G4VSolid(rhs), UnplacedVolume_t(rhs) 284 { 285 kHalfTolerance = 0.5*kCarTolerance; 286 } 287 288 template <class UnplacedVolume_t> 289 G4UAdapter<UnplacedVolume_t>& G4UAdapter<UnplacedVolume_t>:: 290 operator=(const G4UAdapter& rhs) 291 { 292 // Check assignment to self 293 // 294 if (this == &rhs) 295 { 296 return *this; 297 } 298 299 // Copy base class data 300 // 301 G4VSolid::operator=(rhs); 302 UnplacedVolume_t::operator=(rhs); 303 304 // Copy data 305 // 306 fRebuildPolyhedron = false; 307 delete fPolyhedron; fPolyhedron = nullptr; 308 kHalfTolerance = 0.5*kCarTolerance; 309 310 return *this; 311 } 312 313 template <class UnplacedVolume_t> 314 EInside G4UAdapter<UnplacedVolume_t>:: 315 Inside(const G4ThreeVector& p) const 316 { 317 U3Vector pt(p.x(), p.y(), p.z()); 318 vecgeom::EnumInside in_temp; 319 EInside in = kOutside; 320 321 in_temp = UnplacedVolume_t::Inside(pt); 322 323 if (in_temp == vecgeom::EnumInside::eInside) in = kInside; 324 else if (in_temp == vecgeom::EnumInside::eSurface) in = kSurface; 325 326 return in; 327 } 328 329 template <class UnplacedVolume_t> 330 G4ThreeVector G4UAdapter<UnplacedVolume_t>:: 331 SurfaceNormal(const G4ThreeVector& pt) const 332 { 333 U3Vector p(pt.x(), pt.y(), pt.z()); 334 U3Vector n; 335 UnplacedVolume_t::Normal(p, n); 336 return G4ThreeVector(n.x(), n.y(), n.z()); 337 } 338 339 template <class UnplacedVolume_t> 340 G4double G4UAdapter<UnplacedVolume_t>:: 341 DistanceToIn(const G4ThreeVector& pt, const G4ThreeVector& d) const 342 { 343 U3Vector p(pt.x(), pt.y(), pt.z()); 344 U3Vector v(d.x(), d.y(), d.z()); 345 G4double dist = UnplacedVolume_t::DistanceToIn(p, v, kInfinity); 346 347 // apply Geant4 distance conventions 348 // 349 if (dist < kHalfTolerance) return 0.0; 350 return (dist > kInfinity) ? kInfinity : dist; 351 } 352 353 template <class UnplacedVolume_t> 354 G4double G4UAdapter<UnplacedVolume_t>:: 355 DistanceToIn(const G4ThreeVector& pt) const 356 { 357 U3Vector p(pt.x(), pt.y(), pt.z()); 358 G4double dist = UnplacedVolume_t::SafetyToIn(p); 359 360 // Apply Geant4 convention: convert negative values to zero 361 // 362 if (dist < kHalfTolerance) return 0.0; 363 return (dist > kInfinity) ? kInfinity : dist; 364 } 365 366 template <class UnplacedVolume_t> 367 G4double G4UAdapter<UnplacedVolume_t>:: 368 DistanceToOut(const G4ThreeVector& pt, const G4ThreeVector& d, 369 const G4bool calcNorm, G4bool* validNorm, 370 G4ThreeVector* norm) const 371 { 372 U3Vector p(pt.x(), pt.y(), pt.z()); 373 U3Vector v(d.x(), d.y(), d.z()); 374 375 G4double dist = UnplacedVolume_t::DistanceToOut(p, v, kInfinity); 376 if(calcNorm) 377 { 378 *validNorm = UnplacedVolume_t::IsConvex(); 379 U3Vector n, hitpoint = p + dist * v; 380 UnplacedVolume_t::Normal(hitpoint, n); 381 norm->set(n.x(), n.y(), n.z()); 382 } 383 384 // Apply Geant4 distance conventions 385 // 386 if (dist < kHalfTolerance) return 0.0; 387 return (dist > kInfinity) ? kInfinity : dist; 388 } 389 390 template <class UnplacedVolume_t> 391 G4double G4UAdapter<UnplacedVolume_t>:: 392 DistanceToOut(const G4ThreeVector& pt) const 393 { 394 U3Vector p(pt.x(), pt.y(), pt.z()); 395 G4double dist = UnplacedVolume_t::SafetyToOut(p); 396 397 // Apply Geant4 convention: convert negative values to zero 398 // 399 if (dist < kHalfTolerance) return 0.0; 400 return (dist > kInfinity) ? kInfinity : dist; 401 } 402 403 template <class UnplacedVolume_t> 404 G4double G4UAdapter<UnplacedVolume_t>::GetCubicVolume() 405 { 406 return UnplacedVolume_t::Capacity(); 407 } 408 409 template <class UnplacedVolume_t> 410 G4double G4UAdapter<UnplacedVolume_t>::GetSurfaceArea() 411 { 412 return UnplacedVolume_t::SurfaceArea(); 413 } 414 415 template <class UnplacedVolume_t> 416 G4ThreeVector G4UAdapter<UnplacedVolume_t>::GetPointOnSurface() const 417 { 418 U3Vector p = UnplacedVolume_t::SamplePointOnSurface(); 419 return G4ThreeVector(p.x(), p.y(), p.z()); 420 } 421 422 template <class UnplacedVolume_t> 423 G4int G4UAdapter<UnplacedVolume_t>::GetNumOfConstituents() const 424 { 425 return 1; 426 } 427 428 template <class UnplacedVolume_t> 429 G4bool G4UAdapter<UnplacedVolume_t>::IsFaceted() const 430 { 431 return false; 432 } 433 434 // Inline visualization adapters 435 436 namespace 437 { 438 G4Mutex pMutex = G4MUTEX_INITIALIZER; 439 } 440 441 // Free function to enable ostream output 442 template <class UnplacedVolume_t> 443 std::ostream& 444 operator<<(std::ostream& os, const G4UAdapter<UnplacedVolume_t>& uAdapted) 445 { 446 return uAdapted.StreamInfo(os); 447 } 448 449 template <class UnplacedVolume_t> 450 void G4UAdapter<UnplacedVolume_t>:: 451 ComputeDimensions(G4VPVParameterisation*, const G4int, 452 const G4VPhysicalVolume*) 453 { 454 std::ostringstream message; 455 message << "Illegal call to G4UAdapter::ComputeDimensions()" << G4endl 456 << "Method not overloaded by derived class !"; 457 G4Exception("G4UAdapter::ComputeDimensions()", "GeomSolids0003", 458 FatalException, message); 459 } 460 461 template <class UnplacedVolume_t> 462 void G4UAdapter<UnplacedVolume_t>:: 463 DescribeYourselfTo(G4VGraphicsScene& scene) const 464 { 465 scene.AddSolid(*this); 466 } 467 468 template <class UnplacedVolume_t> 469 G4GeometryType G4UAdapter<UnplacedVolume_t>:: 470 GetEntityType() const 471 { 472 473 G4String string = "VSolid"; // UnplacedVolume_t::GetEntityType(); 474 return "G4" + string; 475 } 476 477 template <class UnplacedVolume_t> 478 std::ostream& G4UAdapter<UnplacedVolume_t>:: 479 StreamInfo(std::ostream& os) const 480 { 481 UnplacedVolume_t::Print(os); 482 return os; 483 } 484 485 template <class UnplacedVolume_t> 486 G4VSolid* G4UAdapter<UnplacedVolume_t>::Clone() const 487 { 488 std::ostringstream message; 489 message << "Clone() method not implemented for type: " 490 << GetEntityType() << "!" << G4endl 491 << "Returning NULL pointer!"; 492 G4Exception("G4UAdapter::Clone()", "GeomSolids1001", JustWarning, message); 493 return nullptr; 494 } 495 496 template <class UnplacedVolume_t> 497 G4bool G4UAdapter<UnplacedVolume_t>::CalculateExtent(const EAxis pAxis, 498 const G4VoxelLimits& pVoxelLimit, 499 const G4AffineTransform& pTransform, 500 G4double& pMin, G4double& pMax) const 501 { 502 U3Vector vmin, vmax; 503 UnplacedVolume_t::Extent(vmin,vmax); 504 G4ThreeVector bmin(vmin.x(),vmin.y(),vmin.z()); 505 G4ThreeVector bmax(vmax.x(),vmax.y(),vmax.z()); 506 507 // Check correctness of the bounding box 508 // 509 if (bmin.x() >= bmax.x() || bmin.y() >= bmax.y() || bmin.z() >= bmax.z()) 510 { 511 std::ostringstream message; 512 message << "Bad bounding box (min >= max) for solid: " 513 << GetName() << " - " << GetEntityType() << " !" 514 << "\nmin = " << bmin 515 << "\nmax = " << bmax; 516 G4Exception("G4UAdapter::CalculateExtent()", "GeomMgt0001", 517 JustWarning, message); 518 StreamInfo(G4cout); 519 } 520 521 G4BoundingEnvelope bbox(bmin,bmax); 522 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax); 523 } 524 525 template <class UnplacedVolume_t> 526 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::CreatePolyhedron() const 527 { 528 // Must be implemented in concrete wrappers... 529 530 std::ostringstream message; 531 message << "Visualization not supported for USolid shape " 532 << GetEntityType() << "... Sorry!" << G4endl; 533 G4Exception("G4UAdapter::CreatePolyhedron()", "GeomSolids0003", 534 FatalException, message); 535 return nullptr; 536 } 537 538 template <class UnplacedVolume_t> 539 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::GetPolyhedron() const 540 { 541 if (!fPolyhedron || 542 fRebuildPolyhedron || 543 fPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() != 544 fPolyhedron->GetNumberOfRotationSteps()) 545 { 546 G4AutoLock l(&pMutex); 547 delete fPolyhedron; 548 fPolyhedron = CreatePolyhedron(); 549 fRebuildPolyhedron = false; 550 l.unlock(); 551 } 552 return fPolyhedron; 553 } 554 555 template <class UnplacedVolume_t> 556 G4VisExtent G4UAdapter<UnplacedVolume_t>::GetExtent() const 557 { 558 U3Vector vmin, vmax; 559 UnplacedVolume_t::Extent(vmin,vmax); 560 return G4VisExtent(vmin.x(),vmax.x(), 561 vmin.y(),vmax.y(), 562 vmin.z(),vmax.z()); 563 } 564 565 #endif // G4GEOM_USE_USOLIDS 566 567 #endif // G4UADAPTER_HH 568