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::unique_ptr< const TruthPropagator > propagator_
 
std::shared_ptr< tracking::reco::TrackExtrapolatorTool< TruthPropagator > > trk_extrap_
 
std::string field_map_ {""}
 Path to the magnetic field map.
 
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 43 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 54 of file TruthSeedProcessor.cxx.

54 {
56 parameters.get<std::string>("scoring_hits_coll_name");
57 ecal_sp_coll_name_ = parameters.get<std::string>("ecal_sp_coll_name");
58 sp_pass_name_ = parameters.get<std::string>("sp_pass_name");
60 parameters.get<std::string>("recoil_sim_hits_coll_name");
62 parameters.get<std::string>("tagger_sim_hits_coll_name");
63 input_pass_name_ = parameters.get<std::string>("input_pass_name");
64 sim_particles_coll_name_ =
65 parameters.get<std::string>("sim_particles_coll_name");
66 sim_particles_passname_ =
67 parameters.get<std::string>("sim_particles_passname");
68
69 n_min_hits_tagger_ = parameters.get<int>("n_min_hits_tagger", 11);
70 n_min_hits_recoil_ = parameters.get<int>("n_min_hits_recoil", 7);
71 pdg_ids_ = parameters.get<std::vector<int>>("pdg_ids", {11});
72 z_min_ = parameters.get<double>("z_min", -9999); // mm
73 track_id_ = parameters.get<int>("track_id", -9999);
74 pz_cut_ = parameters.get<double>("pz_cut", -9999); // MeV
75 p_cut_ = parameters.get<double>("p_cut", 0.);
76 p_cut_max_ = parameters.get<double>("p_cut_max", 100000.); // MeV
77 p_cut_ecal_ = parameters.get<double>("p_cut_ecal", -1.); // MeV
78 recoil_sp_ = parameters.get<double>("recoil_sp", true);
79 target_sp_ = parameters.get<double>("tagger_sp", true);
80 seed_smearing_ = parameters.get<bool>("seedSmearing", false);
81 max_track_id_ = parameters.get<int>("max_track_id", 5);
82
83 ldmx_log(info) << "Seed Smearing is set to " << seed_smearing_;
84
85 d0smear_ = parameters.get<std::vector<double>>("d0smear", {0.01, 0.01, 0.01});
86 z0smear_ = parameters.get<std::vector<double>>("z0smear", {0.1, 0.1, 0.1});
87 phismear_ = parameters.get<double>("phismear", 0.001);
88 thetasmear_ = parameters.get<double>("thetasmear", 0.001);
89 relpsmear_ = parameters.get<double>("relpsmear", 0.1);
90
91 // Relative smear factor terms, only used if seed_smearing_ is true.
92 rel_smearfactors_ = parameters.get<std::vector<double>>(
93 "rel_smearfactors", {0.1, 0.1, 0.1, 0.1, 0.1, 0.1});
94 inflate_factors_ = parameters.get<std::vector<double>>(
95 "inflate_factors", {10., 10., 10., 10., 10., 10.});
96
97 // In tracking frame: where do these numbers come from?
98 // These numbers come from approximating the path of the beam up
99 // until it is about to enter the first detector volume (TriggerPad1).
100 // In detector coordinates, (x_,y_,z_) = (-21.7, -883) is
101 // where the beam arrives (if no smearing is applied) and we simply
102 // reorder these values so that they are in tracking coordinates.
103 beam_origin_ = parameters.get<std::vector<double>>("beamOrigin",
104 {-883.0, -21.745876, 0.0});
105
106 // Skip the tagger or recoil trackers if wanted
107 skip_tagger_ = parameters.get<bool>("skip_tagger", false);
108 skip_recoil_ = parameters.get<bool>("skip_recoil", false);
109 particle_hypothesis_ = parameters.get<int>("particle_hypothesis");
110
111 beam_electrons_collection_ =
112 parameters.get<std::string>("beam_electrons_collection");
113 tagger_truth_collection_ =
114 parameters.get<std::string>("tagger_truth_collection");
115 recoil_truth_collection_ =
116 parameters.get<std::string>("recoil_truth_collection");
117 tagger_seeds_collection_ =
118 parameters.get<std::string>("tagger_seeds_collection");
119 recoil_seeds_collection_ =
120 parameters.get<std::string>("recoil_seeds_collection");
121 // Get the field map for the propagator
122 field_map_ = parameters.get<std::string>("field_map");
123}
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 field_map_
Path to the magnetic field map.
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 field_map_, 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 125 of file TruthSeedProcessor.cxx.

127 {
128 std::vector<double> pos{static_cast<double>(hit.getPosition()[0]),
129 static_cast<double>(hit.getPosition()[1]),
130 static_cast<double>(hit.getPosition()[2])};
131 createTruthTrack(pos, hit.getMomentum(), particle.getCharge(), trk,
132 target_surface);
133
134 trk.setTrackID(hit.getTrackID());
135 trk.setPdgID(hit.getPdgID());
136}
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 138 of file TruthSeedProcessor.cxx.

