LDMX Software
tracking::reco::TruthSeedProcessor Class Reference

Create a track seed using truth information extracted from the corresponding SimParticle or SimTrackerHit. More...

#include <TruthSeedProcessor.h>

Public Member Functions

 TruthSeedProcessor (const std::string &name, framework::Process &process)
 Constructor.
 
virtual ~TruthSeedProcessor ()=default
 Destructor.
 
void configure (framework::config::Parameters &parameters) override
 Callback for the EventProcessor to configure itself from the given set of parameters.
 
void onProcessStart () override
 Callback for the EventProcessor to take any necessary action when the processing of events starts.
 
void onNewRun (const ldmx::RunHeader &rh) override
 onNewRun is the first function called for each processor after the conditions are fully configured and accessible.
 
void produce (framework::Event &event) override
 Main loop that creates the seed tracks for both the tagger and recoil tracker.
 
- Public Member Functions inherited from tracking::reco::TrackingGeometryUser
 TrackingGeometryUser (const std::string &name, framework::Process &p)
 
- 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 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.
 
virtual void onProcessEnd ()
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
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 Member Functions

void makeHitCountMap (const std::vector< ldmx::SimTrackerHit > &sim_hits, std::map< int, std::vector< int > > &hit_count_map)
 Create a mapping from the selected scoring plane hit objects to the number of hits they associated particle creates in the tracker.
 
void createTruthTrack (const ldmx::SimParticle &particle, ldmx::Track &trk, const std::shared_ptr< Acts::Surface > &target_surface)
 Use the vertex position of the SimParticle to extract (x_, y_, z_, px, py, pz, q) and create a track seed.
 
void createTruthTrack (const ldmx::SimParticle &particle, const ldmx::SimTrackerHit &hit, ldmx::Track &trk, const std::shared_ptr< Acts::Surface > &target_surface)
 Use the scoring plane hit at the target to extract (x_, y_, z_, px, py, pz) and create a track seed.
 
void createTruthTrack (const std::vector< double > &pos_vec, const std::vector< double > &p_vec, int charge, ldmx::Track &trk, const std::shared_ptr< Acts::Surface > &target_surface)
 Create a seed track from the given position, momentum and charge.
 
bool scoringPlaneHitFilter (const ldmx::SimTrackerHit &hit, const std::vector< ldmx::SimTrackerHit > &ecal_sp_hits)
 Filter that checks if a scoring plane passes specified momentum cuts as well as if the associated SimParticle hits the ECal.
 
ldmx::Track seedFromTruth (const ldmx::Track &tt, bool seed_smearing)
 Create a track seed from a truth track applying a smearing to the truth parameters as well as an inflation to the covariance matrix.
 
ldmx::Track recoilFullSeed (const ldmx::SimParticle &particle, const int trackID, const ldmx::SimTrackerHit &hit, const ldmx::SimTrackerHit &ecal_hit, const std::map< int, std::vector< int > > &hit_count_map, const std::shared_ptr< Acts::Surface > &origin_surface, const std::shared_ptr< Acts::Surface > &target_surface, const std::shared_ptr< Acts::Surface > &ecal_surface)
 
ldmx::Track taggerFullSeed (const ldmx::SimParticle &beam_electron, const int trackID, const ldmx::SimTrackerHit &hit, const std::map< int, std::vector< int > > &hit_count_map, const std::shared_ptr< Acts::Surface > &origin_surface, const std::shared_ptr< Acts::Surface > &target_surface)
 This method retrieves the beam electron and forms a full seed The seed parameters are the truth parameters from the beam electron stored at the beam origin Additionally, the foolowing track states are stored ts_smeared : the truth smeared perigee state at the beam origin ts_truth_target : the truth on-surface state at the target Linear extrapolations are done from the origin of the particle to the reference surfaces This track also contains the list of hits belonging to the beam electron on the sensitive surfaces on the tagger tracker, for acceptance studies.
 

Private Attributes

Acts::GeometryContext gctx_
 The ACTS geometry context properly.
 
std::vector< int > pdg_ids_ {11}
 pdg_ids of the particles we want to select for the seeds
 
std::string scoring_hits_coll_name_ {"TargetScoringPlaneHits"}
 Which scoring plane hits to use for the truth seeds generation.
 
std::string ecal_sp_coll_name_ {"EcalScoringPlaneHits"}
 
std::string sp_pass_name_ {""}
 
std::string tagger_sim_hits_coll_name_ {"TaggerSimHits"}
 Sim hits to check if the truth seed is findable.
 
std::string recoil_sim_hits_coll_name_ {"RecoilSimHits"}
 Sim hits to check if the truth seed is findable.
 
std::string input_pass_name_ {""}
 Pass name for the sim hit collections.
 
std::string sim_particles_coll_name_
 
std::string sim_particles_passname_
 
int n_min_hits_tagger_ {7}
 Minimum number of hits left in the recoil tracker to consider the seed as findable.
 
int n_min_hits_recoil_ {7}
 Minimum number of hits left in the recoil tracker to consider the seed as findable.
 
float z_min_ {-999}
 Min cut on the z_ of the scoring hit.
 
int track_id_ {-999}
 Only select a particular trackID.
 
double pz_cut_ {-9999}
 Ask for a minimum pz for the seeds.
 
double p_cut_ {0.}
 Ask for a minimum p for the seeds.
 
double p_cut_max_ {100000.}
 Ask for a maximum p for the seeds.
 
double p_cut_ecal_ {-1.}
 
bool recoil_sp_ {true}
 
bool target_sp_ {true}
 
bool skip_tagger_ {false}
 
bool skip_recoil_ {false}
 
int max_track_id_ {5}
 
std::shared_ptr< LinPropagator > linpropagator_
 
std::shared_ptr< tracking::reco::TrackExtrapolatorTool< LinPropagator > > trk_extrap_
 
std::default_random_engine generator_
 
std::shared_ptr< std::normal_distribution< float > > normal_
 
bool seed_smearing_ {false}
 
std::vector< double > d0smear_
 
std::vector< double > z0smear_
 
double phismear_
 
double thetasmear_
 
double relpsmear_
 
std::vector< double > rel_smearfactors_
 
std::vector< double > inflate_factors_
 
std::vector< double > beam_origin_ {-880.1, -44., 0.}
 
int particle_hypothesis_
 
std::string beam_electrons_collection_
 
std::string tagger_truth_collection_
 
std::string recoil_truth_collection_
 
std::string tagger_seeds_collection_
 
std::string recoil_seeds_collection_
 

Additional Inherited Members

- Protected Member Functions inherited from tracking::reco::TrackingGeometryUser
const Acts::GeometryContext & geometryContext ()
 
const Acts::MagneticFieldContext & magneticFieldContext ()
 
const Acts::CalibrationContext & calibrationContext ()
 
const geo::TrackersTrackingGeometrygeometry ()
 
- 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

Create a track seed using truth information extracted from the corresponding SimParticle or SimTrackerHit.

When creating seeds in the Tagger tracker, the SimParticle associated with the incident electron (trackID == 1) is used to create the seed from the parameters (x_, y_, z_, px, py, pz, q) at the vertex. For the Recoil tracker, since the electron is produced upstream, the SimParticle can't be used to get any parameters at the target. In this case, the target scoring plane hits are used to extract the parameters above.

Definition at line 42 of file TruthSeedProcessor.h.

Constructor & Destructor Documentation

◆ TruthSeedProcessor()

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

Constructor.

Parameters
nameName for this instance of the class.
processThe Process class associated with EventProcessor, provided by the framework.

Definition at line 7 of file TruthSeedProcessor.cxx.

9 : TrackingGeometryUser(name, process) {}
virtual void process(Event &event) final
Processing an event for a Producer is calling produce.

Member Function Documentation

◆ configure()

void tracking::reco::TruthSeedProcessor::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.

Parameters
parametersParameters for configuration.

Reimplemented from framework::EventProcessor.

