LDMX Software
tracking::reco::Vertexer Class Reference

Public Member Functions

 Vertexer (const std::string &name, framework::Process &process)
 
void onProcessStart () override
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
void onProcessEnd () override
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
void configure (framework::config::Parameters &parameters) override
 Callback for the EventProcessor to configure itself from the given set of parameters.
 
void produce (framework::Event &event) override
 Process the event and put new data products into it.
 
void TaggerRecoilMonitoring (const std::vector< ldmx::Track > &tagger_tracks, const std::vector< ldmx::Track > &recoil_tracks)
 
- 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.
 
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.
 

Private Attributes

Acts::GeometryContext gctx_
 
Acts::MagneticFieldContext bctx_
 
int nevents_ {0}
 
int nvertices_ {0}
 
int nreconstructable_ {0}
 
std::shared_ptr< InterpolatedMagneticField3 > sp_interpolated_bField_
 
std::shared_ptr< Acts::ConstantBField > bField_
 
std::string field_map_
 Path to the magnetic field map.
 
std::string trk_c_name_1 {"TaggerTracks"}
 
std::string trk_c_name_2 {"RecoilTracks"}
 
std::string input_pass_name_ {""}
 
std::shared_ptr< VoidPropagator > propagator_
 
TH1F * h_delta_d0
 
TH1F * h_delta_z0
 
TH1F * h_delta_p
 
TH1F * h_delta_phi
 
TH1F * h_delta_theta
 
TH2F * h_delta_d0_vs_recoil_p
 
TH2F * h_delta_z0_vs_recoil_p
 
TH2F * h_td0_vs_rd0
 
TH2F * h_tz0_vs_rz0
 

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.
 
- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 
- 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.
 

Detailed Description

Definition at line 48 of file Vertexer.h.

Constructor & Destructor Documentation

◆ Vertexer()

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

Definition at line 17 of file Vertexer.cxx.

18 : framework::Producer(name, process) {}
Base class for a module which produces a data product.

Member Function Documentation

◆ configure()

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

Callback for the EventProcessor to configure itself from the given set of parameters.

The parameters a processor has access to are the member variables of the python class in the sequence that has className equal to the EventProcessor class name.

For an example, look at MyProcessor.

Parameters
parametersParameters for configuration.

Reimplemented from framework::EventProcessor.

Definition at line 80 of file Vertexer.cxx.

80 {
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}
std::string field_map_
Path to the magnetic field map.
Definition Vertexer.h:75

◆ onProcessEnd()

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

Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.

Reimplemented from framework::EventProcessor.

Definition at line 202 of file Vertexer.cxx.

202 {
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}
std::string getName() const
Get the processor name.

◆ onProcessStart()

void tracking::reco::Vertexer::onProcessStart ( )
overridevirtual

Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.

Reimplemented from framework::EventProcessor.

Definition at line 20 of file Vertexer.cxx.

20 {
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}

◆ produce()

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

Process the event and put new data products into it.

Parameters
eventThe Event to process.

Implements framework::Producer.

Definition at line 92 of file Vertexer.cxx.

92 {
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}
Implementation of a track object.
Definition Track.h:53

◆ TaggerRecoilMonitoring()

void tracking::reco::Vertexer::TaggerRecoilMonitoring ( const std::vector< ldmx::Track > & tagger_tracks,
const std::vector< ldmx::Track > & recoil_tracks )

Definition at line 225 of file Vertexer.cxx.

227 {
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}

Member Data Documentation

◆ bctx_

Acts::MagneticFieldContext tracking::reco::Vertexer::bctx_
private

Definition at line 66 of file Vertexer.h.

◆ bField_

std::shared_ptr<Acts::ConstantBField> tracking::reco::Vertexer::bField_
private

Definition at line 72 of file Vertexer.h.

◆ field_map_

std::string tracking::reco::Vertexer::field_map_
private

Path to the magnetic field map.

Definition at line 75 of file Vertexer.h.

◆ gctx_

Acts::GeometryContext tracking::reco::Vertexer::gctx_
private

Definition at line 65 of file Vertexer.h.

◆ h_delta_d0

TH1F* tracking::reco::Vertexer::h_delta_d0
private

Definition at line 83 of file Vertexer.h.

◆ h_delta_d0_vs_recoil_p

TH2F* tracking::reco::Vertexer::h_delta_d0_vs_recoil_p
private

Definition at line 90 of file Vertexer.h.

◆ h_delta_p

TH1F* tracking::reco::Vertexer::h_delta_p
private

Definition at line 85 of file Vertexer.h.

◆ h_delta_phi

TH1F* tracking::reco::Vertexer::h_delta_phi
private

Definition at line 87 of file Vertexer.h.

◆ h_delta_theta

TH1F* tracking::reco::Vertexer::h_delta_theta
private

Definition at line 88 of file Vertexer.h.

◆ h_delta_z0

TH1F* tracking::reco::Vertexer::h_delta_z0
private

Definition at line 84 of file Vertexer.h.

◆ h_delta_z0_vs_recoil_p

TH2F* tracking::reco::Vertexer::h_delta_z0_vs_recoil_p
private

Definition at line 91 of file Vertexer.h.

◆ h_td0_vs_rd0

TH2F* tracking::reco::Vertexer::h_td0_vs_rd0
private

Definition at line 93 of file Vertexer.h.

◆ h_tz0_vs_rz0

TH2F* tracking::reco::Vertexer::h_tz0_vs_rz0
private

Definition at line 94 of file Vertexer.h.

◆ input_pass_name_

std::string tracking::reco::Vertexer::input_pass_name_ {""}
private

Definition at line 79 of file Vertexer.h.

79{""};

◆ nevents_

int tracking::reco::Vertexer::nevents_ {0}
private

Definition at line 68 of file Vertexer.h.

68{0};

◆ nreconstructable_

int tracking::reco::Vertexer::nreconstructable_ {0}
private

Definition at line 70 of file Vertexer.h.

70{0};

◆ nvertices_

int tracking::reco::Vertexer::nvertices_ {0}
private

Definition at line 69 of file Vertexer.h.

69{0};

◆ propagator_

std::shared_ptr<VoidPropagator> tracking::reco::Vertexer::propagator_
private

Definition at line 80 of file Vertexer.h.

◆ sp_interpolated_bField_

std::shared_ptr<InterpolatedMagneticField3> tracking::reco::Vertexer::sp_interpolated_bField_
private

Definition at line 71 of file Vertexer.h.

◆ trk_c_name_1

std::string tracking::reco::Vertexer::trk_c_name_1 {"TaggerTracks"}
private

Definition at line 77 of file Vertexer.h.

77{"TaggerTracks"};

◆ trk_c_name_2

std::string tracking::reco::Vertexer::trk_c_name_2 {"RecoilTracks"}
private

Definition at line 78 of file Vertexer.h.

78{"RecoilTracks"};

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