140 {
141 createTruthTrack(particle.getVertex(), particle.getMomentum(),
142 particle.getCharge(), trk, target_surface);
143
144 trk.setPdgID(particle.getPdgID());
145}
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:86
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 147 of file TruthSeedProcessor.cxx.

150 {
151 // Get the position and momentum of the particle at the point of creation.
152 // This only works for the incident electron when creating a tagger tracker
153 // seed. For the recoil tracker, the scoring plane position will need to
154 // be used. For other particles created in the target or tracker planes,
155 // this version of the method can also be used.
156 // These are just Eigen vectors defined as
157 // Eigen::Matrix<double, kSize, 1>;
158 Acts::Vector3 pos{pos_vec[0], pos_vec[1], pos_vec[2]};
159 Acts::Vector3 mom{p_vec[0], p_vec[1], p_vec[2]};
160
161 // Rotate the position and momentum into the ACTS frame.
162 pos = tracking::sim::utils::ldmx2Acts(pos);
163 mom = tracking::sim::utils::ldmx2Acts(mom);
164
165 // Get the charge of the particle.
166 // TODO: Add function that uses the PDG ID to calculate this.
167 double q{charge * Acts::UnitConstants::e};
168
169 // The idea here is:
170 // 1 - Define a bound track state parameters at point P on track. Basically a
171 // curvilinear representation.
172 // 2 - Propagate to target surface to obtain the
173 // BoundTrackState there.
174
175 // Transform the position, momentum and charge to free parameters.
176 auto free_params{tracking::sim::utils::toFreeParameters(pos, mom, q)};
177
178 // Create a line surface at the perigee. The perigee position is extracted
179 // from a particle's vertex or the particle's position at a specific
180 // scoring plane.
181 auto gen_surface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
182 Acts::Vector3(free_params[Acts::eFreePos0], free_params[Acts::eFreePos1],
183 free_params[Acts::eFreePos2]))};
184
185 // ldmx_log(trace)<<"PF:: gen_surface"<<free_params[Acts::eFreePos0]<<"
186 // " <<free_params[Acts::eFreePos1]<<" " <<free_params[Acts::eFreePos2];
187
188 // Transform the parameters to local positions on the perigee surface.
189 auto bound_params{
190 Acts::transformFreeToBoundParameters(free_params, *gen_surface, gctx_)
191 .value()};
192 // Create a particle hypothesis
193 auto part{Acts::GenericParticleHypothesis(
194 Acts::ParticleHypothesis(Acts::PdgParticle(particle_hypothesis_)))};
195 Acts::BoundTrackParameters bound_trk_pars(gen_surface, bound_params,
196 std::nullopt, part);
197
198 auto prop_bound_state =
199 trk_extrap_->extrapolate(bound_trk_pars, target_surface);
200
201 if (!prop_bound_state) {
202 ldmx_log(warn) << "Propagation to target surface failed — "
203 << "track may have exited the B-field map.";
204 return;
205 }
206
207 // Create the seed track object.
208 Acts::Vector3 ref = target_surface->center(geometryContext());
209 Acts::Vector3 ref_ldmx = tracking::sim::utils::acts2Ldmx(ref);
210 trk.setPerigeeLocation(ref_ldmx(0), ref_ldmx(1), ref_ldmx(2));
211
212 auto prop_bound_vec = prop_bound_state->parameters();
213
214 trk.setPerigeeParameters(
215 tracking::sim::utils::convertActsToLdmxPars(prop_bound_vec));
216}
Acts::GeometryContext gctx_
The ACTS geometry context properly.

References gctx_.

◆ 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 472 of file TruthSeedProcessor.cxx.

