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

SD plugins: adapt to new const G4step interface in dd4hep 1.21 #255

Merged
merged 2 commits into from
Jul 27, 2022
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
9 changes: 8 additions & 1 deletion plugins/CaloPreShowerSDAction.cpp
Original file line number Diff line number Diff line change
@@ -1,9 +1,16 @@
#include "DD4hep/Version.h"
#include "DDG4/Geant4SensDetAction.inl"
#include "DDG4/Geant4EventAction.h"
#include "DDG4/Geant4Mapping.h"
#include "G4OpticalPhoton.hh"
#include "G4VProcess.hh"

#if DD4HEP_VERSION_GE(1, 21)
#define GEANT4_CONST_STEP const
#else
#define GEANT4_CONST_STEP
#endif

/// Namespace for the AIDA detector description toolkit
namespace dd4hep {

Expand Down Expand Up @@ -58,7 +65,7 @@ namespace dd4hep {
}

/// Method for generating hit(s) using the information of G4Step object.
template <> bool Geant4SensitiveAction<CalorimeterWithPreShowerLayer>::process(G4Step* step,G4TouchableHistory*) {
template <> bool Geant4SensitiveAction<CalorimeterWithPreShowerLayer>::process(G4Step GEANT4_CONST_STEP * step,G4TouchableHistory*) {
typedef CalorimeterWithPreShowerLayer::Hit Hit;
Geant4StepHandler h(step);
HitContribution contrib = Hit::extractContribution(step);
Expand Down
21 changes: 14 additions & 7 deletions plugins/TPCSDAction.cpp
Original file line number Diff line number Diff line change
@@ -1,9 +1,16 @@
#include "DD4hep/Version.h"
#include "DDG4/Geant4SensDetAction.inl"
#include "DDG4/Geant4EventAction.h"
#include "DDG4/Geant4Mapping.h"
#include "G4OpticalPhoton.hh"
#include "G4VProcess.hh"

#if DD4HEP_VERSION_GE(1, 21)
#define GEANT4_CONST_STEP const
#else
#define GEANT4_CONST_STEP
#endif

/// Namespace for the AIDA detector description toolkit
namespace dd4hep {

Expand Down Expand Up @@ -84,18 +91,18 @@ namespace dd4hep {


/// return the layer number of the volume (either pre or post-position )
int getCopyNumber(G4Step* s, bool usePostPos ){
int getCopyNumber(G4Step GEANT4_CONST_STEP * step, bool usePostPos ){

int cellID = this->volID( s , usePostPos) ;
int cellID = this->volID( step , usePostPos) ;

return this->layerField->value( cellID ) ;
}


/// Returns the volumeID of sensitive volume corresponding to the step (either pre or post-position )
long long int volID( G4Step* s, bool usePostPos=false ) {
long long int volID( G4Step GEANT4_CONST_STEP * step, bool usePostPos=false ) {

Geant4StepHandler h(s);
Geant4StepHandler h(step);

Geant4VolumeManager volMgr = Geant4Mapping::instance().volumeManager();

Expand Down Expand Up @@ -129,7 +136,7 @@ namespace dd4hep {


/// Method for generating hit(s) using the information of G4Step object.
G4bool process(G4Step* step, G4TouchableHistory* ) {
G4bool process(G4Step GEANT4_CONST_STEP * step, G4TouchableHistory* ) {


fHitCollection = sensitive->collection(0) ;
Expand Down Expand Up @@ -450,7 +457,7 @@ namespace dd4hep {
ResetCumulativeVariables();
}

void CumulateLowPtStep(G4Step *step)
void CumulateLowPtStep(G4Step GEANT4_CONST_STEP * step)
{

const G4ThreeVector meanPosition = (step->GetPreStepPoint()->GetPosition() + step->GetPostStepPoint()->GetPosition()) / 2;
Expand Down Expand Up @@ -501,7 +508,7 @@ namespace dd4hep {

/// Method for generating hit(s) using the information of G4Step object.
template <> G4bool
Geant4SensitiveAction<TPCSDData>::process(G4Step* step, G4TouchableHistory* history) {
Geant4SensitiveAction<TPCSDData>::process(G4Step GEANT4_CONST_STEP * step, G4TouchableHistory* history) {
return m_userData.process(step, history);
}

Expand Down