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}
27
28ldmx::HcalID HcalSD::decodeCopyNumber(const std::uint32_t copyNumber,
29 const G4ThreeVector& localPosition,
30 const G4Box* scint) {
31 const unsigned int version{copyNumber / 0x01000000};
32 if (version != 0) {
34 return ldmx::HcalID{Index(copyNumber).field2(), Index(copyNumber).field1(),
35 Index(copyNumber).field0()};
36 }
37 const auto& geometry = getCondition<ldmx::HcalGeometry>(
39 unsigned int strip_id = 0;
40 const unsigned int section = copyNumber / 1000;
41 const unsigned int layer = copyNumber % 1000;
42
43 // 5cm wide bars are HARD-CODED
44 if (section == ldmx::HcalID::BACK) {
45 if (geometry.backLayerIsHorizontal(layer)) {
46 strip_id = int((localPosition.y() + scint->GetYHalfLength()) / 50.0);
47 } else {
48 strip_id = int((localPosition.x() + scint->GetXHalfLength()) / 50.0);
49 }
50 } else {
51 strip_id = int((localPosition.z() + scint->GetZHalfLength()) / 50.0);
52 }
53 return ldmx::HcalID{section, layer, strip_id};
54}
55
56G4bool HcalSD::ProcessHits(G4Step* aStep, G4TouchableHistory* ROhist) {
57 // Get the edep from the step.
58 G4double edep = aStep->GetTotalEnergyDeposit();
59
60 // Skip steps with no energy dep which come from non-Geantino particles.
61 if (edep == 0.0 and not isGeantino(aStep)) {
62 if (verboseLevel > 2) {
63 std::cout << "CalorimeterSD skipping step with zero edep." << std::endl
64 << std::endl;
65 }
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 // Create a new cal hit.
116 ldmx::SimCalorimeterHit& hit{hits_.emplace_back()};
117
118 // Get the scintillator solid box
119 G4Box* scint = nullptr;
120
121 if (aStep) {
122 const auto* pre_step_point = aStep->GetPreStepPoint();
123 if (pre_step_point) {
124 const auto& touchable_handle = pre_step_point->GetTouchableHandle();
125 if (touchable_handle) {
126 const auto* volume = touchable_handle->GetVolume();
127
128 if (volume) {
129 const auto* logical_volume = volume->GetLogicalVolume();
130 if (logical_volume) {
131 auto* solid = logical_volume->GetSolid();
132 if (solid) {
133 scint = static_cast<G4Box*>(solid);
134 }
135 }
136 }
137 }
138 }
139 }
140
141 // Set the step mid-point as the hit position.
142 G4StepPoint* pre_point = aStep->GetPreStepPoint();
143 G4StepPoint* post_point = aStep->GetPostStepPoint();
144 // A Geant4 "touchable" is a way to uniquely identify a particular volume,
145 // short for touchable detector element. See the detector definition and
146 // response section of the Geant4 application developers manual for details.
147 //
148 // The TouchableHandle is just a reference counted pointer to a
149 // G4TouchableHistory object, which is a concrete implementation of a
150 // G4Touchable interface.
151 //
152 auto touchable_history{pre_point->GetTouchableHandle()->GetHistory()};
153 // Affine transform for converting between local and global coordinates
154 auto top_transform{touchable_history->GetTopTransform()};
155 G4ThreeVector position =
156 0.5 * (pre_point->GetPosition() + post_point->GetPosition());
157 G4ThreeVector local_position = top_transform.TransformPoint(position);
158 hit.setPosition(position[0], position[1], position[2]);
159
160 // Create the ID for the hit. Note 2 here corresponds to the "depth" of the
161 // geometry tree. If this changes in the GDML, this would have to be updated
162 // here. Currently, 0 corresponds to the world volume, 1 corresponds to the
163 // Hcal, and 2 to the bars/absorbers
164 int copy_num = touchable_history->GetVolume(2)->GetCopyNo();
165 ldmx::HcalID id = decodeCopyNumber(copy_num, local_position, scint);
166 hit.setID(id.raw());
167
168 // add one contributor for this hit with
169 // ID of ancestor incident on Cal-Region
170 // ID of this track
171 // PDG of this track
172 // EDEP (including birks factor)
173 // time of this hit
174 const G4Track* track = aStep->GetTrack();
175 int track_id = track->GetTrackID();
176 hit.addContrib(getTrackMap().findIncident(track_id), track_id,
177 track->GetParticleDefinition()->GetPDGEncoding(), edep,
178 track->GetGlobalTime());
179 //
180 // Pre/post step details for scintillator response simulation
181
182 // Convert back to mm
183 hit.setPathLength(step_length * CLHEP::cm / CLHEP::mm);
184 hit.setVelocity(track->GetVelocity());
185 const auto& geometry = getCondition<ldmx::HcalGeometry>(
187 // Convert pre/post step position from global coordinates to coordinates
188 // within the scintillator bar
189 const auto local_pre_step_point{
190 top_transform.TransformPoint(pre_point->GetPosition())};
191 const auto local_post_step_point{
192 top_transform.TransformPoint(post_point->GetPosition())};
193
194 // And rotate them to a local coordinate system for the bar that always has
195 // the same x/y/z definitions (see HcalGeometry for details)
196 auto local_pre_position_rotated{geometry.rotateGlobalToLocalBarPosition(
197 {local_pre_step_point[0], local_pre_step_point[1],
198 local_pre_step_point[2]},
199 id)};
200
201 auto local_post_position_rotated{geometry.rotateGlobalToLocalBarPosition(
202 {local_post_step_point[0], local_post_step_point[1],
203 local_post_step_point[2]},
204 id)};
205 hit.setPreStepPosition(local_pre_position_rotated[0],
206 local_pre_position_rotated[1],
207 local_pre_position_rotated[2]);
208 hit.setPostStepPosition(local_post_position_rotated[0],
209 local_post_position_rotated[1],
210 local_post_position_rotated[2]);
211 hit.setPreStepTime(pre_point->GetGlobalTime());
212 hit.setPostStepTime(post_point->GetGlobalTime());
213
214 ldmx_log(trace) << hit;
215
216 return true;
217}
218
219} // namespace simcore
220
221DECLARE_SENSITIVEDETECTOR(simcore::HcalSD)
Class that translates HCal ID into positions of strip hits.
Class that defines an HCal sensitive detector.
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
Stores simulated calorimeter hit information.
void setPosition(const float x, const float y, const float z)
Set the XYZ position of the hit [mm].
Handle to the conditions system, provided at construction to classes which require it.
Class defining a sensitive detector of type HCal.
Definition HcalSD.h:23
ldmx::HcalID decodeCopyNumber(const std::uint32_t copyNumber, const G4ThreeVector &localPosition, const G4Box *scint)
Decode copy number of scintillator bar.
Definition HcalSD.cxx:28
static const std::string COLLECTION_NAME
name of collection to be added to event bus
Definition HcalSD.h:26
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override
Create a hit out of the energy deposition deposited during a step.
Definition HcalSD.cxx:56
HcalSD(const std::string &name, simcore::ConditionsInterface &ci, const framework::config::Parameters &params)
Constructor.
Definition HcalSD.cxx:22
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 ...