474 {
475 for (int i_sim_hit = 0; i_sim_hit < sim_hits.size(); i_sim_hit++) {
476 auto& sim_hit = sim_hits.at(i_sim_hit);
477 // This track never left a hit before
478 if (!hit_count_map.count(sim_hit.getTrackID())) {
479 hit_count_map[sim_hit.getTrackID()].push_back(i_sim_hit);
480 }
481
482 // This track left a hit before.
483 // Check if it's on a different sensor than the others
484
485 else {
486 int sensor_id = tracking::sim::utils::getSensorID(sim_hit);
487 bool found_hit = false;
488
489 for (auto& i_rhit : hit_count_map[sim_hit.getTrackID()]) {
490 int tmp_sensor_id =
491 tracking::sim::utils::getSensorID(sim_hits.at(i_rhit));
492
493 if (sensor_id == tmp_sensor_id) {
494 found_hit = true;
495 break;
496 }
497 } // loop on the already recorded hits
498
499 if (!found_hit) {
500 hit_count_map[sim_hit.getTrackID()].push_back(i_sim_hit);
501 }
502 }
503 } // loop on sim hits
504}

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 // Custom transformation of the interpolated bfield map
16 auto transform_pos = [](const Acts::Vector3& pos_) {
17 Acts::Vector3 rot_pos;
18 rot_pos(0) = pos_(1);
19 rot_pos(1) = pos_(2);
20 rot_pos(2) = pos_(0) + DIPOLE_OFFSET;
21 return rot_pos;
22 };
23
24 auto transform_b_field = [](const Acts::Vector3& field,
25 const Acts::Vector3& /*pos_*/) {
26 Acts::Vector3 rot_field;
27 rot_field(0) = field(2);
28 rot_field(1) = field(0);
29 rot_field(2) = field(1);
30 return rot_field;
31 };
32
33 // Setup the interpolated bfield map
34 const auto map = std::make_shared<InterpolatedMagneticField3>(
35 loadDefaultBField(field_map_, transform_pos, transform_b_field));
36
37 // Setup the stepper and navigator
38 const auto stepper = Acts::EigenStepper<>{map};
39 Acts::Navigator::Config nav_cfg{geometry().getTG()};
40 nav_cfg.resolveMaterial = true;
41 nav_cfg.resolvePassive = true;
42 nav_cfg.resolveSensitive = true;
43 const Acts::Navigator navigator(nav_cfg);
44
45 propagator_ = std::make_unique<TruthPropagator>(
46 stepper, navigator,
47 Acts::getDefaultLogger("TruthPropagator", Acts::Logging::FATAL));
48 trk_extrap_ = std::make_shared<std::decay_t<decltype(*trk_extrap_)>>(
49 *propagator_, geometryContext(), magneticFieldContext());
50 trk_extrap_->setMaxStepSize(200); // mm
51 trk_extrap_->setPathLimit(3000); // mm
52}

References field_map_, and 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 74 of file TruthSeedProcessor.h.

74{};

◆ 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 560 of file TruthSeedProcessor.cxx.

