LDMX Software
tracking::reco::LinearTruthTracking Class Reference

Public Member Functions

 LinearTruthTracking (const std::string &name, framework::Process &process)
 Constructor.
 
virtual ~LinearTruthTracking ()=default
 Destructor.
 
void onProcessEnd () override
 Output event statistics.
 
void configure (framework::config::Parameters &parameters) override
 Configure the processor using the given user specified parameters.
 
void produce (framework::Event &event) override
 Run the processor and create a collection of results which indicate if a charge particle can be found by the recoil tracker.
 
- Public Member Functions inherited from tracking::reco::TrackingGeometryUser
 TrackingGeometryUser (const std::string &name, framework::Process &p)
 
- Public Member Functions inherited from framework::Producer
 Producer (const std::string &name, Process &process)
 Class constructor.
 
virtual void beforeNewRun (ldmx::RunHeader &header)
 Handle allowing producers to modify run headers before the run begins.
 
- Public Member Functions inherited from framework::EventProcessor
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()
 Class destructor.
 
virtual void onNewRun (const ldmx::RunHeader &runHeader)
 Callback for the EventProcessor to take any necessary action when the run being processed changes.
 
virtual void onFileOpen (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a event input ROOT file is closed.
 
virtual void onProcessStart ()
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
template<class T >
const T & getCondition (const std::string &condition_name)
 Access a conditions object for the current event.
 
TDirectory * getHistoDirectory ()
 Access/create a directory in the histogram file for this event processor to create histograms and analysis tuples.
 
void setStorageHint (framework::StorageControl::Hint hint)
 Mark the current event as having the given storage control hint from this module.
 
void setStorageHint (framework::StorageControl::Hint hint, const std::string &purposeString)
 Mark the current event as having the given storage control hint from this module and the given purpose string.
 
int getLogFrequency () const
 Get the current logging frequency from the process.
 
int getRunNumber () const
 Get the run number from the process.
 
std::string getName () const
 Get the processor name.
 
void createHistograms (const std::vector< framework::config::Parameters > &histos)
 Internal function which is used to create histograms passed from the python configuration @parma histos vector of Parameters that configure histograms to create.
 

Protected Member Functions

ldmx::StraightTrack truthTracker (const std::vector< ldmx::SimTrackerHit > &points, std::vector< std::array< double, 3 > > &ecal_points)
 
double calculateDistance (const std::array< double, 3 > &point1, const std::array< double, 3 > &point2)
 
std::tuple< double, double, double, double > fit3DLine (const std::vector< ldmx::SimTrackerHit > &points)
 
double globalChiSquare (const std::vector< ldmx::SimTrackerHit > &points, double a_x, double a_y, double b_x, double b_y)
 
- Protected Member Functions inherited from tracking::reco::TrackingGeometryUser
const Acts::GeometryContext & geometry_context ()
 
const Acts::MagneticFieldContext & magnetic_field_context ()
 
const Acts::CalibrationContext & calibration_context ()
 
const geo::TrackersTrackingGeometrygeometry ()
 
- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 

Protected Attributes

double processing_time_ {0.}
 
long n_events_ {0}
 
unsigned int n_truth_ {0}
 
long n_missing_ {0}
 
long n_empty_ {0}
 
double ecal_first_layer_z_threshold_ {250.0}
 
std::string out_trk_collection_ {"LinearRecoilTruthTracks"}
 The name of the output collection of seeds to be stored.
 
std::string input_hits_collection_ {"RecoilSimHits"}
 The name of the input hits collection to use in finding seeds..
 
std::string input_rec_hits_collection_ {"EcalRecHits"}
 The name of the tagger Tracks (only for Recoil Seeding)
 
std::string input_pass_name_ {""}
 
- Protected Attributes inherited from framework::EventProcessor
HistogramHelper histograms_
 Interface class for making and filling histograms.
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger theLog_
 The logger for this EventProcessor.
 

Additional Inherited Members

- Static Public Member Functions inherited from framework::EventProcessor
static void declare (const std::string &classname, int classtype, EventProcessorMaker *)
 Internal function which is part of the PluginFactory machinery.
 
- Static Public Attributes inherited from framework::Producer
static const int CLASSTYPE {1}
 Constant used to track EventProcessor types by the PluginFactory.
 

Detailed Description

Definition at line 26 of file LinearTruthTracking.h.

Constructor & Destructor Documentation

◆ LinearTruthTracking()

tracking::reco::LinearTruthTracking::LinearTruthTracking ( const std::string & name,
framework::Process & process )

Constructor.

Parameters
nameThe name of the instance of this object.
processThe process running this producer.

Definition at line 8 of file LinearTruthTracking.cxx.

10 : TrackingGeometryUser(name, process) {}

Member Function Documentation

◆ calculateDistance()

double tracking::reco::LinearTruthTracking::calculateDistance ( const std::array< double, 3 > & point1,
const std::array< double, 3 > & point2 )
protected

Definition at line 162 of file LinearTruthTracking.cxx.

163 {
164 return sqrt(pow(point1[1] - point2[1], 2) + pow(point1[2] - point2[2], 2));
165} // calculateDistance

◆ configure()

void tracking::reco::LinearTruthTracking::configure ( framework::config::Parameters & parameters)
overridevirtual

Configure the processor using the given user specified parameters.

Parameters
parametersSet of parameters used to configure this processor.

Reimplemented from framework::EventProcessor.

Definition at line 12 of file LinearTruthTracking.cxx.

12 {
13 // Output seed name
14 out_trk_collection_ = parameters.getParameter<std::string>(
15 "out_trk_collection", "LinearRecoilTruthTracks");
16
17 // Input strip hits
18 input_hits_collection_ = parameters.getParameter<std::string>(
19 "input_hits_collection", "RecoilSimHits");
20 input_rec_hits_collection_ = parameters.getParameter<std::string>(
21 "input_recHits_collection", "EcalRecHits");
22
23 input_pass_name_ =
24 parameters.getParameter<std::string>("input_pass_name", "");
25
26 ecal_first_layer_z_threshold_ =
27 parameters.getParameter<double>("ecal_first_layer_z_threshold");
28} // configure
std::string input_rec_hits_collection_
The name of the tagger Tracks (only for Recoil Seeding)
std::string input_hits_collection_
The name of the input hits collection to use in finding seeds..
std::string out_trk_collection_
The name of the output collection of seeds to be stored.

References input_hits_collection_, input_rec_hits_collection_, and out_trk_collection_.

◆ fit3DLine()

std::tuple< double, double, double, double > tracking::reco::LinearTruthTracking::fit3DLine ( const std::vector< ldmx::SimTrackerHit > & points)
protected

Definition at line 167 of file LinearTruthTracking.cxx.

168 {
169 std::vector<double> z_vals, x_vals, y_vals;
170
171 for (const auto& point : points) {
172 const auto& position = point.getPosition();
173 x_vals.push_back(position[0]);
174 y_vals.push_back(position[1]);
175 z_vals.push_back(position[2]);
176 }
177
178 size_t num_points = points.size();
179 Eigen::MatrixXd a_mat(num_points, 2);
180 Eigen::VectorXd x_vec(num_points), y_vec(num_points);
181
182 for (size_t i = 0; i < num_points; ++i) {
183 a_mat(i, 0) = z_vals[i];
184 a_mat(i, 1) = 1;
185 x_vec(i) = x_vals[i];
186 y_vec(i) = y_vals[i];
187 }
188
189 Eigen::Matrix2d AtA = a_mat.transpose() * a_mat;
190 Eigen::Vector2d At_x = a_mat.transpose() * x_vec;
191 Eigen::Vector2d At_y = a_mat.transpose() * y_vec;
192
193 Eigen::Vector2d x_params = AtA.ldlt().solve(At_x);
194 Eigen::Vector2d y_params = AtA.ldlt().solve(At_y);
195
196 // Return (mx, bx, my, by)
197 return {x_params(0), x_params(1), y_params(0), y_params(1)};
198}

◆ globalChiSquare()

double tracking::reco::LinearTruthTracking::globalChiSquare ( const std::vector< ldmx::SimTrackerHit > & points,
double a_x,
double a_y,
double b_x,
double b_y )
protected

Definition at line 200 of file LinearTruthTracking.cxx.

202 {
203 double chi2_x = 0, chi2_y = 0;
204
205 for (const auto& point : points) {
206 const auto& position = point.getPosition();
207 double x_pos = position[0];
208 double y_pos = position[1];
209 double z_pos = position[2];
210
211 double x_fit = m_x * z_pos + b_x;
212 double y_fit = m_y * z_pos + b_y;
213
214 chi2_x += std::pow((x_pos - x_fit), 2);
215 chi2_y += std::pow((y_pos - y_fit), 2);
216 }
217
218 return chi2_x + chi2_y;
219} // globalChiSquare

◆ onProcessEnd()

void tracking::reco::LinearTruthTracking::onProcessEnd ( )
overridevirtual

Output event statistics.

Reimplemented from framework::EventProcessor.

Definition at line 151 of file LinearTruthTracking.cxx.

151 {
152 ldmx_log(info) << "AVG Time / Event: " << std::fixed << std::setprecision(1)
153 << processing_time_ / n_events_ << " ms";
154 ldmx_log(info) << "Number of Events without enough seed points: "
155 << n_missing_;
156 ldmx_log(info) << "Number of Events with only 1 trackID == 1 point: "
157 << n_empty_;
158 ldmx_log(info) << "Total Truth Tracks / Event: " << n_truth_ << "/"
159 << n_events_;
160} // onProcessEnd

◆ produce()

void tracking::reco::LinearTruthTracking::produce ( framework::Event & event)
overridevirtual

Run the processor and create a collection of results which indicate if a charge particle can be found by the recoil tracker.

Parameters
eventThe event to process.

Implements framework::Producer.

Definition at line 30 of file LinearTruthTracking.cxx.

30 {
31 auto start = std::chrono::high_resolution_clock::now();
32 std::vector<ldmx::StraightTrack> straight_truth_tracks;
33 n_events_++;
34
35 const std::vector<ldmx::SimTrackerHit> recoil_hits =
37 input_pass_name_);
38 const std::vector<ldmx::EcalHit> ecal_rec_hit =
39 event.getCollection<ldmx::EcalHit>(input_rec_hits_collection_,
40 input_pass_name_);
41
42 std::vector<std::array<double, 3>> first_layer_ecal_rec_hits;
43
44 for (const auto& x_ecal : ecal_rec_hit) {
45 if (x_ecal.getZPos() < ecal_first_layer_z_threshold_) {
46 first_layer_ecal_rec_hits.push_back(
47 {x_ecal.getZPos(), x_ecal.getXPos(), x_ecal.getYPos()});
48 } // if first layer of Ecal
49 } // for positions in ecalRecHit
50
51 // Check if we would fit empty seeds.
52 if (recoil_hits.size() < 2) {
53 n_missing_++;
54 n_truth_ += straight_truth_tracks.size();
55 event.add(out_trk_collection_, straight_truth_tracks);
56 return;
57 } // if not enough hits
58
59 std::vector<ldmx::SimTrackerHit> truth_points;
60
61 // Find only recoil electron out of the hit collection
62 for (const auto& point : recoil_hits) {
63 if (point.getTrackID() == 1) {
64 truth_points.push_back(point);
65 } // if trackID == 1
66 } // for point in hit collection
67
68 // should only find 1 track (recoil e-) or 0 track (no recoil e-)
69 if (truth_points.size() > 1) {
70 straight_truth_tracks.push_back(
71 truthTracker(truth_points, first_layer_ecal_rec_hits));
72 } // if we have truth points
73 else {
74 n_empty_++;
75 n_truth_ += straight_truth_tracks.size();
76 return;
77 } // else no truth points
78
79 n_truth_ += straight_truth_tracks.size();
80 event.add(out_trk_collection_, straight_truth_tracks);
81
82 auto end = std::chrono::high_resolution_clock::now();
83 auto diff = end - start;
84 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
85
86 first_layer_ecal_rec_hits.clear();
87 straight_truth_tracks.clear();
88
89} // produce
Stores reconstructed hit information from the ECAL.
Definition EcalHit.h:19
Represents a simulated tracker hit in the simulation.

