LDMX Software
Vertexer.cxx
1#include "Tracking/Reco/Vertexer.h"
2
3#include <chrono>
4
5#include "TFile.h"
6using namespace framework;
7
8// This producer takes in input two track collections and forms all possible
9// vertices from those It can be used to match the tagger and recoil tracks at
10// the target and form the beamspot
11
12// TODO Move all the performance monitor to another processor
13
14namespace tracking {
15namespace reco {
16
17Vertexer::Vertexer(const std::string& name, framework::Process& process)
18 : framework::Producer(name, process) {}
19
20void Vertexer::onProcessStart() {
21 // Monitoring plots
22
23 double d0min = -2;
24 double d0max = 2;
25 double z0min = -2;
26 double z0max = 2;
27
28 h_delta_d0 = new TH1F("h_delta_d0", "h_delta_d0", 400, d0min, d0max);
29 h_delta_z0 = new TH1F("h_delta_z0", "h_delta_z0", 200, z0min, z0max);
30 h_delta_p = new TH1F("h_delta_p", "h_delta_p", 200, -1, 4);
31 // h_delta_pT_vsP = new TH2D("h_delta_pT_vs_p","h_delta_pT_v_p",200,)
32 h_delta_phi = new TH1F("h_delta_phi", "h_delta_phi", 400, -0.2, 0.2);
33 h_delta_theta = new TH1F("h_delta_theta", "h_delta_theta", 200, -0.1, 0.1);
34
35 h_delta_d0_vs_recoil_p =
36 new TH2F("h_delta_d0_vs_recoil_p", "h_delta_d0_vs_recoil_p", 200, 0, 5,
37 400, -1, 1);
38 h_delta_z0_vs_recoil_p =
39 new TH2F("h_delta_z0_vs_recoil_p", "h_delta_z0_vs_recoil_p", 200, 0, 5,
40 400, -1, 1);
41
42 h_td0_vs_rd0 =
43 new TH2F("h_td0_vs_rd0", "h_td0_vs_rd0", 100, -40, 40, 100, -40, 40);
44 h_tz0_vs_rz0 =
45 new TH2F("h_tz0_vs_rz0", "h_tz0_vs_rz0", 100, -40, 40, 100, -40, 40);
46
47 gctx_ = Acts::GeometryContext();
48 bctx_ = Acts::MagneticFieldContext();
49
50 /*
51 * This is unused now, should it be?
52 auto localToGlobalBin_xyz = [](std::array<size_t, 3> bins,
53 std::array<size_t, 3> sizes) {
54 return (bins[0] * (sizes[1] * sizes[2]) + bins[1] * sizes[2] +
55 bins[2]); // xyz - field space
56 // return (bins[1] * (sizes[2] * sizes[0]) + bins[2] * sizes[0] + bins[0]);
57 // //zxy
58 };
59 */
60
61 sp_interpolated_bField_ =
62 std::make_shared<InterpolatedMagneticField3>(loadDefaultBField(
63 field_map_, default_transformPos, default_transformBField));
64
65 // There is a sign issue between the vertexing and the perigee representation
66 Acts::Vector3 b_field(0., 0., -1.5 * Acts::UnitConstants::T);
67 bField_ = std::make_shared<Acts::ConstantBField>(b_field);
68
69 std::cout << "Check if nullptr::" << sp_interpolated_bField_.get()
70 << std::endl;
71
72 // Set up propagator with void navigator
73 // auto&& stepper = Acts::EigenStepper<>{sp_interpolated_bField_};
74 // propagator_ = std::make_shared<VoidPropagator>(stepper);
75
76 auto&& stepper_const = Acts::EigenStepper<>{bField_};
77 propagator_ = std::make_shared<VoidPropagator>(stepper_const);
78}
79
80void Vertexer::configure(framework::config::Parameters& parameters) {
81 // TODO:: the bfield map should be taken automatically
82 field_map_ = parameters.getParameter<std::string>("field_map");
83
84 trk_c_name_1 =
85 parameters.getParameter<std::string>("trk_c_name_1", "TaggerTracks");
86 trk_c_name_2 =
87 parameters.getParameter<std::string>("trk_c_name_2", "RecoilTracks");
88 input_pass_name_ =
89 parameters.getParameter<std::string>("input_pass_name", "");
90}
91
92void Vertexer::produce(framework::Event& event) {
93 nevents_++;
94 // auto start = std::chrono::high_resolution_clock::now();
95
96 // Track linearizer in the proximity of the vertex location
97 using Linearizer = Acts::HelicalTrackLinearizer;
98 Linearizer::Config linearizerConfig;
99 linearizerConfig.bField = bField_;
100 linearizerConfig.propagator = propagator_;
101 Linearizer linearizer(linearizerConfig);
102
103 // Set up Billoir Vertex Fitter
104 using VertexFitter = Acts::FullBilloirVertexFitter;
105
106 // Alternatively one can use
107 // using VertexFitter =
108 // Acts::FullBilloirVertexFitter<tracking::sim::utils::boundTrackParameters,Linearizer>;
109
110 VertexFitter::Config vertexFitterCfg;
111 VertexFitter billoirFitter(vertexFitterCfg);
112 // mg Aug 2024 .. State doesn't exist in v36 and isn't used here anyway
113 // VertexFitter::State state(sp_interpolated_bField_->makeCache(bctx_));
114
115 // Unconstrained fit
116 // See
117 // https://github.com/acts-project/acts/blob/main/Tests/UnitTests/Core/Vertexing/FullBilloirVertexFitterTests.cpp#L149
118 // For constraint implementation
119
120 // Acts::VertexingOptions<Acts::BoundTrackParameters> vfOptions(gctx_,
121 // bctx_);
122 // mg Aug 2024 ... VertexingOptions template change in v36
123 Acts::VertexingOptions vfOptions(gctx_, bctx_);
124
125 // Retrive the two track collections
126
127 const std::vector<ldmx::Track> tracks_1 =
128 event.getCollection<ldmx::Track>(trk_c_name_1, input_pass_name_);
129 const std::vector<ldmx::Track> tracks_2 =
130 event.getCollection<ldmx::Track>(trk_c_name_2, input_pass_name_);
131
132 ldmx_log(debug) << "Retrieved track collections" << std::endl
133 << "Track 1 size:" << tracks_1.size() << std::endl
134 << "Track 2 size:" << tracks_2.size() << std::endl;
135
136 if (tracks_1.size() < 1 || tracks_2.size() < 1) return;
137
138 std::vector<Acts::BoundTrackParameters> billoir_tracks_1, billoir_tracks_2;
139
140 // TODO:: The perigee surface should be common between all tracks.
141
142 std::shared_ptr<Acts::PerigeeSurface> perigeeSurface =
143 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
144 tracks_1.front().getPerigeeX(), tracks_1.front().getPerigeeY(),
145 tracks_1.front().getPerigeeZ()));
146
147 // Monitoring of tagger and recoil tracks
148 TaggerRecoilMonitoring(tracks_1, tracks_2);
149
150 // Start the vertex formation
151 // Form a vertex for each combination of tracks found in the same event
152 // between the two track collections
153
154 for (auto& trk : tracks_1) {
155 billoir_tracks_1.push_back(
156 tracking::sim::utils::boundTrackParameters(trk, perigeeSurface));
157 }
158
159 for (auto& trk : tracks_2) {
160 billoir_tracks_2.push_back(
161 tracking::sim::utils::boundTrackParameters(trk, perigeeSurface));
162 }
163
164 // std::vector<Acts::Vertex<Acts::BoundTrackParameters> > fit_vertices;
165 std::vector<Acts::Vertex> fit_vertices;
166
167 for (auto& b_trk_1 : billoir_tracks_1) {
168 std::vector<const Acts::BoundTrackParameters*> fit_tracks_ptr;
169
170 for (auto& b_trk_2 : billoir_tracks_2) {
171 fit_tracks_ptr.push_back(&b_trk_1);
172 fit_tracks_ptr.push_back(&b_trk_2);
173
174 ldmx_log(debug) << "Calling vertex fitter" << std::endl
175 << "Track 1 parameters" << std::endl
176 << b_trk_1 << std::endl
177 << "Track 2 parameters" << std::endl
178 << b_trk_2 << std::endl;
179
180 // std::cout << "Perigee Surface" << std::endl;
181 // perigeeSurface->toStream(gctx_, std::cout);
182 // std::cout << std::endl;
183
184 } // loop on second set of tracks
185
186 nreconstructable_++;
187 try {
188 // Acts::Vertex<Acts::BoundTrackParameters> fitVtx =
189 // billoirFitter.fit(fit_tracks_ptr, linearizer, vfOptions,
190 // state).value(); fit_vertices.push_back(fitVtx);
191 nvertices_++;
192
193 } catch (...) {
194 ldmx_log(warn) << "Vertex fit failed" << std::endl;
195 }
196
197 } // loop on first set
198
199 // Convert the vertices in the ldmx EDM and store them
200}
201
202void Vertexer::onProcessEnd() {
203 ldmx_log(info) << "Reconstructed " << nvertices_ << " vertices over "
204 << nreconstructable_ << " reconstructable" << std::endl;
205
206 TFile* outfile_ = new TFile((getName() + ".root").c_str(), "RECREATE");
207 outfile_->cd();
208
209 h_delta_d0->Write();
210 h_delta_z0->Write();
211 h_delta_p->Write();
212 h_delta_phi->Write();
213 h_delta_theta->Write();
214
215 h_delta_d0_vs_recoil_p->Write();
216 h_delta_z0_vs_recoil_p->Write();
217
218 h_td0_vs_rd0->Write();
219 h_tz0_vs_rz0->Write();
220
221 outfile_->Close();
222 delete outfile_;
223}
224
225void Vertexer::TaggerRecoilMonitoring(
226 const std::vector<ldmx::Track>& tagger_tracks,
227 const std::vector<ldmx::Track>& recoil_tracks) {
228 // For the moment only check that I have 1 tagger track and one recoil track
229 // To avoid trying to match them
230 // TODO update this logic
231
232 if (tagger_tracks.size() != 1 || recoil_tracks.size() != 1) return;
233
234 ldmx::Track t_trk = tagger_tracks.at(0);
235 ldmx::Track r_trk = recoil_tracks.at(0);
236
237 double t_p, r_p;
238 // these are unsed, should they be? FIXME
239 // double t_d0, r_d0;
240 // double tt_p_phi, r_phi;
241 // double t_theta, r_theta;
242 // double t_z0, r_z0;
243
244 t_p = t_trk.q() / t_trk.getQoP();
245 r_p = r_trk.q() / r_trk.getQoP();
246
247 h_delta_d0->Fill(t_trk.getD0() - r_trk.getD0());
248 h_delta_z0->Fill(t_trk.getZ0() - r_trk.getZ0());
249 h_delta_p->Fill(t_p - r_p);
250 h_delta_phi->Fill(t_trk.getPhi() - r_trk.getPhi());
251 h_delta_theta->Fill(t_trk.getTheta() - r_trk.getTheta());
252
253 // differential plots
254
255 h_delta_d0_vs_recoil_p->Fill(r_p, t_trk.getD0() - r_trk.getD0());
256 h_delta_z0_vs_recoil_p->Fill(r_p, t_trk.getZ0() - r_trk.getZ0());
257
258 //"beamspot"
259 h_td0_vs_rd0->Fill(r_trk.getD0(), t_trk.getD0());
260 h_tz0_vs_rz0->Fill(r_trk.getZ0(), t_trk.getZ0());
261
262 //"pT"
263 // TODO Transverse momentum should obtained orthogonal to the B-Field
264 // direction This assumes to be along Z (which is not very accurate)
265
266 // std::vector<double> r_mom = r_trk.getMomentum();
267 // std::vector<double> t_mom = t_trk.getMomentum();
268
269 // I assume to have a single photon being emitted in the target: I use
270 // momentum conservation p_photon = p_beam - p_recoil
271
272 // h_gamma_px->Fill(t_mom[0] - r_mom[0]);
273 // h_gamma_py->Fill(t_mom[1] - r_mom[1]);
274 // h_gamma_pz->Fill(t_mom[2] - r_mom[2]);
275}
276
277} // namespace reco
278} // namespace tracking
279
280DECLARE_PRODUCER_NS(tracking::reco, Vertexer)
#define DECLARE_PRODUCER_NS(NS, 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
Implementation of a track object.
Definition Track.h:53
All classes in the ldmx-sw project use this namespace.
Definition PerfDict.cxx:45
The measurement calibrator can be a function or a class/struct able to retrieve the sim hits containe...