LDMX Software
VertexProcessor.cxx
1#include "Tracking/Reco/VertexProcessor.h"
2
3#include <chrono>
4
5#include "Acts/MagneticField/ConstantBField.hpp"
6
7using namespace framework;
8
9namespace tracking {
10namespace reco {
11
12VertexProcessor::VertexProcessor(const std::string &name,
13 framework::Process &process)
14 : framework::Producer(name, process) {}
15
17 gctx_ = Acts::GeometryContext();
18 bctx_ = Acts::MagneticFieldContext();
19
20 h_m_ = new TH1F("m", "m", 100, 0., 1.);
21 h_m_truth_filter_ = new TH1F("m_filter", "m", 100, 0., 1.);
22 h_m_truth_ = new TH1F("m_truth", "m_truth", 100, 0., 1.);
23
24 /*
25 * this is unused, should it be? FIXME
26 auto localToGlobalBin_xyz = [](std::array<size_t, 3> bins,
27 std::array<size_t, 3> sizes) {
28 return (bins[0] * (sizes[1] * sizes[2]) + bins[1] * sizes[2] +
29 bins[2]); // xyz - field space
30 // return (bins[1] * (sizes[2] * sizes[0]) + bins[2] * sizes[0] + bins[0]);
31 // //zxy
32 };
33 */
34
35 // Setup a interpolated bfield map
36 sp_interpolated_b_field_ =
37 std::make_shared<InterpolatedMagneticField3>(loadDefaultBField(
38 field_map_, defaultTransformPos, defaultTransformBField));
39
40 ldmx_log(info) << "Check if nullptr::" << sp_interpolated_b_field_.get();
41}
42
44 // TODO:: the bfield map should be taken automatically
45 field_map_ = parameters.get<std::string>("field_map");
46
47 trk_coll_name_ = parameters.get<std::string>("trk_coll_name", "Tracks");
48
49 input_pass_name_ = parameters.get<std::string>("input_pass_name");
50}
51
53 // TODO:: Move this to an external file
54 // And move all this to a single time per processor not for each event!!
55
56 nevents_++;
57 auto start = std::chrono::high_resolution_clock::now();
58 auto &&stepper = Acts::EigenStepper<>{sp_interpolated_b_field_};
59
60 // Set up propagator with void navigator
61 propagator_ = std::make_shared<VoidPropagator>(stepper);
62
63 // Track linearizer in the proximity of the vertex location
64 using Linearizer = Acts::HelicalTrackLinearizer;
65 Linearizer::Config linearizer_config;
66 linearizer_config.bField = sp_interpolated_b_field_;
67 linearizer_config.propagator = propagator_;
68 Linearizer linearizer(linearizer_config);
69
70 // Set up Billoir Vertex Fitter
71 using VertexFitter = Acts::FullBilloirVertexFitter;
72
73 VertexFitter::Config vertex_fitter_cfg;
74
75 VertexFitter billoir_fitter(vertex_fitter_cfg);
76
77 // VertexFitter::State state(sp_interpolated_bField_->makeCache(bctx_));
78
79 // Unconstrained fit
80 // See
81 // https://github.com/acts-project/acts/blob/main/Tests/UnitTests/Core/Vertexing/FullBilloirVertexFitterTests.cpp#L149
82 // For constraint implementation
83
84 Acts::VertexingOptions vf_options(gctx_, bctx_);
85
86 // Retrieve the track collection
87 const std::vector<ldmx::Track> tracks =
88 event.getCollection<ldmx::Track>(trk_coll_name_, input_pass_name_);
89
90 // Retrieve the truth seeds
91 const std::vector<ldmx::Track> seeds =
92 event.getCollection<ldmx::Track>("RecoilTruthSeeds", input_pass_name_);
93
94 if (tracks.size() < 1) return;
95
96 // Transform the EDM ldmx::tracks to the format needed by ACTS
97 std::vector<Acts::BoundTrackParameters> billoir_tracks;
98
99 // TODO:: The perigee surface should be common between all tracks.
100 // So should only be created once in principle.
101 // There should be no perigeeSurface2
102
103 std::shared_ptr<Acts::PerigeeSurface> perigee_surface =
104 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
105 tracks.front().getPerigeeX(), tracks.front().getPerigeeY(),
106 tracks.front().getPerigeeZ()));
107
108 for (unsigned int i_track = 0; i_track < tracks.size(); i_track++) {
109 Acts::BoundVector param_vec;
110 param_vec << tracks.at(i_track).getD0(), tracks.at(i_track).getZ0(),
111 tracks.at(i_track).getPhi(), tracks.at(i_track).getTheta(),
112 tracks.at(i_track).getQoP(), tracks.at(i_track).getT();
113
114 Acts::BoundSquareMatrix cov_mat =
115 tracking::sim::utils::unpackCov(tracks.at(i_track).getPerigeeCov());
116 auto part{Acts::GenericParticleHypothesis(Acts::ParticleHypothesis(
117 Acts::PdgParticle(tracks.at(i_track).getPdgID())))};
118 billoir_tracks.push_back(Acts::BoundTrackParameters(
119 perigee_surface, param_vec, std::move(cov_mat), part));
120 }
121
122 // Select exactly 2 tracks
123 if (billoir_tracks.size() != 2) {
124 return;
125 }
126
127 if (billoir_tracks.at(0).charge() * billoir_tracks.at(1).charge() > 0) return;
128
129 // Pion mass hypothesis
130 double pion_mass = 139.570 * Acts::UnitConstants::MeV;
131
132 TLorentzVector p1, p2;
133 p1.SetXYZM(billoir_tracks.at(0).momentum()(0),
134 billoir_tracks.at(0).momentum()(1),
135 billoir_tracks.at(0).momentum()(2), pion_mass);
136
137 p2.SetXYZM(billoir_tracks.at(1).momentum()(0),
138 billoir_tracks.at(1).momentum()(1),
139 billoir_tracks.at(1).momentum()(2), pion_mass);
140
141 std::vector<TLorentzVector> pion_seeds;
142
143 if (seeds.size() == 2) {
144 for (int i_seed = 0; i_seed < seeds.size(); i_seed++) {
145 std::shared_ptr<Acts::PerigeeSurface> perigee_surface2 =
146 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
147 seeds.at(i_seed).getPerigeeX(), seeds.at(i_seed).getPerigeeY(),
148 seeds.at(i_seed).getPerigeeZ()));
149
150 Acts::BoundVector param_vec;
151 param_vec << seeds.at(i_seed).getD0(), seeds.at(i_seed).getZ0(),
152 seeds.at(i_seed).getPhi(), seeds.at(i_seed).getTheta(),
153 seeds.at(i_seed).getQoP(), seeds.at(i_seed).getT();
154
155 Acts::BoundSquareMatrix cov_mat =
156 tracking::sim::utils::unpackCov(seeds.at(i_seed).getPerigeeCov());
157 int pion_pdg_id = 211; // pi+
158 if (seeds.at(i_seed).q() < 0) pion_pdg_id = -211;
159 // BoundTrackParameters needs the particle hypothesis
160 auto part{Acts::GenericParticleHypothesis(
161 Acts::ParticleHypothesis(Acts::PdgParticle(pion_pdg_id)))};
162 auto bound_seed_params = Acts::BoundTrackParameters(
163 perigee_surface, param_vec, std::move(cov_mat), part);
164
165 TLorentzVector pion4v;
166 pion4v.SetXYZM(bound_seed_params.momentum()(0),
167 bound_seed_params.momentum()(1),
168 bound_seed_params.momentum()(2), pion_mass);
169
170 pion_seeds.push_back(pion4v);
171 } // loops on seeds
172
173 h_m_truth_->Fill((pion_seeds.at(0) + pion_seeds.at(1)).M());
174 }
175
176 if ((pion_seeds.size() == 2) &&
177 (pion_seeds.at(0) + pion_seeds.at(1)).M() > 0.490 &&
178 (pion_seeds.at(0) + pion_seeds.at(1)).M() < 0.510) {
179 // Check if the tracks have opposite charge
180 h_m_truth_filter_->Fill((p1 + p2).M());
181 }
182
183 h_m_->Fill((p1 + p2).M());
184
185 auto end = std::chrono::high_resolution_clock::now();
186 // long long microseconds =
187 // std::chrono::duration_cast<std::chrono::microseconds>(end-start).count();
188 auto diff = end - start;
189 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
190}
191
193 TFile *outfile = new TFile("VertexingResults.root", "RECREATE");
194 outfile->cd();
195
196 h_m_->Write();
197 h_m_truth_->Write();
198 h_m_truth_filter_->Write();
199 outfile->Close();
200 delete outfile;
201
202 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(3)
203 << processing_time_ / nevents_ << " ms";
204}
205
206} // namespace reco
207} // namespace tracking
208
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Implements an event buffer system for storing event data.
Definition Event.h:42
Class which represents the process under execution.
Definition Process.h:36
Base class for a module which produces a data product.
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
Implementation of a track object.
Definition Track.h:53
void onProcessEnd() override
Callback for the EventProcessor to take any necessary action when the processing of events finishes,...
Acts::GeometryContext gctx_
The contexts - TODO: they should move to some global location, I guess.
VertexProcessor(const std::string &name, framework::Process &process)
Constructor.
void configure(framework::config::Parameters &parameters) override
Configure the processor using the given user specified parameters.
std::string field_map_
Path to the magnetic field map.
void produce(framework::Event &event) override
Run the processor.
void onProcessStart() override
Callback for the EventProcessor to take any necessary action when the processing of events starts,...
All classes in the ldmx-sw project use this namespace.
The measurement calibrator can be a function or a class/struct able to retrieve the sim hits containe...