560 {
561 // Retrieve the particleMap
562 auto particle_map{event.getMap<int, ldmx::SimParticle>(
563 sim_particles_coll_name_, sim_particles_passname_)};
564
565 // Retrieve the target scoring hits
566 // Information is extracted using the
567 // scoring plane hit left by the particle at the target.
568
569 const std::vector<ldmx::SimTrackerHit> scoring_hits{
570 event.getCollection<ldmx::SimTrackerHit>(scoring_hits_coll_name_,
571 sp_pass_name_)};
572
573 // Retrieve the scoring plane hits at the ECAL
574 const std::vector<ldmx::SimTrackerHit> scoring_hits_ecal{
575 event.getCollection<ldmx::SimTrackerHit>(ecal_sp_coll_name_,
576 sp_pass_name_)};
577
578 // Retrieve the sim hits in the tagger tracker
579 const std::vector<ldmx::SimTrackerHit> tagger_sim_hits =
582
583 // Retrieve the sim hits in the recoil tracker
584 const std::vector<ldmx::SimTrackerHit> recoil_sim_hits =
587
588 // If sim hit collections are empty throw a warning
589 if (tagger_sim_hits.size() == 0 && !skip_tagger_) {
590 ldmx_log(error) << "Tagger sim hits collection empty for event ";
591 }
592 if (recoil_sim_hits.size() == 0 && !skip_recoil_) {
593 ldmx_log(error) << "Recoil sim hits collection empty for event ";
594 }
595
596 // The map stores which track leaves which sim hits
597 std::map<int, std::vector<int>> hit_count_map_recoil;
598 makeHitCountMap(recoil_sim_hits, hit_count_map_recoil);
599
600 std::map<int, std::vector<int>> hit_count_map_tagger;
601 makeHitCountMap(tagger_sim_hits, hit_count_map_tagger);
602
603 // to keep track of how many sim particles leave hits on the scoring plane
604 std::vector<int> recoil_sh_idxs;
605 std::unordered_map<int, std::vector<int>> recoil_sh_count_map;
606
607 std::vector<int> tagger_sh_idxs;
608 std::unordered_map<int, std::vector<int>> tagger_sh_count_map;
609
610 // Target scoring hits for Tagger will have Z<0, Recoil scoring hits will
611 // have Z>0
612 for (unsigned int i_sh = 0; i_sh < scoring_hits.size(); i_sh++) {
613 const ldmx::SimTrackerHit& hit = scoring_hits.at(i_sh);
614 double zhit = hit.getPosition()[2];
615
616 Acts::Vector3 p_vec{hit.getMomentum()[0], hit.getMomentum()[1],
617 hit.getMomentum()[2]};
618 double tagger_p_max = 0.;
619
620 // Check if it is a tagger track going fwd that passes basic cuts
621 if (zhit < 0.) {
622 // Tagger selection cuts
623 // Negative scoring plane hit, with momentum > p_cut
624 if (p_vec(2) < 0. || p_vec.norm() < p_cut_) continue;
625
626 // Check that the hit was left by a charged particle
627 if (abs(particle_map[hit.getTrackID()].getCharge()) < 1e-8) continue;
628
629 if (p_vec.norm() > tagger_p_max) {
630 tagger_sh_count_map[hit.getTrackID()].push_back(i_sh);
631 }
632 } // Tagger loop
633
634 // Check the recoil hits
635 else {
636 // Recoil selection cuts
637 // Positive scoring plane hit, forward direction with momentum > p_cut
638 if (p_vec(2) < 0. || p_vec.norm() < p_cut_) continue;
639
640 // Check that the hit was left by a charged particle
641 if (abs(particle_map[hit.getTrackID()].getCharge()) < 1e-8) continue;
642
643 recoil_sh_count_map[hit.getTrackID()].push_back(i_sh);
644
645 } // Recoil
646 } // loop on Target scoring plane hits
647
648 for (std::pair<int, std::vector<int>> element : recoil_sh_count_map) {
649 std::sort(
650 element.second.begin(), element.second.end(),
651 [&](const int idx1, int idx2) -> bool {
652 const ldmx::SimTrackerHit& hit1 = scoring_hits.at(idx1);
653 const ldmx::SimTrackerHit& hit2 = scoring_hits.at(idx2);
654
655 Acts::Vector3 phit1{hit1.getMomentum()[0], hit1.getMomentum()[1],
656 hit1.getMomentum()[2]};
657 Acts::Vector3 phit2{hit2.getMomentum()[0], hit2.getMomentum()[1],
658 hit2.getMomentum()[2]};
659
660 return phit1.norm() > phit2.norm();
661 });
662 }
663
664 // Sort tagger hits.
665 for (auto& [_track_id, hit_indices] : tagger_sh_count_map) {
666 std::sort(
667 hit_indices.begin(), hit_indices.end(),
668 [&](const int idx1, int idx2) -> bool {
669 const ldmx::SimTrackerHit& hit1 = scoring_hits.at(idx1);
670 const ldmx::SimTrackerHit& hit2 = scoring_hits.at(idx2);
671
672 Acts::Vector3 phit1{hit1.getMomentum()[0], hit1.getMomentum()[1],
673 hit1.getMomentum()[2]};
674 Acts::Vector3 phit2{hit2.getMomentum()[0], hit2.getMomentum()[1],
675 hit2.getMomentum()[2]};
676
677 return phit1.norm() > phit2.norm();
678 });
679 }
680
681 // Building of the event truth information and the truth seeds
682 // TODO remove the truthtracks in the future as the truth seeds are enough
683
684 std::vector<ldmx::Track> tagger_truth_tracks;
685 std::vector<ldmx::Track> tagger_truth_seeds;
686 std::vector<ldmx::Track> recoil_truth_tracks;
687 std::vector<ldmx::Track> recoil_truth_seeds;
688 std::vector<ldmx::Track> beam_electrons;
689
690 // TODO:: The target should be taken from some conditions DB in the future.
691 // Define the perigee_surface at 0.0.0
692 auto target_surface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
693 Acts::Vector3(0., 0., 0.))};
694
695 // Define the target_surface
696 auto target_unbound_surface = tracking::sim::utils::unboundSurface(0.);
697
698 // ecal
699 auto ecal_surface = tracking::sim::utils::unboundSurface(240.5);
700
701 auto beam_origin_surface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
702 Acts::Vector3(beam_origin_[0], beam_origin_[1], beam_origin_[2]))};
703
704 if (!skip_tagger_) {
705 for (const auto& [track_id, hit_indices] : tagger_sh_count_map) {
706 const ldmx::SimTrackerHit& hit = scoring_hits.at(hit_indices.at(0));
707 const ldmx::SimParticle& phit = particle_map[hit.getTrackID()];
708
709 if (hit_count_map_tagger[hit.getTrackID()].size() > n_min_hits_tagger_) {
710 ldmx::Track truth_tagger_track;
711 createTruthTrack(phit, hit, truth_tagger_track, target_surface);
712 truth_tagger_track.setNhits(
713 hit_count_map_tagger[hit.getTrackID()].size());
714 // Add AtTarget state from the scoring plane hit (pos in mm, mom in MeV,
715 // both already in LDMX global frame)
716 ldmx::Track::TrackState ts_target;
717 ts_target.pos_ = {hit.getPosition()[0], hit.getPosition()[1],
718 hit.getPosition()[2]};
719 ts_target.mom_ = {hit.getMomentum()[0], hit.getMomentum()[1],
720 hit.getMomentum()[2]};
721 ts_target.ts_type_ = ldmx::AtTarget;
722 truth_tagger_track.addTrackState(ts_target);
723 tagger_truth_tracks.push_back(truth_tagger_track);
724
725 if (hit.getPdgID() == 11 && hit.getTrackID() < max_track_id_) {
726 ldmx::Track beam_e_truth_seed =
727 taggerFullSeed(particle_map[hit.getTrackID()], hit.getTrackID(),
728 hit, hit_count_map_tagger, beam_origin_surface,
729 target_unbound_surface);
730 beam_electrons.push_back(beam_e_truth_seed);
731 }
732 }
733 }
734 }
735
736 // Recover the EcalScoring hits
737 std::vector<ldmx::SimTrackerHit> ecal_sp_hits =
738 event.getCollection<ldmx::SimTrackerHit>(ecal_sp_coll_name_,
739 sp_pass_name_);
740 // Select ECAL hits
741 std::vector<ldmx::SimTrackerHit> sel_ecal_sp_hits;
742
743 for (auto sp_hit : ecal_sp_hits) {
744 if (sp_hit.getMomentum()[2] > 0 && ((sp_hit.getID() & 0xfff) == 31)) {
745 sel_ecal_sp_hits.push_back(sp_hit);
746 }
747 }
748
749 // Recoil target surface for truth and seed tracks is the target
750
751 for (std::pair<int, std::vector<int>> element : recoil_sh_count_map) {
752 // Only take the first entry of the vector: it should be the scoring plane
753 // hit with the highest momentum.
754 const ldmx::SimTrackerHit& hit = scoring_hits.at(element.second.at(0));
755 const ldmx::SimParticle& phit = particle_map[hit.getTrackID()];
756 ldmx::SimTrackerHit ecal_hit;
757
758 bool found_ecal_hit = false;
759 for (auto ecal_sp_hit : sel_ecal_sp_hits) {
760 if (ecal_sp_hit.getTrackID() == hit.getTrackID()) {
761 ecal_hit = ecal_sp_hit;
762 found_ecal_hit = true;
763 break;
764 }
765 }
766
767 // Findable particle selection
768 // ldmx_log(trace) << "!!! n_recoil_sim_hits found: " <<
769 // hit_count_map_recoil[hit.getTrackID()].size();
770 if (hit_count_map_recoil[hit.getTrackID()].size() > n_min_hits_recoil_ &&
771 found_ecal_hit && !skip_recoil_) {
772 ldmx::Track truth_recoil_track;
773 createTruthTrack(phit, hit, truth_recoil_track, target_surface);
774 truth_recoil_track.setTrackID(hit.getTrackID());
775
776 // AtTarget state from target scoring plane hit (pos mm, mom MeV, LDMX
777 // frame)
778 ldmx::Track::TrackState ts_target;
779 ts_target.pos_ = {hit.getPosition()[0], hit.getPosition()[1],
780 hit.getPosition()[2]};
781 ts_target.mom_ = {hit.getMomentum()[0], hit.getMomentum()[1],
782 hit.getMomentum()[2]};
783 ts_target.ts_type_ = ldmx::AtTarget;
784 truth_recoil_track.addTrackState(ts_target);
785
786 // Express truth ECAL state in the bound parametrization of ecal_surface
787 // (same surface definition used by CKFProcessor).
788 {
789 Acts::Vector3 ep{ecal_hit.getPosition()[0], ecal_hit.getPosition()[1],
790 ecal_hit.getPosition()[2]};
791 Acts::Vector3 em{ecal_hit.getMomentum()[0], ecal_hit.getMomentum()[1],
792 ecal_hit.getMomentum()[2]};
793 ep = tracking::sim::utils::ldmx2Acts(ep);
794 em = tracking::sim::utils::ldmx2Acts(em);
795 // Linearly extrapolate transverse coordinates to ACTS x = 240.5 mm
796 // (= LDMX z = 240.5 mm), correcting for the track slope over the small
797 // z-offset between the scoring plane and the ECAL surface.
798 if (std::abs(em[0]) > 0) {
799 double delta = 240.5 - ep[0];
800 ep[1] += delta * em[1] / em[0];
801 ep[2] += delta * em[2] / em[0];
802 ep[0] = 240.5;
803 }
804 double q_ecal = phit.getCharge() * Acts::UnitConstants::e;
805 auto ecal_free = tracking::sim::utils::toFreeParameters(ep, em, q_ecal);
806 auto ecal_bound = Acts::transformFreeToBoundParameters(
807 ecal_free, *ecal_surface, gctx_);
808 if (ecal_bound.ok()) {
809 auto part{Acts::GenericParticleHypothesis(Acts::ParticleHypothesis(
810 Acts::PdgParticle(particle_hypothesis_)))};
811 Acts::BoundTrackParameters ecal_pars(
812 ecal_surface, ecal_bound.value(),
813 Acts::BoundSquareMatrix::Identity(), part);
814 truth_recoil_track.addTrackState(tracking::sim::utils::makeTrackState(
815 geometryContext(), ecal_pars, ldmx::AtECAL));
816 }
817 }
818
819 // Attach sim hit indices
820 int nhits = 0;
821 for (auto sim_hit_idx : hit_count_map_recoil.at(hit.getTrackID())) {
822 truth_recoil_track.addMeasurementIndex(sim_hit_idx);
823 nhits++;
824 }
825 truth_recoil_track.setNhits(nhits);
826
827 recoil_truth_tracks.push_back(truth_recoil_track);
828 }
829 }
830
831 /*
832 for (std::pair<int,std::vector<int>> element : recoil_sh_count_map) {
833
834 const ldmx::SimTrackerHit& hit = scoring_hits.at(element.second.at(0));
835 const ldmx::SimParticle& phit = particleMap[hit.getTrackID()];
836
837 if (hit_count_map_recoil[hit.getTrackID()].size() > n_min_hits_recoil_) {
838 ldmx::Track truth_recoil_track;
839 createTruthTrack(phit,hit,truth_recoil_track,targetSurface);
840 truth_recoil_track.setNhits(hit_count_map_recoil[hit.getTrackID()].size());
841 recoil_truth_tracks.push_back(truth_recoil_track);
842 }
843 }
844 */
845
846 // Form a truth seed from a truth track
847
848 for (auto& tt : tagger_truth_tracks) {
849 ldmx::Track seed = seedFromTruth(tt, seed_smearing_);
850
851 tagger_truth_seeds.push_back(seed);
852 }
853
854 ldmx_log(debug) << "Forming seeds from truth";
855 for (auto& tt : recoil_truth_tracks) {
856 ldmx_log(debug) << "Smearing truth track";
857
858 ldmx::Track seed = seedFromTruth(tt, seed_smearing_);
859
860 recoil_truth_seeds.push_back(seed);
861 }
862
863 // even if skip_tagger/recoil_ is true, still make the collections in the
864 // event
865 event.add(beam_electrons_collection_, beam_electrons);
866 event.add(tagger_truth_collection_, tagger_truth_tracks);
867 event.add(recoil_truth_collection_, recoil_truth_tracks);
868 event.add(tagger_seeds_collection_, tagger_truth_seeds);
869 event.add(recoil_seeds_collection_, recoil_truth_seeds);
870}
Class representing a simulated particle.
Definition SimParticle.h:24
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 220 of file TruthSeedProcessor.cxx.

