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 process (Event &event) final
 Processing an event for a Producer is calling produce.
 
- Public Member Functions inherited from framework::EventProcessor
 DECLARE_FACTORY (EventProcessor, EventProcessor *, const std::string &, Process &)
 declare that we have a factory for this class
 
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()=default
 Class destructor.
 
virtual void beforeNewRun (ldmx::RunHeader &run_header)
 Callback for Producers to add parameters to the run header before conditions are initialized.
 
virtual void onNewRun (const ldmx::RunHeader &run_header)
 Callback for the EventProcessor to take any necessary action when the run being processed changes.
 
virtual void onFileOpen (EventFile &event_file)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &event_file)
 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 & geometryContext ()
 
const Acts::MagneticFieldContext & magneticFieldContext ()
 
const Acts::CalibrationContext & calibrationContext ()
 
const geo::TrackersTrackingGeometrygeometry ()
 
void loadBField (const std::string &path, const std::vector< double > &map_offset={0., 0., 0.})
 Load the interpolated B-field map from path and cache it.
 
void loadBField (const std::vector< double > &map_offset={0., 0., 0.})
 Load B-field from the path recorded in the detector GDML.
 
std::shared_ptr< Acts::MagneticFieldProvider > bField () const
 Return the loaded B-field provider.
 
- 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
HistogramPool histograms_
 helper object for making and filling histograms
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger the_log_
 The logger for this EventProcessor.
 

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) {}
virtual void process(Event &event) final
Processing an event for a Producer is calling produce.

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 159 of file LinearTruthTracking.cxx.

160 {
161 return sqrt(pow(point1[1] - point2[1], 2) + pow(point1[2] - point2[2], 2));
162} // 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.get<std::string>("out_trk_collection",
15 "LinearRecoilTruthTracks");
16
17 // Input strip hits_
19 parameters.get<std::string>("input_hits_collection", "RecoilSimHits");
21 parameters.get<std::string>("input_rec_hits_collection", "EcalRecHits");
22
23 input_pass_name_ = parameters.get<std::string>("input_pass_name");
24
25 ecal_first_layer_z_threshold_ =
26 parameters.get<double>("ecal_first_layer_z_threshold");
27} // configure
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
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 framework::config::Parameters::get(), 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 164 of file LinearTruthTracking.cxx.

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

◆ 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 197 of file LinearTruthTracking.cxx.

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

◆ onProcessEnd()

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

Output event statistics.

Reimplemented from framework::EventProcessor.

Definition at line 148 of file LinearTruthTracking.cxx.

148 {
149 ldmx_log(info) << "AVG Time / Event: " << std::fixed << std::setprecision(1)
150 << processing_time_ / n_events_ << " ms";
151 ldmx_log(info) << "Number of Events without enough seed points: "
152 << n_missing_;
153 ldmx_log(info) << "Number of Events with only 1 trackID == 1 point: "
154 << n_empty_;
155 ldmx_log(info) << "Total Truth Tracks / Event: " << n_truth_ << "/"
156 << n_events_;
157} // 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 29 of file LinearTruthTracking.cxx.

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

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