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 beforeNewRun (ldmx::RunHeader &header)
 Handle allowing producers to modify run headers before the run begins.
 
- Public Member Functions inherited from framework::EventProcessor
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()
 Class destructor.
 
virtual void onFileOpen (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a event input ROOT file is closed.
 
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 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.
 
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 seedSmearing_ {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 > beamOrigin_ {-880.1, -44., 0.}
 

Additional Inherited Members

- Static Public Member Functions inherited from framework::EventProcessor
static void declare (const std::string &classname, int classtype, EventProcessorMaker *)
 Internal function which is part of the PluginFactory machinery.
 
- Static Public Attributes inherited from framework::Producer
static const int CLASSTYPE {1}
 Constant used to track EventProcessor types by the PluginFactory.
 
- Protected Member Functions inherited from tracking::reco::TrackingGeometryUser
const Acts::GeometryContext & geometry_context ()
 
const Acts::MagneticFieldContext & magnetic_field_context ()
 
const Acts::CalibrationContext & calibration_context ()
 
const geo::TrackersTrackingGeometrygeometry ()
 
- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 
- Protected Attributes inherited from framework::EventProcessor
HistogramHelper histograms_
 Interface class for making and filling histograms.
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger theLog_
 The logger for this EventProcessor.
 

Detailed Description

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

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 className 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.getParameter<std::string>("scoring_hits_coll_name");
27 sp_pass_name_ = parameters.getParameter<std::string>("sp_pass_name", "");
29 parameters.getParameter<std::string>("recoil_sim_hits_coll_name");
31 parameters.getParameter<std::string>("tagger_sim_hits_coll_name");
33 parameters.getParameter<std::string>("input_pass_name", "");
34
35 n_min_hits_tagger_ = parameters.getParameter<int>("n_min_hits_tagger", 11);
36 n_min_hits_recoil_ = parameters.getParameter<int>("n_min_hits_recoil", 7);
37 pdg_ids_ = parameters.getParameter<std::vector<int>>("pdg_ids", {11});
38 z_min_ = parameters.getParameter<double>("z_min", -9999); // mm
39 track_id_ = parameters.getParameter<int>("track_id", -9999);
40 pz_cut_ = parameters.getParameter<double>("pz_cut", -9999); // MeV
41 p_cut_ = parameters.getParameter<double>("p_cut", 0.);
42 p_cut_max_ = parameters.getParameter<double>("p_cut_max", 100000.); // MeV
43 p_cut_ecal_ = parameters.getParameter<double>("p_cut_ecal", -1.); // MeV
44 recoil_sp_ = parameters.getParameter<double>("recoil_sp", true);
45 target_sp_ = parameters.getParameter<double>("tagger_sp", true);
46 seedSmearing_ = parameters.getParameter<bool>("seedSmearing", false);
47 max_track_id_ = parameters.getParameter<int>("max_track_id", 5);
48
49 ldmx_log(info) << "Seed Smearing is set to " << seedSmearing_;
50
51 d0smear_ = parameters.getParameter<std::vector<double>>("d0smear",
52 {0.01, 0.01, 0.01});
53 z0smear_ =
54 parameters.getParameter<std::vector<double>>("z0smear", {0.1, 0.1, 0.1});
55 phismear_ = parameters.getParameter<double>("phismear", 0.001);
56 thetasmear_ = parameters.getParameter<double>("thetasmear", 0.001);
57 relpsmear_ = parameters.getParameter<double>("relpsmear", 0.1);
58
59 // Relative smear factor terms, only used if seedSmearing_ is true.
60 rel_smearfactors_ = parameters.getParameter<std::vector<double>>(
61 "rel_smearfactors", {0.1, 0.1, 0.1, 0.1, 0.1, 0.1});
62 inflate_factors_ = parameters.getParameter<std::vector<double>>(
63 "inflate_factors", {10., 10., 10., 10., 10., 10.});
64
65 // In tracking frame: where do these numbers come from?
66 // These numbers come from approximating the path of the beam up
67 // until it is about to enter the first detector volume (TriggerPad1).
68 // In detector coordinates, (x,y,z) = (-21.7, -883) is
69 // where the beam arrives (if no smearing is applied) and we simply
70 // reorder these values so that they are in tracking coordinates.
71 beamOrigin_ = parameters.getParameter<std::vector<double>>(
72 "beamOrigin", {-883.0, -21.745876, 0.0});
73
74 // Skip the tagger or recoil trackers if wanted
75 skip_tagger_ = parameters.getParameter<bool>("skip_tagger", false);
76 skip_recoil_ = parameters.getParameter<bool>("skip_recoil", false);
77}
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 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 79 of file TruthSeedProcessor.cxx.

81 {
82 std::vector<double> pos{static_cast<double>(hit.getPosition()[0]),
83 static_cast<double>(hit.getPosition()[1]),
84 static_cast<double>(hit.getPosition()[2])};
85 createTruthTrack(pos, hit.getMomentum(), particle.getCharge(), trk,
86 target_surface);
87
88 trk.setTrackID(hit.getTrackID());
89 trk.setPdgID(hit.getPdgID());
90}
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 see...

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

94 {
95 createTruthTrack(particle.getVertex(), particle.getMomentum(),
96 particle.getCharge(), trk, target_surface);
97
98 trk.setPdgID(particle.getPdgID());
99}
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
posThe 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 101 of file TruthSeedProcessor.cxx.

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

444 {
445 for (int i_sim_hit = 0; i_sim_hit < sim_hits.size(); i_sim_hit++) {
446 auto& sim_hit = sim_hits.at(i_sim_hit);
447 // This track never left a hit before
448 if (!hit_count_map.count(sim_hit.getTrackID())) {
449 hit_count_map[sim_hit.getTrackID()].push_back(i_sim_hit);
450 }
451
452 // This track left a hit before.
453 // Check if it's on a different sensor than the others
454
455 else {
456 int sensorID = tracking::sim::utils::getSensorID(sim_hit);
457 bool foundHit = false;
458
459 for (auto& i_rhit : hit_count_map[sim_hit.getTrackID()]) {
460 int tmp_sensorID =
461 tracking::sim::utils::getSensorID(sim_hits.at(i_rhit));
462
463 if (sensorID == tmp_sensorID) {
464 foundHit = true;
465 break;
466 }
467 } // loop on the already recorded hits
468
469 if (!foundHit) {
470 hit_count_map[sim_hit.getTrackID()].push_back(i_sim_hit);
471 }
472 }
473 } // loop on sim hits
474}

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 navCfg{geometry().getTG()};
17 const Acts::Navigator navigator(navCfg);
18
19 linpropagator_ = std::make_shared<LinPropagator>(stepper, navigator);
20 trk_extrap_ = std::make_shared<std::decay_t<decltype(*trk_extrap_)>>(
21 *linpropagator_, geometry_context(), magnetic_field_context());
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 530 of file TruthSeedProcessor.cxx.

530 {
531 // Retrieve the particleMap
532 auto particleMap{event.getMap<int, ldmx::SimParticle>("SimParticles")};
533
534 // Retrieve the target scoring hits
535 // Information is extracted using the
536 // scoring plane hit left by the particle at the target.
537
538 const std::vector<ldmx::SimTrackerHit> scoring_hits{
539 event.getCollection<ldmx::SimTrackerHit>(scoring_hits_coll_name_,
540 sp_pass_name_)};
541
542 // Retrieve the scoring plane hits at the ECAL
543 const std::vector<ldmx::SimTrackerHit> scoring_hits_ecal{
544 event.getCollection<ldmx::SimTrackerHit>("EcalScoringPlaneHits",
545 sp_pass_name_)};
546
547 // Retrieve the sim hits in the tagger tracker
548 const std::vector<ldmx::SimTrackerHit> tagger_sim_hits =
551
552 // Retrieve the sim hits in the recoil tracker
553 const std::vector<ldmx::SimTrackerHit> recoil_sim_hits =
556
557 // If sim hit collections are empty throw a warning
558 if (tagger_sim_hits.size() == 0 && !skip_tagger_) {
559 ldmx_log(error) << "Tagger sim hits collection empty for event ";
560 }
561 if (recoil_sim_hits.size() == 0 && !skip_recoil_) {
562 ldmx_log(error) << "Recoil sim hits collection empty for event ";
563 }
564
565 // The map stores which track leaves which sim hits
566 std::map<int, std::vector<int>> hit_count_map_recoil;
567 makeHitCountMap(recoil_sim_hits, hit_count_map_recoil);
568
569 std::map<int, std::vector<int>> hit_count_map_tagger;
570 makeHitCountMap(tagger_sim_hits, hit_count_map_tagger);
571
572 // to keep track of how many sim particles leave hits on the scoring plane
573 std::vector<int> recoil_sh_idxs;
574 std::unordered_map<int, std::vector<int>> recoil_sh_count_map;
575
576 std::vector<int> tagger_sh_idxs;
577 std::unordered_map<int, std::vector<int>> tagger_sh_count_map;
578
579 // Target scoring hits for Tagger will have Z<0, Recoil scoring hits will have
580 // Z>0
581 for (unsigned int i_sh = 0; i_sh < scoring_hits.size(); i_sh++) {
582 const ldmx::SimTrackerHit& hit = scoring_hits.at(i_sh);
583 double zhit = hit.getPosition()[2];
584
585 Acts::Vector3 p_vec{hit.getMomentum()[0], hit.getMomentum()[1],
586 hit.getMomentum()[2]};
587 double tagger_p_max = 0.;
588
589 // Check if it is a tagger track going fwd that passes basic cuts
590 if (zhit < 0.) {
591 // Tagger selection cuts
592 // Negative scoring plane hit, with momentum > p_cut
593 if (p_vec(2) < 0. || p_vec.norm() < p_cut_) continue;
594
595 // Check that the hit was left by a charged particle
596 if (abs(particleMap[hit.getTrackID()].getCharge()) < 1e-8) continue;
597
598 if (p_vec.norm() > tagger_p_max) {
599 tagger_sh_count_map[hit.getTrackID()].push_back(i_sh);
600 }
601 } // Tagger loop
602
603 // Check the recoil hits
604 else {
605 // Recoil selection cuts
606 // Positive scoring plane hit, forward direction with momentum > p_cut
607 if (p_vec(2) < 0. || p_vec.norm() < p_cut_) continue;
608
609 // Check that the hit was left by a charged particle
610 if (abs(particleMap[hit.getTrackID()].getCharge()) < 1e-8) continue;
611
612 recoil_sh_count_map[hit.getTrackID()].push_back(i_sh);
613
614 } // Recoil
615 } // loop on Target scoring plane hits
616
617 for (std::pair<int, std::vector<int>> element : recoil_sh_count_map) {
618 std::sort(
619 element.second.begin(), element.second.end(),
620 [&](const int idx1, int idx2) -> bool {
621 const ldmx::SimTrackerHit& hit1 = scoring_hits.at(idx1);
622 const ldmx::SimTrackerHit& hit2 = scoring_hits.at(idx2);
623
624 Acts::Vector3 phit1{hit1.getMomentum()[0], hit1.getMomentum()[1],
625 hit1.getMomentum()[2]};
626 Acts::Vector3 phit2{hit2.getMomentum()[0], hit2.getMomentum()[1],
627 hit2.getMomentum()[2]};
628
629 return phit1.norm() > phit2.norm();
630 });
631 }
632
633 // Sort tagger hits.
634 for (auto& [_track_id, hit_indices] : tagger_sh_count_map) {
635 std::sort(
636 hit_indices.begin(), hit_indices.end(),
637 [&](const int idx1, int idx2) -> bool {
638 const ldmx::SimTrackerHit& hit1 = scoring_hits.at(idx1);
639 const ldmx::SimTrackerHit& hit2 = scoring_hits.at(idx2);
640
641 Acts::Vector3 phit1{hit1.getMomentum()[0], hit1.getMomentum()[1],
642 hit1.getMomentum()[2]};
643 Acts::Vector3 phit2{hit2.getMomentum()[0], hit2.getMomentum()[1],
644 hit2.getMomentum()[2]};
645
646 return phit1.norm() > phit2.norm();
647 });
648 }
649
650 // Building of the event truth information and the truth seeds
651 // TODO remove the truthtracks in the future as the truth seeds are enough
652
653 std::vector<ldmx::Track> tagger_truth_tracks;
654 std::vector<ldmx::Track> tagger_truth_seeds;
655 std::vector<ldmx::Track> recoil_truth_tracks;
656 std::vector<ldmx::Track> recoil_truth_seeds;
657 ldmx::Tracks beam_electrons;
658
659 // TODO:: The target should be taken from some conditions DB in the future.
660 // Define the perigee_surface at 0.0.0
661 auto targetSurface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
662 Acts::Vector3(0., 0., 0.))};
663
664 // Define the target_surface
665 auto targetUnboundSurface = tracking::sim::utils::unboundSurface(0.);
666
667 // ecal
668 auto ecalSurface = tracking::sim::utils::unboundSurface(240.5);
669
670 auto beamOriginSurface{Acts::Surface::makeShared<Acts::PerigeeSurface>(
671 Acts::Vector3(beamOrigin_[0], beamOrigin_[1], beamOrigin_[2]))};
672
673 if (!skip_tagger_) {
674 for (const auto& [track_id, hit_indices] : tagger_sh_count_map) {
675 const ldmx::SimTrackerHit& hit = scoring_hits.at(hit_indices.at(0));
676 const ldmx::SimParticle& phit = particleMap[hit.getTrackID()];
677
678 if (hit_count_map_tagger[hit.getTrackID()].size() > n_min_hits_tagger_) {
679 ldmx::Track truth_tagger_track;
680 createTruthTrack(phit, hit, truth_tagger_track, targetSurface);
681 truth_tagger_track.setNhits(
682 hit_count_map_tagger[hit.getTrackID()].size());
683 tagger_truth_tracks.push_back(truth_tagger_track);
684
685 if (hit.getPdgID() == 11 && hit.getTrackID() < max_track_id_) {
686 ldmx::Track beamETruthSeed = TaggerFullSeed(
687 particleMap[hit.getTrackID()], hit.getTrackID(), hit,
688 hit_count_map_tagger, beamOriginSurface, targetUnboundSurface);
689 beam_electrons.push_back(beamETruthSeed);
690 }
691 }
692 }
693 }
694
695 // Recover the EcalScoring hits
696 std::vector<ldmx::SimTrackerHit> ecal_spHits =
697 event.getCollection<ldmx::SimTrackerHit>("EcalScoringPlaneHits",
698 sp_pass_name_);
699 // Select ECAL hits
700 std::vector<ldmx::SimTrackerHit> sel_ecal_spHits;
701
702 for (auto sp_hit : ecal_spHits) {
703 if (sp_hit.getMomentum()[2] > 0 && ((sp_hit.getID() & 0xfff) == 31)) {
704 sel_ecal_spHits.push_back(sp_hit);
705 }
706 }
707
708 // Recoil target surface for truth and seed tracks is the target
709
710 for (std::pair<int, std::vector<int>> element : recoil_sh_count_map) {
711 // Only take the first entry of the vector: it should be the scoring plane
712 // hit with the highest momentum.
713 const ldmx::SimTrackerHit& hit = scoring_hits.at(element.second.at(0));
714 [[maybe_unused]] const ldmx::SimParticle& phit =
715 particleMap[hit.getTrackID()];
716 ldmx::SimTrackerHit ecal_hit;
717
718 bool foundEcalHit = false;
719 for (auto ecal_sp_hit : sel_ecal_spHits) {
720 if (ecal_sp_hit.getTrackID() == hit.getTrackID()) {
721 ecal_hit = ecal_sp_hit;
722 foundEcalHit = true;
723 break;
724 }
725 }
726
727 // Findable particle selection
728 // std::cout << "!!! n_recoil_sim_hits found: " <<
729 // hit_count_map_recoil[hit.getTrackID()].size() << std::endl;
730 if (hit_count_map_recoil[hit.getTrackID()].size() > n_min_hits_recoil_ &&
731 foundEcalHit && !skip_recoil_) {
732 ldmx::Track truth_recoil_track =
733 RecoilFullSeed(particleMap[hit.getTrackID()], hit.getTrackID(), hit,
734 ecal_hit, hit_count_map_recoil, targetSurface,
735 targetUnboundSurface, ecalSurface);
736 // std::cout << "!!! Recoil track created" << std::endl;
737 // std::cout << "!!! Recoil track momentum: " <<
738 // truth_recoil_track.getMomentum()[0] << " " <<
739 // truth_recoil_track.getMomentum()[1] << " " <<
740 // truth_recoil_track.getMomentum()[2] << std::endl; std::cout << "!!! Hit
741 // momentum: " << hit.getMomentum()[0] << " " << hit.getMomentum()[1] << "
742 // " << hit.getMomentum()[2] << std::endl;
743 recoil_truth_tracks.push_back(truth_recoil_track);
744 }
745 }
746
747 /*
748 for (std::pair<int,std::vector<int>> element : recoil_sh_count_map) {
749
750 const ldmx::SimTrackerHit& hit = scoring_hits.at(element.second.at(0));
751 const ldmx::SimParticle& phit = particleMap[hit.getTrackID()];
752
753 if (hit_count_map_recoil[hit.getTrackID()].size() > n_min_hits_recoil_) {
754 ldmx::Track truth_recoil_track;
755 createTruthTrack(phit,hit,truth_recoil_track,targetSurface);
756 truth_recoil_track.setNhits(hit_count_map_recoil[hit.getTrackID()].size());
757 recoil_truth_tracks.push_back(truth_recoil_track);
758 }
759 }
760 */
761
762 // Form a truth seed from a truth track
763
764 for (auto& tt : tagger_truth_tracks) {
765 ldmx::Track seed = seedFromTruth(tt, seedSmearing_);
766
767 tagger_truth_seeds.push_back(seed);
768 }
769
770 ldmx_log(debug) << "Forming seeds from truth" << std::endl;
771 for (auto& tt : recoil_truth_tracks) {
772 ldmx_log(debug) << "Smearing truth track" << std::endl;
773
774 ldmx::Track seed = seedFromTruth(tt, seedSmearing_);
775
776 recoil_truth_seeds.push_back(seed);
777 }
778
779 // even if skip_tagger/recoil_ is true, still make the collections in the
780 // event
781 event.add("beamElectrons", beam_electrons);
782 event.add("TaggerTruthTracks", tagger_truth_tracks);
783 event.add("RecoilTruthTracks", recoil_truth_tracks);
784 event.add("TaggerTruthSeeds", tagger_truth_seeds);
785 event.add("RecoilTruthSeeds", recoil_truth_seeds);
786}
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
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...
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 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 187 of file TruthSeedProcessor.cxx.