226 {
227 ldmx::Track truth_recoil_track;
228 createTruthTrack(particle, hit, truth_recoil_track, origin_surface);
229 truth_recoil_track.setTrackID(trackID);
230
231 // Seed at the target location
232 ldmx::Track smeared_truth_track = seedFromTruth(truth_recoil_track, false);
233
234 // Add truth track state at the target: use scoring plane hit (already LDMX
235 // frame)
236 ldmx::Track::TrackState ts_truth_target;
237 ts_truth_target.pos_ = {hit.getPosition()[0], hit.getPosition()[1],
238 hit.getPosition()[2]};
239 ts_truth_target.mom_ = {hit.getMomentum()[0], hit.getMomentum()[1],
240 hit.getMomentum()[2]};
241 ts_truth_target.ts_type_ = ldmx::AtTarget;
242 smeared_truth_track.addTrackState(ts_truth_target);
243
244 // Express truth ECAL state in the bound parametrization of ecal_surface
245 // (same surface definition used by CKFProcessor) rather than storing raw
246 // scoring plane hit coordinates.
247 {
248 Acts::Vector3 ep{ecal_hit.getPosition()[0], ecal_hit.getPosition()[1],
249 ecal_hit.getPosition()[2]};
250 Acts::Vector3 em{ecal_hit.getMomentum()[0], ecal_hit.getMomentum()[1],
251 ecal_hit.getMomentum()[2]};
252 ep = tracking::sim::utils::ldmx2Acts(ep);
253 em = tracking::sim::utils::ldmx2Acts(em);
254 // Linearly extrapolate transverse coordinates to ACTS x = 240.5 mm
255 // (= LDMX z = 240.5 mm), correcting for the track slope over the small
256 // z-offset between the scoring plane and the ECAL surface.
257 if (std::abs(em[0]) > 0) {
258 double delta = 240.5 - ep[0];
259 ep[1] += delta * em[1] / em[0];
260 ep[2] += delta * em[2] / em[0];
261 ep[0] = 240.5;
262 }
263 double q_ecal = particle.getCharge() * Acts::UnitConstants::e;
264 auto ecal_free = tracking::sim::utils::toFreeParameters(ep, em, q_ecal);
265 auto ecal_bound =
266 Acts::transformFreeToBoundParameters(ecal_free, *ecal_surface, gctx_);
267 if (ecal_bound.ok()) {
268 auto part{Acts::GenericParticleHypothesis(
269 Acts::ParticleHypothesis(Acts::PdgParticle(particle_hypothesis_)))};
270 Acts::BoundTrackParameters ecal_pars(ecal_surface, ecal_bound.value(),
271 Acts::BoundSquareMatrix::Identity(),
272 part);
273 smeared_truth_track.addTrackState(tracking::sim::utils::makeTrackState(
274 geometryContext(), ecal_pars, ldmx::AtECAL));
275 }
276 }
277
278 // Add the hits
279 int nhits = 0;
280
281 for (auto sim_hit_idx : hit_count_map.at(smeared_truth_track.getTrackID())) {
282 smeared_truth_track.addMeasurementIndex(sim_hit_idx);
283 nhits += 1;
284 }
285
286 smeared_truth_track.setNhits(nhits);
287
288 return smeared_truth_track;
289}