Definition at line 24 of file TruthSeedProcessor.cxx.

24 {
26 parameters.get<std::string>("scoring_hits_coll_name");
27 ecal_sp_coll_name_ = parameters.get<std::string>("ecal_sp_coll_name");
28 sp_pass_name_ = parameters.get<std::string>("sp_pass_name");
30 parameters.get<std::string>("recoil_sim_hits_coll_name");
32 parameters.get<std::string>("tagger_sim_hits_coll_name");
33 input_pass_name_ = parameters.get<std::string>("input_pass_name");
34 sim_particles_coll_name_ =
35 parameters.get<std::string>("sim_particles_coll_name");
36 sim_particles_passname_ =
37 parameters.get<std::string>("sim_particles_passname");
38
39 n_min_hits_tagger_ = parameters.get<int>("n_min_hits_tagger", 11);
40 n_min_hits_recoil_ = parameters.get<int>("n_min_hits_recoil", 7);
41 pdg_ids_ = parameters.get<std::vector<int>>("pdg_ids", {11});
42 z_min_ = parameters.get<double>("z_min", -9999); // mm
43 track_id_ = parameters.get<int>("track_id", -9999);
44 pz_cut_ = parameters.get<double>("pz_cut", -9999); // MeV
45 p_cut_ = parameters.get<double>("p_cut", 0.);
46 p_cut_max_ = parameters.get<double>("p_cut_max", 100000.); // MeV
47 p_cut_ecal_ = parameters.get<double>("p_cut_ecal", -1.); // MeV
48 recoil_sp_ = parameters.get<double>("recoil_sp", true);
49 target_sp_ = parameters.get<double>("tagger_sp", true);
50 seed_smearing_ = parameters.get<bool>("seedSmearing", false);
51 max_track_id_ = parameters.get<int>("max_track_id", 5);
52
53 ldmx_log(info) << "Seed Smearing is set to " << seed_smearing_;
54
55 d0smear_ = parameters.get<std::vector<double>>("d0smear", {0.01, 0.01, 0.01});
56 z0smear_ = parameters.get<std::vector<double>>("z0smear", {0.1, 0.1, 0.1});
57 phismear_ = parameters.get<double>("phismear", 0.001);
58 thetasmear_ = parameters.get<double>("thetasmear", 0.001);
59 relpsmear_ = parameters.get<double>("relpsmear", 0.1);
60
61 // Relative smear factor terms, only used if seed_smearing_ is true.
62 rel_smearfactors_ = parameters.get<std::vector<double>>(
63 "rel_smearfactors", {0.1, 0.1, 0.1, 0.1, 0.1, 0.1});
64 inflate_factors_ = parameters.get<std::vector<double>>(
65 "inflate_factors", {10., 10., 10., 10., 10., 10.});
66
67 // In tracking frame: where do these numbers come from?
68 // These numbers come from approximating the path of the beam up
69 // until it is about to enter the first detector volume (TriggerPad1).
70 // In detector coordinates, (x_,y_,z_) = (-21.7, -883) is
71 // where the beam arrives (if no smearing is applied) and we simply
72 // reorder these values so that they are in tracking coordinates.
73 beam_origin_ = parameters.get<std::vector<double>>("beamOrigin",
74 {-883.0, -21.745876, 0.0});
75
76 // Skip the tagger or recoil trackers if wanted
77 skip_tagger_ = parameters.get<bool>("skip_tagger", false);
78 skip_recoil_ = parameters.get<bool>("skip_recoil", false);
79 particle_hypothesis_ = parameters.get<int>("particle_hypothesis");
80
81 beam_electrons_collection_ =
82 parameters.get<std::string>("beam_electrons_collection");
83 tagger_truth_collection_ =
84 parameters.get<std::string>("tagger_truth_collection");
85 recoil_truth_collection_ =
86 parameters.get<std::string>("recoil_truth_collection");
87 tagger_seeds_collection_ =
88 parameters.get<std::string>("tagger_seeds_collection");
89 recoil_seeds_collection_ =
90 parameters.get<std::string>("recoil_seeds_collection");
91}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
float z_min_
Min cut on the z_ of the scoring hit.
std::string input_pass_name_
Pass name for the sim hit collections.
std::string recoil_sim_hits_coll_name_
Sim hits to check if the truth seed is findable.
std::vector< int > pdg_ids_
pdg_ids of the particles we want to select for the seeds
int track_id_
Only select a particular trackID.
double pz_cut_
Ask for a minimum pz for the seeds.
double p_cut_
Ask for a minimum p for the seeds.
std::string scoring_hits_coll_name_
Which scoring plane hits to use for the truth seeds generation.
std::string tagger_sim_hits_coll_name_
Sim hits to check if the truth seed is findable.
double p_cut_max_
Ask for a maximum p for the seeds.
int n_min_hits_tagger_
Minimum number of hits left in the recoil tracker to consider the seed as findable.
int n_min_hits_recoil_
Minimum number of hits left in the recoil tracker to consider the seed as findable.

References framework::config::Parameters::get(), input_pass_name_, n_min_hits_recoil_, n_min_hits_tagger_, p_cut_, p_cut_max_, pdg_ids_, pz_cut_, recoil_sim_hits_coll_name_, scoring_hits_coll_name_, tagger_sim_hits_coll_name_, track_id_, and z_min_.

◆ createTruthTrack() [1/3]

void tracking::reco::TruthSeedProcessor::createTruthTrack ( const ldmx::SimParticle & particle,
const ldmx::SimTrackerHit & hit,
ldmx::Track & trk,
const std::shared_ptr< Acts::Surface > & target_surface )
private

Use the scoring plane hit at the target to extract (x_, y_, z_, px, py, pz) and create a track seed.

In this case, the SimParticle is used to extract the charge of the particle.

Parameters
particleThe SimParticle to extract the charge from.
hitThe SimTrackerHit used to create the seed.

Definition at line 93 of file TruthSeedProcessor.cxx.

95 {
96 std::vector<double> pos{static_cast<double>(hit.getPosition()[0]),
97 static_cast<double>(hit.getPosition()[1]),
98 static_cast<double>(hit.getPosition()[2])};
99 createTruthTrack(pos, hit.getMomentum(), particle.getCharge(), trk,
100 target_surface);
101
102 trk.setTrackID(hit.getTrackID());
103 trk.setPdgID(hit.getPdgID());
104}
double getCharge() const
Get the charge of this particle.
int getPdgID() const
Get the Sim particle track ID of the hit.
std::vector< float > getPosition() const
Get the XYZ position of the hit [mm].
std::vector< double > getMomentum() const
Get the XYZ momentum of the particle at the position at which the hit took place [MeV].
int getTrackID() const
Get the Sim particle track ID of the hit.
void createTruthTrack(const ldmx::SimParticle &particle, ldmx::Track &trk, const std::shared_ptr< Acts::Surface > &target_surface)
Use the vertex position of the SimParticle to extract (x_, y_, z_, px, py, pz, q) and create a track ...

References createTruthTrack(), ldmx::SimParticle::getCharge(), ldmx::SimTrackerHit::getMomentum(), ldmx::SimTrackerHit::getPdgID(), ldmx::SimTrackerHit::getPosition(), and ldmx::SimTrackerHit::getTrackID().

◆ createTruthTrack() [2/3]

void tracking::reco::TruthSeedProcessor::createTruthTrack ( const ldmx::SimParticle & particle,
ldmx::Track & trk,
const std::shared_ptr< Acts::Surface > & target_surface )
private

Use the vertex position of the SimParticle to extract (x_, y_, z_, px, py, pz, q) and create a track seed.

Parameters
particleThe SimParticle to make a seed from.

Definition at line 106 of file TruthSeedProcessor.cxx.

108 {
109 createTruthTrack(particle.getVertex(), particle.getMomentum(),
110 particle.getCharge(), trk, target_surface);
111
112 trk.setPdgID(particle.getPdgID());
113}
std::vector< double > getVertex() const
Get a vector containing the vertex of this particle in mm.
int getPdgID() const
Get the PDG ID of this particle.
Definition SimParticle.h:85
std::vector< double > getMomentum() const
Get a vector containing the momentum of this particle [MeV].

References createTruthTrack(), ldmx::SimParticle::getCharge(), ldmx::SimParticle::getMomentum(), ldmx::SimParticle::getPdgID(), and ldmx::SimParticle::getVertex().

Referenced by createTruthTrack(), createTruthTrack(), and taggerFullSeed().

◆ createTruthTrack() [3/3]

void tracking::reco::TruthSeedProcessor::createTruthTrack ( const std::vector< double > & pos_vec,
const std::vector< double > & p_vec,
int charge,
ldmx::Track & trk,
const std::shared_ptr< Acts::Surface > & target_surface )
private

Create a seed track from the given position, momentum and charge.

Parameters
pos_The position at which the particle was created.
pThe momentum of the particle at the point of creation.
chargeThe charge of the particle.
target_surfacethe surface to where to express the truth track

Definition at line 115 of file TruthSeedProcessor.cxx.

118 {
119 // Get the position and momentum of the particle at the point of creation.
120 // This only works for the incident electron when creating a tagger tracker
121 // seed. For the recoil tracker, the scoring plane position will need to
122 // be used. For other particles created in the target or tracker planes,
123 // this version of the method can also be used.
124 // These are just Eigen vectors defined as
125 // Eigen::Matrix<double, kSize, 1>;
126 Acts::Vector3 pos{pos_vec[0], pos_vec[1], pos_vec[2]};
127 Acts::Vector3 mom{p_vec[0], p_vec[1], p_vec[2]};
128
129 // Rotate the position and momentum into the ACTS frame.
130 pos = tracking::sim::utils::ldmx2Acts(pos);
131 mom = tracking::sim::utils::ldmx2Acts(mom);
132
133 // Get the charge of the particle.
134 // TODO: Add function that uses the PDG ID to calculate this.
135 double q{charge * Acts::UnitConstants::e};
136
137 // The idea here is:
138 // 1 - Define a bound track state parameters at point P on track. Basically a
139 // curvilinear representation. 2 - Propagate to target surface to obtain the
140 // BoundTrackState there.
141
142 // Transform the position, momentum and charge to free parameters.
143 auto free_params{tracking::sim::utils::toFreeParameters(pos, mom, q)};
144
145 // Create a line surface at the perigee. The perigee position is extracted
146 // from a particle's vertex or the particle's position at a specific
147 // scoring plane.
148 auto gen_surface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
149 Acts::Vector3(free_params[Acts::eFreePos0], free_params[Acts::eFreePos1],
150 free_params[Acts::eFreePos2]))};
151
152 // ldmx_log(trace)<<"PF:: gen_surface"<<free_params[Acts::eFreePos0]<<"
153 // " <<free_params[Acts::eFreePos1]<<" " <<free_params[Acts::eFreePos2];
154
155 // Transform the parameters to local positions on the perigee surface.
156 auto bound_params{
157 Acts::transformFreeToBoundParameters(free_params, *gen_surface, gctx_)
158 .value()};
159 // Create a particle hypothesis
160 auto part{Acts::GenericParticleHypothesis(
161 Acts::ParticleHypothesis(Acts::PdgParticle(particle_hypothesis_)))};
162 Acts::BoundTrackParameters bound_trk_pars(gen_surface, bound_params,
163 std::nullopt, part);
164
165 // CAUTION:: The target surface should be close to the gen surface
166 // Linear propagation to the target surface. I assume 1mm of tolerance
167 Acts::Vector3 tgt_surf_center = target_surface->center(geometryContext());
168 Acts::Vector3 gen_surf_center = gen_surface->center(geometryContext());
169 // Tolerance
170 double tol = 1; // mm
171
172 if (abs(tgt_surf_center(0) - gen_surf_center(0)) > tol)
173 ldmx_log(error) << "Linear extrapolation to a far away surface in B field."
174 << " This will cause inaccuracies in track parameters"
175 << " Distance extrapolated = "
176 << (tgt_surf_center(0) - gen_surf_center(0));
177
178 auto prop_bound_state =
179 trk_extrap_->extrapolate(bound_trk_pars, target_surface);
180
181 // Create the seed track object.
182 Acts::Vector3 ref = target_surface->center(geometryContext());
183 // trk.setPerigeeLocation(free_params[Acts::eFreePos0],
184 // free_params[Acts::eFreePos1],
185 // free_params[Acts::eFreePos2]);
186
187 trk.setPerigeeLocation(ref(0), ref(1), ref(2));
188
189 auto prop_bound_vec = (prop_bound_state.value()).parameters();
190
192 tracking::sim::utils::convertActsToLdmxPars(prop_bound_vec));
193
194 trk.setPosition(pos(0), pos(1), pos(2));
195 trk.setMomentum(mom(0), mom(1), mom(2));
196}
void setPerigeeParameters(const std::vector< double > &par)
d_0 z_0 phi_0 theta q/p t
Definition Track.h:161
Acts::GeometryContext gctx_
The ACTS geometry context properly.

References gctx_, and ldmx::Track::setPerigeeParameters().

◆ makeHitCountMap()

void tracking::reco::TruthSeedProcessor::makeHitCountMap ( const std::vector< ldmx::SimTrackerHit > & sim_hits,
std::map< int, std::vector< int > > & hit_count_map )
private

Create a mapping from the selected scoring plane hit objects to the number of hits they associated particle creates in the tracker.

Parameters
sim_hitsvector
hit_count_mapfilled with the hits lefts by each track

Definition at line 454 of file TruthSeedProcessor.cxx.

456 {
457 for (int i_sim_hit = 0; i_sim_hit < sim_hits.size(); i_sim_hit++) {
458 auto& sim_hit = sim_hits.at(i_sim_hit);
459 // This track never left a hit before
460 if (!hit_count_map.count(sim_hit.getTrackID())) {
461 hit_count_map[sim_hit.getTrackID()].push_back(i_sim_hit);
462 }
463
464 // This track left a hit before.
465 // Check if it's on a different sensor than the others
466
467 else {
468 int sensor_id = tracking::sim::utils::getSensorID(sim_hit);
469 bool found_hit = false;
470
471 for (auto& i_rhit : hit_count_map[sim_hit.getTrackID()]) {
472 int tmp_sensor_id =
473 tracking::sim::utils::getSensorID(sim_hits.at(i_rhit));
474
475 if (sensor_id == tmp_sensor_id) {
476 found_hit = true;
477 break;
478 }
479 } // loop on the already recorded hits
480
481 if (!found_hit) {
482 hit_count_map[sim_hit.getTrackID()].push_back(i_sim_hit);
483 }
484 }
485 } // loop on sim hits
486}

Referenced by produce().

◆ onNewRun()

void tracking::reco::TruthSeedProcessor::onNewRun ( const ldmx::RunHeader & rh)
overridevirtual

onNewRun is the first function called for each processor after the conditions are fully configured and accessible.

This is where you could create single-processors, multi-event calculation objects.

Reimplemented from framework::EventProcessor.

Definition at line 11 of file TruthSeedProcessor.cxx.

11 {
12 gctx_ = Acts::GeometryContext();
13 normal_ = std::make_shared<std::normal_distribution<float>>(0., 1.);
14
15 Acts::StraightLineStepper stepper;
16 Acts::Navigator::Config nav_cfg{geometry().getTG()};
17 const Acts::Navigator navigator(nav_cfg);
18
19 linpropagator_ = std::make_shared<LinPropagator>(stepper, navigator);
20 trk_extrap_ = std::make_shared<std::decay_t<decltype(*trk_extrap_)>>(
21 *linpropagator_, geometryContext(), magneticFieldContext());
22}

References gctx_.

◆ onProcessStart()

void tracking::reco::TruthSeedProcessor::onProcessStart ( )
inlineoverridevirtual

