ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UAdapter.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4UAdapter.hh
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 <base/Global.h>
61 #include <base/Vector3D.h>
62 
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 
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 G4GeometryType GetEntityType() const override;
169  // Provide identification of the class of an object.
170  // (required for persistency)
171 
172  virtual G4VSolid* Clone() const override;
173  // Returns a pointer of a dynamically allocated copy of the solid.
174  // Returns NULL pointer with warning in case the concrete solid does not
175  // implement this method. The caller has responsibility for ownership.
176 
177  virtual std::ostream& StreamInfo(std::ostream& os) const override;
178  // Dumps contents of the solid to a stream.
179 
180  virtual void DescribeYourselfTo(G4VGraphicsScene& scene) const override;
181  // A "double dispatch" function which identifies the solid
182  // to the graphics scene for visualization.
183 
184  virtual G4VisExtent GetExtent() const override;
185  // Provide extent (bounding box) as possible hint to the graphics view.
186  virtual G4Polyhedron* CreatePolyhedron() const override;
187  // Create Polyhedron used for Visualisation
188  virtual G4Polyhedron* GetPolyhedron() const override;
189  // Smart access function - creates on request and stores for future
190  // access. A null pointer means "not available".
191 
192  public: // without description
193 
194  G4UAdapter(__void__&);
195  // Fake default constructor for usage restricted to direct object
196  // persistency for clients requiring preallocation of memory for
197  // persistifiable objects.
198 
199  G4UAdapter(const G4UAdapter& rhs);
200  G4UAdapter& operator=(const G4UAdapter& rhs);
201  // Copy constructor and assignment operator.
202 
203  public: // VecGeom overridden methods
204 
205  vecgeom::Precision
206  DistanceToOut(U3Vector const& position, U3Vector const& direction,
207  vecgeom::Precision stepMax = kInfinity) const override
208  {
209  return UnplacedVolume_t::DistanceToOut(position, direction, stepMax);
210  }
211 
212  vecgeom::EnumInside
213  Inside(U3Vector const& aPoint) const override
214  {
215  return UnplacedVolume_t::Inside(aPoint);
216  }
217 
218  vecgeom::Precision
219  DistanceToIn(U3Vector const& position, U3Vector const& direction,
220  const vecgeom::Precision step_max = kInfinity) const override
221  {
222  return UnplacedVolume_t::DistanceToIn(position, direction, step_max);
223  }
224 
225  G4bool Normal(U3Vector const& aPoint, U3Vector& aNormal) const override
226  {
227  return UnplacedVolume_t::Normal(aPoint, aNormal);
228  }
229 
230  void Extent(U3Vector& aMin, U3Vector& aMax) const override
231  {
232  return UnplacedVolume_t::Extent(aMin, aMax);
233  }
234 
235  U3Vector SamplePointOnSurface() const override
236  {
237  return UnplacedVolume_t::SamplePointOnSurface();
238  }
239 
240  protected: // data
241 
242  mutable G4bool fRebuildPolyhedron = false;
243  mutable G4Polyhedron* fPolyhedron = nullptr;
244 
245  G4double kHalfTolerance; // Cached geometrical tolerance
246 
247  using UnplacedVolume_t::DistanceToOut;
248  using UnplacedVolume_t::DistanceToIn;
249 };
250 
251 // Inline implementations
252 
253 template <class UnplacedVolume_t>
254 G4UAdapter<UnplacedVolume_t>::G4UAdapter(__void__& a)
255  : G4VSolid(a), UnplacedVolume_t(*this),
256  kHalfTolerance(0.5*kCarTolerance)
257 {
258 }
259 
260 template <class UnplacedVolume_t>
261 G4UAdapter<UnplacedVolume_t>::~G4UAdapter()
262 {
263  delete fPolyhedron; fPolyhedron = nullptr;
264 }
265 
266 template <class UnplacedVolume_t>
268 operator==(const G4UAdapter& rhs) const
269 {
270  return (this == &rhs) ? true : false;
271 }
272 
273 template <class UnplacedVolume_t>
274 G4UAdapter<UnplacedVolume_t>::
275 G4UAdapter(const G4UAdapter& rhs)
276  : G4VSolid(rhs), UnplacedVolume_t(rhs)
277 {
278  kHalfTolerance = 0.5*kCarTolerance;
279 }
280 
281 template <class UnplacedVolume_t>
282 G4UAdapter<UnplacedVolume_t>& G4UAdapter<UnplacedVolume_t>::
283 operator=(const G4UAdapter& rhs)
284 {
285  // Check assignment to self
286  //
287  if (this == &rhs)
288  {
289  return *this;
290  }
291 
292  // Copy base class data
293  //
294  G4VSolid::operator=(rhs);
295  UnplacedVolume_t::operator=(rhs);
296 
297  // Copy data
298  //
299  fRebuildPolyhedron = false;
300  delete fPolyhedron; fPolyhedron = nullptr;
301  kHalfTolerance = 0.5*kCarTolerance;
302 
303  return *this;
304 }
305 
306 template <class UnplacedVolume_t>
307 EInside G4UAdapter<UnplacedVolume_t>::
308 Inside(const G4ThreeVector& p) const
309 {
310  U3Vector pt(p.x(), p.y(), p.z());
311  vecgeom::EnumInside in_temp;
312  EInside in = kOutside;
313 
314  in_temp = UnplacedVolume_t::Inside(pt);
315 
316  if (in_temp == vecgeom::EnumInside::eInside) in = kInside;
317  else if (in_temp == vecgeom::EnumInside::eSurface) in = kSurface;
318 
319  return in;
320 }
321 
322 template <class UnplacedVolume_t>
323 G4ThreeVector G4UAdapter<UnplacedVolume_t>::
324 SurfaceNormal(const G4ThreeVector& pt) const
325 {
326  U3Vector p(pt.x(), pt.y(), pt.z());
327  U3Vector n;
328  UnplacedVolume_t::Normal(p, n);
329  return G4ThreeVector(n.x(), n.y(), n.z());
330 }
331 
332 template <class UnplacedVolume_t>
333 G4double G4UAdapter<UnplacedVolume_t>::
334 DistanceToIn(const G4ThreeVector& pt, const G4ThreeVector& d) const
335 {
336  U3Vector p(pt.x(), pt.y(), pt.z());
337  U3Vector v(d.x(), d.y(), d.z());
338  G4double dist = UnplacedVolume_t::DistanceToIn(p, v, kInfinity);
339 
340  // apply Geant4 distance conventions
341  //
342  if (dist < kHalfTolerance) return 0.0;
343  return (dist > kInfinity) ? kInfinity : dist;
344 }
345 
346 template <class UnplacedVolume_t>
347 G4double G4UAdapter<UnplacedVolume_t>::
348 DistanceToIn(const G4ThreeVector& pt) const
349 {
350  U3Vector p(pt.x(), pt.y(), pt.z());
351  G4double dist = UnplacedVolume_t::SafetyToIn(p);
352 
353  // Apply Geant4 convention: convert negative values to zero
354  //
355  if (dist < kHalfTolerance) return 0.0;
356  return (dist > kInfinity) ? kInfinity : dist;
357 }
358 
359 template <class UnplacedVolume_t>
360 G4double G4UAdapter<UnplacedVolume_t>::
361 DistanceToOut(const G4ThreeVector& pt, const G4ThreeVector& d,
362  const G4bool calcNorm, G4bool* validNorm,
363  G4ThreeVector* norm) const
364 {
365  U3Vector p(pt.x(), pt.y(), pt.z());
366  U3Vector v(d.x(), d.y(), d.z());
367 
368  G4double dist = UnplacedVolume_t::DistanceToOut(p, v, kInfinity);
369  if(calcNorm) // *norm=n, but only after calcNorm check and if convex volume
370  {
371  if (UnplacedVolume_t::IsConvex())
372  {
373  U3Vector n, hitpoint = p + dist * v;
374  UnplacedVolume_t::Normal(hitpoint, n);
375  *validNorm = true;
376  norm->set(n.x(), n.y(), n.z());
377  }
378  else
379  {
380  *validNorm = false;
381  }
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 // Inline visualization adapters
423 
424 namespace
425 {
427 }
428 
429 template <class UnplacedVolume_t>
430 void G4UAdapter<UnplacedVolume_t>::
431 ComputeDimensions(G4VPVParameterisation*, const G4int,
432  const G4VPhysicalVolume*)
433 {
434  std::ostringstream message;
435  message << "Illegal call to G4UAdapter::ComputeDimensions()" << G4endl
436  << "Method not overloaded by derived class !";
437  G4Exception("G4UAdapter::ComputeDimensions()", "GeomSolids0003",
438  FatalException, message);
439 }
440 
441 template <class UnplacedVolume_t>
442 void G4UAdapter<UnplacedVolume_t>::
443 DescribeYourselfTo(G4VGraphicsScene& scene) const
444 {
445  scene.AddSolid(*this);
446 }
447 
448 template <class UnplacedVolume_t>
449 G4GeometryType G4UAdapter<UnplacedVolume_t>::
450 GetEntityType() const
451 {
452 
453  G4String string = "VSolid"; // UnplacedVolume_t::GetEntityType();
454  return "G4" + string;
455 }
456 
457 template <class UnplacedVolume_t>
458 std::ostream& G4UAdapter<UnplacedVolume_t>::
459 StreamInfo(std::ostream& os) const
460 {
462  return os;
463 }
464 
465 template <class UnplacedVolume_t>
466 G4VSolid* G4UAdapter<UnplacedVolume_t>::Clone() const
467 {
468  std::ostringstream message;
469  message << "Clone() method not implemented for type: "
470  << GetEntityType() << "!" << G4endl
471  << "Returning NULL pointer!";
472  G4Exception("G4UAdapter::Clone()", "GeomSolids1001", JustWarning, message);
473  return nullptr;
474 }
475 
476 template <class UnplacedVolume_t>
477 G4bool G4UAdapter<UnplacedVolume_t>::CalculateExtent(const EAxis pAxis,
478  const G4VoxelLimits& pVoxelLimit,
479  const G4AffineTransform& pTransform,
480  G4double& pMin, G4double& pMax) const
481 {
482  U3Vector vmin, vmax;
483  UnplacedVolume_t::Extent(vmin,vmax);
484  G4ThreeVector bmin(vmin.x(),vmin.y(),vmin.z());
485  G4ThreeVector bmax(vmax.x(),vmax.y(),vmax.z());
486 
487  // Check correctness of the bounding box
488  //
489  if (bmin.x() >= bmax.x() || bmin.y() >= bmax.y() || bmin.z() >= bmax.z())
490  {
491  std::ostringstream message;
492  message << "Bad bounding box (min >= max) for solid: "
493  << GetName() << " - " << GetEntityType() << " !"
494  << "\nmin = " << bmin
495  << "\nmax = " << bmax;
496  G4Exception("G4UAdapter::CalculateExtent()", "GeomMgt0001",
497  JustWarning, message);
498  StreamInfo(G4cout);
499  }
500 
501  G4BoundingEnvelope bbox(bmin,bmax);
502  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
503 }
504 
505 template <class UnplacedVolume_t>
506 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::CreatePolyhedron() const
507 {
508  // Must be implemented in concrete wrappers...
509 
510  std::ostringstream message;
511  message << "Visualization not supported for USolid shape "
512  << GetEntityType() << "... Sorry!" << G4endl;
513  G4Exception("G4UAdapter::CreatePolyhedron()", "GeomSolids0003",
514  FatalException, message);
515  return nullptr;
516 }
517 
518 template <class UnplacedVolume_t>
519 G4Polyhedron* G4UAdapter<UnplacedVolume_t>::GetPolyhedron() const
520 {
521  if (!fPolyhedron ||
522  fRebuildPolyhedron ||
523  fPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
524  fPolyhedron->GetNumberOfRotationSteps())
525  {
526  G4AutoLock l(&pMutex);
527  delete fPolyhedron;
528  fPolyhedron = CreatePolyhedron();
529  fRebuildPolyhedron = false;
530  l.unlock();
531  }
532  return fPolyhedron;
533 }
534 
535 template <class UnplacedVolume_t>
536 G4VisExtent G4UAdapter<UnplacedVolume_t>::GetExtent() const
537 {
538  U3Vector vmin, vmax;
539  UnplacedVolume_t::Extent(vmin,vmax);
540  return G4VisExtent(vmin.x(),vmax.x(),
541  vmin.y(),vmax.y(),
542  vmin.z(),vmax.z());
543 }
544 
545 #endif // G4GEOM_USE_USOLIDS
546 
547 #endif // G4UADAPTER_HH