◆ 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 506 of file TruthSeedProcessor.cxx.

508 {
509 // Clean some of the hits we don't want
510 if (hit.getPosition()[2] < z_min_) return false;
511
512 // Check if the track_id was requested
513 if (track_id_ > 0 && hit.getTrackID() != track_id_) return false;
514
515 // Check if we are requesting particular particles
516 if (std::find(pdg_ids_.begin(), pdg_ids_.end(), hit.getPdgID()) ==
517 pdg_ids_.end())
518 return false;
519
520 Acts::Vector3 p_vec{hit.getMomentum()[0], hit.getMomentum()[1],
521 hit.getMomentum()[2]};
522
523 // p cut
524 if (p_cut_ >= 0. && p_vec.norm() < p_cut_) return false;
525
526 // p cut Max
527 if (p_cut_ < 100000. && p_vec.norm() > p_cut_max_) return false;
528
529 // pz cut
530 if (pz_cut_ > -9999 && p_vec(2) < pz_cut_) return false;
531
532 // Check the ecal scoring plane
533 bool pass_ecal_scoring_plane = true;
534
535 if (p_cut_ecal_ > 0) { // only check if we care about it.
536
537 for (auto& e_sp_hit : ecal_sp_hits) {
538 if (e_sp_hit.getTrackID() == hit.getTrackID() &&
539 e_sp_hit.getPdgID() == hit.getPdgID()) {
540 Acts::Vector3 e_sp_p{e_sp_hit.getMomentum()[0],
541 e_sp_hit.getMomentum()[1],
542 e_sp_hit.getMomentum()[2]};
543
544 if (e_sp_p.norm() < p_cut_ecal_) pass_ecal_scoring_plane = false;
545
546 // Skip the rest of the scoring plane hits since we already found the
547 // track we care about
548 break;
549
550 } // check that the hit belongs to the inital particle from the target
551 // scoring plane hit
552 } // loop on Ecal scoring plane hits
553 } // pcutEcal
554
555 if (!pass_ecal_scoring_plane) return false;
556
557 return true;
558}

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 356 of file TruthSeedProcessor.cxx.

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

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 291 of file TruthSeedProcessor.cxx.

296 {
297 ldmx::Track truth_track;
298
299 createTruthTrack(beam_electron, truth_track, origin_surface);
300 truth_track.setTrackID(trackID);
301
302 // Smeared track at the beam origin
303 ldmx::Track smeared_truth_track = seedFromTruth(truth_track, true);
304
305 ldmx_log(debug) << "Truth parameters at beam origin";
306 for (auto par : truth_track.getPerigeeParameters())
307 ldmx_log(debug) << par << " ";
308 ldmx_log(debug);
309
310 // Add the truth track state at the target using the scoring plane hit
311 // (position and momentum are already in LDMX global frame)
312 ldmx::Track::TrackState ts_truth_target;
313 ts_truth_target.pos_ = {hit.getPosition()[0], hit.getPosition()[1],
314 hit.getPosition()[2]};
315 ts_truth_target.mom_ = {hit.getMomentum()[0], hit.getMomentum()[1],
316 hit.getMomentum()[2]};
317 ts_truth_target.ts_type_ = ldmx::AtTarget;
318 smeared_truth_track.addTrackState(ts_truth_target);
319
320 ldmx_log(debug) << "Truth position at target: " << hit.getPosition()[0] << " "
321 << hit.getPosition()[1] << " " << hit.getPosition()[2];
322
323 // Add the truth track state at the beam origin using particle vertex/momentum
324 // (SimParticle vertex and momentum are in LDMX global frame)
325 ldmx::Track::TrackState ts_truth_beam_origin;
326 ts_truth_beam_origin.pos_ = {beam_electron.getVertex()[0],
327 beam_electron.getVertex()[1],
328 beam_electron.getVertex()[2]};
329 ts_truth_beam_origin.mom_ = {beam_electron.getMomentum()[0],
330 beam_electron.getMomentum()[1],
331 beam_electron.getMomentum()[2]};
332 ts_truth_beam_origin.ts_type_ = ldmx::AtBeamOrigin;
333 smeared_truth_track.addTrackState(ts_truth_beam_origin);
334
335 ldmx_log(debug) << "Smeared parameters at origin";
336 for (auto par : smeared_truth_track.getPerigeeParameters())
337 ldmx_log(debug) << par << " ";
338
339 // assign the sim hit indices
340 // TODO this is not fully correct as the sim hits
341 // might be duplicated on sensors
342 // and should be merged if that is the case
343
344 int nhits = 0;
345
346 for (auto sim_hit_idx : hit_count_map.at(smeared_truth_track.getTrackID())) {
347 smeared_truth_track.addMeasurementIndex(sim_hit_idx);
348 nhits += 1;
349 }
350
351 smeared_truth_track.setNhits(nhits);
352
353 return smeared_truth_track;
354}

