LDMX Software
Public Member Functions | Static Public Attributes | Private Attributes | List of all members
simcore::EcalSD Class Reference

ECal sensitive detector that uses an EcalHexReadout to create the hits. More...

#include <EcalSD.h>

Public Member Functions

 EcalSD (const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &p)
 Class constructor.
 
virtual ~EcalSD ()=default
 Class destructor.
 
virtual bool isSensDet (G4LogicalVolume *vol) const override
 Should the input volume be consider apart of this sensitive detector?
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override
 Process steps to create hits.
 
virtual void saveHits (framework::Event &event) override
 Add our hits to the event bus.
 
virtual void OnFinishedEvent () override
 Clear the map of hits we have accumulated.
 
- Public Member Functions inherited from simcore::SensitiveDetector
 SensitiveDetector (const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &parameters)
 Constructor.
 
virtual ~SensitiveDetector ()=default
 Destructor.
 
virtual void EndOfEvent (G4HCofThisEvent *) override
 This is Geant4's handle to tell us the event is ending.
 

Static Public Attributes

static const std::string COLLECTION_NAME = "EcalSimHits"
 Name of output collection of hits.
 

Private Attributes

std::map< ldmx::EcalID, ldmx::SimCalorimeterHithits_
 map of hits to add to the event (will be squashed)
 
bool enableHitContribs_
 enable hit contribs
 
bool compressHitContribs_
 compress hit contribs
 

Additional Inherited Members

- Public Types inherited from simcore::SensitiveDetector
using Factory = ::simcore::Factory< SensitiveDetector, SensitiveDetector *, const std::string &, simcore::ConditionsInterface &, const framework::config::Parameters & >
 The SD Factory.
 
- Protected Member Functions inherited from simcore::SensitiveDetector
template<class T >
const T & getCondition (const std::string &condition_name)
 Record the configuration of this detector into the run header.
 
bool isGeantino (const G4Step *step) const
 Check if the passed step is a step of a geantino.
 
const TrackMapgetTrackMap () const
 Get a handle to the current track map.
 

Detailed Description

ECal sensitive detector that uses an EcalHexReadout to create the hits.

Definition at line 30 of file EcalSD.h.

Constructor & Destructor Documentation

◆ EcalSD()

simcore::EcalSD::EcalSD ( const std::string &  name,
simcore::ConditionsInterface ci,
const framework::config::Parameters p 
)

Class constructor.

Parameters
nameThe name of the sensitive detector.
theCollectionNameThe name of the hits collection.
subDetIDThe subdetector ID.

Definition at line 19 of file EcalSD.cxx.

21 : SensitiveDetector(name, ci, p) {
22 enableHitContribs_ = p.getParameter<bool>("enableHitContribs");
23 compressHitContribs_ = p.getParameter<bool>("compressHitContribs");
24}
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
bool compressHitContribs_
compress hit contribs
Definition EcalSD.h:85
bool enableHitContribs_
enable hit contribs
Definition EcalSD.h:83
SensitiveDetector(const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &parameters)
Constructor.

References compressHitContribs_, enableHitContribs_, and framework::config::Parameters::getParameter().

Member Function Documentation

◆ isSensDet()

virtual bool simcore::EcalSD::isSensDet ( G4LogicalVolume *  vol) const
inlineoverridevirtual

Should the input volume be consider apart of this sensitive detector?

Note
Dependent on names defined in GDML!

Implements simcore::SensitiveDetector.

Definition at line 54 of file EcalSD.h.

54 {
55 auto region = vol->GetRegion();
56 if (region and region->GetName().contains("CalorimeterRegion")) {
57 return vol->GetName().contains("Si");
58 }
59 return false;
60 }

◆ OnFinishedEvent()

virtual void simcore::EcalSD::OnFinishedEvent ( )
inlineoverridevirtual

Clear the map of hits we have accumulated.

Implements simcore::SensitiveDetector.

Definition at line 77 of file EcalSD.h.

77{ hits_.clear(); }
std::map< ldmx::EcalID, ldmx::SimCalorimeterHit > hits_
map of hits to add to the event (will be squashed)
Definition EcalSD.h:81

References hits_.

◆ ProcessHits()

G4bool simcore::EcalSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)
overridevirtual

Process steps to create hits.

Parameters
aStepThe step information.
ROhistThe readout history.

DEBUG this printout is helpful when developing the GDML and/or EcalGeometry since Geant4 will probe where exactly the GDML sensitive volumes are std::cout << "(" << position[0] << ", " << position[1] << ", " << position[2] << ") " << cpynum << " -> layer " << layerNumber << " module " << module_position << std::endl;

convert position to center of cell position

This is the behavior that has been done in the past, although it is completely redundant with the ID information already deduced. It would probably help us more if we persisted the actual simulated position of the hit rather than the cell center; however, that is up for more discussion.

Implements simcore::SensitiveDetector.

Definition at line 26 of file EcalSD.cxx.

