C++ API

This reference is generated from Doxygen XML and rendered in Sphinx via Breathe.

class AdePTConfiguration

Create and configure instances of an AdePT transport implementation.

Public Functions

AdePTConfiguration()
~AdePTConfiguration()
inline void SetNumThreads(int numThreads)
inline void SetTrackInAllRegions(bool trackInAllRegions)
inline void SetCallUserSteppingAction(bool callUserSteppingAction)
inline void SetCallUserTrackingAction(bool callUserTrackingAction)
inline void AddGPURegionName(std::string name)
inline void RemoveGPURegionName(std::string name)
inline void AddWDTRegionName(std::string name)
inline void AddDeadRegionName(std::string name)
inline void SetVerbosity(int verbosity)
inline void SetMillionsOfTrackSlots(double millionSlots)
inline void SetMillionsOfHitSlots(double millionSlots)
inline void SetHitBufferFlushThreshold(float threshold)
inline void SetCPUCapacityFactor(float CPUCapacityFactor)
inline void SetHitBufferSafetyFactor(double HitBufferSafetyFactor)
inline void SetAdePTSeed(int seed)
inline void SetCUDAStackLimit(int limit)
inline void SetCUDAHeapLimit(int limit)
inline void SetLastNParticlesOnCPU(int Nparticles)
inline void SetMaxWDTIter(int maxIter)
inline void SetWDTKineticEnergyLimit(double ekin)
inline void SetSpeedOfLight(bool speedOfLight)
inline void SetMultipleStepsInMSCWithTransportation(bool setMultipleSteps)
inline void SetEnergyLossFluctuation(bool setELossFluct)
inline void SetVecGeomGDML(std::string filename)
inline void SetCovfieBfieldFile(std::string filename)
inline bool GetTrackInAllRegions() const
inline bool GetCallUserSteppingAction() const
inline bool GetCallUserTrackingAction() const
inline bool GetSpeedOfLight() const
inline bool GetMultipleStepsInMSCWithTransportation() const
inline bool GetEnergyLossFluctuation() const
inline int GetNumThreads() const
inline int GetVerbosity() const
inline int GetCUDAStackLimit() const
inline int GetCUDAHeapLimit() const
inline uint64_t GetAdePTSeed() const
inline unsigned short GetLastNParticlesOnCPU() const
inline unsigned short GetMaxWDTIter() const
inline double GetWDTKineticEnergyLimit() const
inline float GetHitBufferFlushThreshold() const
inline float GetCPUCapacityFactor() const
inline double GetHitBufferSafetyFactor() const
inline double GetMillionsOfTrackSlots() const
inline double GetMillionsOfHitSlots() const
inline const std::vector<std::string> *GetGPURegionNames() const
inline const std::vector<std::string> *GetCPURegionNames() const
inline const std::vector<std::string> &GetWDTRegionNames() const
inline const std::vector<std::string> &GetDeadRegionNames() const
inline std::string GetVecGeomGDML() const
inline std::string GetCovfieBfieldFile() const

Private Members

bool fTrackInAllRegions = {false}
bool fCallUserSteppingAction = {false}
bool fCallUserTrackingAction = {false}
bool fSpeedOfLight = {false}
bool fSetMultipleStepsInMSCWithTransportation = {false}
bool fSetEnergyLossFluctuation = {false}
int fNumThreads = {-1}
int fVerbosity = {0}
int fCUDAStackLimit = {0}
int fCUDAHeapLimit = {0}
uint64_t fAdePTSeed = {1234567}
float fHitBufferFlushThreshold = {0.8}
float fCPUCapacityFactor = {2.5}
double fHitBufferSafetyFactor = {1.5}
double fMillionsOfTrackSlots = {1}
double fMillionsOfHitSlots = {1}
unsigned short fLastNParticlesOnCPU = {0}
unsigned short fMaxWDTIter = {5}
double fWDTKineticEnergyLimit = {0.2}
std::vector<std::string> fGPURegionNames = {}
std::vector<std::string> fCPURegionNames = {}
std::vector<std::string> fWDTRegionNames = {}
std::vector<std::string> fDeadRegionNames = {}
std::string fVecGeomGDML = {""}
std::string fCovfieBfieldFile = {""}
std::unique_ptr<AdePTConfigurationMessenger> fAdePTConfigurationMessenger
class AdePTConfigurationMessenger : public G4UImessenger

Public Functions

AdePTConfigurationMessenger(AdePTConfiguration*)

Implementation of the AdePTConfigurationMessenger class.

~AdePTConfigurationMessenger() override
void SetNewValue(G4UIcommand*, G4String) override

Private Members

AdePTConfiguration *fAdePTConfiguration
std::unique_ptr<G4UIdirectory> fDir
std::unique_ptr<G4UIcmdWithAnInteger> fSetCUDAStackLimitCmd
std::unique_ptr<G4UIcmdWithAnInteger> fSetCUDAHeapLimitCmd
std::unique_ptr<G4UIcmdWithAnInteger> fSetAdePTSeedCmd
std::unique_ptr<G4UIcmdWithAnInteger> fSetFinishOnCpuCmd
std::unique_ptr<G4UIcmdWithAnInteger> fSetMaxWDTIterCmd
std::unique_ptr<G4UIcmdWithADouble> fSetWDTKineticEnergyLimitCmd
std::unique_ptr<G4UIcmdWithABool> fSetTrackInAllRegionsCmd
std::unique_ptr<G4UIcmdWithABool> fSetCallUserSteppingActionCmd
std::unique_ptr<G4UIcmdWithABool> fSetCallUserTrackingActionCmd
std::unique_ptr<G4UIcmdWithABool> fSetSpeedOfLightCmd
std::unique_ptr<G4UIcmdWithABool> fSetMultipleStepsInMSCWithTransportationCmd
std::unique_ptr<G4UIcmdWithABool> fSetEnergyLossFluctuationCmd
std::unique_ptr<G4UIcmdWithAString> fAddRegionCmd
std::unique_ptr<G4UIcmdWithAString> fRemoveRegionCmd
std::unique_ptr<G4UIcmdWithAString> fAddWDTRegionCmd
std::unique_ptr<G4UIcmdWithAnInteger> fSetVerbosityCmd
std::unique_ptr<G4UIcmdWithADouble> fSetMillionsOfTrackSlotsCmd
std::unique_ptr<G4UIcmdWithADouble> fSetMillionsOfHitSlotsCmd
std::unique_ptr<G4UIcmdWithADouble> fSetHitBufferFlushThresholdCmd
std::unique_ptr<G4UIcmdWithADouble> fSetCPUCapacityFactorCmd
std::unique_ptr<G4UIcmdWithADouble> fSetHitBufferSafetyFactorCmd
std::unique_ptr<G4UIcmdWithAString> fSetGDMLCmd
std::unique_ptr<G4UIcmdWithAString> fSetCovfieFileCmd
class AdePTG4HepEmState

Owns the prepared host-side G4HepEm inputs used by transport.

The Geant4 integration side prepares one of these objects before the shared transport is created. This wrapper owns both:

  • the rebuilt G4HepEmData

  • a deep copy of the G4HepEmParameters taken from the provided config

Cleanup is intentionally split:

  • DataDeleter performs the deep cleanup of the owned G4HepEmData and then deletes the outer G4HepEmData allocation.

  • ParametersDeleter performs the deep cleanup of the owned G4HepEmParameters, including the GPU mirror created by AdePT, and then deletes the outer G4HepEmParameters allocation.

Public Functions

explicit AdePTG4HepEmState(G4HepEmConfig *hepEmConfig)

Build the AdePT-owned G4HepEmData and G4HepEmParameters copies from the supplied config.

Rebuild a complete AdePT-owned set of host-side G4HepEm inputs from the supplied config.

AdePTG4HepEmState owns two different G4HepEm objects:

  • a deep copy of the G4HepEmParameters stored in the supplied G4HepEmConfig

  • a freshly rebuilt G4HepEmData derived from that copied parameter block

We must copy G4HepEmParameters because the original object remains owned by the worker-local G4HepEmConfig, while the shared AdePT transport can outlive the worker that first created it. G4HepEmData is rebuilt here directly, so it is already fully owned by AdePT and does not need a second copy step.

~AdePTG4HepEmState()

Destroy the owned G4HepEmData and G4HepEmParameters copies.

AdePTG4HepEmState(const AdePTG4HepEmState&) = delete
AdePTG4HepEmState &operator=(const AdePTG4HepEmState&) = delete
AdePTG4HepEmState(AdePTG4HepEmState&&) noexcept = default
AdePTG4HepEmState &operator=(AdePTG4HepEmState&&) noexcept = default
inline G4HepEmData *GetData() const

Access the owned host-side HepEm data tables.

inline G4HepEmParameters *GetParameters() const

Access the owned HepEm parameter copy.

Private Members

std::unique_ptr<G4HepEmData, DataDeleter> fData

Owned G4HepEmData rebuilt for AdePT.

std::unique_ptr<G4HepEmParameters, ParametersDeleter> fParameters

Owned deep copy of G4HepEmParameters used to build and upload transport data.

class AdePTGeant4Integration

Public Types

enum class DeferredStepType : unsigned char

Values:

enumerator ReplayStep
enumerator ReturnTrack

Public Functions

inline explicit AdePTGeant4Integration()
~AdePTGeant4Integration()
AdePTGeant4Integration(const AdePTGeant4Integration&) = delete
AdePTGeant4Integration &operator=(const AdePTGeant4Integration&) = delete
AdePTGeant4Integration(AdePTGeant4Integration&&) = default
AdePTGeant4Integration &operator=(AdePTGeant4Integration&&) = default
void ProcessGPUStep(std::span<const GPUStep> gpuSteps, bool const callUserSteppingAction = false, bool const callUserTrackingaction = false)

Reconstructs GPU steps on host and calls the user-defined sensitive detector code.

void ReturnDeferredTrack(std::span<const GPUStep> gpuSteps, bool const callUserActions = false)

Return a deferred parent track to Geant4 without rebuilding the visible G4 step.

This is only used for returned gamma handoff steps with one parent step, no secondaries, and zero deposited energy. In that case there is no GPU step to process on the host, so only the parent G4Track is rebuilt from the post-step state and pushed back to the Geant4 stack.

std::vector<float> GetUniformField() const

Returns the Z value of the user-defined uniform magnetic field.

This function can only be called when the user-defined field is a G4UniformMagField

inline int GetEventID() const
inline int GetThreadID() const
inline HostTrackDataMapper &GetHostTrackDataMapper()
void QueueDeferredStep(std::span<const GPUStep> gpuSteps, DeferredStepType type = DeferredStepType::ReplayStep)

Defer a returned step for later sorted replay on the host.

DeferredStepStore TakeDeferredSteps()

Transfer ownership of the currently queued deferred steps.

This drains the integration-local deferred-step storage without copying it.

inline void SetHepEmTrackingManager(G4HepEmTrackingManagerSpecialized *hepEmTrackingManager)

Private Functions

void FillG4NavigationHistory(const vecgeom::NavigationState &aNavState, G4NavigationHistory &aG4NavigationHistory) const

Reconstruct G4TouchableHistory from a VecGeom Navigation index.

G4TouchableHandle MakeTouchableFromNavState(vecgeom::NavigationState const &navState) const
G4Track *ConstructSecondaryTrackInPlace(GPUStep const *secStep) const

Construct the temporary secondary track that is attached to the secondary vector of the parent step.

void InitSecondaryHostTrackDataFromParent(GPUStep const *secStep, HostTrackData &secTData, int g4ParentID, G4TouchableHandle &preTouchable) const
void FillG4Track(GPUStep const *aGPUStep, G4Track *aG4Track, const HostTrackData &hostTData, G4TouchableHandle &aPreG4TouchableHandle, G4TouchableHandle &aPostG4TouchableHandle) const
void FillG4Step(GPUStep const *aGPUStep, G4Step *aG4Step, const HostTrackData &hostTData, G4TouchableHandle &aPreG4TouchableHandle, G4TouchableHandle &aPostG4TouchableHandle, G4StepStatus aPreStepStatus, G4StepStatus aPostStepStatus, bool callUserTrackingAction, bool callUserSteppingAction) const
G4Track *MakeTrackForCPUStacking(const G4Track &track) const

Create a heap-owned track that can be pushed onto the Geant4 stack.

This is only used as a fallback for gamma/lepton nuclear when no Geant4 nuclear process is attached. In that case there is no temporary nuclear replay track to continue on the CPU, and the visible reconstructed track cannot be handed to the stack manager because it is reused integration storage.

G4Track *MakeReturnedTrackFromStep(GPUStep const &parentStep, const HostTrackData &hostTData, bool setStopButAlive) const

Recreate a track to from a returned parent step to be continued on CPU.

Out-of-GPU-region and finish-on-CPU steps used to hand Geant4 a returned track built from the post-step state. The visible reconstructed step keeps the transported GPU-step data, but the continued CPU track must still match that old current-state handoff.

Private Members

G4HepEmTrackingManagerSpecialized *fHepEmTrackingManager = {nullptr}
std::unique_ptr<HostTrackDataMapper> fHostTrackDataMapper
std::unique_ptr<AdePTGeant4Integration_detail::StepReconstructionObjects, AdePTGeant4Integration_detail::Deleter> fStepReconstructionObjects = {nullptr}
std::vector<GPUStep> fDeferredGPUSteps
std::vector<DeferredStep> fDeferredSteps
class AdePTGeometryBridge

Global bridge between Geant4 host geometry and the VecGeom world used by AdePT.

The geometry world and the VecGeom-to-Geant4 lookup tables are process-global, so this bridge provides the corresponding global services needed during setup and touchable reconstruction.

Public Static Functions

static void CreateVecGeomWorld(G4VPhysicalVolume const *physvol)

Converts the Geant4 world volume into the process-global VecGeom world.

Convert the Geant4 world into the process-global VecGeom world via G4VG.

Throws:

std::runtime_error – if the input Geant4 world or resulting VecGeom world is null.

static void CheckGeometry(G4HepEmData const *hepEmData)

Verifies that the Geant4 and VecGeom geometries match.

Compare the Geant4 and VecGeom host geometries for consistency.

static void InitVolAuxData(adeptint::VolAuxData *volAuxData, G4HepEmData const *hepEmData, G4HepEmTrackingManagerSpecialized *hepEmTM, bool trackInAllRegions, std::vector<std::string> const *gpuRegionNames, std::vector<std::string> const &deadRegionNames, adeptint::WDTHostRaw &wdtRaw)

Fills the auxiliary per-volume data needed by AdePT.

Fill the auxiliary per-volume transport metadata used by AdePT.

static adeptint::WDTHostPacked PackWDT(adeptint::WDTHostRaw const &wdtRaw)

Pack the Woodcock tracking data from the sparse host-side map into arrays that can be copied to the GPU.

Parameters:

wdtRaw – Raw WDT data collected during geometry traversal.

Returns:

Packed, dense WDT data ready to be handed to the transport for device upload.

static G4VPhysicalVolume const *GetG4PhysicalVolume(vecgeom::VPlacedVolume const *placedVolume)

Returns the Geant4 placed volume matching a VecGeom placed volume.

Resolve the Geant4 placed volume associated with a VecGeom placed volume.

Throws:

std::runtime_error – if the VecGeom placed volume is not present in the global lookup table.

Private Static Functions

static void MapVecGeomToG4(std::vector<G4VPhysicalVolume const*> &vecgeomPvToG4Map, std::vector<G4LogicalVolume const*> &vecgeomLvToG4Map)

Builds the lookup tables from VecGeom placed/logical volume ids to Geant4 volumes.

Private Static Attributes

static std::vector<G4VPhysicalVolume const*> fGlobalVecGeomPvToG4Map
static std::vector<G4LogicalVolume const*> fGlobalVecGeomLvToG4Map
class AdePTTrackingManager : public G4VTrackingManager

Public Types

using AdePTTransport = adept::transport::AdePTTransport

Public Functions

explicit AdePTTrackingManager(AdePTConfiguration *config, int verbosity = 0)
~AdePTTrackingManager()
void BuildPhysicsTable(const G4ParticleDefinition&) override
void PreparePhysicsTable(const G4ParticleDefinition&) override
void HandOverOneTrack(G4Track *aTrack) override
void FlushEvent() override
void InitializeAdePT()
inline G4HepEmConfig *GetG4HepEmConfig()

Private Functions

void ProcessTrack(G4Track *aTrack)

Steps a particle using the generic G4 tracking, until it dies or enters a user-defined GPU region, in which case tracking is delegated to AdePT.

In order to be able to send tracks back and forth between the custom tracking and the generic G4 tracking, this function has to provide the functionality of G4TrackingManager::ProcessOneTrack, in addition to checking whether the track is in a GPU Region

const vecgeom::NavigationState GetVecGeomFromG4State(const G4NavigationHistory &aG4NavigationHistory)

Get the corresponding VecGeom NavigationState from the G4NavigationHistory. The boundary status will be set to false by default.

Parameters:

aG4NavigationHistory – the given G4NavigationHistory

Returns:

the corresponding vecgeom::NavigationState

const vecgeom::NavigationState GetVecGeomFromG4State(const G4Track &aG4Track, const G4NavigationHistory *aG4NavigationHistory = nullptr)

Get the corresponding VecGeom NavigationState from the track’s G4NavigationHistory, and set the boundary status based on the track’s state.

Parameters:
  • aG4Track – The G4Track from which to extract the NavigationState

  • aG4NavigationHistory – Navigation history that is used to define the navState. If not provided, it will default to aG4Track.GetNextTouchableHandle()->GetHistory()

Returns:

The corresponding vecgeom::NavigationState

void InitializeSharedAdePTTransport()

Perform the one-time shared AdePT transport initialization on the first Geant4 worker.

The first worker prepares all host-side inputs needed by transport:

  • the uniform magnetic-field values

  • the AdePT-owned AdePTG4HepEmState

  • geometry consistency checks

  • VolAuxData

  • packed WDT metadata

Once that host-side preparation is complete, the worker creates the shared AdePTTransport. The transport constructor then performs the corresponding one-time device initialization and upload.

void ProcessReturnedGPUSteps(int threadId, int eventId)

Drain returned GPU-step batches from transport and reconstruct the corresponding Geant4 steps on the CPU.

Transport still owns the batch lifetime and iterates over the available batches. This helper provides the Geant4-side reconstruction logic for each batch.

Private Members

std::unique_ptr<G4HepEmTrackingManagerSpecialized> fHepEmTrackingManager
AdePTGeant4Integration fGeant4Integration
std::set<G4Region const*> fGPURegions = {}
std::shared_ptr<AdePTTransport> fAdeptTransport
AdePTConfiguration *const fAdePTConfiguration
int fVerbosity = {0}
unsigned int fTrackCounter = {0}
int fCurrentEventID = {0}
bool fAdePTInitialized = {false}
bool fSpeedOfLight = {false}

Private Static Attributes

static int fNumThreads = {0}
class AdePTTransport
#include <AdePTTransport.hh>

Public Functions

AdePTTransport(const AdePTTransportConfig &configuration, std::unique_ptr<AdePTG4HepEmState> adeptG4HepEmState, adeptint::VolAuxData *auxData, const adeptint::WDTHostPacked &wdtPacked, const std::vector<float> &uniformFieldValues)
AdePTTransport(const AdePTTransport &other) = delete
~AdePTTransport()
void AddTrack(int pdg, uint64_t trackId, uint64_t parentId, double energy, double x, double y, double z, double dirx, double diry, double dirz, double globalTime, double localTime, double properTime, float weight, unsigned short stepCounter, int threadId, unsigned int eventId, vecgeom::NavigationState &&state)

Adds a track to the buffer.

inline int GetDebugLevel() const
template<typename Callback>
void HandleGPUStepBatchesWith(int threadId, int eventId, Callback &&callback)

Handle the currently available returned GPU-step batches for one thread and event.

Transport retains ownership of the step-buffer lifetime. For each available batch, callback is invoked with a std::span<const GPUStep> view and the batch is released again when the callback returns.

In this code path, the callback is the AdePTTrackingManager logic that reconstructs Geant4 steps from the returned GPU steps.

void RequestFlush(int threadId)

Request that the device flush all pending work for the given worker.

void WaitForFlushProgress()

Wait until the transport threads make further flush progress.

bool AreReturnedStepsFlushed(int threadId) const

Check whether all returned GPU-step batches have been made available on the host.

void MarkHostFlushed(int threadId)

Mark the worker event as fully handled on the host side after replaying the returned steps.

Public Members

uint64_t fAdePTSeed = 1234567

Private Functions

void Initialize(adeptint::VolAuxData *auxData, const adeptint::WDTHostPacked &wdtPacked, const std::vector<float> &uniformFieldValues)
void InitBVH()
bool InitializeGeometry(const vecgeom::cxx::VPlacedVolume *world)
bool InitializePhysics()

Private Members

unsigned int fNThread = {0}

Number of G4 workers.

unsigned int fTrackCapacity = {0}

Number of track slots to allocate on device.

unsigned int fStepCapacity = {0}

Number of step slots to allocate on device.

int fDebugLevel = {0}

Debug level.

int fCUDAStackLimit = {0}

CUDA device stack limit.

int fCUDAHeapLimit = {0}

CUDA device heap limit.

unsigned short fLastNParticlesOnCPU = {0}

Number N of last N particles that are finished on CPU.

unsigned short fMaxWDTIter = {5}

Maximum number of Woodcock tracking iterations per step.

std::unique_ptr<GPUstate, GPUstateDeleter> fGPUstate = {nullptr}

CUDA state placeholder.

std::unique_ptr<TrackBuffer> fBuffer = {nullptr}

Buffers for transferring tracks between host and device.

std::unique_ptr<AdePTG4HepEmState> fAdePTG4HepEmState

Transport-owned wrapper around G4HepEmData and copied G4HepEmParameters

adeptint::WDTDeviceBuffers fWDTDev = {}

device buffers for Woodcock tracking data

std::thread fGPUWorker

Thread to manage GPU.

std::condition_variable fCV_G4Workers

Communicate with G4 workers.

std::mutex fMutex_G4Workers

Mutex associated to the condition variable.

std::vector<std::atomic<EventState>> fEventStates

State machine for each G4 worker.

bool fHasWDTRegions = false
bool fReturnAllSteps = false
bool fReturnFirstAndLastStep = false
std::string fBfieldFile = {""}

Path to magnetic field file (in the covfie format)

double fCPUCapacityFactor{2.5}

Factor by which the host step buffer is larger than the device step buffer. Must be at least 2 Filling fraction of the host step buffer when the steps are copied out and not taken directly by the G4 workers.

double fCPUCopyFraction = {0.5}

Needed to stall the GPU, in case the nPartInFlight * fStepBufferSafetyFactor > available StepSlots.

double fStepBufferSafetyFactor = {1.5}
struct AdePTTransportConfig

Plain configuration values required to initialize the AdePT transport engine.

Public Members

uint64_t adeptSeed = {1234567}
unsigned int numThreads = {0}
unsigned int trackCapacity = {0}
unsigned int stepCapacity = {0}
int debugLevel = {0}
int cudaStackLimit = {0}
int cudaHeapLimit = {0}
unsigned short lastNParticlesOnCPU = {0}
unsigned short maxWDTIter = {5}
bool returnAllSteps = {false}
bool returnFirstAndLastStep = {false}
std::string bfieldFile = {}
double cpuCapacityFactor = {2.5}
double cpuCopyFraction = {0.5}
double stepBufferSafetyFactor = {1.5}
struct AllowFinishOffEventArray
#include <TransportStats.hh>

Array of flags whether the event can be finished off.

Public Functions

inline unsigned short operator[](unsigned int idx) const

Public Members

const unsigned short *flags = {nullptr}
template<typename Type, typename Enable = void>
struct Atomic_t

Atomic_t generic implementation specialized using SFINAE mechanism.

template<typename Type>
struct Atomic_t<Type, typename std::enable_if<std::is_integral<Type>::value>::type> : public adept::AtomicBase_t<Type>
#include <Atomic.h>

Specialization for integral types.

Public Functions

inline Type operator=(Type desired)

Standard library atomic operations for integral types.

Atomically assigns the desired value to the atomic variable.

inline Type fetch_add(Type arg)

Atomically replaces the current value with the result of arithmetic addition of the value and arg.

inline Type fetch_sub(Type arg)

Atomically replaces the current value with the result of arithmetic subtraction of the value and arg.

inline Type fetch_and(Type arg)

Atomically replaces the current value with the result of bitwise AND of the value and arg.

inline Type fetch_or(Type arg)

Atomically replaces the current value with the result of bitwise OR of the value and arg.

inline Type fetch_xor(Type arg)

Atomically replaces the current value with the result of bitwise XOR of the value and arg.

inline Type operator+=(Type arg)

Performs atomic pre-increment.

inline Type operator-=(Type arg)

Performs atomic subtract.

inline Type operator++()

Performs atomic pre-increment.

inline Type operator++(int)

Performs atomic post-increment.

inline Type operator--()

Performs atomic pre-decrement.

inline Type operator--(int)

Performs atomic post-decrement.

template<typename Type>
struct Atomic_t<Type, typename std::enable_if<!std::is_integral<Type>::value>::type> : public adept::AtomicBase_t<Type>
#include <Atomic.h>

Specialization for non-integral types.

Public Types

using Base_t = AtomicBase_t<Type>

Public Functions

inline Type operator=(Type desired)

Atomically assigns the desired value to the atomic variable.

inline Type fetch_add(Type arg)

Atomically replaces the current value with the result of arithmetic addition of the value and arg.

inline Type fetch_sub(Type arg)

Atomically replaces the current value with the result of arithmetic subtraction of the value and arg.

inline Type operator+=(Type arg)

Performs atomic pre-increment.

inline Type operator-=(Type arg)

Performs atomic subtract.

inline Type operator++()

Performs atomic pre-increment.

inline Type operator++(int)

Performs atomic post-increment.

inline Type operator--()

Performs atomic pre-decrement.

inline Type operator--(int)

Performs atomic post-decrement.

Public Members

AtomicType_t fData

Atomic data.

template<typename Type>
struct AtomicBase_t
#include <Atomic.h>

A portable atomic type.

Not all base types are supported by CUDA.

Subclassed by adept::Atomic_t< Type, typename std::enable_if< std::is_integral< Type >::value >::type >, adept::Atomic_t< Type, typename std::enable_if<!std::is_integral< Type >::value >::type >

Public Types

using AtomicType_t = std::atomic<Type>

Standard library atomic behaviour.

Public Functions

inline AtomicBase_t()

Constructor taking an address.

inline AtomicBase_t(AtomicBase_t const &other)

Copy constructor.

inline AtomicBase_t &operator=(AtomicBase_t const &other)

Assignment.

inline void store(Type desired)

Implementation-dependent.

Atomically replaces the current value with desired.

inline Type load() const

Atomically loads and returns the current value of the atomic variable.

inline Type exchange(Type desired)

Atomically replaces the underlying value with desired.

inline bool compare_exchange_strong(Type &expected, Type desired)

Atomically compares the stored value with the expected one, and if equal stores desired.

If not loads the old value into expected. Returns true if swap was successful.

Public Members

AtomicType_t fData = {0}

Atomic data.

Public Static Functions

static inline AtomicBase_t *MakeInstanceAt(void *addr)

Emplace the data at a given address.

template<typename Type>
class BlockData : protected copcore::VariableSizeObjectInterface<BlockData<Type>, Type>
#include <BlockData.h>

A container of data that adopts memory.

It is caller responsibility to allocate at least SizeOfInstance bytes for a single object, and SizeOfAlignAware for multiple objects. Write access to data elements in the block is given atomically, up to the capacity.

Public Types

using AtomicInt_t = adept::Atomic_t<int>
using Queue_t = adept::mpmc_bounded_queue<int>
using Value_t = Type
using Base_t = copcore::VariableSizeObjectInterface<BlockData<Value_t>, Value_t>
using ArrayData_t = copcore::VariableSizeObj<Value_t>

Public Functions

inline int Capacity() const

Maximum number of elements.

inline void Clear()

Clear the content.

inline Type const &operator[](const int index) const

Read-only index operator.

inline Type &operator[](const int index)

Read/write index operator.

inline Type *NextElement()

Dispatch next free element, nullptr if none left.

inline void ReleaseElement(int index)

Release an element.

inline int GetNused()

Number of elements currently distributed.

inline int GetNholes()

Number of holes in the block.

inline bool IsFull() const

Check if container is fully distributed.

Public Static Functions

static inline size_t SizeOfInstance(int capacity)

Returns the size in bytes of a BlockData object with given capacity.

static inline Cont *MakeCopy(const Cont &other)

< Enumerate the part of the private interface, we want to expose.

static inline Cont *MakeCopy(size_t new_size, const Cont &other)

< Enumerate the part of the private interface, we want to expose.

