LDMX Software
SimCalorimeterHit.h
Go to the documentation of this file.
1
7#ifndef SIMCORE_EVENT_SIMCALORIMETERHIT_H_
8#define SIMCORE_EVENT_SIMCALORIMETERHIT_H_
9
10// ROOT
11#include "TObject.h" //For ClassDef
12
13// LDMX
14#include "SimCore/Event/SimParticle.h"
15
16namespace ldmx {
17
32 public:
34 static const std::string ECAL_COLLECTION;
35
37 static const std::string HCAL_COLLECTION;
38
43 struct Contrib {
55 int incidentID{-1};
56
58 int trackID{-1};
59
61 int pdgCode{0};
62
64 float edep{0};
65
67 float time{0};
68 };
69
73 SimCalorimeterHit() = default;
74
78 virtual ~SimCalorimeterHit() = default;
79
83 void Clear();
84
88 void Print() const;
89
94 int getID() const { return id_; }
95
100 void setID(const int id) { this->id_ = id; }
101
106 float getEdep() const { return edep_; }
107
112 void setEdep(const float edep) { this->edep_ = edep; }
113
118 std::vector<float> getPosition() const { return {x_, y_, z_}; }
119
125 std::vector<float> getPreStepPosition() const {
126 return {preStepX_, preStepY_, preStepZ_};
127 }
133 std::vector<float> getPostStepPosition() const {
134 return {postStepX_, postStepY_, postStepZ_};
135 }
142 void setPosition(const float x, const float y, const float z) {
143 this->x_ = x;
144 this->y_ = y;
145 this->z_ = z;
146 }
154 void setPreStepPosition(const float x, const float y, const float z) {
155 preStepX_ = x;
156 preStepY_ = y;
157 preStepZ_ = z;
158 }
166 void setPostStepPosition(const float x, const float y, const float z) {
167 postStepX_ = x;
168 postStepY_ = y;
169 postStepZ_ = z;
170 }
171
176 void setPathLength(const float length) { pathLength_ = length; }
181 float getPathLength() const { return pathLength_; }
186 void setPreStepTime(const float time) { preStepTime_ = time; }
191 void setPostStepTime(const float time) { postStepTime_ = time; }
192
197 void setVelocity(float velocity) { velocity_ = velocity; }
198
203 float getTime() const { return time_; }
208 float getPreStepTime() const { return preStepTime_; }
213 float getPostStepTime() const { return postStepTime_; }
214
219 float getVelocity() const { return velocity_; }
220
225 void setTime(const float time) { this->time_ = time; }
226
231 unsigned getNumberOfContribs() const { return nContribs_; }
232
242 void addContrib(int incidentID, int trackID, int pdgCode, float edep,
243 float time);
244
250 Contrib getContrib(int i) const;
251
258 int findContribIndex(int trackID, int pdgCode) const;
259
267 void updateContrib(int i, float edep, float time);
268
272 bool operator<(const SimCalorimeterHit &rhs) const {
273 return this->getTime() < rhs.getTime();
274 }
275
276 private:
283 int id_{0};
284
288 float edep_{0};
289
293 float x_{0};
294
298 float y_{0};
299
303 float z_{0};
304
308 float time_{0};
309
313 std::vector<int> trackIDContribs_;
314
318 std::vector<int> incidentIDContribs_;
319
323 std::vector<int> pdgCodeContribs_;
324
328 std::vector<float> edepContribs_;
329
333 std::vector<float> timeContribs_;
334
338 unsigned nContribs_{0};
339
340 /*
341 * Parameters used only for hits corresponding to a single interactions
342 * (currently Hcal and TS).
343 */
344
349 float pathLength_{-1};
350
355 float preStepX_{0};
356 float preStepY_{0};
357 float preStepZ_{0};
361 float preStepTime_{0};
366 float postStepX_{0};
367 float postStepY_{0};
368 float postStepZ_{0};
373
377 float velocity_{-1};
378
382 ClassDef(SimCalorimeterHit, 4)
383};
384} // namespace ldmx
385
386#endif
Stores simulated calorimeter hit information.
float edep_
The energy deposition.
virtual ~SimCalorimeterHit()=default
Class destructor.
void setPreStepTime(const float time)
Set global pre-step time of the hit [ns].
float getVelocity() const
Get the track velocity of the hit [mm/ns].
float getPathLength() const
Get the physical path length for the interaction [mm].
void Clear()
Clear the data in the object.
void setPathLength(const float length)
Set the physical path length for the interaction [mm].
std::vector< float > getPostStepPosition() const
Get the XYZ post-step position of the hit in the coordinate frame of the sensitive volume [mm].
std::vector< int > incidentIDContribs_
The list of incident IDs contributing to the hit.
float getPostStepTime() const
Get the post-step time of the hit [ns].
float preStepX_
The X, Y, and Z positions [mm] before the interaction in the coordinate frame of the sensitive volume...
float velocity_
The track velocity [mm/ns].
int findContribIndex(int trackID, int pdgCode) const
Find the index of a hit contribution from a SimParticle and PDG code.
float z_
The Z position.
unsigned getNumberOfContribs() const
Get the number of hit contributions.
float time_
The global time of the hit.
float y_
The Y position.
float preStepTime_
The global time before the interaction [ns].
float getTime() const
Get the global time of the hit [ns].
void setEdep(const float edep)
Set the energy deposition of the hit [MeV].
void setPosition(const float x, const float y, const float z)
Set the XYZ position of the hit [mm].
SimCalorimeterHit()=default
Class constructor.
Contrib getContrib(int i) const
Get a hit contribution by index.
unsigned nContribs_
The number of hit contributions.
void addContrib(int incidentID, int trackID, int pdgCode, float edep, float time)
Add a hit contribution from a SimParticle.
std::vector< float > timeContribs_
The list of times contributing to the hit.
void setPostStepPosition(const float x, const float y, const float z)
Set the XYZ post-step position of the hit in the coordinate frame of the sensitive volume [mm].
static const std::string HCAL_COLLECTION
name of the hcal sim collection, should match gdml
float x_
The X position.
void setID(const int id)
Set the detector ID.
float pathLength_
The true path length [mm].
std::vector< int > trackIDContribs_
The list of track IDs contributing to the hit.
std::vector< float > getPosition() const
Get the XYZ position of the hit [mm].
int getID() const
Get the detector ID.
void setPostStepTime(const float time)
Set global post-step time of the hit [ns].
bool operator<(const SimCalorimeterHit &rhs) const
Sort by time of hit.
void setTime(const float time)
Set the time of the hit [ns].
std::vector< int > pdgCodeContribs_
The list of PDG codes contributing to the hit.
std::vector< float > edepContribs_
The list of energy depositions contributing to the hit.
std::vector< float > getPreStepPosition() const
Get the XYZ pre-step position of the hit in the coordinate frame of the sensitive volume [mm].
float getEdep() const
Get the energy deposition of the hit [MeV].
float postStepX_
The X, Y, and Z positions [mm] after the interaction in the coordinate frame of the sensitive volume.
static const std::string ECAL_COLLECTION
name of the ecal sim collection, should match gdml
void updateContrib(int i, float edep, float time)
Update an existing hit contribution by incrementing its edep and setting the time if the new time is ...
float postStepTime_
The global time after the interaction [ns].
int id_
Member variables used in all calorimeter types.
void setPreStepPosition(const float x, const float y, const float z)
Set the XYZ pre-step position of the hit in the coordinate frame of the sensitive volume [mm].
void Print() const
Print out the object.
void setVelocity(float velocity)
Set the velocity of the track [mm/ns].
float getPreStepTime() const
Get the pre-step time of the hit [ns].
Information about a contribution to the hit in the associated cell.
int trackID
track ID of this contributor
int pdgCode
PDG ID of this contributor.
float time
Time this contributor made the hit (global Geant4 time)
int incidentID
trackID of incident particle that is an ancestor of the contributor
float edep
Energy depostied by this contributor.