26 {
27 static const int layer_depth = 2; // index depends on GDML implementation
28 const auto& geometry = getCondition<ldmx::EcalGeometry>(
29 ldmx::EcalGeometry::CONDITIONS_OBJECT_NAME);
30
31 // Get the edep from the step.
32 G4double edep = aStep->GetTotalEnergyDeposit();
33
34 // Skip steps with no energy dep which come from non-Geantino particles.
35 if (edep == 0.0 and not isGeantino(aStep)) {
36 if (verboseLevel > 2) {
37 G4cout << "CalorimeterSD skipping step with zero edep." << G4endl
38 << G4endl;
39 }
40 return false;
41 }
42
43 // Compute the hit position
44 G4StepPoint* prePoint = aStep->GetPreStepPoint();
45 G4StepPoint* postPoint = aStep->GetPostStepPoint();
46 G4ThreeVector position =
47 0.5 * (prePoint->GetPosition() + postPoint->GetPosition());
48
49 // Create the ID for the hit.
50 int cpynum{0}; // Initialize cpynum to 0
51
52 auto preStepPoint = aStep->GetPreStepPoint();
53 if (preStepPoint) {
54 auto touchableHandle = preStepPoint->GetTouchableHandle();
55 if (touchableHandle) {
56 auto history = touchableHandle->GetHistory();
57 if (history) {
58 auto volume = history->GetVolume(layer_depth);
59 if (volume) {
60 cpynum = volume->GetCopyNo();
61 }
62 }
63 }
64 }
65 int layerNumber;
66 layerNumber = cpynum / 7;
67 int module_position = cpynum % 7;
79 // fastest, but need to trust module number between GDML and EcalGeometry
80 // match
81 ldmx::EcalID id =
82 geometry.getID(position[0], position[1], layerNumber, module_position);
83
84 // medium, only need to trust z-layer positions in GDML and EcalGeometry match
85 // helpful for debugging any issues where transverse position is not
86 // matching between the GDML and EcalGeometry
87 // ldmx::EcalID id = geometry.getID(position[0], position[1], layerNumber);
88
89 // slowest, completely rely on EcalGeometry
90 // this is helpful for validating the EcalGeometry implementation and
91 // configuration since this will be called with any hit position that
92 // is inside of the configured SD volumes from Geant4's point of view
93 // ldmx::EcalID id = geometry.getID(position[0], position[1], position[2]);
94
95 if (hits_.find(id) == hits_.end()) {
96 // hit in empty cell
97 auto& hit = hits_[id];
98 hit.setID(id.raw());
108 auto [x, y, z] = geometry.getPosition(id);
109 hit.setPosition(x, y, z);
110 }
111
112 auto& hit = hits_[id];
113
114 // hit variables
115 auto track = aStep->GetTrack();
116 auto time = track->GetGlobalTime();
117 auto track_id = track->GetTrackID();
118 auto pdg = track->GetParticleDefinition()->GetPDGEncoding();
119
120 if (enableHitContribs_) {
121 int contrib_i = hit.findContribIndex(track_id, pdg);
122 if (compressHitContribs_ and contrib_i != -1) {
123 hit.updateContrib(contrib_i, edep, time);
124 } else {
125 hit.addContrib(getTrackMap().findIncident(track_id), track_id, pdg, edep,
126 time);
127 }
128 } else {
129 // no hit contribs and hit already exists
130 hit.setEdep(hit.getEdep() + edep);
131 if (time < hit.getTime() or hit.getTime() == 0) {
132 hit.setTime(time);
133 }
134 }
135
136 return true;
137}
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
Definition EcalID.h:20
bool isGeantino(const G4Step *step) const
Check if the passed step is a step of a geantino.
const TrackMap & getTrackMap() const
Get a handle to the current track map.
int findIncident(int trackID) const
Find a trajectory's nearest parent that is incident on the calorimeter region.
Definition TrackMap.cxx:35

References compressHitContribs_, enableHitContribs_, simcore::TrackMap::findIncident(), simcore::SensitiveDetector::getTrackMap(), hits_, and simcore::SensitiveDetector::isGeantino().

◆ saveHits()

void simcore::EcalSD::saveHits ( framework::Event event)
overridevirtual

Add our hits to the event bus.

Implements simcore::SensitiveDetector.

Definition at line 139 of file EcalSD.cxx.

139 {
140 // squash hits into list
141 std::vector<ldmx::SimCalorimeterHit> hits;
142 hits.reserve(hits_.size());
143 for (const auto& [id, hit] : hits_) hits.push_back(hit);
144 event.add(COLLECTION_NAME, hits);
145}
static const std::string COLLECTION_NAME
Name of output collection of hits.
Definition EcalSD.h:33

References COLLECTION_NAME, and hits_.

Member Data Documentation

◆ COLLECTION_NAME

const std::string simcore::EcalSD::COLLECTION_NAME = "EcalSimHits"
static

Name of output collection of hits.

Definition at line 33 of file EcalSD.h.

Referenced by saveHits().

◆ compressHitContribs_

bool simcore::EcalSD::compressHitContribs_
private

compress hit contribs

Definition at line 85 of file EcalSD.h.

Referenced by EcalSD(), and ProcessHits().

◆ enableHitContribs_

bool simcore::EcalSD::enableHitContribs_
private

enable hit contribs

Definition at line 83 of file EcalSD.h.

Referenced by EcalSD(), and ProcessHits().

◆ hits_

std::map<ldmx::EcalID, ldmx::SimCalorimeterHit> simcore::EcalSD::hits_
private

map of hits to add to the event (will be squashed)

Definition at line 81 of file EcalSD.h.

Referenced by OnFinishedEvent(), ProcessHits(), and saveHits().


The documentation for this class was generated from the following files: