Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Unify common part of RecordHit #341

Merged
merged 2 commits into from
Feb 3, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
14 changes: 6 additions & 8 deletions include/AdePT/core/AdePTScoringTemplate.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -18,14 +18,12 @@ void FreeGPU(Scoring *scoring, Scoring *scoring_dev);

template <typename Scoring>
__device__ void RecordHit(Scoring *scoring_dev, int aParentID, char aParticleType, double aStepLength,
double aTotalEnergyDeposit, vecgeom::NavigationState const *aPreState,
vecgeom::Vector3D<Precision> const *aPrePosition,
vecgeom::Vector3D<Precision> const *aPreMomentumDirection,
vecgeom::Vector3D<Precision> const *aPrePolarization, double aPreEKin, double aPreCharge,
vecgeom::NavigationState const *aPostState, vecgeom::Vector3D<Precision> const *aPostPosition,
vecgeom::Vector3D<Precision> const *aPostMomentumDirection,
vecgeom::Vector3D<Precision> const *aPostPolarization, double aPostEKin, double aPostCharge,
unsigned int eventId, short threadId);
double aTotalEnergyDeposit, vecgeom::NavigationState const &aPreState,
vecgeom::Vector3D<Precision> const &aPrePosition,
vecgeom::Vector3D<Precision> const &aPreMomentumDirection, double aPreEKin, double aPreCharge,
vecgeom::NavigationState const &aPostState, vecgeom::Vector3D<Precision> const &aPostPosition,
vecgeom::Vector3D<Precision> const &aPostMomentumDirection, double aPostEKin,
double aPostCharge, unsigned int eventId, short threadId);

template <typename Scoring>
__device__ void AccountProduced(Scoring *scoring_dev, int num_ele, int num_pos, int num_gam);
Expand Down
3 changes: 0 additions & 3 deletions include/AdePT/core/AsyncAdePTTransport.icc
Original file line number Diff line number Diff line change
Expand Up @@ -322,9 +322,6 @@ void AsyncAdePTTransport<IntegrationLayer>::Flush(G4int threadId, G4int eventId)

if (tracks.empty()) {
async_adept_impl::FlushScoring(fScoring[threadId]);

// TODO: Most likely this needs to go to AsyncAdePTTransport in a function that
// we will forward-declare here (Otherwise AdePTScoring is undefined)
fEventStates[threadId].store(EventState::ScoringRetrieved, std::memory_order_release);
}

Expand Down
46 changes: 10 additions & 36 deletions include/AdePT/core/HostScoringImpl.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -53,15 +53,6 @@ __device__ __forceinline__ GPUHit *GetNextFreeHit(HostScoring *hostScoring_dev)
return &hostScoring_dev->fGPUHitsBuffer_dev[aHitIndex];
}

/// @brief Utility function to copy a 3D vector, used for filling the Step Points
__device__ __forceinline__ void Copy3DVector(vecgeom::Vector3D<Precision> const *source,
vecgeom::Vector3D<Precision> *destination)
{
destination->x() = source->x();
destination->y() = source->y();
destination->z() = source->z();
}

/// @brief Copy the hits buffer to the host
void CopyHitsToHost(HostScoring &hostScoring, HostScoring *hostScoring_dev, cudaStream_t &stream)
{
Expand Down Expand Up @@ -159,37 +150,20 @@ void FreeGPU(HostScoring *hostScoring, HostScoring *hostScoring_dev)
/// @brief Record a hit
template <>
__device__ void RecordHit(HostScoring *hostScoring_dev, int aParentID, char aParticleType, double aStepLength,
double aTotalEnergyDeposit, vecgeom::NavigationState const *aPreState,
vecgeom::Vector3D<Precision> const *aPrePosition,
vecgeom::Vector3D<Precision> const *aPreMomentumDirection,
vecgeom::Vector3D<Precision> const *aPrePolarization, double aPreEKin, double aPreCharge,
vecgeom::NavigationState const *aPostState, vecgeom::Vector3D<Precision> const *aPostPosition,
vecgeom::Vector3D<Precision> const *aPostMomentumDirection,
vecgeom::Vector3D<Precision> const *aPostPolarization, double aPostEKin, double aPostCharge,
unsigned int, short)
double aTotalEnergyDeposit, vecgeom::NavigationState const &aPreState,
vecgeom::Vector3D<Precision> const &aPrePosition,
vecgeom::Vector3D<Precision> const &aPreMomentumDirection, double aPreEKin, double aPreCharge,
vecgeom::NavigationState const &aPostState, vecgeom::Vector3D<Precision> const &aPostPosition,
vecgeom::Vector3D<Precision> const &aPostMomentumDirection, double aPostEKin,
double aPostCharge, unsigned int, short)
{
// Acquire a hit slot
GPUHit *aGPUHit = GetNextFreeHit(hostScoring_dev);
GPUHit &aGPUHit = *GetNextFreeHit(hostScoring_dev);

// Fill the required data
aGPUHit->fParentID = aParentID;
aGPUHit->fParticleType = aParticleType;
aGPUHit->fStepLength = aStepLength;
aGPUHit->fTotalEnergyDeposit = aTotalEnergyDeposit;
// Pre step point
aGPUHit->fPreStepPoint.fNavigationState = *aPreState;
Copy3DVector(aPrePosition, &(aGPUHit->fPreStepPoint.fPosition));
Copy3DVector(aPreMomentumDirection, &(aGPUHit->fPreStepPoint.fMomentumDirection));
// Copy3DVector(aPrePolarization, aGPUHit.fPreStepPoint.fPolarization);
aGPUHit->fPreStepPoint.fEKin = aPreEKin;
aGPUHit->fPreStepPoint.fCharge = aPreCharge;
// Post step point
aGPUHit->fPostStepPoint.fNavigationState = *aPostState;
Copy3DVector(aPostPosition, &(aGPUHit->fPostStepPoint.fPosition));
Copy3DVector(aPostMomentumDirection, &(aGPUHit->fPostStepPoint.fMomentumDirection));
// Copy3DVector(aPostPolarization, aGPUHit.fPostStepPoint.fPolarization);
aGPUHit->fPostStepPoint.fEKin = aPostEKin;
aGPUHit->fPostStepPoint.fCharge = aPostCharge;
FillHit(aGPUHit, aParentID, aParticleType, aStepLength, aTotalEnergyDeposit, aPreState, aPrePosition,
aPreMomentumDirection, aPreEKin, aPreCharge, aPostState, aPostPosition, aPostMomentumDirection, aPostEKin,
aPostCharge, 0, 0);
}

/// @brief Account for the number of produced secondaries
Expand Down
46 changes: 9 additions & 37 deletions include/AdePT/core/PerEventScoringImpl.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -238,51 +238,23 @@ void PerEventScoring::CopyToHost(cudaStream_t cudaStream)

namespace adept_scoring {

/// @brief Utility function to copy a 3D vector, used for filling the Step Points
__device__ __forceinline__ void Copy3DVector(vecgeom::Vector3D<Precision> const *source,
vecgeom::Vector3D<Precision> *destination)
{
destination->x() = source->x();
destination->y() = source->y();
destination->z() = source->z();
}

/// @brief Record a hit
template <>
__device__ void RecordHit(AsyncAdePT::PerEventScoring * /*scoring*/, int aParentID, char aParticleType,
double aStepLength, double aTotalEnergyDeposit, vecgeom::NavigationState const *aPreState,
vecgeom::Vector3D<Precision> const *aPrePosition,
vecgeom::Vector3D<Precision> const *aPreMomentumDirection,
vecgeom::Vector3D<Precision> const * /*aPrePolarization*/, double aPreEKin, double aPreCharge,
vecgeom::NavigationState const *aPostState, vecgeom::Vector3D<Precision> const *aPostPosition,
vecgeom::Vector3D<Precision> const *aPostMomentumDirection,
vecgeom::Vector3D<Precision> const * /*aPostPolarization*/, double aPostEKin,
double aStepLength, double aTotalEnergyDeposit, vecgeom::NavigationState const &aPreState,
vecgeom::Vector3D<Precision> const &aPrePosition,
vecgeom::Vector3D<Precision> const &aPreMomentumDirection, double aPreEKin, double aPreCharge,
vecgeom::NavigationState const &aPostState, vecgeom::Vector3D<Precision> const &aPostPosition,
vecgeom::Vector3D<Precision> const &aPostMomentumDirection, double aPostEKin,
double aPostCharge, unsigned int eventID, short threadID)
{
// Acquire a hit slot
GPUHit &aGPUHit = AsyncAdePT::gHitScoringBuffer_dev.GetNextSlot();
aGPUHit.fEventId = eventID;
aGPUHit.threadId = threadID;
GPUHit &aGPUHit = AsyncAdePT::gHitScoringBuffer_dev.GetNextSlot();

// Fill the required data
aGPUHit.fParentID = aParentID;
aGPUHit.fParticleType = aParticleType;
aGPUHit.fStepLength = aStepLength;
aGPUHit.fTotalEnergyDeposit = aTotalEnergyDeposit;
// Pre step point
aGPUHit.fPreStepPoint.fNavigationState = *aPreState;
Copy3DVector(aPrePosition, &(aGPUHit.fPreStepPoint.fPosition));
Copy3DVector(aPreMomentumDirection, &(aGPUHit.fPreStepPoint.fMomentumDirection));
// Copy3DVector(aPrePolarization, aGPUHit.fPreStepPoint.fPolarization);
aGPUHit.fPreStepPoint.fEKin = aPreEKin;
aGPUHit.fPreStepPoint.fCharge = aPreCharge;
// Post step point
aGPUHit.fPostStepPoint.fNavigationState = *aPostState;
Copy3DVector(aPostPosition, &(aGPUHit.fPostStepPoint.fPosition));
Copy3DVector(aPostMomentumDirection, &(aGPUHit.fPostStepPoint.fMomentumDirection));
// Copy3DVector(aPostPolarization, aGPUHit.fPostStepPoint.fPolarization);
aGPUHit.fPostStepPoint.fEKin = aPostEKin;
aGPUHit.fPostStepPoint.fCharge = aPostCharge;
FillHit(aGPUHit, aParentID, aParticleType, aStepLength, aTotalEnergyDeposit, aPreState, aPrePosition,
aPreMomentumDirection, aPreEKin, aPreCharge, aPostState, aPostPosition, aPostMomentumDirection, aPostEKin,
aPostCharge, eventID, threadID);
}

/// @brief Account for the number of produced secondaries
Expand Down
43 changes: 43 additions & 0 deletions include/AdePT/core/ScoringCommons.hh
Original file line number Diff line number Diff line change
@@ -1,6 +1,8 @@
// SPDX-FileCopyrightText: 2024 CERN
// SPDX-License-Identifier: Apache-2.0

// This file contains elements that can be shared between several scoring implementations

#ifndef SCORING_COMMONS_HH
#define SCORING_COMMONS_HH

Expand Down Expand Up @@ -55,4 +57,45 @@ struct GlobalCounters {
}
};