References input_hits_collection_, input_rec_hits_collection_, and out_trk_collection_.

◆ truthTracker()

ldmx::StraightTrack tracking::reco::LinearTruthTracking::truthTracker ( const std::vector< ldmx::SimTrackerHit > & points,
std::vector< std::array< double, 3 > > & ecal_points )
protected

Definition at line 91 of file LinearTruthTracking.cxx.

93 {
95
96 // We have pre-selected hits corresponding to the Recoil e-, so we only fit
97 // these (no RecHit)
98 auto [m_x, b_x, m_y, b_y] = fit3DLine(points);
99
100 trk.setSlopeX(m_x);
101 trk.setInterceptX(b_x);
102 trk.setSlopeY(m_y);
103 trk.setInterceptY(b_y);
104 trk.setTheta(std::atan2(m_y, std::sqrt(1 + m_x * m_x)));
105 trk.setPhi(std::atan2(m_x, 1.0));
106 trk.setTargetLocation(0.0, b_x, b_y);
107 trk.setNhits(points.size());
108 trk.setNdf(points.size() - 2);
109
110 trk.setTrackID(points[0].getTrackID());
111 trk.setPdgID(points[0].getPdgID());
112 trk.setTruthProb(1.0);
113
114 trk.setAllTruthSensorPoints(points);
115
116 trk.setChi2(globalChiSquare(points, m_x, m_y, b_x, b_y));
117
118 // Z position from the first point in ecalPoints
119 if (ecal_points.size() > 0) {
120 double ecal_first_layer_z = ecal_points[0][0];
121
122 // Extrapolate track to first layer of Ecal
123 std::array<double, 3> extrapolated_point = {ecal_first_layer_z,
124 m_x * ecal_first_layer_z + b_x,
125 m_y * ecal_first_layer_z + b_y};
126
127 const std::array<double, 3>* closest_rec_hit = nullptr;
128 double min_distance = std::numeric_limits<double>::max();
129
130 // Loop through ecalPoints to find the closest recHit to our track
131 for (const auto& ecal_rec_hit : ecal_points) {
132 double temp_distance =
133 calculateDistance(extrapolated_point, ecal_rec_hit);
134
135 if (temp_distance < min_distance) {
136 min_distance = temp_distance;
137 closest_rec_hit = &ecal_rec_hit;
138 } // if compare min distance
139 } // for recHits
140
141 if (closest_rec_hit != nullptr) {
142 trk.setFirstLayerEcalRecHit(*closest_rec_hit);
143 trk.setDistancetoRecHit(min_distance);
144 trk.setEcalLayer1Location(extrapolated_point);
145 } // if set trk info
146 } // make sure there are ecalPoints to work with
147
148 return trk;
149} // truthTracker