193 {
194 ldmx::Track truth_recoil_track;
195 createTruthTrack(particle, hit, truth_recoil_track, origin_surface);
196 truth_recoil_track.setTrackID(trackID);
197
198 // Seed at the target location
199 ldmx::Track smearedTruthTrack = seedFromTruth(truth_recoil_track, false);
200
201 // Add the track state at the target
202 ldmx::Track truth_track_target;
203 createTruthTrack(particle, hit, truth_track_target, target_surface);
204
205 // Store the truth track state on the seed track
206 ldmx::Track::TrackState ts_truth_target;
207 Acts::Vector3 ref = target_surface->center(geometry_context());
208 ts_truth_target.refX = ref(0);
209 ts_truth_target.refY = ref(1);
210 ts_truth_target.refZ = ref(2);
211 ts_truth_target.params = truth_track_target.getPerigeeParameters();
212 // empty cov
213 ts_truth_target.ts_type = ldmx::TrackStateType::AtTarget;
214 smearedTruthTrack.addTrackState(ts_truth_target);
215
216 // Add the track state at the ecal
217 ldmx::Track truth_track_ecal;
218 createTruthTrack(particle, ecal_hit, truth_track_ecal, ecal_surface);
219
220 ldmx::Track::TrackState ts_truth_ecal;
221 Acts::Vector3 ref_ecal = ecal_surface->center(geometry_context());
222 ts_truth_ecal.refX = ref_ecal(0);
223 ts_truth_ecal.refY = ref_ecal(1);
224 ts_truth_ecal.refZ = ref_ecal(2);
225 ts_truth_ecal.params = truth_track_ecal.getPerigeeParameters();
226 // empty cov
227 ts_truth_ecal.ts_type = ldmx::TrackStateType::AtECAL;
228 smearedTruthTrack.addTrackState(ts_truth_ecal);
229
230 // Add the hits
231 int nhits = 0;
232
233 for (auto sim_hit_idx : hit_count_map.at(smearedTruthTrack.getTrackID())) {
234 smearedTruthTrack.addMeasurementIndex(sim_hit_idx);
235 nhits += 1;
236 }
237
238 smearedTruthTrack.setNhits(nhits);
239
240 return smearedTruthTrack;
241}

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

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

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

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

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

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