Callback for the EventProcessor to take any necessary action when the processing of events starts.

For this class, the callback is used to retrieve the GeometryContext from ACTS.

Reimplemented from framework::EventProcessor.

Definition at line 73 of file TruthSeedProcessor.h.

73{};

◆ produce()

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

Main loop that creates the seed tracks for both the tagger and recoil tracker.

Parameters
eventThe event containing the collections to process.

Implements framework::Producer.

Definition at line 542 of file TruthSeedProcessor.cxx.

542 {
543 // Retrieve the particleMap
544 auto particle_map{event.getMap<int, ldmx::SimParticle>(
545 sim_particles_coll_name_, sim_particles_passname_)};
546
547 // Retrieve the target scoring hits
548 // Information is extracted using the
549 // scoring plane hit left by the particle at the target.
550
551 const std::vector<ldmx::SimTrackerHit> scoring_hits{
552 event.getCollection<ldmx::SimTrackerHit>(scoring_hits_coll_name_,
553 sp_pass_name_)};
554
555 // Retrieve the scoring plane hits at the ECAL
556 const std::vector<ldmx::SimTrackerHit> scoring_hits_ecal{
557 event.getCollection<ldmx::SimTrackerHit>(ecal_sp_coll_name_,
558 sp_pass_name_)};
559
560 // Retrieve the sim hits in the tagger tracker
561 const std::vector<ldmx::SimTrackerHit> tagger_sim_hits =
564
565 // Retrieve the sim hits in the recoil tracker
566 const std::vector<ldmx::SimTrackerHit> recoil_sim_hits =
569
570 // If sim hit collections are empty throw a warning
571 if (tagger_sim_hits.size() == 0 && !skip_tagger_) {
572 ldmx_log(error) << "Tagger sim hits collection empty for event ";
573 }
574 if (recoil_sim_hits.size() == 0 && !skip_recoil_) {
575 ldmx_log(error) << "Recoil sim hits collection empty for event ";
576 }
577
578 // The map stores which track leaves which sim hits
579 std::map<int, std::vector<int>> hit_count_map_recoil;
580 makeHitCountMap(recoil_sim_hits, hit_count_map_recoil);
581
582 std::map<int, std::vector<int>> hit_count_map_tagger;
583 makeHitCountMap(tagger_sim_hits, hit_count_map_tagger);
584
585 // to keep track of how many sim particles leave hits on the scoring plane
586 std::vector<int> recoil_sh_idxs;
587 std::unordered_map<int, std::vector<int>> recoil_sh_count_map;
588
589 std::vector<int> tagger_sh_idxs;
590 std::unordered_map<int, std::vector<int>> tagger_sh_count_map;
591
592 // Target scoring hits for Tagger will have Z<0, Recoil scoring hits will
593 // have Z>0
594 for (unsigned int i_sh = 0; i_sh < scoring_hits.size(); i_sh++) {
595 const ldmx::SimTrackerHit& hit = scoring_hits.at(i_sh);
596 double zhit = hit.getPosition()[2];
597
598 Acts::Vector3 p_vec{hit.getMomentum()[0], hit.getMomentum()[1],
599 hit.getMomentum()[2]};
600 double tagger_p_max = 0.;
601
602 // Check if it is a tagger track going fwd that passes basic cuts
603 if (zhit < 0.) {
604 // Tagger selection cuts
605 // Negative scoring plane hit, with momentum > p_cut
606 if (p_vec(2) < 0. || p_vec.norm() < p_cut_) continue;
607
608 // Check that the hit was left by a charged particle
609 if (abs(particle_map[hit.getTrackID()].getCharge()) < 1e-8) continue;
610
611 if (p_vec.norm() > tagger_p_max) {
612 tagger_sh_count_map[hit.getTrackID()].push_back(i_sh);
613 }
614 } // Tagger loop
615
616 // Check the recoil hits
617 else {
618 // Recoil selection cuts
619 // Positive scoring plane hit, forward direction with momentum > p_cut
620 if (p_vec(2) < 0. || p_vec.norm() < p_cut_) continue;
621
622 // Check that the hit was left by a charged particle
623 if (abs(particle_map[hit.getTrackID()].getCharge()) < 1e-8) continue;
624
625 recoil_sh_count_map[hit.getTrackID()].push_back(i_sh);
626
627 } // Recoil
628 } // loop on Target scoring plane hits
629
630 for (std::pair<int, std::vector<int>> element : recoil_sh_count_map) {
631 std::sort(
632 element.second.begin(), element.second.end(),
633 [&](const int idx1, int idx2) -> bool {
634 const ldmx::SimTrackerHit& hit1 = scoring_hits.at(idx1);
635 const ldmx::SimTrackerHit& hit2 = scoring_hits.at(idx2);
636
637 Acts::Vector3 phit1{hit1.getMomentum()[0], hit1.getMomentum()[1],
638 hit1.getMomentum()[2]};
639 Acts::Vector3 phit2{hit2.getMomentum()[0], hit2.getMomentum()[1],
640 hit2.getMomentum()[2]};
641
642 return phit1.norm() > phit2.norm();
643 });
644 }
645
646 // Sort tagger hits.
647 for (auto& [_track_id, hit_indices] : tagger_sh_count_map) {
648 std::sort(
649 hit_indices.begin(), hit_indices.end(),
650 [&](const int idx1, int idx2) -> bool {
651 const ldmx::SimTrackerHit& hit1 = scoring_hits.at(idx1);
652 const ldmx::SimTrackerHit& hit2 = scoring_hits.at(idx2);
653
654 Acts::Vector3 phit1{hit1.getMomentum()[0], hit1.getMomentum()[1],
655 hit1.getMomentum()[2]};
656 Acts::Vector3 phit2{hit2.getMomentum()[0], hit2.getMomentum()[1],
657 hit2.getMomentum()[2]};
658
659 return phit1.norm() > phit2.norm();
660 });
661 }
662
663 // Building of the event truth information and the truth seeds
664 // TODO remove the truthtracks in the future as the truth seeds are enough
665
666 std::vector<ldmx::Track> tagger_truth_tracks;
667 std::vector<ldmx::Track> tagger_truth_seeds;
668 std::vector<ldmx::Track> recoil_truth_tracks;
669 std::vector<ldmx::Track> recoil_truth_seeds;
670 ldmx::Tracks beam_electrons;
671
672 // TODO:: The target should be taken from some conditions DB in the future.
673 // Define the perigee_surface at 0.0.0
674 auto target_surface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
675 Acts::Vector3(0., 0., 0.))};
676
677 // Define the target_surface
678 auto target_unbound_surface = tracking::sim::utils::unboundSurface(0.);
679
680 // ecal
681 auto ecal_surface = tracking::sim::utils::unboundSurface(240.5);
682
683 auto beam_origin_surface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
684 Acts::Vector3(beam_origin_[0], beam_origin_[1], beam_origin_[2]))};
685
686 if (!skip_tagger_) {
687 for (const auto& [track_id, hit_indices] : tagger_sh_count_map) {
688 const ldmx::SimTrackerHit& hit = scoring_hits.at(hit_indices.at(0));
689 const ldmx::SimParticle& phit = particle_map[hit.getTrackID()];
690
691 if (hit_count_map_tagger[hit.getTrackID()].size() > n_min_hits_tagger_) {
692 ldmx::Track truth_tagger_track;
693 createTruthTrack(phit, hit, truth_tagger_track, target_surface);
694 truth_tagger_track.setNhits(
695 hit_count_map_tagger[hit.getTrackID()].size());
696 tagger_truth_tracks.push_back(truth_tagger_track);
697
698 if (hit.getPdgID() == 11 && hit.getTrackID() < max_track_id_) {
699 ldmx::Track beam_e_truth_seed =
700 taggerFullSeed(particle_map[hit.getTrackID()], hit.getTrackID(),
701 hit, hit_count_map_tagger, beam_origin_surface,
702 target_unbound_surface);
703 beam_electrons.push_back(beam_e_truth_seed);
704 }
705 }
706 }
707 }
708
709 // Recover the EcalScoring hits
710 std::vector<ldmx::SimTrackerHit> ecal_sp_hits =
711 event.getCollection<ldmx::SimTrackerHit>(ecal_sp_coll_name_,
712 sp_pass_name_);
713 // Select ECAL hits
714 std::vector<ldmx::SimTrackerHit> sel_ecal_sp_hits;
715
716 for (auto sp_hit : ecal_sp_hits) {
717 if (sp_hit.getMomentum()[2] > 0 && ((sp_hit.getID() & 0xfff) == 31)) {
718 sel_ecal_sp_hits.push_back(sp_hit);
719 }
720 }
721
722 // Recoil target surface for truth and seed tracks is the target
723
724 for (std::pair<int, std::vector<int>> element : recoil_sh_count_map) {
725 // Only take the first entry of the vector: it should be the scoring plane
726 // hit with the highest momentum.
727 const ldmx::SimTrackerHit& hit = scoring_hits.at(element.second.at(0));
728 [[maybe_unused]] const ldmx::SimParticle& phit =
729 particle_map[hit.getTrackID()];
730 ldmx::SimTrackerHit ecal_hit;
731
732 bool found_ecal_hit = false;
733 for (auto ecal_sp_hit : sel_ecal_sp_hits) {
734 if (ecal_sp_hit.getTrackID() == hit.getTrackID()) {
735 ecal_hit = ecal_sp_hit;
736 found_ecal_hit = true;
737 break;
738 }
739 }
740
741 // Findable particle selection
742 // ldmx_log(trace) << "!!! n_recoil_sim_hits found: " <<
743 // hit_count_map_recoil[hit.getTrackID()].size();
744 if (hit_count_map_recoil[hit.getTrackID()].size() > n_min_hits_recoil_ &&
745 found_ecal_hit && !skip_recoil_) {
746 ldmx::Track truth_recoil_track =
747 recoilFullSeed(particle_map[hit.getTrackID()], hit.getTrackID(), hit,
748 ecal_hit, hit_count_map_recoil, target_surface,
749 target_unbound_surface, ecal_surface);
750 // ldmx_log(trace) << "!!! Recoil track created";
751 // ldmx_log(trace) << "!!! Recoil track momentum: " <<
752 // truth_recoil_track.getMomentum()[0] << " " <<
753 // truth_recoil_track.getMomentum()[1] << " " <<
754 // truth_recoil_track.getMomentum()[2]; ldmx_log(trace) << "!!! Hit
755 // momentum: " << hit.getMomentum()[0] << " " << hit.getMomentum()[1] << "
756 // " << hit.getMomentum()[2];
757 recoil_truth_tracks.push_back(truth_recoil_track);
758 }
759 }
760
761 /*
762 for (std::pair<int,std::vector<int>> element : recoil_sh_count_map) {
763
764 const ldmx::SimTrackerHit& hit = scoring_hits.at(element.second.at(0));
765 const ldmx::SimParticle& phit = particleMap[hit.getTrackID()];
766
767 if (hit_count_map_recoil[hit.getTrackID()].size() > n_min_hits_recoil_) {
768 ldmx::Track truth_recoil_track;
769 createTruthTrack(phit,hit,truth_recoil_track,targetSurface);
770 truth_recoil_track.setNhits(hit_count_map_recoil[hit.getTrackID()].size());
771 recoil_truth_tracks.push_back(truth_recoil_track);
772 }
773 }
774 */
775
776 // Form a truth seed from a truth track
777
778 for (auto& tt : tagger_truth_tracks) {
779 ldmx::Track seed = seedFromTruth(tt, seed_smearing_);
780
781 tagger_truth_seeds.push_back(seed);
782 }
783
784 ldmx_log(debug) << "Forming seeds from truth";
785 for (auto& tt : recoil_truth_tracks) {
786 ldmx_log(debug) << "Smearing truth track";
787
788 ldmx::Track seed = seedFromTruth(tt, seed_smearing_);
789
790 recoil_truth_seeds.push_back(seed);
791 }
792
793 // even if skip_tagger/recoil_ is true, still make the collections in the
794 // event
795 event.add(beam_electrons_collection_, beam_electrons);
796 event.add(tagger_truth_collection_, tagger_truth_tracks);
797 event.add(recoil_truth_collection_, recoil_truth_tracks);
798 event.add(tagger_seeds_collection_, tagger_truth_seeds);
799 event.add(recoil_seeds_collection_, recoil_truth_seeds);
800}
Class representing a simulated particle.
Definition SimParticle.h:23
Represents a simulated tracker hit in the simulation.
Implementation of a track object.
Definition Track.h:53
void makeHitCountMap(const std::vector< ldmx::SimTrackerHit > &sim_hits, std::map< int, std::vector< int > > &hit_count_map)
Create a mapping from the selected scoring plane hit objects to the number of hits they associated pa...
ldmx::Track taggerFullSeed(const ldmx::SimParticle &beam_electron, const int trackID, const ldmx::SimTrackerHit &hit, const std::map< int, std::vector< int > > &hit_count_map, const std::shared_ptr< Acts::Surface > &origin_surface, const std::shared_ptr< Acts::Surface > &target_surface)
This method retrieves the beam electron and forms a full seed The seed parameters are the truth param...
ldmx::Track seedFromTruth(const ldmx::Track &tt, bool seed_smearing)
Create a track seed from a truth track applying a smearing to the truth parameters as well as an infl...

References ldmx::SimTrackerHit::getMomentum(), ldmx::SimTrackerHit::getPosition(), ldmx::SimTrackerHit::getTrackID(), input_pass_name_, makeHitCountMap(), p_cut_, recoil_sim_hits_coll_name_, scoring_hits_coll_name_, and tagger_sim_hits_coll_name_.

◆ recoilFullSeed()

ldmx::Track tracking::reco::TruthSeedProcessor::recoilFullSeed ( const ldmx::SimParticle & particle,
const int trackID,
const ldmx::SimTrackerHit & hit,
const ldmx::SimTrackerHit & ecal_hit,
const std::map< int, std::vector< int > > & hit_count_map,
const std::shared_ptr< Acts::Surface > & origin_surface,
const std::shared_ptr< Acts::Surface > & target_surface,
const std::shared_ptr< Acts::Surface > & ecal_surface )
private

Definition at line 200 of file TruthSeedProcessor.cxx.

206 {
207 ldmx::Track truth_recoil_track;
208 createTruthTrack(particle, hit, truth_recoil_track, origin_surface);
209 truth_recoil_track.setTrackID(trackID);
210
211 // Seed at the target location
212 ldmx::Track smeared_truth_track = seedFromTruth(truth_recoil_track, false);
213
214 // Add the track state at the target
215 ldmx::Track truth_track_target;
216 createTruthTrack(particle, hit, truth_track_target, target_surface);
217
218 // Store the truth track state on the seed track
219 ldmx::Track::TrackState ts_truth_target;
220 Acts::Vector3 ref = target_surface->center(geometryContext());
221 ts_truth_target.ref_x_ = ref(0);
222 ts_truth_target.ref_y_ = ref(1);
223 ts_truth_target.ref_z_ = ref(2);
224 ts_truth_target.params_ = truth_track_target.getPerigeeParameters();
225 // empty cov
226 ts_truth_target.ts_type_ = ldmx::TrackStateType::AtTarget;
227 smeared_truth_track.addTrackState(ts_truth_target);
228
229 // Add the track state at the ecal
230 ldmx::Track truth_track_ecal;
231 createTruthTrack(particle, ecal_hit, truth_track_ecal, ecal_surface);
232
233 ldmx::Track::TrackState ts_truth_ecal;
234 Acts::Vector3 ref_ecal = ecal_surface->center(geometryContext());
235 ts_truth_ecal.ref_x_ = ref_ecal(0);
236 ts_truth_ecal.ref_y_ = ref_ecal(1);
237 ts_truth_ecal.ref_z_ = ref_ecal(2);
238 ts_truth_ecal.params_ = truth_track_ecal.getPerigeeParameters();
239 // empty cov
240 ts_truth_ecal.ts_type_ = ldmx::TrackStateType::AtECAL;
241 smeared_truth_track.addTrackState(ts_truth_ecal);
242
243 // Add the hits
244 int nhits = 0;
245
246 for (auto sim_hit_idx : hit_count_map.at(smeared_truth_track.getTrackID())) {
247 smeared_truth_track.addMeasurementIndex(sim_hit_idx);
248 nhits += 1;
249 }
250
251 smeared_truth_track.setNhits(nhits);
252
253 return smeared_truth_track;
254}

◆ scoringPlaneHitFilter()

bool tracking::reco::TruthSeedProcessor::scoringPlaneHitFilter ( const ldmx::SimTrackerHit & hit,
const std::vector< ldmx::SimTrackerHit > & ecal_sp_hits )
private

Filter that checks if a scoring plane passes specified momentum cuts as well as if the associated SimParticle hits the ECal.

Parameters
hitThe target scoring plane hit to check.
ecal_sp_hitsThe ECal scoring plane hit used to check if the associated particle hits the ECal.

Definition at line 488 of file TruthSeedProcessor.cxx.

490 {
491 // Clean some of the hits we don't want
492 if (hit.getPosition()[2] < z_min_) return false;
493
494 // Check if the track_id was requested
495 if (track_id_ > 0 && hit.getTrackID() != track_id_) return false;
496
497 // Check if we are requesting particular particles
498 if (std::find(pdg_ids_.begin(), pdg_ids_.end(), hit.getPdgID()) ==
499 pdg_ids_.end())
500 return false;
501
502 Acts::Vector3 p_vec{hit.getMomentum()[0], hit.getMomentum()[1],
503 hit.getMomentum()[2]};
504
505 // p cut
506 if (p_cut_ >= 0. && p_vec.norm() < p_cut_) return false;
507
508 // p cut Max
509 if (p_cut_ < 100000. && p_vec.norm() > p_cut_max_) return false;
510
511 // pz cut
512 if (pz_cut_ > -9999 && p_vec(2) < pz_cut_) return false;
513
514 // Check the ecal scoring plane
515 bool pass_ecal_scoring_plane = true;
516
517 if (p_cut_ecal_ > 0) { // only check if we care about it.
518
519 for (auto& e_sp_hit : ecal_sp_hits) {
520 if (e_sp_hit.getTrackID() == hit.getTrackID() &&
521 e_sp_hit.getPdgID() == hit.getPdgID()) {
522 Acts::Vector3 e_sp_p{e_sp_hit.getMomentum()[0],
523 e_sp_hit.getMomentum()[1],
524 e_sp_hit.getMomentum()[2]};
525
526 if (e_sp_p.norm() < p_cut_ecal_) pass_ecal_scoring_plane = false;
527
528 // Skip the rest of the scoring plane hits since we already found the
529 // track we care about
530 break;
531
532 } // check that the hit belongs to the inital particle from the target
533 // scoring plane hit
534 } // loop on Ecal scoring plane hits
535 } // pcutEcal
536
537 if (!pass_ecal_scoring_plane) return false;
538
539 return true;
540}

References ldmx::SimTrackerHit::getMomentum(), ldmx::SimTrackerHit::getPdgID(), ldmx::SimTrackerHit::getPosition(), ldmx::SimTrackerHit::getTrackID(), p_cut_, p_cut_max_, pdg_ids_, pz_cut_, track_id_, and z_min_.

◆ seedFromTruth()

ldmx::Track tracking::reco::TruthSeedProcessor::seedFromTruth ( const ldmx::Track & tt,
bool seed_smearing )
private

Create a track seed from a truth track applying a smearing to the truth parameters as well as an inflation to the covariance matrix.

Parameters
ttTruthTrack to be used to form a seed
Returns
seed The seed track

Definition at line 338 of file TruthSeedProcessor.cxx.

339 {
340 ldmx::Track seed = ldmx::Track();
341 seed.setPerigeeLocation(tt.getPerigeeLocation()[0],
342 tt.getPerigeeLocation()[1],
343 tt.getPerigeeLocation()[2]);
344 seed.setChi2(0.);
345 seed.setNhits(tt.getNhits());
346 seed.setNdf(0);
347 seed.setNsharedHits(0);
348 seed.setTrackID(tt.getTrackID());
349 seed.setPdgID(tt.getPdgID());
350 seed.setTruthProb(1.);
351
352 Acts::BoundVector bound_params;
353 Acts::BoundVector stddev;
354
355 if (seed_smearing) {
356 ldmx_log(debug) << "Smear track and inflate covariance";
357
358 /*
359 double sigma_d0 = rel_smearfactors_[Acts::eBoundLoc0] * tt.getD0();
360 double sigma_z0 = rel_smearfactors_[Acts::eBoundLoc1] * tt.getZ0();
361 double sigma_phi = rel_smearfactors_[Acts::eBoundPhi] * tt.getPhi();
362 double sigma_theta = rel_smearfactors_[Acts::eBoundTheta] *
363 tt.getTheta(); double sigma_p = rel_smearfactors_[Acts::eBoundQOverP]
364 * abs(1/tt.getQoP()); double sigma_t =
365 rel_smearfactors_[Acts::eBoundTime] * tt.getT();
366 */
367
368 double sigma_d0 = d0smear_[0];
369 double sigma_z0 = z0smear_[0];
370 double sigma_phi = phismear_;
371 double sigma_theta = thetasmear_;
372 double sigma_p = relpsmear_ * abs(1 / tt.getQoP());
373 double sigma_t = 1. * Acts::UnitConstants::ns;
374
375 double smear = (*normal_)(generator_);
376 double d0smear = tt.getD0() + smear * sigma_d0;
377
378 smear = (*normal_)(generator_);
379 double z0smear = tt.getZ0() + smear * sigma_z0;
380
381 smear = (*normal_)(generator_);
382 double phismear = tt.getPhi() + smear * sigma_phi;
383
384 smear = (*normal_)(generator_);
385 double thetasmear = tt.getTheta() + smear * sigma_theta;
386
387 double p = std::abs(1. / tt.getQoP());
388 smear = (*normal_)(generator_);
389 double psmear = p + smear * sigma_p;
390
391 double q = tt.getQoP() < 0 ? -1. : 1.;
392 double qo_psmear = q / psmear;
393
394 smear = (*normal_)(generator_);
395 double tsmear = tt.getT() + smear * sigma_t;
396
397 bound_params << d0smear, z0smear, phismear, thetasmear, qo_psmear, tsmear;
398
399 stddev[Acts::eBoundLoc0] =
400 inflate_factors_[Acts::eBoundLoc0] * sigma_d0 * Acts::UnitConstants::mm;
401 stddev[Acts::eBoundLoc1] =
402 inflate_factors_[Acts::eBoundLoc1] * sigma_z0 * Acts::UnitConstants::mm;
403 stddev[Acts::eBoundPhi] = inflate_factors_[Acts::eBoundPhi] * sigma_phi;
404 stddev[Acts::eBoundTheta] =
405 inflate_factors_[Acts::eBoundTheta] * sigma_theta;
406 stddev[Acts::eBoundQOverP] =
407 inflate_factors_[Acts::eBoundQOverP] * (1. / p) * (1. / p) * sigma_p;
408 stddev[Acts::eBoundTime] =
409 inflate_factors_[Acts::eBoundTime] * sigma_t * Acts::UnitConstants::ns;
410
411 ldmx_log(debug) << stddev;
412
413 std::vector<double> v_seed_params(
414 (bound_params).data(),
415 bound_params.data() + bound_params.rows() * bound_params.cols());
416
417 Acts::BoundSquareMatrix bound_cov =
418 stddev.cwiseProduct(stddev).asDiagonal();
419 std::vector<double> v_seed_cov;
420 tracking::sim::utils::flatCov(bound_cov, v_seed_cov);
421 seed.setPerigeeParameters(v_seed_params);
422 seed.setPerigeeCov(v_seed_cov);
423
424 } else {
425 // Do not smear the seed
426
427 bound_params << tt.getD0(), tt.getZ0(), tt.getPhi(), tt.getTheta(),
428 tt.getQoP(), tt.getT();
429
430 std::vector<double> v_seed_params(
431 (bound_params).data(),
432 bound_params.data() + bound_params.rows() * bound_params.cols());
433
434 double p = std::abs(1. / tt.getQoP());
435 double sigma_p = 0.75 * p * Acts::UnitConstants::GeV;
436 stddev[Acts::eBoundLoc0] = 2 * Acts::UnitConstants::mm;
437 stddev[Acts::eBoundLoc1] = 5 * Acts::UnitConstants::mm;
438 stddev[Acts::eBoundTime] = 1000 * Acts::UnitConstants::ns;
439 stddev[Acts::eBoundPhi] = 5 * Acts::UnitConstants::degree;
440 stddev[Acts::eBoundTheta] = 5 * Acts::UnitConstants::degree;
441 stddev[Acts::eBoundQOverP] = (1. / p) * (1. / p) * sigma_p;
442
443 Acts::BoundSquareMatrix bound_cov =
444 stddev.cwiseProduct(stddev).asDiagonal();
445 std::vector<double> v_seed_cov;
446 tracking::sim::utils::flatCov(bound_cov, v_seed_cov);
447 seed.setPerigeeParameters(v_seed_params);
448 seed.setPerigeeCov(v_seed_cov);
449 }
450
451 return seed;
452}

References ldmx::Track::setPerigeeParameters().

Referenced by taggerFullSeed().

◆ taggerFullSeed()

ldmx::Track tracking::reco::TruthSeedProcessor::taggerFullSeed ( const ldmx::SimParticle & beam_electron,
const int trackID,
const ldmx::SimTrackerHit & hit,
const std::map< int, std::vector< int > > & hit_count_map,
const std::shared_ptr< Acts::Surface > & origin_surface,
const std::shared_ptr< Acts::Surface > & target_surface )
private

This method retrieves the beam electron and forms a full seed The seed parameters are the truth parameters from the beam electron stored at the beam origin Additionally, the foolowing track states are stored ts_smeared : the truth smeared perigee state at the beam origin ts_truth_target : the truth on-surface state at the target Linear extrapolations are done from the origin of the particle to the reference surfaces This track also contains the list of hits belonging to the beam electron on the sensitive surfaces on the tagger tracker, for acceptance studies.

Parameters
beam_electron: the beam electron particle
hit: the scoring hit at the target from the beam electron particle survived
hit_count_map: the sim hit on track map
origin_surface: where to express the track origin parameters. Can be perigee, plane...
target_surface: the target surface for the truth target state

Definition at line 256 of file TruthSeedProcessor.cxx.

261 {
262 ldmx::Track truth_track;
263
264 createTruthTrack(beam_electron, truth_track, origin_surface);
265 truth_track.setTrackID(trackID);
266
267 // Smeared track at the beam origin
268 ldmx::Track smeared_truth_track = seedFromTruth(truth_track, true);
269 // ldmx_log(trace) << "!!!!! Smeared truth track momentums: "
270 // << smearedTruthTrack.getMomentum()[0] << " "
271 // << smearedTruthTrack.getMomentum()[1] << " "
272 // << smearedTruthTrack.getMomentum()[2];
273 // ldmx_log(trace) << "!!!!! Smeared truth track Q/P: " <<
274 // smearedTruthTrack.getQoP()
275 // ;
276 ldmx_log(debug) << "Truth parameters at beam origin";
277 for (auto par : truth_track.getPerigeeParameters())
278 ldmx_log(debug) << par << " ";
279 ldmx_log(debug);
280
281 // Add the truth track state at the target
282 // Truth track target will be obtained from the scoring plane hit then
283 // extrapolated linearly to the target surface
284
285 ldmx::Track truth_track_target;
286 createTruthTrack(beam_electron, hit, truth_track_target, target_surface);
287
288 // Store the truth track state on the seed track
289 ldmx::Track::TrackState ts_truth_target;
290 Acts::Vector3 ref = target_surface->center(geometryContext());
291 ts_truth_target.ref_x_ = ref(0);
292 ts_truth_target.ref_y_ = ref(1);
293 ts_truth_target.ref_z_ = ref(2);
294 ts_truth_target.params_ = truth_track_target.getPerigeeParameters();
295 // empty cov
296 ts_truth_target.ts_type_ = ldmx::TrackStateType::AtTarget;
297 smeared_truth_track.addTrackState(ts_truth_target);
298
299 ldmx_log(debug) << "Truth parameters at target";
300 for (auto par : truth_track_target.getPerigeeParameters())
301 ldmx_log(debug) << par << " ";
302 ldmx_log(debug);
303
304 // This is the un-smeared truth track that can be used for pulls and residuals
305 ldmx::Track seed_truth_track = seedFromTruth(truth_track, false);
306
307 ldmx::Track::TrackState ts_truth_beam_origin;
308 Acts::Vector3 ref_origin = origin_surface->center(geometryContext());
309 ts_truth_beam_origin.ref_x_ = ref_origin(0);
310 ts_truth_beam_origin.ref_y_ = ref_origin(1);
311 ts_truth_beam_origin.ref_z_ = ref_origin(2);
312 ts_truth_beam_origin.params_ = seed_truth_track.getPerigeeParameters();
313 // ts_truth_beam_origin.cov = seedTruthTrack.getPerigeeCov();
314 ts_truth_beam_origin.ts_type_ = ldmx::TrackStateType::AtBeamOrigin;
315 smeared_truth_track.addTrackState(ts_truth_beam_origin);
316
317 ldmx_log(debug) << "Smeared parameters at origin";
318 for (auto par : smeared_truth_track.getPerigeeParameters())
319 ldmx_log(debug) << par << " ";
320
321 // assign the sim hit indices
322 // TODO this is not fully correct as the sim hits
323 // might be duplicated on sensors
324 // and should be merged if that is the case
325
326 int nhits = 0;
327
328 for (auto sim_hit_idx : hit_count_map.at(smeared_truth_track.getTrackID())) {
329 smeared_truth_track.addMeasurementIndex(sim_hit_idx);
330 nhits += 1;
331 }
332
333 smeared_truth_track.setNhits(nhits);
334
335 return smeared_truth_track;
336}

References createTruthTrack(), and seedFromTruth().

Member Data Documentation

◆ beam_electrons_collection_

std::string tracking::reco::TruthSeedProcessor::beam_electrons_collection_
private

Definition at line 285 of file TruthSeedProcessor.h.

◆ beam_origin_

std::vector<double> tracking::reco::TruthSeedProcessor::beam_origin_ {-880.1, -44., 0.}
private

Definition at line 282 of file TruthSeedProcessor.h.

282{-880.1, -44., 0.};

◆ d0smear_

std::vector<double> tracking::reco::TruthSeedProcessor::d0smear_
private

Definition at line 275 of file TruthSeedProcessor.h.

◆ ecal_sp_coll_name_

std::string tracking::reco::TruthSeedProcessor::ecal_sp_coll_name_ {"EcalScoringPlaneHits"}
private

Definition at line 199 of file TruthSeedProcessor.h.

199{"EcalScoringPlaneHits"};

◆ gctx_

Acts::GeometryContext tracking::reco::TruthSeedProcessor::gctx_
private

The ACTS geometry context properly.

Definition at line 192 of file TruthSeedProcessor.h.

Referenced by createTruthTrack(), and onNewRun().

◆ generator_

std::default_random_engine tracking::reco::TruthSeedProcessor::generator_
private

Definition at line 270 of file TruthSeedProcessor.h.

◆ inflate_factors_

std::vector<double> tracking::reco::TruthSeedProcessor::inflate_factors_
private

Definition at line 281 of file TruthSeedProcessor.h.

◆ input_pass_name_

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

Pass name for the sim hit collections.

Definition at line 209 of file TruthSeedProcessor.h.

209{""};

Referenced by configure(), and produce().

◆ linpropagator_

std::shared_ptr<LinPropagator> tracking::reco::TruthSeedProcessor::linpropagator_
private

Definition at line 262 of file TruthSeedProcessor.h.

◆ max_track_id_

int tracking::reco::TruthSeedProcessor::max_track_id_ {5}
private

Definition at line 260 of file TruthSeedProcessor.h.

260{5};

◆ n_min_hits_recoil_

int tracking::reco::TruthSeedProcessor::n_min_hits_recoil_ {7}
private

Minimum number of hits left in the recoil tracker to consider the seed as findable.

Definition at line 224 of file TruthSeedProcessor.h.

224{7};

Referenced by configure().

◆ n_min_hits_tagger_

int tracking::reco::TruthSeedProcessor::n_min_hits_tagger_ {7}
private

Minimum number of hits left in the recoil tracker to consider the seed as findable.

Definition at line 218 of file TruthSeedProcessor.h.

218{7};

Referenced by configure().

◆ normal_

std::shared_ptr<std::normal_distribution<float> > tracking::reco::TruthSeedProcessor::normal_
private

Definition at line 271 of file TruthSeedProcessor.h.

◆ p_cut_

double tracking::reco::TruthSeedProcessor::p_cut_ {0.}
private

Ask for a minimum p for the seeds.

Definition at line 239 of file TruthSeedProcessor.h.

239{0.};

Referenced by configure(), produce(), and scoringPlaneHitFilter().

◆ p_cut_ecal_

double tracking::reco::TruthSeedProcessor::p_cut_ecal_ {-1.}
private

Definition at line 245 of file TruthSeedProcessor.h.

245{-1.};

◆ p_cut_max_

double tracking::reco::TruthSeedProcessor::p_cut_max_ {100000.}
private

Ask for a maximum p for the seeds.

Definition at line 242 of file TruthSeedProcessor.h.

242{100000.};

Referenced by configure(), and scoringPlaneHitFilter().

◆ particle_hypothesis_

int tracking::reco::TruthSeedProcessor::particle_hypothesis_
private

Definition at line 283 of file TruthSeedProcessor.h.

◆ pdg_ids_

std::vector<int> tracking::reco::TruthSeedProcessor::pdg_ids_ {11}
private

pdg_ids of the particles we want to select for the seeds

Definition at line 195 of file TruthSeedProcessor.h.

195{11};

Referenced by configure(), and scoringPlaneHitFilter().

◆ phismear_

double tracking::reco::TruthSeedProcessor::phismear_
private

Definition at line 277 of file TruthSeedProcessor.h.

◆ pz_cut_

double tracking::reco::TruthSeedProcessor::pz_cut_ {-9999}
private

Ask for a minimum pz for the seeds.

Definition at line 236 of file TruthSeedProcessor.h.

236{-9999};

Referenced by configure(), and scoringPlaneHitFilter().

◆ recoil_seeds_collection_

std::string tracking::reco::TruthSeedProcessor::recoil_seeds_collection_
private

Definition at line 289 of file TruthSeedProcessor.h.

◆ recoil_sim_hits_coll_name_

std::string tracking::reco::TruthSeedProcessor::recoil_sim_hits_coll_name_ {"RecoilSimHits"}
private

Sim hits to check if the truth seed is findable.

Definition at line 206 of file TruthSeedProcessor.h.

206{"RecoilSimHits"};

Referenced by configure(), and produce().

◆ recoil_sp_

bool tracking::reco::TruthSeedProcessor::recoil_sp_ {true}
private

Definition at line 248 of file TruthSeedProcessor.h.

248{true};

◆ recoil_truth_collection_

std::string tracking::reco::TruthSeedProcessor::recoil_truth_collection_
private

Definition at line 287 of file TruthSeedProcessor.h.

◆ rel_smearfactors_

std::vector<double> tracking::reco::TruthSeedProcessor::rel_smearfactors_
private

Definition at line 280 of file TruthSeedProcessor.h.

◆ relpsmear_

double tracking::reco::TruthSeedProcessor::relpsmear_
private

Definition at line 279 of file TruthSeedProcessor.h.

◆ scoring_hits_coll_name_

std::string tracking::reco::TruthSeedProcessor::scoring_hits_coll_name_ {"TargetScoringPlaneHits"}
private

Which scoring plane hits to use for the truth seeds generation.

Definition at line 198 of file TruthSeedProcessor.h.

198{"TargetScoringPlaneHits"};

Referenced by configure(), and produce().

◆ seed_smearing_

bool tracking::reco::TruthSeedProcessor::seed_smearing_ {false}
private

Definition at line 273 of file TruthSeedProcessor.h.

273{false};

◆ sim_particles_coll_name_

std::string tracking::reco::TruthSeedProcessor::sim_particles_coll_name_
private

Definition at line 211 of file TruthSeedProcessor.h.

◆ sim_particles_passname_

std::string tracking::reco::TruthSeedProcessor::sim_particles_passname_
private

Definition at line 212 of file TruthSeedProcessor.h.

◆ skip_recoil_

bool tracking::reco::TruthSeedProcessor::skip_recoil_ {false}
private

Definition at line 257 of file TruthSeedProcessor.h.

257{false};

◆ skip_tagger_

bool tracking::reco::TruthSeedProcessor::skip_tagger_ {false}
private

Definition at line 254 of file TruthSeedProcessor.h.

254{false};

◆ sp_pass_name_

std::string tracking::reco::TruthSeedProcessor::sp_pass_name_ {""}
private

Definition at line 200 of file TruthSeedProcessor.h.

200{""};

◆ tagger_seeds_collection_

std::string tracking::reco::TruthSeedProcessor::tagger_seeds_collection_
private

Definition at line 288 of file TruthSeedProcessor.h.

◆ tagger_sim_hits_coll_name_

std::string tracking::reco::TruthSeedProcessor::tagger_sim_hits_coll_name_ {"TaggerSimHits"}
private

Sim hits to check if the truth seed is findable.

Definition at line 203 of file TruthSeedProcessor.h.

203{"TaggerSimHits"};

Referenced by configure(), and produce().

◆ tagger_truth_collection_

std::string tracking::reco::TruthSeedProcessor::tagger_truth_collection_
private

Definition at line 286 of file TruthSeedProcessor.h.

◆ target_sp_

bool tracking::reco::TruthSeedProcessor::target_sp_ {true}
private

Definition at line 251 of file TruthSeedProcessor.h.

251{true};

◆ thetasmear_

double tracking::reco::TruthSeedProcessor::thetasmear_
private

Definition at line 278 of file TruthSeedProcessor.h.

◆ track_id_

int tracking::reco::TruthSeedProcessor::track_id_ {-999}
private

Only select a particular trackID.

Definition at line 233 of file TruthSeedProcessor.h.

233{-999};

Referenced by configure(), and scoringPlaneHitFilter().

◆ trk_extrap_

std::shared_ptr<tracking::reco::TrackExtrapolatorTool<LinPropagator> > tracking::reco::TruthSeedProcessor::trk_extrap_
private

Definition at line 266 of file TruthSeedProcessor.h.

◆ z0smear_

std::vector<double> tracking::reco::TruthSeedProcessor::z0smear_
private

Definition at line 276 of file TruthSeedProcessor.h.

◆ z_min_

float tracking::reco::TruthSeedProcessor::z_min_ {-999}
private

Min cut on the z_ of the scoring hit.

It could be used to clean the scoring hits if desired.

Definition at line 230 of file TruthSeedProcessor.h.

230{-999};

Referenced by configure(), and scoringPlaneHitFilter().


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