Member Data Documentation

◆ ecal_first_layer_z_threshold_

double tracking::reco::LinearTruthTracking::ecal_first_layer_z_threshold_ {250.0}
protected

Definition at line 83 of file LinearTruthTracking.h.

83{250.0};

◆ input_hits_collection_

std::string tracking::reco::LinearTruthTracking::input_hits_collection_ {"RecoilSimHits"}
protected

The name of the input hits collection to use in finding seeds..

Definition at line 88 of file LinearTruthTracking.h.

88{"RecoilSimHits"};

Referenced by configure(), and produce().

◆ input_pass_name_

std::string tracking::reco::LinearTruthTracking::input_pass_name_ {""}
protected

Definition at line 91 of file LinearTruthTracking.h.

91{""};

◆ input_rec_hits_collection_

std::string tracking::reco::LinearTruthTracking::input_rec_hits_collection_ {"EcalRecHits"}
protected

The name of the tagger Tracks (only for Recoil Seeding)

Definition at line 90 of file LinearTruthTracking.h.

90{"EcalRecHits"};

Referenced by configure(), and produce().

◆ n_empty_

long tracking::reco::LinearTruthTracking::n_empty_ {0}
protected

Definition at line 80 of file LinearTruthTracking.h.

80{0};

◆ n_events_

long tracking::reco::LinearTruthTracking::n_events_ {0}
protected

Definition at line 77 of file LinearTruthTracking.h.

77{0};

◆ n_missing_

long tracking::reco::LinearTruthTracking::n_missing_ {0}
protected

Definition at line 79 of file LinearTruthTracking.h.

79{0};

◆ n_truth_

unsigned int tracking::reco::LinearTruthTracking::n_truth_ {0}
protected

Definition at line 78 of file LinearTruthTracking.h.

78{0};

◆ out_trk_collection_

std::string tracking::reco::LinearTruthTracking::out_trk_collection_ {"LinearRecoilTruthTracks"}
protected

The name of the output collection of seeds to be stored.

Definition at line 86 of file LinearTruthTracking.h.

86{"LinearRecoilTruthTracks"};

Referenced by configure(), and produce().

◆ processing_time_

double tracking::reco::LinearTruthTracking::processing_time_ {0.}
protected

Definition at line 76 of file LinearTruthTracking.h.

76{0.};

The documentation for this class was generated from the following files: