LDMX Software
SimCalorimeterHit.h
Go to the documentation of this file.
1
8#ifndef SIMCORE_EVENT_SIMCALORIMETERHIT_H_
9#define SIMCORE_EVENT_SIMCALORIMETERHIT_H_
10
11// ROOT
12#include "TObject.h" //For ClassDef
13
14// LDMX
15#include "SimCore/Event/SimParticle.h"
16
17namespace ldmx {
18
33 public:
35 static const std::string ECAL_COLLECTION;
36
38 static const std::string HCAL_COLLECTION;
39
44 struct Contrib {
56 int incidentID{-1};
57
59 int trackID{-1};
60
62 int pdgCode{0};
63
65 float edep{0};
66
68 float time{0};
69 };
70
74 SimCalorimeterHit() = default;
75
79 virtual ~SimCalorimeterHit() = default;
80
84 void Clear();
85
89 void Print() const;
90
95 int getID() const { return id_; }
96
101 void setID(const int id) { this->id_ = id; }
102
107 float getEdep() const { return edep_; }
108
113 void setEdep(const float edep) { this->edep_ = edep; }
114
119 std::vector<float> getPosition() const { return {x_, y_, z_}; }
120
126 std::vector<float> getPreStepPosition() const {
127 return {preStepX_, preStepY_, preStepZ_};
128 }
134 std::vector<float> getPostStepPosition() const {
135 return {postStepX_, postStepY_, postStepZ_};
136 }
143 void setPosition(const float x, const float y, const float z) {
144 this->x_ = x;
145 this->y_ = y;
146 this->z_ = z;
147 }
155 void setPreStepPosition(const float x, const float y, const float z) {
156 preStepX_ = x;
157 preStepY_ = y;
158 preStepZ_ = z;
159 }
167 void setPostStepPosition(const float x, const float y, const float z) {
168 postStepX_ = x;
169 postStepY_ = y;
170 postStepZ_ = z;
171 }
172
177 void setPathLength(const float length) { pathLength_ = length; }
182 float getPathLength() const { return pathLength_; }
187 void setPreStepTime(const float time) { preStepTime_ = time; }
192 void setPostStepTime(const float time) { postStepTime_ = time; }
193
198 void setVelocity(float velocity) { velocity_ = velocity; }
199
204 float getTime() const { return time_; }
209 float getPreStepTime() const { return preStepTime_; }
214 float getPostStepTime() const { return postStepTime_; }
215
220 float getVelocity() const { return velocity_; }
221
226 void setTime(const float time) { this->time_ = time; }
227
232 unsigned getNumberOfContribs() const { return nContribs_; }
233
243 void addContrib(int incidentID, int trackID, int pdgCode, float edep,
244 float time);
245
251 Contrib getContrib(int i) const;
252
259 int findContribIndex(int trackID, int pdgCode) const;
260
268 void updateContrib(int i, float edep, float time);
269
273 bool operator<(const SimCalorimeterHit &rhs) const {
274 return this->getTime() < rhs.getTime();
275 }
276
280 std::vector<int> getTrackIds() const { return trackIDContribs_; }
284 std::vector<int> getIncidentIds() const { return incidentIDContribs_; }
285
289 std::vector<int> getPdgIds() const { return pdgCodeContribs_; }
290
294 std::vector<float> getEdeps() const { return edepContribs_; }
295
299 std::vector<float> getTimes() const { return timeContribs_; }
300
301 private:
308 int id_{0};
309
313 float edep_{0};
314
318 float x_{0};
319
323 float y_{0};
324
328 float z_{0};
329
333 float time_{0};
334
338 std::vector<int> trackIDContribs_;
339
343 std::vector<int> incidentIDContribs_;
344
348 std::vector<int> pdgCodeContribs_;
349
353 std::vector<float> edepContribs_;
354
358 std::vector<float> timeContribs_;
359
363 unsigned nContribs_{0};
364
365 /*
366 * Parameters used only for hits corresponding to a single interactions
367 * (currently Hcal and TS).
368 */
369
374 float pathLength_{-1};
375
380 float preStepX_{0};
381 float preStepY_{0};
382 float preStepZ_{0};
386 float preStepTime_{0};
391 float postStepX_{0};
392 float postStepY_{0};
393 float postStepZ_{0};
398
402 float velocity_{-1};
403
407 ClassDef(SimCalorimeterHit, 5)
408};
409} // namespace ldmx
410
411#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.
std::vector< float > getEdeps() const
Get the list of energy depositions contributing to the hit.
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].
std::vector< int > getPdgIds() const
Get the list of PDG codes contributing to the hit.
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< int > getIncidentIds() const
Get the list of incident 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 ...
std::vector< float > getTimes() const
Get the list of times contributing to the hit.
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.
std::vector< int > getTrackIds() const
Get the list of track IDs contributing to the hit.
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.