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 incident_id_{-1};
57
59 int track_id_{-1};
60
62 int pdg_code_{0};
63
65 float edep_{0};
66
68 float time_{0};
69
71 int origin_id_{-1};
72 };
73
77 SimCalorimeterHit() = default;
78
82 virtual ~SimCalorimeterHit() = default;
83
87 void clear();
88
92 friend std::ostream &operator<<(std::ostream &o, const SimCalorimeterHit &d);
93
98 int getID() const { return id_; }
99
104 void setID(const int id) { this->id_ = id; }
105
110 float getEdep() const { return edep_; }
111
116 void setEdep(const float edep) { this->edep_ = edep; }
117
122 std::vector<float> getPosition() const { return {x_, y_, z_}; }
123
129 std::vector<float> getPreStepPosition() const {
130 return {pre_step_x_, pre_step_y_, pre_step_z_};
131 }
137 std::vector<float> getPostStepPosition() const {
138 return {post_step_x_, post_step_y_, post_step_z_};
139 }
146 void setPosition(const float x, const float y, const float z) {
147 this->x_ = x;
148 this->y_ = y;
149 this->z_ = z;
150 }
158 void setPreStepPosition(const float x, const float y, const float z) {
159 pre_step_x_ = x;
160 pre_step_y_ = y;
161 pre_step_z_ = z;
162 }
170 void setPostStepPosition(const float x, const float y, const float z) {
171 post_step_x_ = x;
172 post_step_y_ = y;
173 post_step_z_ = z;
174 }
175
180 void setPathLength(const float length) { path_length_ = length; }
185 float getPathLength() const { return path_length_; }
190 void setPreStepTime(const float time) { pre_step_time_ = time; }
195 void setPostStepTime(const float time) { post_step_time_ = time; }
196
201 void setVelocity(float velocity) { velocity_ = velocity; }
202
207 float getTime() const { return time_; }
212 float getPreStepTime() const { return pre_step_time_; }
217 float getPostStepTime() const { return post_step_time_; }
218
223 float getVelocity() const { return velocity_; }
224
229 void setTime(const float time) { this->time_ = time; }
230
235 unsigned getNumberOfContribs() const { return n_contribs_; }
236
246 void addContrib(int incidentID, int trackID, int pdgCode, float edep,
247 float time, int originID = -1);
248
254 Contrib getContrib(int i) const;
255
262 int findContribIndex(int trackID, int pdgCode) const;
263
271 void updateContrib(int i, float edep, float time);
272
276 bool operator<(const SimCalorimeterHit &rhs) const {
277 return this->getTime() < rhs.getTime();
278 }
279
283 std::vector<int> getTrackIds() const { return track_id_contribs_; }
287 std::vector<int> getIncidentIds() const { return incident_id_contribs_; }
288
292 std::vector<int> getPdgIds() const { return pdg_code_contribs_; }
293
297 std::vector<float> getEdeps() const { return edep_contribs_; }
298
302 std::vector<float> getTimes() const { return time_contribs_; }
303
304 private:
311 int id_{0};
312
316 float edep_{0};
317
321 float x_{0};
322
326 float y_{0};
327
331 float z_{0};
332
336 float time_{0};
337
341 std::vector<int> track_id_contribs_;
342
346 std::vector<int> incident_id_contribs_;
347
351 std::vector<int> pdg_code_contribs_;
352
356 std::vector<float> edep_contribs_;
357
361 std::vector<float> time_contribs_;
362
366 std::vector<int> origin_contribs_;
367
371 unsigned n_contribs_{0};
372
373 /*
374 * Parameters used only for hits corresponding to a single interactions
375 * (currently Hcal and TS).
376 */
377
382 float path_length_{-1};
383
388 float pre_step_x_{0};
389 float pre_step_y_{0};
390 float pre_step_z_{0};
399 float post_step_x_{0};
400 float post_step_y_{0};
401 float post_step_z_{0};
406
410 float velocity_{-1};
411
415 ClassDef(SimCalorimeterHit, 7)
416};
417} // namespace ldmx
418
419#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< float > time_contribs_
The list of times contributing to the hit.
float getPostStepTime() const
Get the post-step time of the hit [ns].
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.
float pre_step_time_
The global time before the interaction [ns].
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 getTime() const
Get the global time of the hit [ns].
void setEdep(const float edep)
Set the energy deposition of the hit [MeV].
friend std::ostream & operator<<(std::ostream &o, const SimCalorimeterHit &d)
Print out the object.
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.
std::vector< float > edep_contribs_
The list of energy depositions contributing to the hit.
std::vector< int > track_id_contribs_
The list of track IDs contributing to the hit.
SimCalorimeterHit()=default
Class constructor.
float pre_step_x_
The X, Y, and Z positions [mm] before the interaction in the coordinate frame of the sensitive volume...
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.
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].
std::vector< int > incident_id_contribs_
The list of incident IDs contributing to the hit.
float path_length_
The true path length [mm].
Contrib getContrib(int i) const
Get a hit contribution by index_.
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 addContrib(int incidentID, int trackID, int pdgCode, float edep, float time, int originID=-1)
Add a hit contribution from a SimParticle.
void setTime(const float time)
Set the time of the hit [ns].
std::vector< int > pdg_code_contribs_
The list of PDG codes 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].
unsigned n_contribs_
The number of hit contributions.
float getEdep() const
Get the energy deposition of the hit [MeV].
std::vector< int > origin_contribs_
The list of origin IDs contributing to the hit.
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 post_step_x_
The X, Y, and Z positions [mm] after the interaction in the coordinate frame of the sensitive volume.
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].
std::vector< int > getTrackIds() const
Get the list of track IDs contributing to the hit.
float post_step_time_
The global time after the interaction [ns].
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.
float time_
Time this contributor made the hit (global Geant4 time)
float edep_
Energy depostied by this contributor.
int track_id_
track ID of this contributor
int incident_id_
trackID of incident particle that is an ancestor of the contributor
int origin_id_
Saves ID of electron incident particle came from.
int pdg_code_
PDG ID of this contributor.