References createTruthTrack(), ldmx::SimParticle::getMomentum(), ldmx::SimTrackerHit::getMomentum(), ldmx::SimTrackerHit::getPosition(), ldmx::SimParticle::getVertex(), and seedFromTruth().

Member Data Documentation

◆ beam_electrons_collection_

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

Definition at line 289 of file TruthSeedProcessor.h.

◆ beam_origin_

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

Definition at line 286 of file TruthSeedProcessor.h.

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

◆ d0smear_

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

Definition at line 279 of file TruthSeedProcessor.h.

◆ ecal_sp_coll_name_

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

Definition at line 200 of file TruthSeedProcessor.h.

200{"EcalScoringPlaneHits"};

◆ field_map_

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

Path to the magnetic field map.

Definition at line 270 of file TruthSeedProcessor.h.

270{""};

Referenced by configure(), and onNewRun().

◆ gctx_

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

The ACTS geometry context properly.

Definition at line 193 of file TruthSeedProcessor.h.

Referenced by createTruthTrack(), and onNewRun().

◆ generator_

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

Definition at line 274 of file TruthSeedProcessor.h.

◆ inflate_factors_

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

Definition at line 285 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 210 of file TruthSeedProcessor.h.

210{""};

Referenced by configure(), and produce().

◆ max_track_id_

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

Definition at line 261 of file TruthSeedProcessor.h.

261{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 225 of file TruthSeedProcessor.h.

225{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 219 of file TruthSeedProcessor.h.

219{7};

Referenced by configure().

◆ normal_

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

Definition at line 275 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 240 of file TruthSeedProcessor.h.

240{0.};

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

◆ p_cut_ecal_

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

Definition at line 246 of file TruthSeedProcessor.h.

246{-1.};

◆ p_cut_max_

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

Ask for a maximum p for the seeds.

Definition at line 243 of file TruthSeedProcessor.h.

243{100000.};

Referenced by configure(), and scoringPlaneHitFilter().

◆ particle_hypothesis_

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

Definition at line 287 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 196 of file TruthSeedProcessor.h.

196{11};

Referenced by configure(), and scoringPlaneHitFilter().

◆ phismear_

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

Definition at line 281 of file TruthSeedProcessor.h.

◆ propagator_

std::unique_ptr<const TruthPropagator> tracking::reco::TruthSeedProcessor::propagator_
private

Definition at line 263 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 237 of file TruthSeedProcessor.h.

237{-9999};

Referenced by configure(), and scoringPlaneHitFilter().

◆ recoil_seeds_collection_

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

Definition at line 293 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 207 of file TruthSeedProcessor.h.

207{"RecoilSimHits"};

Referenced by configure(), and produce().

◆ recoil_sp_

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

Definition at line 249 of file TruthSeedProcessor.h.

249{true};

◆ recoil_truth_collection_

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

Definition at line 291 of file TruthSeedProcessor.h.

◆ rel_smearfactors_

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

Definition at line 284 of file TruthSeedProcessor.h.

◆ relpsmear_

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

Definition at line 283 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 199 of file TruthSeedProcessor.h.

199{"TargetScoringPlaneHits"};

Referenced by configure(), and produce().

◆ seed_smearing_

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

Definition at line 277 of file TruthSeedProcessor.h.

277{false};

◆ sim_particles_coll_name_

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

Definition at line 212 of file TruthSeedProcessor.h.

◆ sim_particles_passname_

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

Definition at line 213 of file TruthSeedProcessor.h.

◆ skip_recoil_

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

Definition at line 258 of file TruthSeedProcessor.h.

258{false};

◆ skip_tagger_

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

Definition at line 255 of file TruthSeedProcessor.h.

255{false};

◆ sp_pass_name_

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

Definition at line 201 of file TruthSeedProcessor.h.

201{""};

◆ tagger_seeds_collection_

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

Definition at line 292 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 204 of file TruthSeedProcessor.h.

204{"TaggerSimHits"};

Referenced by configure(), and produce().

◆ tagger_truth_collection_

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

Definition at line 290 of file TruthSeedProcessor.h.

◆ target_sp_

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

Definition at line 252 of file TruthSeedProcessor.h.

252{true};

◆ thetasmear_

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

Definition at line 282 of file TruthSeedProcessor.h.

◆ track_id_

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

Only select a particular trackID.

Definition at line 234 of file TruthSeedProcessor.h.

234{-999};

Referenced by configure(), and scoringPlaneHitFilter().

◆ trk_extrap_

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

Definition at line 267 of file TruthSeedProcessor.h.

◆ z0smear_

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

Definition at line 280 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 231 of file TruthSeedProcessor.h.

231{-999};

Referenced by configure(), and scoringPlaneHitFilter().


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