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 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.
 
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_b_field_
 
std::shared_ptr< Acts::ConstantBField > b_field_
 
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

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

Base class for a module which produces a data product.
virtual void process(Event &event) final
Processing an event for a Producer is calling produce.

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 class_name 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.get<std::string>("field_map");
83
84 trk_c_name_1_ = parameters.get<std::string>("trk_c_name_1", "TaggerTracks");
85 trk_c_name_2_ = parameters.get<std::string>("trk_c_name_2", "RecoilTracks");
86 input_pass_name_ = parameters.get<std::string>("input_pass_name");
87}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
std::string field_map_
Path to the magnetic field map.
Definition Vertexer.h:75

References framework::config::Parameters::get().

◆ 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 200 of file Vertexer.cxx.

200 {
201 ldmx_log(info) << "Reconstructed " << nvertices_ << " vertices over "
202 << nreconstructable_ << " reconstructable" << std::endl;
203
204 TFile* outfile = new TFile((getName() + ".root").c_str(), "RECREATE");
205 outfile->cd();
206
207 h_delta_d0_->Write();
208 h_delta_z0_->Write();
209 h_delta_p_->Write();
210 h_delta_phi_->Write();
211 h_delta_theta_->Write();
212
213 h_delta_d0_vs_recoil_p_->Write();
214 h_delta_z0_vs_recoil_p_->Write();
215
216 h_td0_vs_rd0_->Write();
217 h_tz0_vs_rz0_->Write();
218
219 outfile->Close();
220 delete outfile;
221}
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_b_field_ =
62 std::make_shared<InterpolatedMagneticField3>(loadDefaultBField(
63 field_map_, defaultTransformPos, defaultTransformBField));
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 b_field_ = std::make_shared<Acts::ConstantBField>(b_field);
68
69 std::cout << "Check if nullptr::" << sp_interpolated_b_field_.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<>{b_field_};
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 89 of file Vertexer.cxx.

89 {
90 nevents_++;
91 // auto start = std::chrono::high_resolution_clock::now();
92
93 // Track linearizer in the proximity of the vertex location
94 using Linearizer = Acts::HelicalTrackLinearizer;
95 Linearizer::Config linearizer_config;
96 linearizer_config.bField = b_field_;
97 linearizer_config.propagator = propagator_;
98 Linearizer linearizer(linearizer_config);
99
100 // Set up Billoir Vertex Fitter
101 using VertexFitter = Acts::FullBilloirVertexFitter;
102
103 // Alternatively one can use
104 // using VertexFitter =
105 // Acts::FullBilloirVertexFitter<tracking::sim::utils::boundTrackParameters,Linearizer>;
106
107 VertexFitter::Config vertex_fitter_cfg;
108 VertexFitter billoir_fitter(vertex_fitter_cfg);
109 // mg Aug 2024 .. State doesn't exist in v36 and isn't used here anyway
110 // VertexFitter::State state(sp_interpolated_bField_->makeCache(bctx_));
111
112 // Unconstrained fit
113 // See
114 // https://github.com/acts-project/acts/blob/main/Tests/UnitTests/Core/Vertexing/FullBilloirVertexFitterTests.cpp#L149
115 // For constraint implementation
116
117 // Acts::VertexingOptions<Acts::BoundTrackParameters> vfOptions(gctx_,
118 // bctx_);
119 // mg Aug 2024 ... VertexingOptions template change in v36
120 Acts::VertexingOptions vf_options(gctx_, bctx_);
121
122 // Retrive the two track collections
123
124 const std::vector<ldmx::Track> tracks_1 =
125 event.getCollection<ldmx::Track>(trk_c_name_1_, input_pass_name_);
126 const std::vector<ldmx::Track> tracks_2 =
127 event.getCollection<ldmx::Track>(trk_c_name_2_, input_pass_name_);
128
129 ldmx_log(debug) << "Retrieved track collections" << std::endl
130 << "Track 1 size:" << tracks_1.size() << std::endl
131 << "Track 2 size:" << tracks_2.size() << std::endl;
132
133 if (tracks_1.size() < 1 || tracks_2.size() < 1) return;
134
135 std::vector<Acts::BoundTrackParameters> billoir_tracks_1, billoir_tracks_2;
136
137 // TODO:: The perigee surface should be common between all tracks.
138
139 Acts::Vector3 perigee_acts = tracking::sim::utils::ldmx2Acts(Acts::Vector3(
140 tracks_1.front().getPerigeeX(), tracks_1.front().getPerigeeY(),
141 tracks_1.front().getPerigeeZ()));
142 std::shared_ptr<Acts::PerigeeSurface> perigee_surface =
143 Acts::Surface::makeShared<Acts::PerigeeSurface>(perigee_acts);
144
145 // Monitoring of tagger and recoil tracks
146 taggerRecoilMonitoring(tracks_1, tracks_2);
147
148 // Start the vertex formation
149 // Form a vertex for each combination of tracks found in the same event
150 // between the two track collections
151
152 for (auto& trk : tracks_1) {
153 billoir_tracks_1.push_back(
154 tracking::sim::utils::boundTrackParameters(trk, perigee_surface));
155 }
156
157 for (auto& trk : tracks_2) {
158 billoir_tracks_2.push_back(
159 tracking::sim::utils::boundTrackParameters(trk, perigee_surface));
160 }
161
162 // std::vector<Acts::Vertex<Acts::BoundTrackParameters> > fit_vertices;
163 std::vector<Acts::Vertex> fit_vertices;
164
165 for (auto& b_trk_1 : billoir_tracks_1) {
166 std::vector<const Acts::BoundTrackParameters*> fit_tracks_ptr;
167
168 for (auto& b_trk_2 : billoir_tracks_2) {
169 fit_tracks_ptr.push_back(&b_trk_1);
170 fit_tracks_ptr.push_back(&b_trk_2);
171
172 ldmx_log(debug) << "Calling vertex fitter" << std::endl
173 << "Track 1 parameters" << std::endl
174 << b_trk_1 << std::endl
175 << "Track 2 parameters" << std::endl
176 << b_trk_2 << std::endl;
177
178 // std::cout << "Perigee Surface" << std::endl;
179 // perigeeSurface->toStream(gctx_, std::cout);
180 // std::cout << std::endl;
181
182 } // loop on second set of tracks
183
184 nreconstructable_++;
185 try {
186 // Acts::Vertex<Acts::BoundTrackParameters> fitVtx =
187 // billoirFitter.fit(fit_tracks_ptr, linearizer, vfOptions,
188 // state).value(); fit_vertices.push_back(fitVtx);
189 nvertices_++;
190
191 } catch (...) {
192 ldmx_log(warn) << "Vertex fit failed" << std::endl;
193 }
194
195 } // loop on first set
196
197 // Convert the vertices in the ldmx EDM and store them
198}
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 223 of file Vertexer.cxx.

