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_truthFilter_ = 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_bField_ =
37 std::make_shared<InterpolatedMagneticField3>(loadDefaultBField(
38 field_map_, default_transformPos, default_transformBField));
39
40 ldmx_log(info) << "Check if nullptr::" << sp_interpolated_bField_.get();
41}
42
44 // TODO:: the bfield map should be taken automatically
45 field_map_ = parameters.getParameter<std::string>("field_map");
46
47 trk_coll_name_ =
48 parameters.getParameter<std::string>("trk_coll_name", "Tracks");
49
50 input_pass_name_ =
51 parameters.getParameter<std::string>("input_pass_name", "");
52}
53
55 // TODO:: Move this to an external file
56 // And move all this to a single time per processor not for each event!!
57
58 nevents_++;
59 auto start = std::chrono::high_resolution_clock::now();
60 auto &&stepper = Acts::EigenStepper<>{sp_interpolated_bField_};
61
62 // Set up propagator with void navigator
63 propagator_ = std::make_shared<VoidPropagator>(stepper);
64
65 // Track linearizer in the proximity of the vertex location
66 using Linearizer = Acts::HelicalTrackLinearizer;
67 Linearizer::Config linearizerConfig;
68 linearizerConfig.bField = sp_interpolated_bField_;
69 linearizerConfig.propagator = propagator_;
70 Linearizer linearizer(linearizerConfig);
71
72 // Set up Billoir Vertex Fitter
73 using VertexFitter = Acts::FullBilloirVertexFitter;
74
75 VertexFitter::Config vertexFitterCfg;
76
77 VertexFitter billoirFitter(vertexFitterCfg);
78
79 // VertexFitter::State state(sp_interpolated_bField_->makeCache(bctx_));
80
81 // Unconstrained fit
82 // See
83 // https://github.com/acts-project/acts/blob/main/Tests/UnitTests/Core/Vertexing/FullBilloirVertexFitterTests.cpp#L149
84 // For constraint implementation
85
86 Acts::VertexingOptions vfOptions(gctx_, bctx_);
87
88 // Retrieve the track collection
89 const std::vector<ldmx::Track> tracks =
90 event.getCollection<ldmx::Track>(trk_coll_name_, input_pass_name_);
91
92 // Retrieve the truth seeds
93 const std::vector<ldmx::Track> seeds =
94 event.getCollection<ldmx::Track>("RecoilTruthSeeds", input_pass_name_);
95
96 if (tracks.size() < 1) return;
97
98 // Transform the EDM ldmx::tracks to the format needed by ACTS
99 std::vector<Acts::BoundTrackParameters> billoir_tracks;
100
101 // TODO:: The perigee surface should be common between all tracks.
102 // So should only be created once in principle.
103 // There should be no perigeeSurface2
104
105 std::shared_ptr<Acts::PerigeeSurface> perigeeSurface =
106 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
107 tracks.front().getPerigeeX(), tracks.front().getPerigeeY(),
108 tracks.front().getPerigeeZ()));
109
110 for (unsigned int iTrack = 0; iTrack < tracks.size(); iTrack++) {
111 Acts::BoundVector paramVec;
112 paramVec << tracks.at(iTrack).getD0(), tracks.at(iTrack).getZ0(),
113 tracks.at(iTrack).getPhi(), tracks.at(iTrack).getTheta(),
114 tracks.at(iTrack).getQoP(), tracks.at(iTrack).getT();
115
116 Acts::BoundSquareMatrix covMat =
117 tracking::sim::utils::unpackCov(tracks.at(iTrack).getPerigeeCov());
118 auto part{Acts::GenericParticleHypothesis(Acts::ParticleHypothesis(
119 Acts::PdgParticle(tracks.at(iTrack).getPdgID())))};
120 billoir_tracks.push_back(Acts::BoundTrackParameters(
121 perigeeSurface, paramVec, std::move(covMat), part));
122 }
123
124 // Select exactly 2 tracks
125 if (billoir_tracks.size() != 2) {
126 return;
127 }
128
129 if (billoir_tracks.at(0).charge() * billoir_tracks.at(1).charge() > 0) return;
130
131 // Pion mass hypothesis
132 double pion_mass = 139.570 * Acts::UnitConstants::MeV;
133
134 TLorentzVector p1, p2;
135 p1.SetXYZM(billoir_tracks.at(0).momentum()(0),
136 billoir_tracks.at(0).momentum()(1),
137 billoir_tracks.at(0).momentum()(2), pion_mass);
138
139 p2.SetXYZM(billoir_tracks.at(1).momentum()(0),
140 billoir_tracks.at(1).momentum()(1),
141 billoir_tracks.at(1).momentum()(2), pion_mass);
142
143 std::vector<TLorentzVector> pion_seeds;
144
145 if (seeds.size() == 2) {
146 for (int iSeed = 0; iSeed < seeds.size(); iSeed++) {
147 std::shared_ptr<Acts::PerigeeSurface> perigeeSurface2 =
148 Acts::Surface::makeShared<Acts::PerigeeSurface>(Acts::Vector3(
149 seeds.at(iSeed).getPerigeeX(), seeds.at(iSeed).getPerigeeY(),
150 seeds.at(iSeed).getPerigeeZ()));
151
152 Acts::BoundVector paramVec;
153 paramVec << seeds.at(iSeed).getD0(), seeds.at(iSeed).getZ0(),
154 seeds.at(iSeed).getPhi(), seeds.at(iSeed).getTheta(),
155 seeds.at(iSeed).getQoP(), seeds.at(iSeed).getT();
156
157 Acts::BoundSquareMatrix covMat =
158 tracking::sim::utils::unpackCov(seeds.at(iSeed).getPerigeeCov());
159 int pionPdgId = 211; // pi+
160 if (seeds.at(iSeed).q() < 0) pionPdgId = -211;
161 // BoundTrackParameters needs the particle hypothesis
162 auto part{Acts::GenericParticleHypothesis(
163 Acts::ParticleHypothesis(Acts::PdgParticle(pionPdgId)))};
164 auto boundSeedParams = Acts::BoundTrackParameters(
165 perigeeSurface, paramVec, std::move(covMat), part);
166
167 TLorentzVector pion4v;
168 pion4v.SetXYZM(boundSeedParams.momentum()(0),
169 boundSeedParams.momentum()(1),
170 boundSeedParams.momentum()(2), pion_mass);
171
172 pion_seeds.push_back(pion4v);
173 } // loops on seeds
174
175 h_m_truth_->Fill((pion_seeds.at(0) + pion_seeds.at(1)).M());
176 }
177
178 if ((pion_seeds.size() == 2) &&
179 (pion_seeds.at(0) + pion_seeds.at(1)).M() > 0.490 &&
180 (pion_seeds.at(0) + pion_seeds.at(1)).M() < 0.510) {
181 // Check if the tracks have opposite charge
182 h_m_truthFilter_->Fill((p1 + p2).M());
183 }
184
185 h_m_->Fill((p1 + p2).M());
186
187 auto end = std::chrono::high_resolution_clock::now();
188 // long long microseconds =
189 // std::chrono::duration_cast<std::chrono::microseconds>(end-start).count();
190 auto diff = end - start;
191 processing_time_ += std::chrono::duration<double, std::milli>(diff).count();
192}
193
195 TFile *outfile = new TFile("VertexingResults.root", "RECREATE");
196 outfile->cd();
197
198 h_m_->Write();
199 h_m_truth_->Write();
200 h_m_truthFilter_->Write();
201 outfile->Close();
202 delete outfile;
203
204 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(3)
205 << processing_time_ / nevents_ << " ms";
206}
207
208} // namespace reco
209} // namespace tracking
210
211DECLARE_PRODUCER_NS(tracking::reco, VertexProcessor)
#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
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.
Definition PerfDict.cxx:45
The measurement calibrator can be a function or a class/struct able to retrieve the sim hits containe...