static inline Cont *MakeCopyAt(const Cont &other, void *addr)
static inline Cont *MakeCopyAt(size_t new_size, const Cont &other, void *addr)
template<typename ...T>
static inline Cont *MakeInstance(size_t nvalues, const T&... params)
template<typename ...T>
static inline Cont *MakeInstanceAt(size_t nvalues, void *addr, const T&... params)
static inline void ReleaseInstance(Cont *obj)
static inline constexpr size_t SizeOf(size_t nvalues)
static inline constexpr size_t SizeOfAlignAware(size_t nvalues)

Private Functions

inline ArrayData_t &GetVariableData()

Functions required by VariableSizeObjectInterface.

inline const ArrayData_t &GetVariableData() const
inline BlockData(size_t nvalues)
inline BlockData(BlockData const &other)
inline BlockData(size_t new_size, BlockData const &other)
inline ~BlockData()

Private Members

int fCapacity = {0}

Maximum number of elements.

AtomicInt_t fNbooked

Number of booked elements.

AtomicInt_t fNused

Number of used elements.

Queue_t *fHoles = {nullptr}

Queue of holes.

ArrayData_t fData

Data follows, has to be last.

friend Base_t

Private Static Functions

static inline size_t SizeOfExtra(int capacity)

Returns the size in bytes of extra queue data needed by BlockData object with given capacity.

template<typename Type>
struct Cell_t

Internal data structure to handle the data sequence.

Public Functions

inline Cell_t()

Cell constructor.

Public Members

adept::Atomic_t<int> fSequence

Atomic sequence counter.

Type fData

Data stored in the cell.

class ConstBzFieldStepper

A very simple stepper treating the propagation of particles in a constant Bz magnetic field ( neglecting energy loss of particle ) This class is roughly equivalent to TGeoHelix in ROOT.

Public Functions