225 {
226 // For the moment only check that I have 1 tagger track and one recoil track
227 // To avoid trying to match them
228 // TODO update this logic
229
230 if (tagger_tracks.size() != 1 || recoil_tracks.size() != 1) return;
231
232 ldmx::Track t_trk = tagger_tracks.at(0);
233 ldmx::Track r_trk = recoil_tracks.at(0);
234
235 double t_p, r_p;
236 // these are unsed, should they be? FIXME
237 // double t_d0, r_d0;
238 // double tt_p_phi, r_phi;
239 // double t_theta, r_theta;
240 // double t_z0, r_z0;
241
242 t_p = t_trk.getCharge() / t_trk.getQoP();
243 r_p = r_trk.getCharge() / r_trk.getQoP();
244
245 h_delta_d0_->Fill(t_trk.getD0() - r_trk.getD0());
246 h_delta_z0_->Fill(t_trk.getZ0() - r_trk.getZ0());
247 h_delta_p_->Fill(t_p - r_p);
248 h_delta_phi_->Fill(t_trk.getPhi() - r_trk.getPhi());
249 h_delta_theta_->Fill(t_trk.getTheta() - r_trk.getTheta());
250
251 // differential plots
252
253 h_delta_d0_vs_recoil_p_->Fill(r_p, t_trk.getD0() - r_trk.getD0());
254 h_delta_z0_vs_recoil_p_->Fill(r_p, t_trk.getZ0() - r_trk.getZ0());
255
256 //"beamspot"
257 h_td0_vs_rd0_->Fill(r_trk.getD0(), t_trk.getD0());
258 h_tz0_vs_rz0_->Fill(r_trk.getZ0(), t_trk.getZ0());
259
260 //"pT"
261 // TODO Transverse momentum should obtained orthogonal to the B-Field
262 // direction This assumes to be along Z (which is not very accurate)
263
264 // std::vector<double> r_mom = r_trk.getMomentum();
265 // std::vector<double> t_mom = t_trk.getMomentum();
266
267 // I assume to have a single photon being emitted in the target: I use
268 // momentum conservation p_photon = p_beam - p_recoil
269
270 // h_gamma_px->Fill(t_mom[0] - r_mom[0]);
271 // h_gamma_py->Fill(t_mom[1] - r_mom[1]);
272 // h_gamma_pz->Fill(t_mom[2] - r_mom[2]);
273}

Member Data Documentation

◆ b_field_

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

Definition at line 72 of file Vertexer.h.

◆ bctx_

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

Definition at line 66 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_b_field_

std::shared_ptr<InterpolatedMagneticField3> tracking::reco::Vertexer::sp_interpolated_b_field_
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: