Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-30 10:26:05

0001 // This file is part of VecGeom and is distributed under the
0002 // conditions in the file LICENSE.txt in the top directory.
0003 // For the full list of authors see CONTRIBUTORS.txt and `git log`.
0004 
0005 /// \brief Class storing the global navigation index lookup table.
0006 /// \file management/NavIndexTable.h
0007 /// \author Andrei Gheata (andrei.gheata@cern.ch)
0008 
0009 #ifndef VECGEOM_MANAGEMENT_NAVINDEXTABLE_H_
0010 #define VECGEOM_MANAGEMENT_NAVINDEXTABLE_H_
0011 
0012 #include "VecGeom/management/GeoVisitor.h"
0013 #include "VecGeom/navigation/NavigationState.h"
0014 
0015 #include <new>
0016 
0017 namespace vecgeom {
0018 inline namespace VECGEOM_IMPL_NAMESPACE {
0019 
0020 class NavIndexTable;
0021 
0022 class BuildNavIndexVisitor : public GeoVisitorNavIndex {
0023 private:
0024   size_t fTableSize = sizeof(
0025       NavIndex_t); ///< The iteration needs to be done in two passes, the first one is just to compute the array size
0026   NavIndex_t *fNavInd = nullptr; ///< Array storing the navigation info related to a given state
0027   int fLimitDepth     = 0;       ///< limit depth to scache the transformations, 0 means unlimited
0028   // NavIndex_t fWorld   = 1;       ///< index for the world volume
0029   NavIndex_t fCurrent = 1;     ///< Current navigation index being filled.
0030   bool fDoCount       = true;  ///< First pass to compute the table size
0031   bool fValidate      = false; ///< If set, this flag will force a thorough validation upon visiting
0032 
0033 public:
0034   BuildNavIndexVisitor(int depth_limit, bool do_count)
0035       : GeoVisitorNavIndex(), fLimitDepth(depth_limit), fDoCount(do_count)
0036   {
0037   }
0038 
0039   size_t GetTableSize() const { return fTableSize; }
0040   void SetTable(NavIndex_t *table) { fNavInd = table; }
0041   void SetDoCount(bool flag) { fDoCount = flag; }
0042   void SetValidate(bool flag) { fValidate = flag; }
0043 
0044   NavIndex_t apply(NavStatePath *state, int level, NavIndex_t mother, int dind);
0045 };
0046 
0047 class NavIndexTable {
0048 private:
0049   NavIndex_t *fNavInd       = nullptr; ///< address of the table
0050   VPlacedVolume *fVolBuffer = nullptr; ///< placed volume buffer
0051   NavIndex_t fWorld         = 1;       ///< index for the world volume
0052   size_t fTableSize         = 0;       ///< table size in bytes
0053   size_t fDepthLimit        = 0;       ///< depth limnit to which transformations will be cached
0054 
0055   NavIndexTable(NavIndex_t *table, size_t table_size) : fNavInd(table), fTableSize(table_size) {}
0056 
0057 public:
0058   ~NavIndexTable() { CleanTable(); }
0059   static NavIndexTable *Instance(NavIndex_t *nav_table = nullptr, size_t table_size = 0)
0060   {
0061     static NavIndexTable instance(nav_table, table_size);
0062     return &instance;
0063   }
0064 
0065   void CleanTable()
0066   {
0067     delete[] fNavInd;
0068     fNavInd = nullptr;
0069   }
0070 
0071   VECCORE_ATT_HOST_DEVICE
0072   VECGEOM_FORCE_INLINE
0073   size_t GetTableSize() const { return fTableSize; }
0074 
0075   VECCORE_ATT_HOST_DEVICE
0076   VECGEOM_FORCE_INLINE
0077   NavIndex_t *GetTable() const { return fNavInd; }
0078 
0079   VECCORE_ATT_HOST_DEVICE
0080   VECGEOM_FORCE_INLINE
0081   void SetVolumeBuffer(VPlacedVolume *buffer) { fVolBuffer = buffer; }
0082 
0083   VECCORE_ATT_HOST_DEVICE
0084   VECGEOM_FORCE_INLINE
0085   NavIndex_t GetWorld() const { return fWorld; }
0086 
0087   bool AllocateTable(size_t bytes)
0088   {
0089     size_t nelem       = bytes / sizeof(NavIndex_t) + 1;
0090     NavIndex_t *buffer = new (std::nothrow) NavIndex_t[nelem];
0091     if (buffer) {
0092       assert((reinterpret_cast<uintptr_t>(buffer) % sizeof(Precision) == 0) &&
0093              "Buffer must be aligned to at least Precision, because we store Transformation3D's translation and "
0094              "rotation in it.");
0095       fNavInd    = buffer;
0096       fNavInd[0] = 0;
0097       fTableSize = bytes;
0098     } else {
0099       std::cout << "=== EEE === AlocateTable bad_alloc intercepted while trying to allocate " << bytes << " bytes\n";
0100     }
0101     return buffer != nullptr;
0102   }
0103 
0104   static size_t ComputeTableSize(VPlacedVolume const *top, int maxdepth, int depth_limit)
0105   {
0106     NavStatePath *state = NavStatePath::MakeInstance(maxdepth);
0107     BuildNavIndexVisitor visitor(depth_limit, true); // just count table size
0108     visitAllPlacedVolumesNavIndex(top, &visitor, state);
0109     size_t table_size = visitor.GetTableSize();
0110     NavStatePath::ReleaseInstance(state);
0111     return table_size;
0112   }
0113 
0114   bool CreateTable(VPlacedVolume const *top, int maxdepth, int depth_limit)
0115   {
0116     fDepthLimit         = depth_limit;
0117     NavStatePath *state = NavStatePath::MakeInstance(maxdepth);
0118     BuildNavIndexVisitor visitor(depth_limit, true); // just count table size
0119 
0120     visitAllPlacedVolumesNavIndex(top, &visitor, state);
0121     bool hasTable = AllocateTable(visitor.GetTableSize());
0122     if (!hasTable) return false;
0123 
0124     visitor.SetTable(fNavInd);
0125     visitor.SetDoCount(false);
0126 
0127     state->Clear();
0128     visitAllPlacedVolumesNavIndex(top, &visitor, state);
0129     NavStatePath::ReleaseInstance(state);
0130 
0131     return true;
0132   }
0133 
0134   bool Validate(VPlacedVolume const *top, int maxdepth) const
0135   {
0136     NavStatePath *state = NavStatePath::MakeInstance(maxdepth);
0137     state->Clear();
0138     auto visitor = new BuildNavIndexVisitor(0, false);
0139     visitor->SetValidate(true);
0140     visitAllPlacedVolumesNavIndex(top, visitor, state);
0141     NavStatePath::ReleaseInstance(state);
0142     return true;
0143   }
0144 
0145   NavIndex_t ValidateState(NavStatePath *state);
0146 
0147   // vecgeom::cuda::NavIndexTable *CopyToGPU() const
0148 
0149   /// Traverses the geometry tree keeping track of the state context (volume path or navigation state)
0150   /// and applies the injected Visitor
0151   template <typename Visitor>
0152   static void visitAllPlacedVolumesNavIndex(VPlacedVolume const *currentvolume, Visitor *visitor, NavStatePath *state,
0153                                             int level = 0, NavIndex_t mother = 0, int dind = 0)
0154   {
0155     if (currentvolume != NULL) {
0156       state->Push(currentvolume);
0157       NavIndex_t new_mother = visitor->apply(state, level, mother, dind);
0158       int size              = currentvolume->GetDaughters().size();
0159       for (int i = 0; i < size; ++i) {
0160         visitAllPlacedVolumesNavIndex(currentvolume->GetDaughters().operator[](i), visitor, state, level + 1,
0161                                       new_mother, i);
0162       }
0163       state->Pop();
0164     }
0165   }
0166 };
0167 
0168 } // namespace VECGEOM_IMPL_NAMESPACE
0169 } // namespace vecgeom
0170 
0171 #endif