References createTruthTrack(), and seedFromTruth().

Member Data Documentation

◆ beamOrigin_

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

Definition at line 278 of file TruthSeedProcessor.h.

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

◆ d0smear_

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

Definition at line 271 of file TruthSeedProcessor.h.

◆ 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 266 of file TruthSeedProcessor.h.

◆ inflate_factors_

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

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

208{""};

Referenced by configure(), and produce().

◆ linpropagator_

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

Definition at line 258 of file TruthSeedProcessor.h.

◆ max_track_id_

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

Definition at line 256 of file TruthSeedProcessor.h.

256{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 220 of file TruthSeedProcessor.h.

220{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 214 of file TruthSeedProcessor.h.

214{7};

Referenced by configure().

◆ normal_

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

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

235{0.};

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

◆ p_cut_ecal_

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

Definition at line 241 of file TruthSeedProcessor.h.

241{-1.};

◆ p_cut_max_

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

Ask for a maximum p for the seeds.

Definition at line 238 of file TruthSeedProcessor.h.

238{100000.};

Referenced by configure(), and scoringPlaneHitFilter().

◆ 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 273 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 232 of file TruthSeedProcessor.h.

232{-9999};

Referenced by configure(), and scoringPlaneHitFilter().

◆ 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 205 of file TruthSeedProcessor.h.

205{"RecoilSimHits"};

Referenced by configure(), and produce().

◆ recoil_sp_

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

Definition at line 244 of file TruthSeedProcessor.h.

244{true};

◆ rel_smearfactors_

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

Definition at line 276 of file TruthSeedProcessor.h.

◆ relpsmear_

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

Definition at line 275 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().

◆ seedSmearing_

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

Definition at line 269 of file TruthSeedProcessor.h.

269{false};

◆ skip_recoil_

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

Definition at line 253 of file TruthSeedProcessor.h.

253{false};

◆ skip_tagger_

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

Definition at line 250 of file TruthSeedProcessor.h.

250{false};

◆ sp_pass_name_

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

Definition at line 199 of file TruthSeedProcessor.h.

199{""};

◆ 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 202 of file TruthSeedProcessor.h.

202{"TaggerSimHits"};

Referenced by configure(), and produce().

◆ target_sp_

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

Definition at line 247 of file TruthSeedProcessor.h.

247{true};

◆ thetasmear_

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

Definition at line 274 of file TruthSeedProcessor.h.

◆ track_id_

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

Only select a particular trackID.

Definition at line 229 of file TruthSeedProcessor.h.

229{-999};

Referenced by configure(), and scoringPlaneHitFilter().

◆ trk_extrap_

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

Definition at line 262 of file TruthSeedProcessor.h.

◆ z0smear_

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

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

226{-999};

Referenced by configure(), and scoringPlaneHitFilter().


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