inline ConstBzFieldStepper(float Bz = 0.)
inline void SetBz(double Bz)
inline double GetBz() const
template<typename BaseType, typename BaseIType> inline  __attribute__ ((always_inline)) void DoStep(BaseType const &

this function propagates the track along the helix solution by a step input: current position, current direction, some particle properties output: new position, new direction of particle

template<typename Vector3D, typename BaseType, typename BaseIType>
inline void DoStep(Vector3D const &pos, Vector3D const &dir, BaseIType const &charge, BaseType const &momentum, BaseType const &step, Vector3D &newpos, Vector3D &newdir) const

this function propagates the track along the helix solution by a step input: current position, current direction, some particle properties output: new position, new direction of particle

Public Members

BaseType const BaseType const BaseType const BaseType const BaseType const BaseIType const BaseType const BaseType const BaseType BaseType BaseType BaseType BaseType BaseType & const

Private Members

double fBz
class ConstFieldHelixStepper

A simple stepper treating the propagation of particles in a constant magnetic field ( not just along the z-axis&#8212; for that there is ConstBzFieldHelixStepper )

Public Functions

inline ConstFieldHelixStepper(double Bx, double By, double Bz)
inline ConstFieldHelixStepper(double Bfield[3])
inline ConstFieldHelixStepper(Vector3D<double> const &Bfield)
inline void SetB(double Bx, double By, double Bz)
inline Vector3D<double> const GetFieldVec() const
template<typename Real_t, typename Int_t>
inline void DoStep(Real_t const &posx, Real_t const &posy, Real_t const &posz, Real_t const &dirx, Real_t const &diry, Real_t const &dirz, Int_t const &charge, Real_t const &momentum, Real_t const &step, Real_t &newposx, Real_t &newposy, Real_t &newposz, Real_t &newdirx, Real_t &newdiry, Real_t &newdirz) const

this function propagates the track along the helix solution by a step input: current position, current direction, some particle properties output: new position, new direction of particle

this function propagates the track along the “helix-solution” by a step step input: current position (x0, y0, z0), current direction ( dirX0, dirY0, dirZ0 ), some particle properties output: new position, new direction of particle

template<typename Real_t, typename Int_t>
inline void DoStep(Vector3D<Real_t> const &position, Vector3D<Real_t> const &direction, Int_t const &charge, Real_t const &momentum, Real_t const &step, Vector3D<Real_t> &endPosition, Vector3D<Real_t> &endDirection) const

this function propagates the track along the helix solution by a step input: current position, current direction, some particle properties output: new position, new direction of particle

template<typename Real_t, typename Int_t>
inline void PrintStep(vecgeom::Vector3D<Real_t> const &startPosition, vecgeom::Vector3D<Real_t> const &startDirection, Int_t const &charge, Real_t const &momentum, Real_t const &step, vecgeom::Vector3D<Real_t> &endPosition, vecgeom::Vector3D<Real_t> &endDirection) const
template<typename Real_t, typename Int_t>
inline void DoStep(vecgeom::Vector3D<Real_t> const &startPosition, vecgeom::Vector3D<Real_t> const &startDirection, Int_t const &charge, Real_t const &momentum, Real_t const &step, vecgeom::Vector3D<Real_t> &endPosition, vecgeom::Vector3D<Real_t> &endDirection) const

Protected Functions

inline void CalculateDerived(Vector3D<double> Bvec)
template<typename Real_t>
inline bool CheckModulus(Real_t &newdirX_v, Real_t &newdirY_v, Real_t &newdirZ_v) const

Private Types

template<typename T>
using Vector3D = vecgeom::Vector3D<T>

Private Members

Vector3D<double> fUnit = {1, 0, 0}
double fBmag = 0
template<class T>
struct CudaDeleter

Public Functions

inline void operator()(T *ptr) const
template<class T>
struct CudaHostDeleter

Public Functions

inline void operator()(T *ptr) const
struct DataDeleter

Deletes the outer G4HepEmData object after first freeing all tables it owns.

FreeG4HepEmData releases both the host-side tables and any device-side mirrors embedded in the G4HepEmData object, but it does not delete the outer G4HepEmData allocation itself. This deleter performs both steps for the owned fData member. It does not touch the separately owned G4HepEmParameters copy stored in fParameters.

Public Functions

void operator()(G4HepEmData *data) const

Release the tables owned by G4HepEmData and then delete the outer object.

This is the cleanup path for the fully owned fData member. It is separate from the class destructor because std::unique_ptr needs a deleter for the deep cleanup before the outer G4HepEmData allocation itself can be deleted.

struct DeferredStep

Stored work for a returned step that is replayed later on the host.

The host collects returned GPU-step blocks during transport and replays them later in one fixed order, so the Geant4-side work stays reproducible from run to run.

Public Members

std::size_t firstGPUStep = {0}
std::size_t numGPUSteps = {0}
DeferredStepType type = {DeferredStepType::ReplayStep}
struct DeferredStepStore

Owns the deferred returned-step data drained from the integration.

Public Members

std::vector<GPUStep> gpuSteps = {}
std::vector<DeferredStep> steps = {}
struct Deleter

Public Functions

void operator()(StepReconstructionObjects *ptr)
template<class Equation_t, class T_Field, unsigned int Nvar, typename Real_t>
class DormandPrinceRK45

Public Static Functions

static void EvaluateDerivatives(const T_Field &field, const Real_t y[], int charge, Real_t dydx[])
static inline void StepWithErrorEstimate(const T_Field &field, const Real_t *yIn, const Real_t *dydx, int charge, Real_t Step, Real_t *yOut, Real_t *yerr, Real_t *next_dydx)

Public Static Attributes

static constexpr unsigned int kMethodOrder = 4
class ErrorEstimatorRK
#include <ErrorEstimatorRK.h>

Public Functions

inline ErrorEstimatorRK(float eps_rel_max, int noComponents = 6)
template<class Real_t>
Real_t EstimateSquareError(const Real_t yError[], const Real_t &hStep, const Real_t &magMomentumSq) const
inline float GetMaxRelativeError() const

Public Static Attributes

static constexpr const double tinyValue = 1.0e-30

Private Members

const float fEpsRelMax
const float fInvEpsilonRelSq
class fieldPropagatorConstBany

Public Functions

inline void stepInField (ConstFieldHelixStepper &helixAnyB, double kinE, double mass, int charge, double step, vecgeom::Vector3D< vecgeom::double > &position, vecgeom::Vector3D< vecgeom::double > &direction)
class fieldPropagatorConstBz

Public Functions

inline fieldPropagatorConstBz(double Bz)
inline ~fieldPropagatorConstBz()
double ComputeSafeLength(double momentumMag, int charge, const Vector3D &direction)
template<class Navigator = AdePTNavigator>
double ComputeStepAndNextVolume(double kinE, double mass, int charge, double physicsStep, Vector3D &position, Vector3D &direction, vecgeom::NavigationState const &current_state, vecgeom::NavigationState &new_state, long &hitsurf_index, bool &propagated, const double safety = 0.0, const int max_iteration = 100)

Private Types

using Vector3D = vecgeom::Vector3D<double>

Private Members

double BzValue
template<class Field_t, class RkDriver_t, typename Real_t, class Navigator>
class fieldPropagatorRungeKutta

Runge-Kutta propagator in magnetic field.

Template Parameters:
  • Field_t – Field type

  • RkDriver_t – Driver type

  • Real_t – Precision type

  • Navigator – Navigator type

Public Static Functions

static inline double ComputeStepAndNextVolume(Field_t const &magneticField, double kinE, double mass, int charge, double physicsStep, double safeLength, vecgeom::Vector3D<double> &position, vecgeom::Vector3D<double> &direction, vecgeom::NavigationState const &current_state, vecgeom::NavigationState &next_state, long &hitsurf_index, bool &propagated, const Real_t &safetyIn, const int max_iterations, int &iterDone, int threadId, bool &zero_first_step, bool verbose = false)

Propagating method attempting to push a particle with the step proposed by physics.

Parameters:
  • magneticField – Field map

  • kinE – Kinetc energy

  • mass – Particle mass

  • charge – Particle charge

  • physicsStep – Step to propagate with

  • safeLength – B-field safety depending only on track curvature

  • position – Particle position

  • direction – Particle direction

  • current_state[in] Current geometry state

  • next_state[out] Geometry state after propagation

  • hitsurf_index[out] Index of the hit surface (surface model only)

  • propagated[out] Checks if the step was fully propagated

  • safetyIn – Geometric isotropic safety

  • max_iterations – Maximum allowed iterations

  • iterDone[out] Number of iterations performed

  • threadId – Thread id

  • zero_first_step – Detected zero first step

  • verbose – Verbosity

Returns:

Length of the step made

static inline Real_t ComputeSafeLength(vecgeom::Vector3D<double> &momentumVec, vecgeom::Vector3D<Real_t> &BfieldVec, int charge)

Protected Static Attributes

static constexpr unsigned int fMaxTrials = 100
static constexpr unsigned int Nvar = 6
static constexpr Real_t kPush = 0.
static constexpr Real_t kDistCheckPush = kPush + Navigator::kBoundaryPush
class G4EmStandardPhysics_AdePT : public G4EmStandardPhysics

Public Functions

explicit G4EmStandardPhysics_AdePT(G4int ver = 1, const G4String &name = "")
~G4EmStandardPhysics_AdePT() override
void ConstructProcess() override

Private Members

AdePTTrackingManager *fTrackingManager = {nullptr}
AdePTConfiguration *fAdePTConfiguration = {nullptr}
class G4EmStandardPhysics_HepEm : public G4EmStandardPhysics

Public Functions

explicit G4EmStandardPhysics_HepEm(G4int ver = 1, const G4String &name = "")
~G4EmStandardPhysics_HepEm() override
void ConstructProcess() override

Private Members

G4HepEmTrackingManager *fTrackingManager = {nullptr}
class G4HepEmTrackingManagerSpecialized : public G4HepEmTrackingManager

4HepEmTrackingManagerSpecialized class: The derived class from G4HepEmTrackingManager must implement HandOverOneTrack and CheckEarlyTrackingExit to allow for an early exit of the tracking loop in the G4HepEmTrackingManager

Public Functions

G4HepEmTrackingManagerSpecialized()
~G4HepEmTrackingManagerSpecialized()
inline void SetGPURegions(const std::set<G4Region const*> &gpuRegions)
inline void SetTrackInAllRegions(bool trackInAllRegions)

Set whether AdePT should transport particles across the whole geometry.

inline bool GetTrackInAllRegions() const
void HandOverOneTrack(G4Track *aTrack) override
bool CheckEarlyTrackingExit(G4Track *track, G4EventManager *evtMgr, G4UserTrackingAction *userTrackingAction, G4TrackVector &secondaries) const override

The function checks within the TrackElectron and TrackGamma calls in G4HepEmTracking manager, if a barrier is hit.

Parameters:
  • secondaries – Secondary track container produced during tracking.

  • track – a G4track to be checked

  • evtMgr – the G4EventManager to stack secondaries

  • userTrackingAction – the UserTrackingAction

Returns:

True if tracking should exit early and hand back control

Returns:

inline void ResetFinishEventOnCPUSize(int num_threads)
inline void SetFinishEventOnCPU(int threadid, int eventid)
inline int GetFinishEventOnCPU(int threadid) const

Private Members

bool fTrackInAllRegions = false

Whether the whole geometry is a GPU region.

std::set<G4Region const*> fGPURegions = {}

List of GPU regions.

std::vector<int> fFinishEventOnCPU

Vector over worker threads for events finished on CPU.

struct GPUstateDeleter
#include <EventState.hh>

Public Functions

void operator()(GPUstate *ptr)
struct GPUStep
#include <GPUStep.hh>

Public Functions

inline bool operator<(GPUStep const &other) const

Public Members

GPUStepPoint fPreStepPoint
GPUStepPoint fPostStepPoint
double fStepLength = {0}
double fTotalEnergyDeposit = {0}
double fGlobalTime = {0.}
float fLocalTime = {0.f}
float fProperTime = {0.f}
double fPreGlobalTime = {0.}
float fTrackWeight = {1}
uint64_t fTrackID = {0}
uint64_t fParentID = {0}
short fStepLimProcessId = {-1}
int fEventId = {0}
short threadId = {-1}
unsigned short fStepCounter = {0}
bool fLastStepOfTrack = {false}
ParticleType fParticleType = {ParticleType::Electron}
unsigned char fNumSecondaries = {0}
struct GPUStepPoint
#include <GPUStep.hh>

Public Members

vecgeom::Vector3D<double> fPosition
vecgeom::Vector3D<double> fMomentumDirection
double fEKin
vecgeom::NavigationState fNavigationState = {0}
class HostCircularBuffer

Tracks occupied ranges in the pinned host buffer used for returned GPU steps.

The HostCircularBuffer manages the memory of the pinned HostBuffer in returned-step transfer. It keeps track of the used space in the sorted vector of segments fSegments. The StepProcessingThread adds segments when it submits items to the StepQueue (in fact, the memory is already allocated as soon as the copy from TransferStepsToHost is done). The StepProcessingThread can delete segments when it copies the steps to the holdoutBuffer and the G4Worker finish their work. Since both StepProcessingThread and G4Workers can change the segments, a mutex is used to lock the access to fSegments. In TransferStepsToHost, the StepProcessingThread checks whether there is enough contiguous memory in the HostCircularBuffer before the copy can start.

Public Functions

explicit HostCircularBuffer(std::size_t capacity)
bool AddSegment(std::size_t begin, std::size_t end)

Adds a segment at the provided position. Ensures segments remain sorted.

void RemoveSegment(std::size_t segmentBegin)

Removes a segment based on its starting offset.

std::size_t GetFreeContiguousSlots(std::size_t transferSize)

Returns the contiguous free space in front of the write offset (or at the beginning of the buffer in case of a wraparound).

double GetFillFractionAfterLastRequest() const

Return the fill fraction according to the existing transfer-manager bookkeeping.

This intentionally preserves the historical semantics: the value is based on fFreeContiguousSpace, which is updated by GetFreeContiguousSlots(transferSize) and includes the pending transfer size.

std::size_t GetWriteOffset() const

Private Functions

bool checkForOverlaps() const
void printSegments(const std::string &msg) const

Private Members

std::size_t fCapacity = {0}
std::size_t fWriteOffset = {0}
std::size_t fFreeContiguousSpace
std::vector<Segment> fSegments
mutable std::mutex bufferManagerMutex
struct HostTrackData

A helper struct to store the data that is stored exclusively on the CPU.

Public Functions

HostTrackData() = default
HostTrackData(HostTrackData&&) noexcept = default
HostTrackData &operator=(HostTrackData&&) noexcept = default
HostTrackData(HostTrackData const&) = delete
HostTrackData &operator=(HostTrackData const&) = delete

Public Members

int g4id = 0
int g4parentid = 0
uint64_t gpuId = 0
G4PrimaryParticle *primary = nullptr
G4VProcess *creatorProcess = nullptr
G4VUserTrackInformation *userTrackInfo = nullptr
G4LogicalVolume *logicalVolumeAtVertex = nullptr
G4ThreeVector vertexPosition
G4ThreeVector vertexMomentumDirection
G4double vertexKineticEnergy = 0.0
ParticleType particleType = {ParticleType::Electron}
class HostTrackDataMapper

Public Functions

inline ~HostTrackDataMapper()
inline void beginEvent (int eventID, size_t expectedTracks=1 '000 '000)

Call once at the start of each event, so we can clear and reserve.

inline HostTrackData &get(uint64_t gpuId) noexcept

HOT PATH: 1 hash + bucket probe, returns a reference into the table.

inline HostTrackData &create(uint64_t gpuId, bool useNewId = true)
inline HostTrackData &activateForGPU(uint64_t gpuId, int g4id, bool haveReverse) noexcept

Ensure a live HostTrackData slot exists for this GPU track. If the slot was retired (or never existed), create it and bind it to the given G4 id. Preserves reproducibility by keeping the original g4id.

Returns:

reference to the active HostTrackData.

inline bool getGPUId(int g4id, uint64_t &gpuid)

Sets the gpuid by reference and returns whether the entry already existed.

Parameters:
  • g4id – int G4 id that is checked

  • gpuid – uint64 gpu id that is returned

Returns:

whether the GPU id already existed

inline void removeTrack(uint64_t gpuId)

Call when a track (with given gpuId) is completely done:

inline void retireToCPU(uint64_t gpuId)
inline bool contains(uint64_t gpuId) const

Whether an entry exists in the GPU to Index map for the given GPU id.

Parameters:

gpuId – GPU id to be checked

Returns:

true if the value exists

Private Types

using GpuToIndexMap = std::unordered_map<uint64_t, int>

Private Functions

inline void eraseHostTrackData(GpuToIndexMap::iterator it, bool keepReverseMap, bool deleteUserTrackInfo)

Private Members

std::unordered_map<uint64_t, int> gpuToIndex
std::unordered_map<int, uint64_t> g4idToGpuId
std::vector<HostTrackData> hostDataVec
int currentGpuReturnG4ID = std::numeric_limits<int>::max()
int currentEventID = -1
template<typename MagneticField_t>
class MagneticFieldEquation

Public Functions

template<typename Real_t>
inline void EvaluateDerivativesReturnB(MagneticField_t const &magField, vecgeom::Vector3D<float> const &position, vecgeom::Vector3D<float> const &momentum, int charge, Real_t dy_ds[], vecgeom::Vector3D<float> &BfieldVal)
template<typename Real_t>
inline void EvaluateDerivativesReturnB(MagField_t const &magField, const Real_t y[], int charge, Real_t dy_ds[], vecgeom::Vector3D<float> &BfieldVec)
template<typename Real_t>
inline void EvaluateDerivativesReturnB(MagField_t const &magField, vecgeom::Vector3D<float> const &position, vecgeom::Vector3D<float> const &momentum, int charge, Real_t dy_ds[], vecgeom::Vector3D<float> &BfieldVec)

Public Static Functions

template<typename Real_t>
static inline void EvaluateDerivatives(MagneticField_t const &magField, const Real_t y[], int charge, Real_t dy_ds[])
template<typename Real_t>
static inline void EvaluateDerivativesGivenB(const Real_t y[], const Real_t Bfield[3], int charge, Real_t dy_ds[])
template<typename Real_t>
static inline void EvaluateDerivatives(const Real_t y[], vecgeom::Vector3D<Real_t> const &Bvec, int charge, Real_t dy_ds[])
template<typename Real_t>
static void EvaluateDerivativesReturnB(MagneticField_t const &magField, const Real_t y[], int charge, Real_t dy_ds[], vecgeom::Vector3D<float> &BfieldVal)
template<typename Real_t, typename Int_t>
static inline void EvaluateRhsGivenB(const Real_t y[], const Int_t &charge, const vecgeom::Vector3D<Real_t> &Bvalue, Real_t dy_ds[])

Public Static Attributes

static constexpr double gCof = copcore::units::kCLight
static constexpr int Nvar = 6

Private Static Functions

template<typename Real_t, typename Int_t>
static inline void Force(vecgeom::Vector3D<Real_t> const &momentum, Int_t const &charge, vecgeom::Vector3D<Real_t> const &Bfield, Real_t &dx_ds, Real_t &dy_ds, Real_t &dz_ds, Real_t &dpx_ds, Real_t &dpy_ds, Real_t &dpz_ds)
template<typename T>
class MParrayT : protected copcore::VariableSizeObjectInterface<MParrayT<T>, T>
#include <MParrayT.h>

A variable-size array having elements added in an atomic way.

Public Types

using value_type = T
using pointer = value_type*
using const_pointer = const value_type*
using reference = value_type&
using const_reference = const value_type&
using iterator = value_type*
using const_iterator = const value_type*
using size_t = std::size_t
using AtomicInt_t = adept::Atomic_t<int>
using Base_t = copcore::VariableSizeObjectInterface<MParrayT<T>, T>
using ArrayData_t = copcore::VariableSizeObj<T>

Public Functions

inline size_t size() const

Maximum number of elements.

inline constexpr size_t max_size() const

Maximum number of elements.

inline void clear()

Clear the content.

inline const_reference operator[](size_t index) const

Read-only index operator.

inline bool push_back(const_reference val)

Dispatch next free element, nullptr if none left.

inline bool full() const

Check if container is fully distributed.

inline const_iterator begin() const
inline const_iterator end() const
inline const_reference front() const
inline const_reference back() const
inline const_pointer data() const

Public Static Functions

static inline size_t SizeOfInstance(int capacity)

Returns the size in bytes of a BlockData object with given capacity.

Private Functions

inline ArrayData_t &GetVariableData()

Functions required by VariableSizeObjectInterface.

inline const ArrayData_t &GetVariableData() const
inline MParrayT(size_t nvalues)
inline MParrayT(MParrayT const &other)
inline MParrayT(size_t new_size, MParrayT const &other)
inline ~MParrayT()

Private Members

size_t fCapacity = {0}

Maximum number of elements.

AtomicInt_t fNbooked

Number of booked elements.

AtomicInt_t fNused

Number of used elements.

ArrayData_t fData

Data follows, has to be last.

friend Base_t
template<typename Type>
class mpmc_bounded_queue : protected copcore::VariableSizeObjectInterface<mpmc_bounded_queue<Type>, internal::Cell_t<Type>>

Class MPMC bounded queue.

Public Types

using AtomicInt_t = adept::Atomic_t<int>
using Value_t = internal::Cell_t<Type>
using Base_t = copcore::VariableSizeObjectInterface<mpmc_bounded_queue<Type>, Value_t>
using ArrayData_t = copcore::VariableSizeObj<Value_t>

Public Functions

inline int Capacity() const

Maximum number of elements.

inline void clear()

Size function.

inline int size() const

Size function.

inline bool enqueue(Type const &data)

MPMC enqueue function.

inline bool dequeue(Type &data)

MPMC dequeue function.

inline Type &operator[](const int index)

Usage as array rather than a queue.

This mode only allowed if no pending concurrent enqueue/dequeue operations. The index should be smaller than the number of stored elements. The range is only checked in debug mode.

inline Type const &operator[](const int index) const

Public Static Functions

static inline size_t SizeOfInstance(int capacity)

Returns the size in bytes of a BlockData object with given capacity.

static inline Cont *MakeCopy(const Cont &other)

< Enumerate the part of the private interface, we want to expose.

static inline Cont *MakeCopy(size_t new_size, const Cont &other)

< Enumerate the part of the private interface, we want to expose.

static inline Cont *MakeCopyAt(const Cont &other, void *addr)
static inline Cont *MakeCopyAt(size_t new_size, const Cont &other, void *addr)
template<typename ...T>
static inline Cont *MakeInstance(size_t nvalues, const T&... params)
template<typename ...T>
static inline Cont *MakeInstanceAt(size_t nvalues, void *addr, const T&... params)
static inline void ReleaseInstance(Cont *obj)
static inline constexpr size_t SizeOf(size_t nvalues)
static inline constexpr size_t SizeOfAlignAware(size_t nvalues)

Private Functions

inline ArrayData_t &GetVariableData()

Functions required by VariableSizeObjectInterface.

inline const ArrayData_t &GetVariableData() const
inline mpmc_bounded_queue(int nvalues)

MPMC bounded queue constructor.

Parameters:

nvalues – Maximum number of elements in the queue. Must be a power of two and >= 2.

inline mpmc_bounded_queue(int, mpmc_bounded_queue const&)
inline mpmc_bounded_queue(mpmc_bounded_queue const &other)

MPMC bounded queue copy constructor.

inline mpmc_bounded_queue(size_t new_size, mpmc_bounded_queue const &other)

MPMC bounded queue copy constructor with given size.

void operator=(mpmc_bounded_queue const&) = delete

Operator =.

Private Members

int const fCapacity

Capacity of the queue.

int const fMask

Mask used to navigate fast in the circular buffer.

AtomicInt_t fEnqueue

Enqueueing index.

AtomicInt_t fDequeue

Dequeueing index.

AtomicInt_t fNstored

Number of stored elements.

ArrayData_t fBuffer

Buffer of cells with atomic access to data elements.

friend Base_t
struct ParametersDeleter

Deletes the outer G4HepEmParameters object after first freeing all host/device allocations it owns.

The copied parameter block owns its fParametersPerRegion host array and the GPU mirror pointed to by fParametersPerRegion_gpu after transport upload. FreeG4HepEmParameters releases those nested allocations, while this deleter also deletes the outer G4HepEmParameters allocation.

Public Functions

void operator()(G4HepEmParameters *parameters) const

Release the copied HepEm parameters and then delete the outer object.

The copied G4HepEmParameters block is fully owned by AdePTG4HepEmState. This deep cleanup therefore releases both the host-side per-region array and the device-side mirror created during transport upload before deleting the outer G4HepEmParameters allocation itself.

class RanluxppDouble : public RanluxppEngineImpl<48>
#include <Ranluxpp.h>

Public Functions

inline RanluxppDouble(uint64_t seed = 314159265)
inline double Rndm()

Generate a double-precision random number with 48 bits of randomness.

inline double operator()()

Generate a double-precision random number (non-virtual method)

inline uint64_t IntRndm()

Generate a random integer value with 48 bits.

inline uint64_t IntRndm64()

Generate a uniformly random 64-bit integer by concatenating two 32-bit words.

inline RanluxppDouble BranchNoAdvance()

Branch a new RNG state, also advancing the current one. The caller must Advance() the branched RNG state to decorrelate the produced numbers.

inline RanluxppDouble Branch()

Branch a new RNG state, also advancing the current one.

template<int w>
class RanluxppEngineImpl
#include <Ranluxpp.h>

Public Functions

RanluxppEngineImpl() = default
inline void __attribute__((noinline)) Advance()

Produce next block of random bits.

inline uint64_t NextRandomBits()

Return the next random bits, generate a new block if necessary.

inline double NextRandomFloat()

Return a floating point number, converted from the next random bits.

inline void SetSeed(uint64_t s)

Initialize and seed the state of the generator.

inline void Skip(uint64_t n)

Skip n random numbers without generating them.

Protected Functions

inline void SaveState(uint64_t *state) const
inline void XORstate(const uint64_t *state)

Private Members

uint64_t fState[9]

RANLUX state of the generator.

unsigned fCarry

Carry bit of the RANLUX state.

int fPosition = 0

Current position in bits.

Private Static Attributes

static constexpr const uint64_t *kA = kA_2048
static constexpr int kMaxPos = 9 * 64
template<class Stepper_t, typename Real_t, typename Int_t, class Equation_t, class MagField_t>
class RkIntegrationDriver

A simple stepper treating the propagation of particles in a constant magnetic field ( not just along the z-axis&#8212; for that there is ConstBzFieldHelixStepper )

Public Functions

template<int Verbose>
inline bool Advance(vecgeom::Vector3D<double> &position, vecgeom::Vector3D<double> &momentumVec, Int_t const &chargeInt, Real_t const &length, MagField_t const &magField, Real_t dydx_next[Nvar], Real_t &hgood, unsigned int maxTrials, int cordIters)

Public Static Functions

template<int Verbose = 1>
static inline bool Advance(vecgeom::Vector3D<double> &position, vecgeom::Vector3D<double> &momentumVec, Int_t const &charge, Real_t const &step, MagField_t const &magField, Real_t dydx_next[], Real_t &hgood, unsigned int maxTrials = 5, int cordIters = 0)
static inline void PrintStep(vecgeom::Vector3D<Real_t> const &startPosition, vecgeom::Vector3D<Real_t> const &startDirection, Int_t const &charge, Real_t const &momentum, Real_t const &step, vecgeom::Vector3D<Real_t> &endPosition, vecgeom::Vector3D<Real_t> &endDirection)
static inline bool IntegrateStep(const Real_t yStart[], const Real_t dydx[], int charge, Real_t &xCurrent, Real_t htry, const MagField_t &magField, Real_t yEnd[], Real_t next_dydx[], Real_t &hnext)

Public Static Attributes

static constexpr int Nvar = 6
static constexpr Real_t fEpsilonRelativeMax = 0.001
static constexpr Real_t fMinimumStep = 0.01 * copcore::units::millimeter
static constexpr Real_t kSmall = 1.0e-30

Protected Static Functions

static inline bool CheckModulus(Real_t &newdirX_v, Real_t &newdirY_v, Real_t &newdirZ_v)
struct SecondaryInitData
#include <GPUStep.hh>

Minimal data struct that is needed along with the parent track to provide the initial track information that is sent back to the CPU.

Public Members

uint64_t trackId
vecgeom::Vector3D<double> dir
double eKin
short creatorProcessId = {-1}
ParticleType particleType = {ParticleType::Electron}
struct Segment

Public Functions

inline bool operator<(const Segment &other) const

Public Members

std::size_t begin = {0}
std::size_t end = {0}
struct Stats
#include <TransportStats.hh>

Public Members

int inFlight[GPUQueueIndex::NumSpecies]
int wdtGammasInFlight
float queueFillLevel[GPUQueueIndex::NumParticleQueues]
float slotFillLevel[GPUQueueIndex::NumSpecies]
unsigned int *perEventInFlight = {nullptr}

Updated asynchronously; sized to the G4 worker count.

unsigned int *perEventInFlightPrevious = {nullptr}

Used in transport kernels; sized to the G4 worker count.

unsigned int stepBufferOccupancy
struct StepReconstructionObjects

This struct holds temporary Geant4 objects needed to replay GPU steps on the host. These are allocated as members or using placement new to go around G4’s pool allocators, which cause a destruction order fiasco when the pools are destroyed before the objects, and then the objects’ destructors are called. This also keeps all these objects much closer in memory.

Public Functions

inline void ResetSecondaryTracks()
inline StepReconstructionObjects()

Public Members

G4NavigationHistory fPreG4NavigationHistory
G4NavigationHistory fPostG4NavigationHistory
std::aligned_storage<sizeof(G4StepPoint), alignof(G4StepPoint)>::type stepPointStorage[2]
std::aligned_storage<sizeof(G4Step), alignof(G4Step)>::type stepStorage
std::aligned_storage<sizeof(G4TouchableHistory), alignof(G4TouchableHistory)>::type toucheableHistoryStorage[2]
std::aligned_storage<sizeof(G4TouchableHandle), alignof(G4TouchableHandle)>::type toucheableHandleStorage[2]
G4Step *fG4Step = nullptr
G4TouchableHandle *fPreG4TouchableHistoryHandle = nullptr
G4TouchableHandle *fPostG4TouchableHistoryHandle = nullptr
std::aligned_storage<sizeof(G4DynamicParticle), alignof(G4DynamicParticle)>::type dynParticleStorage[3]
std::aligned_storage<sizeof(G4DynamicParticle), alignof(G4DynamicParticle)>::type secDynStorage[kMaxTotalSecondaries]
std::aligned_storage<sizeof(G4Track), alignof(G4Track)>::type trackStorage[3]
std::aligned_storage<sizeof(G4Track), alignof(G4Track)>::type secTrackStorage[kMaxTotalSecondaries]
G4Track *fElectronTrack = nullptr
G4Track *fPositronTrack = nullptr
G4Track *fGammaTrack = nullptr
G4DynamicParticle *secDyn[kMaxTotalSecondaries] = {nullptr}
G4Track *secTrk[kMaxTotalSecondaries] = {nullptr}
std::aligned_storage<sizeof(G4TrackVector), alignof(G4TrackVector)>::type secondaryStorage
G4TrackVector *fSecondaryVector = nullptr
int secEBase = 0
int secPBase = secEBase + 1
int secGBase = secPBase + 1
int secEUsed = 0
int secPUsed = 0
int secGUsed = 0

Public Static Attributes

static constexpr int kMaxSecElectrons = 1
static constexpr int kMaxSecPositrons = 1
static constexpr int kMaxSecGammas = 3
static constexpr int kMaxTotalSecondaries = kMaxSecElectrons + kMaxSecPositrons + kMaxSecGammas
static constexpr int kMaxSecondaries = 3
struct StreamlikeIdentity
#include <Portability.hh>

Allow passing a single string into a streamlike operator for device-compatible ADEPT_VALIDATE messages.

Public Functions

inline ADEPT_ATT_HOST_DEVICE operator char const*() const
template<BackendType T>
struct StreamType
#include <Global.h>

Public Types

using value_type = int

Public Static Functions

static inline void CreateStream(value_type &stream)
template<typename Real_t>
class SurfNavigator
#include <SurfNavigator.h>

Public Types

using Vector3D = vecgeom::Vector3D<double>
using SurfData = vgbrep::SurfData<Real_t>
using Real_b = typename SurfData::Real_b

Public Static Functions

static inline int LocatePointIn(int pvol_id, Vector3D const &point, vecgeom::NavigationState &path, bool top, int *exclude = nullptr)

Locates the point in the geometry volume tree.

Parameters:
  • pvol_id – Placed volume id to be checked first

  • point – Point to be checked, in the local frame of pvol

  • path – Path to a parent of pvol that must contain the point

  • top – Check first if pvol contains the point

  • exclude – Placed volume id to exclude from the search

Returns:

Index of the placed volume that contains the point

static inline double ComputeSafety(Vector3D const &globalpoint, vecgeom::NavigationState const &state, double limit = vecgeom::InfinityLength<Real_t>())

Computes the isotropic safety from the globalpoint.

Parameters:
  • globalpoint – Point in global coordinates

  • state – Path where to compute safety

  • limit – Limit to which safety should be computed accurately

Returns:

Isotropic safe distance

static inline double ComputeStepAndNextVolume(Vector3D const &globalpoint, Vector3D const &globaldir, double step_limit, vecgeom::NavigationState const &in_state, vecgeom::NavigationState &out_state, long &hitsurf_index, double push = 0)
static inline double ComputeStepAndPropagatedState(Vector3D const &globalpoint, Vector3D const &globaldir, double step_limit, vecgeom::NavigationState const &in_state, vecgeom::NavigationState &out_state, long &hitsurf_index, double push = 0)
static inline void RelocateToNextVolume(Vector3D const &globalpoint, Vector3D const &globaldir, long hitsurf_index, vecgeom::NavigationState &out_state)

Public Static Attributes

static constexpr double kBoundaryPush = 10 * vecgeom::kTolerance
struct ToDeviceBuffer
#include <TrackBuffer.hh>

Public Members

adeptint::TrackData *tracks
unsigned int maxTracks
std::atomic_uint nTrack
mutable std::shared_mutex mutex
struct TrackBuffer
#include <TrackBuffer.hh>

Buffer holding input tracks to be transported on GPU.

Public Functions

TrackBuffer(unsigned int numToDevice)
inline ToDeviceBuffer &getActiveBuffer()
inline void swapToDeviceBuffers()
inline TrackHandle createToDeviceSlot()

Create a handle with lock for tracks that go to the device. Create a shared_lock and a reference to a track.

Returns:

TrackHandle with lock and reference to track slot.

Public Members

std::array<ToDeviceBuffer, 2> toDeviceBuffer
std::atomic_short toDeviceIndex = {0}
unsigned int fNumToDevice = {0}

number of slots in the toDevice buffer

unique_ptr_cuda<adeptint::TrackData, CudaHostDeleter<adeptint::TrackData>> toDevice_host

Tracks to be transported to the device.

unique_ptr_cuda<adeptint::TrackData> toDevice_dev

toDevice buffer of tracks

struct TrackData
#include <TrackData.h>

Track data uploaded from Geant4 to AdePT for GPU transport.

Public Functions

TrackData() = default
inline TrackData(int pdg_id, uint64_t trackId, uint64_t parentId, double ene, double x, double y, double z, double dirx, double diry, double dirz, double gTime, double lTime, double pTime, float weight, unsigned short stepCounter, vecgeom::NavigationState &&state, unsigned int eventId = 0, short threadId = -1)

Public Members

vecgeom::NavigationState navState
double position[3]
double direction[3]
double eKin = {0}
double globalTime = {0}
double localTime = {0}
double properTime = {0}
float weight = {0}
int pdg = {0}
uint64_t trackId = {0}

track id (non-consecutive, reproducible)

uint64_t parentId = {0}
unsigned short stepCounter = {0}
unsigned int eventId = {0}
short threadId = {-1}
struct TrackHandle
#include <TrackBuffer.hh>

A handle to access TrackData vectors while holding a lock.

Public Members

adeptint::TrackData &track
std::shared_lock<std::shared_mutex> lock
struct TransportLoopCounters
#include <TransportStats.hh>

Host-only counters accumulating transport-loop stop/stall/flush action reasons across the full run. These are incremented on the host transport thread and printed at shutdown when verbosity >= 1.

Public Members

unsigned long long totalIterations = {0}

Total transport iterations executed.

unsigned long long leakExtractionByQueuePressure = {0}

Iterations where leak queue exceeded 50% threshold.

unsigned long long leakExtractionByEventFlush = {0}

Iterations where an event flush requested leak extraction.

unsigned long long leakExtractionBlocked = {0}

Times transport stalled waiting for in-progress extraction.

unsigned long long eventDrainedToStepFlush = {0}

Events that transitioned to RequestStepFlush (queues drained)

unsigned long long stepBufferSwaps = {0}

Total step-buffer swaps performed.

unsigned long long stepBufferSwapByOccupancy = {0}

Swaps triggered by occupancy >= half capacity.

unsigned long long stepBufferSwapByOccupancy10k = {0}

Swaps triggered by occupancy >= 10000.

unsigned long long stepBufferSwapByPressure = {0}

Swaps triggered by nextStepMightFail (overflow risk)

unsigned long long stepBufferSwapByEventFlush = {0}

Swaps triggered by event RequestStepFlush.

template<typename V>
class VariableSizeObj
#include <VariableSizeObj.h>

Public Types

using Index_t = unsigned int

Public Functions

inline VariableSizeObj(TRootIOCtor*)

Beginning address of the array values &#8212; real size: [fN].

inline VariableSizeObj(unsigned int nvalues)
inline VariableSizeObj(const VariableSizeObj &other)
inline VariableSizeObj(size_t new_size, const VariableSizeObj &other)
inline V *GetValues()
inline const V *GetValues() const
inline V &operator[](Index_t index)
inline const V &operator[](Index_t index) const
inline VariableSizeObj &operator=(const VariableSizeObj &rhs)

Public Members

bool fSelfAlloc
const Index_t fN

Record whether the memory is externally managed.

V fRealArray[1]
template<class T, BackendType backend>
class VariableSizeObjAllocator
template<class T>
class VariableSizeObjAllocator<T, BackendType::CPU>

Partial variable-size allocator specialization for the CPU backend.

Public Types

using value_type = T

Public Functions

inline VariableSizeObjAllocator(std::size_t capacity, int device = 0)
inline VariableSizeObjAllocator()
VariableSizeObjAllocator(const VariableSizeObjAllocator&) = default
template<class U>
inline VariableSizeObjAllocator(const VariableSizeObjAllocator<U, BackendType::CPU> &other)
inline bool operator==(const VariableSizeObjAllocator &other) const
inline bool operator!=(const VariableSizeObjAllocator &other) const
template<typename ...P>
inline value_type *allocate(std::size_t n, const P&... params) const
inline void deallocate(value_type *ptr, std::size_t n = 0) const
inline int GetDevice() const
inline void SetCapacity(std::size_t capacity)

Private Members

std::size_t fCapacity = {0}

Capacity of each VariableSizeObj container.

int fDeviceId = {0}

Device id.

template<typename Cont, typename V>
class VariableSizeObjectInterface
#include <VariableSizeObj.h>

Public Static Functions

template<typename ...T>
static inline Cont *MakeInstance(size_t nvalues, const T&... params)
template<typename ...T>
static inline Cont *MakeInstanceAt(size_t nvalues, void *addr, const T&... params)
static inline Cont *MakeCopy(const Cont &other)
static inline Cont *MakeCopy(size_t new_size, const Cont &other)
static inline Cont *MakeCopyAt(const Cont &other, void *addr)
static inline Cont *MakeCopyAt(size_t new_size, const Cont &other, void *addr)
static inline void ReleaseInstance(Cont *obj)
static inline constexpr size_t SizeOf(size_t nvalues)
static inline constexpr size_t SizeOfExtra(size_t)
static inline constexpr size_t SizeOfAlignAware(size_t nvalues)

Protected Functions

VariableSizeObjectInterface() = default
~VariableSizeObjectInterface() = default

Private Static Functions

static inline constexpr size_t FillUp(size_t nvalues)
static inline constexpr size_t RealFillUp(size_t nvalues)
struct VolAuxArray
#include <GeometryAuxData.hh>

Structure holding the arrays of auxiliary volume data on host and device.

Public Members

int fNumVolumes = {0}
VolAuxData *fAuxData = {nullptr}

array of auxiliary volume data on host

VolAuxData *fAuxData_dev = {nullptr}

array of auxiliary volume data on device

Public Static Functions

static inline VolAuxArray &GetInstance()
struct VolAuxData
#include <GeometryAuxData.hh>

Auxiliary logical volume data. This stores in the same structure the material-cuts couple index, the sensitive volume handler index and the flag if the region is active for AdePT.

Public Members

int fSensIndex = {-1}

index of handler for sensitive volumes (-1 means non-sensitive)

int fMCIndex = {0}

material-cut couple index in G4HepEm

int fGPUregionId = {-1}

GPU region index, corresponds to G4Region.instanceID if tracked on GPU, -1 otherwise.

struct WDTDeviceBuffers
#include <WoodcockData.hh>

Owned device buffers to manage lifetime of Woodcock tracking data.

Public Members

WDTRoot *d_roots = nullptr
WDTRegion *d_regions = nullptr
int *d_map = nullptr
struct WDTDeviceView
#include <WoodcockData.hh>

Device view pointing to the uploaded Woodcock data.

Public Members

const WDTRoot *roots

[nRoots]

const WDTRegion *regions

[nRegions] (only WDT-enabled regions)

const int *regionToWDT

[regionToWDTLen], regionId -> bucket (index into regions[]) or -1

int nRoots
int nRegions
unsigned short maxIter

maximum number of Woodcock iterations

struct WDTHostPacked
#include <WoodcockData.hh>

Compact, upload-ready representation of Woodcock data.

Public Members

std::vector<WDTRoot> roots

packed per-region contiguous

std::vector<WDTRegion> regions

one per WDT region

std::vector<int> regionToWDT

dense by regionId (size = number of G4 regions)

struct WDTHostRaw
#include <WoodcockData.hh>

Sparse host-side Woodcock collection built during geometry traversal.

Public Members

std::unordered_map<int, std::vector<int>> regionToRootIndices

regionId -> list of indices into roots

std::vector<WDTRoot> roots

one entry per root placed volume

float ekinMin = {0.f}

global Woodcock minimum kinetic energy

struct WDTRegion
#include <WoodcockData.hh>

Compact Woodcock region header.

Public Members

int offset

first index in roots[]

int count

number of roots for this region

float ekinMin

kinetic energy threshold

struct WDTRoot
#include <WoodcockData.hh>

Root Volume data of a Woodcock tracking region: Navigation index + G4HepEm material cut couple index.

Public Members

vecgeom::NavigationState root

NavState of the root placed volume.

int hepemIMC

G4HepEm mat-cut index for this root.

namespace adept

Typedefs

using MParray = MParrayT<int>
namespace internal
namespace transport

The AdePT transport service. This provides the interfaces for:

  • initializing geometry and physics on the AdePT side

  • filling the buffer with tracks to be transported on the GPU

  • driving the GPU transport worker

Typedefs

template<typename T = void, typename Deleter = CudaDeleter<T>>
using unique_ptr_cuda = std::unique_ptr<T, Deleter>

Enums

enum class EventState : unsigned char

Values:

enumerator NewTracksFromG4
enumerator G4RequestsFlush
enumerator Inject
enumerator InjectionCompleted
enumerator Transporting
enumerator WaitingForTransportToFinish
enumerator RequestStepFlush
enumerator SwappingStepBuffers
enumerator FlushingSteps
enumerator StepsFlushed
enumerator DeviceFlushed

Functions

void InitVolAuxArray(adeptint::VolAuxArray &array)
void freeCuda(void *ptr)
void freeCudaHost(void *ptr)
void freeCudaStream(void *stream)
void freeCudaEvent(void *event)
namespace detail

Functions

void setDeviceLimits(int stackLimit = 0, int heapLimit = 0)
void CopySurfaceModelToGPU()
void InitWDTOnDevice(const adeptint::WDTHostPacked&, adeptint::WDTDeviceBuffers&, unsigned short)
void UploadG4HepEmToGPU(G4HepEmData *hepEmData, G4HepEmParameters *hepEmParameters)
std::thread LaunchGPUWorker(int, int, int, adept::transport::TrackBuffer&, adept::transport::GPUstate&, std::vector<std::atomic<adept::transport::EventState>>&, std::condition_variable&, int, int, bool, bool, unsigned short, const double, bool)
std::unique_ptr<adept::transport::GPUstate, adept::transport::GPUstateDeleter> InitializeGPU(int trackCapacity, int stepCapacity, int numThreads, adept::transport::TrackBuffer &trackBuffer, double CPUCapacityFactor, double CPUCopyFraction, std::string &generalBfieldFile, const std::vector<float> &uniformBfieldValues)
void FreeGPU(std::unique_ptr<adept::transport::GPUstate, adept::transport::GPUstateDeleter>&, std::thread&, adeptint::WDTDeviceBuffers&)
namespace adept_integration

Functions

inline int GetThreadId()

Return the non-negative Geant4 worker slot used by AdePT per-thread host arrays.

In sequential Geant4, G4GetThreadId() returns MASTER_ID (-1) or SEQUENTIAL_ID (-2). AdePT uses non-negative slot indices 0..N-1, so map any negative Geant4 thread id to slot 0.

namespace AdePTGeant4Integration_detail

The AdePT-Geant4 integration layer.

  • Initialization and checking of VecGeom geometry

  • G4Step reconstruction from GPU steps

  • Processing of reconstructed steps using the Geant4 sensitive detector

namespace adeptint
namespace copcore

Enums

enum BackendType

Backend types enumeration.

Values:

enumerator CPU
enumerator CUDA
enumerator HIP

Functions

static inline void error_check(int, const char*, int)

CUDA error checking.

static inline int get_num_SMs()

Get number of SMs on the current device.

template<typename Backend>
const char *BackendName(Backend const &backend)

Getting the backend name in templated constructs.

namespace units

Variables

static constexpr double kAvogadro = 6.02214076e+23 / mole
static constexpr double kCLight = 2.99792458e+8 * m / s
static constexpr double kCLightSquare = kCLight * kCLight
static constexpr double kHPlanck = 6.62607015e-34 * joule * s
static constexpr double kHBarPlanck = kHPlanck / kTwoPi
static constexpr double kHBarPlanckCLight = kHBarPlanck * kCLight
static constexpr double kHBarPlanckCLightSquare = kHBarPlanckCLight * kHBarPlanckCLight
static constexpr double kElectronCharge = -eplus
static constexpr double kUnitChargeSquare = eplus * eplus
static constexpr double kElectronMassC2 = 0.510998910 * MeV
static constexpr double kInvElectronMassC2 = 1.0 / kElectronMassC2
static constexpr double kProtonMassC2 = 938.272013 * MeV
static constexpr double kNeutronMassC2 = 939.56536 * MeV
static constexpr double kAtomicMassUnitC2 = 931.494028 * MeV
static constexpr double kAtomicMassUnit = kAtomicMassUnitC2 / kCLightSquare
static constexpr double kMu0 = 4 * kPi * 1.e-7 * henry / m
static constexpr double kEpsilon0 = 1. / (kCLightSquare * kMu0)
static constexpr double kEMCoupling = kUnitChargeSquare / (4 * kPi * kEpsilon0)
static constexpr double kFineStructConst = kEMCoupling / kHBarPlanckCLight
static constexpr double kClassicElectronRadius = kEMCoupling / kElectronMassC2
static constexpr double kRedElectronComptonWLenght = kHBarPlanckCLight / kElectronMassC2
static constexpr double kBohrRadius = kRedElectronComptonWLenght / kFineStructConst
static constexpr double kBoltzmann = 8.617333e-11 * MeV / kelvin
static constexpr double kSTPTemperature = 273.15 * kelvin
static constexpr double kNTPTemperature = 293.15 * kelvin
static constexpr double kSTPPressure = 1. * atmosphere
static constexpr double kGasThreshold = 10. * mg / cm3
static constexpr double kUniverseMeanDensity = 1.e-25 * g / cm3
static constexpr double kPi = 3.14159265358979323846
static constexpr double kTwoPi = 2. * kPi
static constexpr double kHalfPi = kPi / 2.
static constexpr double kPiSquare = kPi * kPi
static constexpr double millimeter = 1.
static constexpr double millimeter2 = millimeter * millimeter
static constexpr double millimeter3 = millimeter * millimeter * millimeter
static constexpr double centimeter = 10. * millimeter
static constexpr double centimeter2 = centimeter * centimeter
static constexpr double centimeter3 = centimeter * centimeter * centimeter
static constexpr double meter = 1000. * millimeter
static constexpr double meter2 = meter * meter
static constexpr double meter3 = meter * meter * meter
static constexpr double kilometer = 1000. * meter
static constexpr double kilometer2 = kilometer * kilometer
static constexpr double kilometer3 = kilometer * kilometer * kilometer
static constexpr double parsec = 3.0856775807e+16 * meter
static constexpr double micrometer = 1.e-6 * meter
static constexpr double nanometer = 1.e-9 * meter
static constexpr double angstrom = 1.e-10 * meter
static constexpr double fermi = 1.e-15 * meter
static constexpr double barn = 1.e-28 * meter2
static constexpr double millibarn = 1.e-3 * barn
static constexpr double microbarn = 1.e-6 * barn
static constexpr double nanobarn = 1.e-9 * barn
static constexpr double picobarn = 1.e-12 * barn
static constexpr double nm = nanometer
static constexpr double um = micrometer
static constexpr double mm = millimeter
static constexpr double mm2 = millimeter2
static constexpr double mm3 = millimeter3
static constexpr double cm = centimeter
static constexpr double cm2 = centimeter2
static constexpr double cm3 = centimeter3
static constexpr double liter = 1.e+3 * cm3
static constexpr double L = liter
static constexpr double dL = 1.e-1 * liter
static constexpr double cL = 1.e-2 * liter
static constexpr double mL = 1.e-3 * liter
static constexpr double m = meter
static constexpr double m2 = meter2
static constexpr double m3 = meter3
static constexpr double km = kilometer
static constexpr double km2 = kilometer2
static constexpr double km3 = kilometer3
static constexpr double pc = parsec
static constexpr double radian = 1.
static constexpr double milliradian = 1.e-3 * radian
static constexpr double degree = (kPi / 180.0) * radian
static constexpr double steradian = 1.
static constexpr double rad = radian
static constexpr double mrad = milliradian
static constexpr double sr = steradian
static constexpr double deg = degree
static constexpr double nanosecond = 1.
static constexpr double second = 1.e+9 * nanosecond
static constexpr double millisecond = 1.e-3 * second
static constexpr double microsecond = 1.e-6 * second
static constexpr double picosecond = 1.e-12 * second
static constexpr double hertz = 1. / second
static constexpr double kilohertz = 1.e+3 * hertz
static constexpr double megahertz = 1.e+6 * hertz
static constexpr double s = second
static constexpr double ms = millisecond
static constexpr double ns = nanosecond
static constexpr double eplus = 1.
static constexpr double e_SI = 1.602176634e-19
static constexpr double coulomb = eplus / e_SI
static constexpr double megaelectronvolt = 1.
static constexpr double electronvolt = 1.e-6 * megaelectronvolt
static constexpr double kiloelectronvolt = 1.e-3 * megaelectronvolt
static constexpr double gigaelectronvolt = 1.e+3 * megaelectronvolt
static constexpr double teraelectronvolt = 1.e+6 * megaelectronvolt
static constexpr double petaelectronvolt = 1.e+9 * megaelectronvolt
static constexpr double joule = electronvolt / e_SI
static constexpr double eV = electronvolt
static constexpr double keV = kiloelectronvolt
static constexpr double MeV = megaelectronvolt
static constexpr double GeV = gigaelectronvolt
static constexpr double TeV = teraelectronvolt
static constexpr double PeV = petaelectronvolt
static constexpr double kilogram = joule * second * second / (meter * meter)
static constexpr double gram = 1.e-3 * kilogram
static constexpr double milligram = 1.e-3 * gram
static constexpr double kg = kilogram
static constexpr double g = gram
static constexpr double mg = milligram
static constexpr double watt = joule / second
static constexpr double newton = joule / meter
static constexpr double pascal = newton / m2
static constexpr double bar = 1.e+5 * pascal
static constexpr double atmosphere = 101325 * pascal
static constexpr double ampere = coulomb / second
static constexpr double milliampere = 1.e-3 * ampere
static constexpr double microampere = 1.e-6 * ampere
static constexpr double nanoampere = 1.e-9 * ampere
static constexpr double megavolt = megaelectronvolt / eplus
static constexpr double gigavolt = 1.e+3 * megavolt
static constexpr double kilovolt = 1.e-3 * megavolt
static constexpr double volt = 1.e-6 * megavolt
static constexpr double ohm = volt / ampere
static constexpr double farad = coulomb / volt
static constexpr double millifarad = 1.e-3 * farad
static constexpr double microfarad = 1.e-6 * farad
static constexpr double nanofarad = 1.e-9 * farad
static constexpr double picofarad = 1.e-12 * farad
static constexpr double weber = volt * second
static constexpr double tesla = volt * second / meter2
static constexpr double gauss = 1.e-4 * tesla
static constexpr double kilogauss = 1.e-1 * tesla
static constexpr double henry = weber / ampere
static constexpr double kelvin = 1.
static constexpr double mole = 1.
static constexpr double becquerel = 1. / second
static constexpr double curie = 3.7e+10 * becquerel
static constexpr double kilobecquerel = 1.e+3 * becquerel
static constexpr double megabecquerel = 1.e+6 * becquerel
static constexpr double gigabecquerel = 1.e+9 * becquerel
static constexpr double millicurie = 1.e-3 * curie
static constexpr double microcurie = 1.e-6 * curie
static constexpr double Bq = becquerel
static constexpr double kBq = kilobecquerel
static constexpr double MBq = megabecquerel
static constexpr double GBq = gigabecquerel
static constexpr double Ci = curie
static constexpr double mCi = millicurie
static constexpr double uCi = microcurie
static constexpr double gray = joule / kilogram
static constexpr double kilogray = 1.e+3 * gray
static constexpr double milligray = 1.e-3 * gray
static constexpr double microgray = 1.e-6 * gray
static constexpr double candela = 1.
static constexpr double lumen = candela * steradian
static constexpr double lux = lumen / meter2
static constexpr double perCent = 0.01
static constexpr double perThousand = 0.001
static constexpr double perMillion = 0.000001
namespace COPCORE_IMPL

Typedefs

using LoopNavigator = vecgeom::LoopNavigator

Functions

double getDoubleOpt(char **begin, char **end, const std::string &option, double defaultval)
double getIntOpt(char **begin, char **end, const std::string &option, int defaultval)
bool getBoolOpt(char **begin, char **end, const std::string &option, bool defaultval)
std::string getStringOpt(char **begin, char **end, const std::string &option, const std::string &defaultval)
namespace fieldConstants

Variables

static constexpr double kB2C = -0.299792458 * (copcore::units::GeV / (copcore::units::tesla * copcore::units::meter))
static constexpr float deltaIntersection = 0.001 * copcore::units::millimeter
static constexpr float deltaChord = 0.25 * copcore::units::millimeter
namespace portability

Functions

std::logic_error make_debug_error(char const *which, char const *condition, char const *file, int line)
std::runtime_error make_runtime_error(char const *which, char const *what, char const *condition, char const *file, int line)
inline ADEPT_ATT_HOST_DEVICE void unreachable ()

Invoke undefined behavior.

namespace detail

Functions

inline ADEPT_ATT_HOST_DEVICE char const * operator<< (StreamlikeIdentity const &, char const *s)
namespace PrintFieldVectors

Functions

template<typename Real_t>
static inline void Print3vec(const Real_t y[3], const char *name)
template<typename Real_t>
static inline void PrintScalar(const Real_t yVal, const char *name[10])
template<typename Real_t>
static inline void PrintSixvec(const Real_t y[NvarPrint], Real_t originalMomentum = -1.0, bool newline = true)
template<typename Real_t>
static inline void PrintSixvec(vecgeom::Vector3D<Real_t> const &position, vecgeom::Vector3D<Real_t> const &momentum, Real_t origMomentumMag = -1.0)
template<typename Real_t>
static inline void PrintLineSixvec(const Real_t y[NvarPrint])
template<typename Real_t>
static inline void PrintSixvecAndDyDx(const Real_t y[NvarPrint], int charge, const Real_t Bfield[3], Real_t const dydx[NvarPrint])
template<typename Real_t>
static inline void PrintLineSixvecDyDx(const Real_t y[NvarPrint], int charge, const Real_t Bfield[3], Real_t const dydx[NvarPrint])

Variables

constexpr int prec = 6
constexpr int nd = prec + 4
file AdePTConfiguration.hh
file AdePTConfiguration.hh
#include <cstdint>
#include <memory>
#include <string>
#include <vector>
file AdePTTrackingManager.hh
#include “G4VTrackingManager.hh”
#include “G4RegionStore.hh”
#include “globals.hh”
#include <memory>
#include <vector>
file AdePTTrackingManager.hh
file AdePTConfigurationMessenger.hh
#include “G4UImessenger.hh”
#include “globals.hh”
#include <memory>
file G4EmStandardPhysics_AdePT.hh
#include “G4EmStandardPhysics.hh”
file G4EmStandardPhysics_AdePT.hh
file G4EmStandardPhysics_HepEm.hh
#include “G4EmStandardPhysics.hh”
file G4EmStandardPhysics_HepEm.hh
file AdePTGeometryBridge.hh
#include <VecGeom/base/Config.h>
#include <string>
#include <vector>
file AdePTGeant4Integration.hh
#include <G4EventManager.hh>
#include <G4Event.hh>
#include <G4Track.hh>
#include <span>
#include <vector>
file HostTrackDataMapper.hh
#include “G4VUserTrackInformation.hh”
#include “G4VProcess.hh”
#include <VecGeom/navigation/NavigationState.h>
#include <unordered_map>
#include <cstdint>
#include <algorithm>
#include <limits>
#include <memory>
file AdePTThreadId.hh
#include “G4Threading.hh”
file G4HepEmTrackingManagerSpecialized.hh
#include “G4HepEmTrackingManager.hh”
file AdePTTransport.cuh
file AdePTTransport.hh
#include <VecGeom/base/Config.h>
#include <VecGeom/management/CudaManager.h>
#include <condition_variable>
#include <mutex>
#include <memory>
#include <span>
#include <thread>
#include <unordered_map>
#include <optional>
#include “AdePTTransport.icc
file AdePTTransport.icc
file AdePTTransportConfig.hh
#include <cstdint>
#include <string>
file Atomic.h
#include <cassert>
#include <atomic>

Portable atomic structure.

Author

Andrei Gheata (andrei.gheata@cern.ch)

file BlockData.h

Templated data structure storing a contiguous block of data.

Author

Andrei Gheata (andrei.gheata@cern.ch)

file MParray.h

Multi-producer array for handling track indices.

Author

Andrei Gheata (andrei.gheata@cern.ch)

file MParrayT.h

Multi-producer array of arbitrary type that can be filled concurrently.

Author

Andrei Gheata (andrei.gheata@cern.ch)

file mpmc_bounded_queue.h
#include <stdint.h>
#include <cassert>

Implementation of CUDA-aware multi-producer, multi-consumer bounded queue.

Author

Andrei Gheata based on http://www.1024cores.net/home/lock-free-algorithms/queues/bounded-mpmc-queue

file SlotManager.cuh
file VariableSizeObj.h
#include <string.h>
#include <cassert>

Array/vector of items with all information stored contiguously in memory.

Resizing operations require explicit new creation and copy. This is decomposed in the content part (VariableSizeObj) and the interface part (VariableSizeObjectInterface) so that the variable size part can be placed at the end of the object. Derived types must decorate the implementation methods with appropriate host device qualifiers (or use the macros pre-defined in VecCore) for usage on GPU.

Author

Philippe Canal (pcanal@fnal.gov). Copied/adapted for VecGeom library (andrei.gheata@cern.ch).

Use in derived classes

The derivation should be done on VariableSizeObjectInterface like:

class MyVariableSizeType : protected copcore::VariableSizeObjectInterface<MyVariableSizeType, V>

where V is the variable data type. The derived class must have the LAST data member (say fData) of type copcore::VariableSizeObj<V>. It has to declare the base class as friend to its private methods, and has to mandatory implement the following private methods:

using ArrayData_t = copcore::VariableSizeObj<Value_t>; ArrayData_t &GetVariableData() { return fData; } const ArrayData_t &GetVariableData() const { return fData; }

but also a private constructor initializing the capacity of the variable size array:

MyVariableSizeType(size_t nvalues, …param) : … , fData(nvalues)

The derived class should have no public constructors, but should enumerate the public base class static construction methods:

using Base_t = copcore::VariableSizeObjectInterface<MyVariableSizeType, V>; using Base_t::MakeInstance; using Base_t::MakeInstanceAt;

These methods internally invoke the private constructor, their invokation being dependent on the extra parameters taken by this constructor and explained in the Usage section.

In case the derived type should be copiable, it has to expose the Base_t copy interfaces:

 using Base_t::MakeCopy;
 using Base_t::MakeCopyAt;
and implement private copy constructors:

MyVariableSizeType(MyVariableSizeType const &other) : MyVariableSizeType(other.fData.fN, other) {…} // Copy and resize MyVariableSizeType(size_t new_size, MyVariableSizeType const &other) : fData(new_size, other.fData) {…}

Derived types must implement a public static method returning the size of objects for a given number of elements, which must be called mandatory by the user before allocating memory for the object. This may be just an alias to the SizeOf helper method from the base class as below, or has to be extended to handle additional custom requirements:

static size_t SizeOfInstance(int capacity) { return Base_t::SizeOf(capacity); }
Important note: the derived type should not contain any pointer data members that have to be dynamically allocated, or if it does then it has to be supported in the constructors and in the SizeOfInstance method. The base helper class is taylored for data member types for which the size can be computed at compilation time (i.e. sizeof will return the actual object size).

Composition using data members of types derived from VariableSizeObject is allowed. In such case it is mandatory to implement the private method:

static size_t SizeOfExtra(int capacity) { …; return extra_size; }

The extra size has to be computed by summing the sizes of all VariableSizeObjectInterface-derived data members as returned by their method SizeOfAlignedAware(capacity).

Implementation examples for derived VariableSizeObjectInterface types can be seen here

Usage of VariableSize objects

The most important use case the helper was designed for is allocation of variable size objects in a pre-allocated memory buffer. To do this, one needs to reteive the container size for the desired array capacity via:

size_t buffer_size = MyVariableSizeType::SizeOfInstance(size_t capacity);

or:

size_t buffer_size = MyVariableSizeType::SizeOfAlignAware(size_t capacity);

The latter version has to be used for allocating multiple objects of type MyVariableSizeType in a contiguous buffer.

char *buffer = nullptr; host_or_device_malloc_method(buffer, [N *] buffer_size);

To allocate a single object, use:

auto obj = MyVariableSizeType::MakeInstanceAt(capacity, buffer, params…);

where params… are the private constructor parameters. In case of multiple contiguous objects, the padding in the buffer must be equal to buffer_size.

In case standard allocation on the heap is needed, one has to use the method:

auto obj = MyVariableSizeType::MakeInstance(capacity, param…);

The cleanup of VariableSize objects should not be done via delete, but rather using:

MyVariableSizeType::ReleaseInstance(obj);

This will call the actual delete in case the object is allocated on the heap, or do nothing for the buffer adoption case, when the buffer has to be de-allocated separately.

Copying VariableSize objects is possible in case the copy constructors are defined as described above. The copying can be done together with resizing the container, using:

auto copy = MyVariableSizeType::MakeCopy([new_size], MyVariableSizeType const &source); // using malloc auto copy = MyVariableSizeType::MakeCopy([new_size], MyVariableSizeType const &source, void *addr); // adopting memory

If omitting new_size above, the copy will have the same size as the source.

file VariableSizeObjAllocator.h
#include <cstddef>
#include <stdexcept>
#include <iostream>

Backend-aware allocator for variable size objects.

A standard allocator providing allocate/deallocate interface for VariableSizeObj objects. Specializations are provided for CPU, CUDA and TODO HIP.

Author

Andrei Gheata (andrei.gheata@cern.ch).

file AdePTG4HepEmState.hh
#include <memory>
file GeometryAuxData.hh
file AdePTSteppingAction.cuh
file AdePTSteppingActionSelector.cuh
file electrons.cuh
file electrons_split.cuh
file gammas.cuh
file gammas_split.cuh
file gammasWDT.cuh
file gammasWDT_split.cuh
file WoodcockHelper.cuh
file CompareResponses.h
#include <VecGeom/base/Vector3D.h>

Functions

static bool CompareResponseVector3D(int id, vecgeom::Vector3D<double> const &originalVec, vecgeom::Vector3D<double> const &baselineVec, vecgeom::Vector3D<double> const &resultVec, const char *vecName, double thresholdRel)
static void ReportSameMoveVector3D(int id, vecgeom::Vector3D<double> const &originalVec, vecgeom::Vector3D<double> const &resultVec, const char *vecName)
file ConstBzFieldStepper.h
#include “VecGeom/base/Global.h”
#include “fieldConstants.h
file ConstFieldHelixStepper.h
#include <VecGeom/base/Vector3D.h>
#include “fieldConstants.h
file DormandPrinceRK45.h
file ErrorEstimatorRK.h
file fieldConstants.h
file fieldPropagatorConstBany.h
#include <VecGeom/base/Vector3D.h>
file fieldPropagatorConstBz.h
#include <VecGeom/base/Vector3D.h>
#include “fieldConstants.h
file fieldPropagatorRungeKutta.h
#include <VecGeom/base/Vector3D.h>
#include “fieldConstants.h
#include <VecGeom/navigation/NavigationState.h>
file GeneralMagneticField.cuh
file MagneticFieldEquation.h
#include <VecGeom/base/Vector3D.h>
#include <iostream>
file PrintFieldVectors.h
#include <iostream>
#include <iomanip>
#include “VecGeom/base/Vector3D.h”

Variables

static constexpr int NvarPrint = 6
file RkIntegrationDriver.h
#include “ErrorEstimatorRK.h
#include <VecCore/VecMath.h>
#include <VecGeom/base/Vector3D.h>
#include <iostream>
#include <iomanip>
file UniformMagneticField.cuh
file AdePTNavigator.h

Top switch for different navigators.

Typedefs

using AdePTNavigator = BVHNavigator
file BVHNavigator.h
#include <VecGeom/base/Global.h>
#include <VecGeom/navigation/BVHNavigator.h>

Navigation methods for geometry.

Typedefs

using BVHNavigator = vecgeom::BVHNavigator
file LoopNavigator.h
#include <VecGeom/base/Global.h>
#include <VecGeom/navigation/LoopNavigator.h>

Navigation methods for geometry.

file SurfNavigator.h
#include <VecGeom/base/Global.h>
#include <VecGeom/base/Vector3D.h>
#include <VecGeom/navigation/NavigationState.h>
#include <VecGeom/surfaces/BVHSurfNavigator.h>

Navigation methods using the surface model.

file ParticleManager.cuh
file ParticleQueues.cuh
file TrackBuffer.hh
#include <array>
#include <atomic>
#include <chrono>
#include <iostream>
#include <mutex>
#include <shared_mutex>
#include <thread>
file G4HepEmRandomEngineDeviceImpl.hh
#include <G4HepEmRandomEngine.hh>
file Ranluxpp.h
#include “ranluxpp/mulmod.h
#include “ranluxpp/ranlux_lcg.h
#include <cassert>
#include <cstdint>
file helpers.h
#include <cstdint>

Functions

static inline uint64_t add_overflow(uint64_t a, uint64_t b, unsigned &overflow)

Compute a + b and set overflow accordingly.

static inline uint64_t add_carry(uint64_t a, uint64_t b, unsigned &carry)

Compute a + b and increment carry if there was an overflow.

static inline uint64_t sub_overflow(uint64_t a, uint64_t b, unsigned &overflow)

Compute a - b and set overflow accordingly.

static inline uint64_t sub_carry(uint64_t a, uint64_t b, unsigned &carry)

Compute a - b and increment carry if there was an overflow.

static inline int64_t compute_r(const uint64_t *upper, uint64_t *r)

Update r = r - (t1 + t2) + (t3 + t2) * b ** 10.

This function also yields cbar = floor(r / m) as its return value (int64_t because the value can be -1). With an initial value of r = t0, this can be used for computing the remainder after division by m (see the function mod_m in mulmod.h). The function to_ranlux passes r = 0 and uses only the return value to obtain the decimal expansion after divison by m.

file mulmod.h
#include “helpers.h
#include <cstdint>

Functions

static void multiply9x9(const uint64_t *in1, const uint64_t *in2, uint64_t *out)

Multiply two 576 bit numbers, stored as 9 numbers of 64 bits each.

Parameters:
  • in1[in] first factor as 9 numbers of 64 bits each

  • in2[in] second factor as 9 numbers of 64 bits each

  • out[out] result with 18 numbers of 64 bits each

static void mod_m(const uint64_t *mul, uint64_t *out)

Compute a value congruent to mul modulo m less than 2 ** 576.

\( m = 2^{576} - 2^{240} + 1 \)

The result in out is guaranteed to be smaller than the modulus.

Parameters:
  • mul[in] product from multiply9x9 with 18 numbers of 64 bits each

  • out[out] result with 9 numbers of 64 bits each

static void mulmod(const uint64_t *in1, uint64_t *inout)

Combine multiply9x9 and mod_m with internal temporary storage.

The result in inout is guaranteed to be smaller than the modulus.

Parameters:
  • in1[in] first factor with 9 numbers of 64 bits each

  • inout[inout] second factor and also the output of the same size

static void powermod(const uint64_t *base, uint64_t *res, uint64_t n)

Compute base to the n modulo m.

The arguments base and res may point to the same location.

Parameters:
  • base[in] with 9 numbers of 64 bits each

  • res[out] output with 9 numbers of 64 bits each

  • n[in] exponent

file ranlux_lcg.h
#include “helpers.h
#include <cstdint>

Functions

static void to_lcg(const uint64_t *ranlux, unsigned c, uint64_t *lcg)

Convert RANLUX numbers to an LCG state.

\( m = 2^{576} - 2^{240} + 1 \)

Parameters:
  • ranlux[in] the RANLUX numbers as 576 bits

  • lcg[out] the 576 bits of the LCG state, smaller than m

  • c[in] the carry bit of the RANLUX state

static void to_ranlux(const uint64_t *lcg, uint64_t *ranlux, unsigned &c_out)

Convert an LCG state to RANLUX numbers.

\( m = 2^{576} - 2^{240} + 1 \)

Parameters:
  • lcg[in] the 576 bits of the LCG state, must be smaller than m

  • ranlux[out] the RANLUX numbers as 576 bits

  • c_out[out] the carry bit of the RANLUX state

file DeviceGlobals.cuh
file EventState.hh
file GPUState.cuh
file TransportStats.hh
file DeviceStepBuffer.cuh
file GPUStep.hh
#include “VecGeom/navigation/NavigationState.h”

Functions

inline void Copy3DVector(vecgeom::Vector3D<double> const &source, vecgeom::Vector3D<double> &destination)

Utility function to copy a 3D vector, used for filling the Step Points.

inline void FillGPUStep(GPUStep &aGPUStep, uint64_t aTrackID, uint64_t aParentID, short aStepLimProcessId, ParticleType aParticleType, double aStepLength, double aTotalEnergyDeposit, float aTrackWeight, vecgeom::NavigationState const &aPreState, vecgeom::Vector3D<double> const &aPrePosition, vecgeom::Vector3D<double> const &aPreMomentumDirection, double aPreEKin, vecgeom::NavigationState const &aPostState, vecgeom::Vector3D<double> const &aPostPosition, vecgeom::Vector3D<double> const &aPostMomentumDirection, double aPostEKin, double aGlobalTime, float aLocalTime, float aProperTime, double aPreGlobalTime, unsigned int eventID, short threadID, bool isLastStep, unsigned short stepCounter, unsigned char aNumSecondaries)

Fill the provided GPU step with the given data.

Variables

constexpr short kAdePTTransportationProcess = 10

AdePT-specific step-limiting process ids stored in GPUStep::fStepLimProcessId.

constexpr short kAdePTOutOfGPURegionProcess = 11

Returned step for a track that leaves the GPU region and continues on the CPU.

constexpr short kAdePTFinishOnCPUProcess = 12

Returned step for a track that is intentionally finished on the CPU.

file GPUStepRecording.cuh
file GPUStepTransferManager.cuh
file HostCircularBuffer.hh
#include <cstddef>
#include <mutex>
#include <string>
#include <vector>
file AdePTPrecision.hh

Typedefs

using rk_integration_t = double
file ArgParser.h
#include <algorithm>
#include <iostream>

Defines

OPTION_INT(name, defaultval)
OPTION_DOUBLE(name, defaultval)
OPTION_BOOL(name, defaultval)
OPTION_STRING(name, defaultval)
file Global.h
#include <cstdio>
#include <type_traits>

CopCore global macros and types.

Defines

COPCORE_IMPL
COPCORE_CUDA_CHECK(err)
COPCORE_EXCEPTION(message)

Trigger a runtime error depending on the backend.

COPCORE_REQUIRES(...)

Macro to template-specialize on a specific compile-time requirement.

COPCORE_CALLABLE_FUNC(FUNC)

macro to declare device callable functions usable in executors

COPCORE_CALLABLE_DECLARE(HVAR, FUNC)

macro to pass callable function to executors

COPCORE_CALLABLE_IN_NAMESPACE_DECLARE(HVAR, NAMESPACE, FUNC)
file Macros.h

CopCore global macros.

Compiling mixtures of C/C++/CUDA, we run into the cases of compiling for the host and the device. This involves both marking functions as accessible from the host/device/both, as well as (potential) different compile paths when compiling for the host vs for the device. The macros in this file assist this markup process. Only C/C++ and NVidia CUDA are currently supported, but other “backends” such as HIP can be added as time goes on.

AdePT follows LHCb’s Allen project in using macros for host/device/inline functions that match the CUDA keywords.

Defines

DISABLE_PEDANTIC_WARNINGS

Macro documentation notes.

COPCORE_CUDA_COMPILER Defined when using a CUDA compiler (effectively a proxy for __CUDACC__).

COPCORE_DEVICE_COMPILATION Defined when compiling CUDA code for device (effectively a proxy for __CUDA_ARCH__).

__host__ Marks a function as callable from host code. For non-CUDA compilation this macro expands to nothing.

__device__ Marks a function as callable from device code. For non-CUDA compilation this macro expands to nothing.

__forceinline__ Marks a function to always inline. If not provided by CUDA headers, it expands to inline __attribute__((always_inline)) on host compilers.

ENABLE_PEDANTIC_WARNINGS
file PhysicalConstants.h

Physical constants in internal units.

This is porting the corresponding file from gitlab.cern.ch/GeantV, adopting CopCore namespace. Original file was taken from CLHEP and changed according to GeantV. Authors of the original version: M. Maire

Author

M Novak, A Ribon

Date

december 2015

The basic units are :

  • centimeter (centimeter)

  • second (second)

  • giga electron volt (GeV)

  • positron charge (eplus)

  • degree Kelvin (kelvin)

  • the amount of substance (mole)

  • luminous intensity (candela)

  • radian (radian)

  • steradian (steradian)

You can add your own constants. // CLHEP 17.07.20 use PDG 2019 values

file Portability.hh
#include <stdexcept>
#include <sstream>

Defines

ADEPT_ATT_HOST
ADEPT_ATT_DEVICE
ADEPT_ATT_HOST_DEVICE
ADEPT_DEBUG
ADEPT_DEVICE_DEBUG
ADEPT_DEVICE_PLATFORM

API prefix token for the device offloading type.

ADEPT_DEVICE_PLATFORM_UPPER_STR
ADEPT_DEVICE_API_SYMBOL(TOK)

Add a prefix “hip” or “cuda” to a code token.

ADEPT_UNLIKELY(COND)

Mark the result of this condition to be “unlikely”.

This asks the compiler to move the section of code to a “cold” part of the instructions, improving instruction locality. It should be used primarily for error checking conditions.

ADEPT_DEBUG_THROW_(MSG, WHICH)
ADEPT_DEBUG_FAIL(MSG, WHICH)

Throw a debug assertion regardless of the ADEPT_DEBUG setting. This is used internally but is also useful for catching subtle programming errors in downstream code.

ADEPT_DEBUG_ASSERT_(COND, WHICH)
ADEPT_NOASSERT_(COND)
ADEPT_ASSERT(COND)

Internal debug assertion macro. This replaces standard assert usage.

ADEPT_ASSERT_UNREACHABLE()

Throw an assertion if the code point is reached. When debug assertions are turned off, this changes to a compiler hint that improves optimization (and may force the code to exit uncermoniously if the point is encountered, rather than continuing on with undefined behavior).

ADEPT_RUNTIME_THROW(WHICH, WHAT, COND)
ADEPT_VALIDATE(COND, MSG)

Always-on runtime assertion macro. This can check user input and input data consistency, and will raise std::runtime_error on failure with a descriptive error message that is streamed as the second argument. If used in __device__ -annotated code, the second argument must be a single C string.

An always-on debug-type assertion without a detailed message can be constructed by omitting the stream (but leaving the comma):

ADEPT_VALIDATE(file_stream,);

ADEPT_NOT_CONFIGURED(WHAT)

Assert if the code point is reached because an optional feature is disabled. This generally should be used for the constructors of dummy class definitions in, e.g., Foo.nocuda.cc:

Foo::Foo()
{
    ADEPT_NOT_CONFIGURED("CUDA");
}

ADEPT_ERROR_CHECK(RESULT, STMT)
ADEPT_DEVICE_API_CALL(STMT)

Safely and portably dispatch a CUDA/HIP API call.

Example:

ADEPT_DEVICE_API_CALL(Malloc(&ptr_gpu, 100 * sizeof(float)));
ADEPT_DEVICE_API_CALL(DeviceSynchronize());
file ResourceManagement.cuh
file ResourceManagement.hh
#include <memory>
file ScopedTransportProfiler.hh
file SystemOfUnits.h

System of units.

This is porting the corresponding file from gitlab.cern.ch/GeantV, adopting CopCore namespace. Original file was taken from CLHEP and changed according to GeantV. Authors of the original version: M. Maire, S. Giani

Author

M Novak, A Ribon

Date

december 2015

The basic units are :

  • centimeter (centimeter)

  • second (second)

  • giga electron volt (GeV)

  • positron charge (eplus)

  • degree Kelvin (kelvin)

  • the amount of substance (mole)

  • luminous intensity (candela)

  • radian (radian)

  • steradian (steradian)

It’s a non exhaustive list of derived and pratical units (i.e. mostly the SI units).

The value of \( \pi \) (kPi) is defined here as it is needed for radian to degree conversion. The other physical constants are defined in the header file : Geant/PhysicalConstants.h

You can add your own units.

file ParticleTypes.hh

Enums

enum class ParticleType : char

Strongly-typed enum for the three EM particle species tracked by AdePT.

Used in GPUStep, SecondaryInitData, HostTrackData, and anywhere a step or track needs to carry its particle species as data.

Note: The numeric values (0, 1, 2) intentionally match GPUQueueIndex (transport/queues/ParticleQueues.cuh), which is used as an array index for GPU queue/state arrays. Both enums represent the same three species, but GPUQueueIndex also carries NumSpecies / GammaWDT / NumParticleQueues sentinels for loop bounds and array sizing. Use ParticleType for physics data (steps, tracks); use GPUQueueIndex::{Electron,Positron,Gamma} for GPU state array indexing. A static_assert in ParticleQueues.cuh guarantees the values stay in sync.

Values:

enumerator Electron
enumerator Positron
enumerator Gamma
file Track.cuh
file TrackData.h
#include “VecGeom/navigation/NavigationState.h”
#include <cstdint>
#include <utility>

Typedefs

using MParrayTracks = adept::MParrayT<adeptint::TrackData>
file TrackDebug.cuh
file WoodcockData.hh
#include <VecGeom/navigation/NavigationState.h>
#include <unordered_map>
#include <vector>
file AdePTConfiguration.cc
#include <memory>
file AdePTConfigurationMessenger.cc
#include “G4UIdirectory.hh”
#include “G4UIcmdWithAnInteger.hh”
#include “G4UIcmdWithAString.hh”
#include “G4UIcmdWithABool.hh”
#include “G4UIcmdWithADouble.hh”
#include “G4Tokenizer.hh”
file G4IntegrationRdcFacade.cu
file AdePTGeometryBridge.cpp
#include <VecGeom/management/GeoManager.h>
#include <VecGeom/navigation/NavigationState.h>
#include <G4HepEmData.hh>
#include <G4Exception.hh>
#include <G4HepEmMatCutData.hh>
#include <G4LogicalVolume.hh>
#include <G4MaterialCutsCouple.hh>
#include <G4PVReplica.hh>
#include <G4RegionStore.hh>
#include <G4ReplicaNavigation.hh>
#include <G4SystemOfUnits.hh>
#include <G4TransportationManager.hh>
#include <G4VG.hh>
#include <G4VSensitiveDetector.hh>
#include <G4ios.hh>
#include <algorithm>
#include <cmath>
#include <functional>
#include <stdexcept>
file G4EmStandardPhysics_AdePT.cc
#include <G4HepEmConfig.hh>
#include “G4Electron.hh”
#include “G4Gamma.hh”
#include “G4Positron.hh”
file G4EmStandardPhysics_HepEm.cc
#include “G4Electron.hh”
#include “G4Gamma.hh”
#include “G4HepEmTrackingManager.hh”
#include “G4Positron.hh”
file AdePTGeant4Integration.cpp
#include <VecGeom/navigation/NavigationState.h>
#include <G4ios.hh>
#include <G4SystemOfUnits.hh>
#include <G4TransportationManager.hh>
#include <G4VSensitiveDetector.hh>
#include <G4UniformMagField.hh>
#include <G4FieldManager.hh>
#include <G4TouchableHandle.hh>
#include <G4StepStatus.hh>
#include <G4HepEmNoProcess.hh>
#include <G4HepEmConfig.hh>
#include <G4HepEmParameters.hh>
#include “G4Electron.hh”
#include “G4Positron.hh”
#include “G4Gamma.hh”
#include <new>
#include <type_traits>
file AdePTTrackingManager.cc
#include “G4Threading.hh”
#include “G4Track.hh”
#include “G4EventManager.hh”
#include “G4Event.hh”
#include “G4RunManager.hh”
#include “G4MTRunManager.hh”
#include “G4TransportationManager.hh”
#include “G4EmParameters.hh”
#include “G4Electron.hh”
#include “G4Gamma.hh”
#include “G4Positron.hh”
#include <VecGeom/management/GeoManager.h>
#include <algorithm>
#include <cassert>
file G4HepEmTrackingManagerSpecialized.cc
#include “G4Electron.hh”
#include “G4Gamma.hh”
#include “G4Positron.hh”
file AdePTTransport.cc
#include <VecGeom/management/BVHManager.h>
#include <VecGeom/management/GeoManager.h>
#include <G4HepEmData.hh>
#include <G4HepEmParameters.hh>
#include <cassert>
#include <chrono>
#include <iostream>
#include <stdexcept>
file AdePTTransport.cu
file AdePTG4HepEmState.cpp
#include <G4HepEmConfig.hh>
#include <G4HepEmData.hh>
#include <G4HepEmElectronInit.hh>
#include <G4HepEmGammaInit.hh>
#include <G4HepEmMaterialInit.hh>
#include <G4HepEmMatCutData.hh>
#include <G4HepEmParameters.hh>
#include <G4ios.hh>
#include <algorithm>
#include <stdexcept>
file HostCircularBuffer.cpp
#include <algorithm>
#include <cassert>
#include <iostream>
#include <iterator>

Defines

RESET
BOLD_RED
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/g4integration/config
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/config
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration/config
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/containers
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/core
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration/cuda
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/g4hepem
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/transport/g4hepem
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/g4integration
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/g4integration/geometry
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/geometry
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration/geometry
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/integration
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/kernels
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/magneticfield
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/navigation
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration/physics_constructors
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/queues
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/random
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/random/ranluxpp
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/g4integration/returned_steps
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration/returned_steps
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/state
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/steps
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/transport/steps
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/support
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/g4integration/tracking_managers
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/g4integration/tracking_managers
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/tracks
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/src/transport
dir /home/docs/checkouts/readthedocs.org/user_builds/adept-project/checkouts/latest/include/AdePT/transport/woodcock