/// @brief Utility function to copy a 3D vector, used for filling the Step Points
__device__ __forceinline__ void Copy3DVector(vecgeom::Vector3D<Precision> const &source,
vecgeom::Vector3D<Precision> &destination)
{
destination.x() = source.x();
destination.y() = source.y();
destination.z() = source.z();
};

/// @brief Fill the provided hit with the given data
__device__ __forceinline__ void FillHit(GPUHit &aGPUHit, int aParentID, char aParticleType, double aStepLength,
double aTotalEnergyDeposit, vecgeom::NavigationState const &aPreState,
vecgeom::Vector3D<Precision> const &aPrePosition,
vecgeom::Vector3D<Precision> const &aPreMomentumDirection, double aPreEKin,
double aPreCharge, vecgeom::NavigationState const &aPostState,
vecgeom::Vector3D<Precision> const &aPostPosition,
vecgeom::Vector3D<Precision> const &aPostMomentumDirection, double aPostEKin,
double aPostCharge, unsigned int eventID, short threadID)
{
aGPUHit.fEventId = eventID;
aGPUHit.threadId = threadID;

// Fill the required data
aGPUHit.fParentID = aParentID;
aGPUHit.fParticleType = aParticleType;
aGPUHit.fStepLength = aStepLength;
aGPUHit.fTotalEnergyDeposit = aTotalEnergyDeposit;
// Pre step point
aGPUHit.fPreStepPoint.fNavigationState = aPreState;
Copy3DVector(aPrePosition, aGPUHit.fPreStepPoint.fPosition);
Copy3DVector(aPreMomentumDirection, aGPUHit.fPreStepPoint.fMomentumDirection);
aGPUHit.fPreStepPoint.fEKin = aPreEKin;
aGPUHit.fPreStepPoint.fCharge = aPreCharge;
// Post step point
aGPUHit.fPostStepPoint.fNavigationState = aPostState;
Copy3DVector(aPostPosition, aGPUHit.fPostStepPoint.fPosition);
Copy3DVector(aPostMomentumDirection, aGPUHit.fPostStepPoint.fMomentumDirection);
aGPUHit.fPostStepPoint.fEKin = aPostEKin;
aGPUHit.fPostStepPoint.fCharge = aPostCharge;
};

#endif
2 changes: 1 addition & 1 deletion include/AdePT/integration/AdePTGeant4Integration.hh
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@
#define ADEPTGEANT4_INTEGRATION_H

#include <AdePT/core/CommonStruct.h>
#include <AdePT/core/HostScoringStruct.cuh>
#include <AdePT/core/ScoringCommons.hh>

#include <G4EventManager.hh>
#include <G4Event.hh>
Expand Down
14 changes: 6 additions & 8 deletions include/AdePT/kernels/electrons.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -600,16 +600,14 @@ static __device__ __forceinline__ void TransportElectrons(adept::TrackManager<Tr
static_cast<char>(IsElectron ? 0 : 1), // Particle type
elTrack.GetPStepLength(), // Step length
energyDeposit, // Total Edep
&navState, // Pre-step point navstate
&preStepPos, // Pre-step point position
&preStepDir, // Pre-step point momentum direction
nullptr, // Pre-step point polarization
navState, // Pre-step point navstate
preStepPos, // Pre-step point position
preStepDir, // Pre-step point momentum direction
preStepEnergy, // Pre-step point kinetic energy
IsElectron ? -1 : 1, // Pre-step point charge
&nextState, // Post-step point navstate
&pos, // Post-step point position
&dir, // Post-step point momentum direction
nullptr, // Post-step point polarization
nextState, // Post-step point navstate
pos, // Post-step point position
dir, // Post-step point momentum direction
eKin, // Post-step point kinetic energy
IsElectron ? -1 : 1, // Post-step point charge
currentTrack.eventId, currentTrack.threadId); // eventID and threadID (not needed here)
Expand Down
Loading