LDMX Software
HcalSD.cxx
1#include "SimCore/SDs/HcalSD.h"
2
3/*~~~~~~~~~~~~~~*/
4/* DetDescr */
5/*~~~~~~~~~~~~~~*/
7#include "DetDescr/HcalID.h"
8
9// STL
10#include <iostream>
11
12// Geant4
13#include "G4Box.hh"
14#include "G4ParticleTypes.hh"
15#include "G4Step.hh"
16#include "G4StepPoint.hh"
17
18namespace simcore {
19
20const std::string HcalSD::COLLECTION_NAME = "HcalSimHits";
21
22HcalSD::HcalSD(const std::string& name, simcore::ConditionsInterface& ci,
24 : SensitiveDetector(name, ci, p), birksc1_(1.29e-2), birksc2_(9.59e-6) {
25 gdml_identifiers_ = {p.get<std::vector<std::string>>("gdml_identifiers")};
26 enable_hit_contribs_ = p.get<bool>("enable_hit_contribs");
27 compress_hit_contribs_ = p.get<bool>("compress_hit_contribs");
28 max_origin_track_id_ = p.get<int>("max_origin_track_id");
29}
30
31ldmx::HcalID HcalSD::decodeCopyNumber(const std::uint32_t copyNumber,
32 const G4ThreeVector& localPosition,
33 const G4Box* scint) {
34 const unsigned int version{copyNumber / 0x01000000};
35 if (version != 0) {
37 return ldmx::HcalID{Index(copyNumber).field2(), Index(copyNumber).field1(),
38 Index(copyNumber).field0()};
39 }
40 const auto& geometry = getCondition<ldmx::HcalGeometry>(
42 unsigned int strip_id = 0;
43 const unsigned int section = copyNumber / 1000;
44 const unsigned int layer = copyNumber % 1000;
45
46 // 5cm wide bars are HARD-CODED
47 if (section == ldmx::HcalID::BACK) {
48 if (geometry.backLayerIsHorizontal(layer)) {
49 strip_id = int((localPosition.y() + scint->GetYHalfLength()) / 50.0);
50 } else {
51 strip_id = int((localPosition.x() + scint->GetXHalfLength()) / 50.0);
52 }
53 } else {
54 strip_id = int((localPosition.z() + scint->GetZHalfLength()) / 50.0);
55 }
56 return ldmx::HcalID{section, layer, strip_id};
57}
58
59G4bool HcalSD::ProcessHits(G4Step* aStep, G4TouchableHistory* ROhist) {
60 // Get the edep from the step.
61 G4double edep = aStep->GetTotalEnergyDeposit();
62
63 // Skip steps with no energy dep which come from non-Geantino particles.
64 if (edep == 0.0 and not isGeantino(aStep)) {
65 ldmx_log(trace) << "CalorimeterSD skipping step with zero edep.";
66 return false;
67 }
68
69 //---------------------------------------------------------------------------------------------------
70 // Birks' Law
71 // ===========
72 //
73 // In the case of Scintillator as active medium, we can
74 // describe the quenching effects with the Birks' law,
75 // using the expression and the coefficients taken from
76 // the paper NIM 80 (1970) 239-244 for the organic
77 // scintillator NE-102:
78 // S*dE/dr
79 // dL/dr = -----------------------------------
80 // 1 + C1*(dE/dr)
81 // with:
82 // S=1
83 // C1 = 1.29 x 10^-2 g*cm^-2*MeV^-1
84 // C2 = 9.59 x 10^-6 g^2*cm^-4*MeV^-2
85 // These are the same values used by ATLAS TileCal
86 // and CMS HCAL (and also the default in Geant3).
87 // You can try different values for these parameters,
88 // to have an idea on the uncertainties due to them,
89 // by uncommenting one of the lines below.
90 // To get the "dE/dr" that appears in the formula,
91 // which has the dimensions
92 // [ dE/dr ] = MeV * cm^2 / g
93 // we have to divide the energy deposit in MeV by the
94 // product of the step length (in cm) and the density
95 // of the scintillator:
96
97 G4double birks_factor(1.0);
98 G4double step_length = aStep->GetStepLength() / CLHEP::cm;
99 // Do not apply Birks for gamma deposits!
100 // Check, cut if necessary.
101 if (step_length > 1.0e-6) {
102 G4double rho = aStep->GetPreStepPoint()->GetMaterial()->GetDensity() /
103 (CLHEP::g / CLHEP::cm3);
104 G4double dedx = edep / (rho * step_length); //[MeV*cm^2/g]
105 birks_factor = 1.0 / (1.0 + birksc1_ * dedx + birksc2_ * dedx * dedx);
106 if (aStep->GetTrack()->GetDefinition() == G4Gamma::GammaDefinition())
107 birks_factor = 1.0;
108 if (aStep->GetTrack()->GetDefinition() == G4Neutron::NeutronDefinition())
109 birks_factor = 1.0;
110 }
111
112 // update edep to include birksFactor
113 edep *= birks_factor;
114
115 // Get the scintillator solid box
116 G4Box* scint = nullptr;
117
118 if (aStep) {
119 const auto* pre_step_point = aStep->GetPreStepPoint();
120 if (pre_step_point) {
121 const auto& touchable_handle = pre_step_point->GetTouchableHandle();
122 if (touchable_handle) {
123 const auto* volume = touchable_handle->GetVolume();
124
125 if (volume) {
126 const auto* logical_volume = volume->GetLogicalVolume();
127 if (logical_volume) {
128 auto* solid = logical_volume->GetSolid();
129 if (solid) {
130 scint = static_cast<G4Box*>(solid);
131 }
132 }
133 }
134 }
135 }
136 }
137
138 // Set the step mid-point as the hit position.
139 G4StepPoint* pre_point = aStep->GetPreStepPoint();
140 G4StepPoint* post_point = aStep->GetPostStepPoint();
141 // A Geant4 "touchable" is a way to uniquely identify a particular volume,
142 // short for touchable detector element. See the detector definition and
143 // response section of the Geant4 application developers manual for details.
144 //
145 // The TouchableHandle is just a reference counted pointer to a
146 // G4TouchableHistory object, which is a concrete implementation of a
147 // G4Touchable interface.
148 //
149 auto touchable_history{pre_point->GetTouchableHandle()->GetHistory()};
150 // Affine transform for converting between local and global coordinates
151 auto top_transform{touchable_history->GetTopTransform()};
152 G4ThreeVector position =
153 0.5 * (pre_point->GetPosition() + post_point->GetPosition());
154 G4ThreeVector local_position = top_transform.TransformPoint(position);
155
156 // Create the ID for the hit. Note 2 here corresponds to the "depth" of the
157 // geometry tree. If this changes in the GDML, this would have to be updated
158 // here. Currently, 0 corresponds to the world volume, 1 corresponds to the
159 // Hcal, and 2 to the bars/absorbers
160 int copy_num = touchable_history->GetVolume(2)->GetCopyNo();
161 ldmx::HcalID id = decodeCopyNumber(copy_num, local_position, scint);
162
163 if (hits_.find(id) == hits_.end()) {
164 // hit in empty cell/bar
165 auto& new_hit = hits_[id];
166 new_hit.setID(id.raw());
167 new_hit.setPosition(position[0], position[1], position[2]);
168 }
169
170 auto& hit = hits_[id];
171
172 // hit variables
173 const G4Track* track = aStep->GetTrack();
174 auto time = track->GetGlobalTime();
175 auto track_id = track->GetTrackID();
176 auto pdg = track->GetParticleDefinition()->GetPDGEncoding();
177
179 int contrib_i = hit.findContribIndex(track_id, pdg);
180 if (compress_hit_contribs_ and contrib_i != -1) {
181 hit.updateContrib(contrib_i, edep, time);
182 } else {
183 auto map{getTrackMap()};
184 auto incident{map.findIncident(track_id)};
185 // default "origin" is just the same as incident
186 // "origin" checks if a hit "originates" from one of the earliest
187 // track IDs (i.e. probably one of the primaries)
188 int origin{incident};
189 for (int i{1}; i < max_origin_track_id_; ++i) {
190 if (map.isDescendant(track_id, i, 100)) {
191 origin = i;
192 break;
193 }
194 }
195 hit.addContrib(incident, track_id, pdg, edep, time, origin);
196 }
197 } else {
198 // no hit contribs and hit already exists
199 hit.setEdep(hit.getEdep() + edep);
200 if (time < hit.getTime() or hit.getTime() == 0) {
201 hit.setTime(time);
202 }
203 }
204 // Pre/post step details for scintillator response simulation
205 // Note: These are set for each step, so the last step's details will be used
206 // if multiple steps occur in the same bar
207
208 // Convert back to mm
209 hit.setPathLength(step_length * CLHEP::cm / CLHEP::mm);
210 hit.setVelocity(track->GetVelocity());
211 const auto& geometry = getCondition<ldmx::HcalGeometry>(
213 // Convert pre/post step position from global coordinates to coordinates
214 // within the scintillator bar
215 const auto local_pre_step_point{
216 top_transform.TransformPoint(pre_point->GetPosition())};
217 const auto local_post_step_point{
218 top_transform.TransformPoint(post_point->GetPosition())};
219
220 // And rotate them to a local coordinate system for the bar that always has
221 // the same x/y/z definitions (see HcalGeometry for details)
222 auto local_pre_position_rotated{geometry.rotateGlobalToLocalBarPosition(
223 {local_pre_step_point[0], local_pre_step_point[1],
224 local_pre_step_point[2]},
225 id)};
226
227 auto local_post_position_rotated{geometry.rotateGlobalToLocalBarPosition(
228 {local_post_step_point[0], local_post_step_point[1],
229 local_post_step_point[2]},
230 id)};
231 hit.setPreStepPosition(local_pre_position_rotated[0],
232 local_pre_position_rotated[1],
233 local_pre_position_rotated[2]);
234 hit.setPostStepPosition(local_post_position_rotated[0],
235 local_post_position_rotated[1],
236 local_post_position_rotated[2]);
237 hit.setPreStepTime(pre_point->GetGlobalTime());
238 hit.setPostStepTime(post_point->GetGlobalTime());
239
240 ldmx_log(trace) << hit;
241
242 return true;
243}
244
246 // squash hits into list
247 std::vector<ldmx::SimCalorimeterHit> hits;
248 hits.reserve(hits_.size());
249 for (const auto& [id, hit] : hits_) hits.push_back(hit);
250 event.add(COLLECTION_NAME, hits);
251}
252
253} // namespace simcore
254
255DECLARE_SENSITIVEDETECTOR(simcore::HcalSD)
Class that translates HCal ID into positions of strip hits.
Class that defines an HCal sensitive detector.
Implements an event buffer system for storing event data.
Definition Event.h:42
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
static constexpr const char * CONDITIONS_OBJECT_NAME
Conditions object: The name of the python configuration calling this class (Hcal/python/HcalGeometry....
Implements detector ids for HCal subdetector.
Definition HcalID.h:19
A maximally-packed index of up to four different fields.
Definition PackedIndex.h:32
Handle to the conditions system, provided at construction to classes which require it.
Class defining a sensitive detector of type HCal.
Definition HcalSD.h:24
ldmx::HcalID decodeCopyNumber(const std::uint32_t copyNumber, const G4ThreeVector &localPosition, const G4Box *scint)
Decode copy number of scintillator bar.
Definition HcalSD.cxx:31
virtual void saveHits(framework::Event &event) override
Add our hits to the event bus and then reset the container.
Definition HcalSD.cxx:245
bool compress_hit_contribs_
compress hit contribs
Definition HcalSD.h:109
static const std::string COLLECTION_NAME
name of collection to be added to event bus
Definition HcalSD.h:27
std::map< ldmx::HcalID, ldmx::SimCalorimeterHit > hits_
map of hits to add to the event (will be squashed)
Definition HcalSD.h:105
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override
Create a hit out of the energy deposition deposited during a step.
Definition HcalSD.cxx:59
HcalSD(const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &params)
Constructor.
Definition HcalSD.cxx:22
bool enable_hit_contribs_
enable hit contribs
Definition HcalSD.h:107
int max_origin_track_id_
maximum track ID to be considered an "origin"
Definition HcalSD.h:111
Dynamically loaded Geant4 SensitiveDetector for saving hits in specific volumes within the simulation...
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.
const T & getCondition(const std::string &condition_name)
Record the configuration of this detector into the run header.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...