Track3DKalmanHit_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Class: Track3DKalmanHit
4 // Module Type: producer
5 // File: Track3DKalmanHit_module.cc
6 //
7 // This module produces RecoBase/Track objects using KalmanFilterAlgorithm.
8 //
9 // Configuration parameters:
10 //
11 // Hist - Histogram flag (generate histograms if true).
12 // UseClusterHits - Use clustered hits if true.
13 // UsePFParticleHits - Use PFParticle hits if true.
14 // UsePFParticleSeeds - If true, use seeds associated with PFParticles.
15 // HitModuleLabel - Module label for unclustered Hits.
16 // ClusterModuleLabel - Module label for Clusters.
17 // PFParticleModuleLabel - Module label for PFParticles.
18 // Track3DKalmanHitAlg - Algorithm class for Track3DKalmanHit module
19 // SpacePointAlg - Parmaeter set for space points.
20 //
21 // Usage Notes.
22 //
23 // 1. It is an error if UseClusterHits and UsePFParticleHits are both
24 // true (throw exception in that case). However, both can be false,
25 // in which case all hits are used as input.
26 //
27 // 2. If clustered hits are used as input, then tracks can span clusters
28 // (cluster structure of hits is ignored).
29 //
30 // 3. If PFParticle hits are used as input, then tracks do not span
31 // PFParticles. However, it is possible for one the hits from
32 // one PFParticle to give rise to multiple Tracks.
33 //
34 // 4. UsePFParticleSeeds has no effect unless UsePFParticleHits is true.
35 //
36 ////////////////////////////////////////////////////////////////////////
37 
38 #include <cmath>
39 #include <deque>
40 
46 #include "art_root_io/TFileDirectory.h"
47 #include "art_root_io/TFileService.h"
49 #include "canvas/Persistency/Common/FindManyP.h"
52 #include "fhiclcpp/ParameterSet.h"
54 
68 
69 #include "TH1F.h"
70 
71 namespace trkf {
72 
74  public:
75  // Copnstructors, destructor.
76  explicit Track3DKalmanHit(fhicl::ParameterSet const& pset);
77 
78  private:
79  void produce(art::Event& e) override;
80  void beginJob() override;
81  void endJob() override;
82 
83  //Member functions that depend on art::event and use art::Assns
84  //note that return types are move aware
85  void prepareForInput();
86  KalmanInputs getInput(const art::Event& evt) const;
87  Hits getClusteredHits(const art::Event& evt) const;
88  KalmanInputs getPFParticleStuff(const art::Event& evt) const;
89  Hits getAllHits(const art::Event& evt) const;
90  void createOutputs(const art::Event& evt,
91  detinfo::DetectorPropertiesData const& detProp,
92  const std::vector<KalmanOutput>& outputs,
93  const KalmanInputs& inputs,
94  std::vector<recob::Track>& tracks,
95  std::vector<recob::SpacePoint>& spts,
101  void fillHistograms(std::vector<KalmanOutput>& outputs);
102 
103  // Fcl parameters.
104  bool fHist; ///< Make histograms.
105  bool fUseClusterHits; ///< Use clustered hits as input.
106  bool fUsePFParticleHits; ///< Use PFParticle hits as input.
107  bool fUsePFParticleSeeds; ///< Use PFParticle seeds.
108  std::string fHitModuleLabel; ///< Unclustered Hits.
109  std::string fClusterModuleLabel; ///< Clustered Hits.
110  std::string fPFParticleModuleLabel; ///< PFParticle label.
111  Track3DKalmanHitAlg fTKHAlg; ///< Track3DKalmanHit algorithm.
112  SpacePointAlg fSpacePointAlg; ///< Space point algorithm.
113 
114  // Histograms.
115  TH1F* fHIncChisq; ///< Incremental chisquare.
116  TH1F* fHPull; ///< Hit pull.
117 
118  // Statistics.
119  int fNumEvent; ///< Number of events seen.
120  };
122 } // namespace trkf
123 
124 //----------------------------------------------------------------------------
125 /// Constructor.
126 ///
127 /// Arguments:
128 ///
129 /// p - Fcl parameters.
130 ///
132  : EDProducer{pset}
133  , fHist(false)
134  , fUseClusterHits(false)
135  , fUsePFParticleHits(false)
136  , fUsePFParticleSeeds(false)
137  , fTKHAlg(pset.get<fhicl::ParameterSet>("Track3DKalmanHitAlg"))
138  , fSpacePointAlg(pset.get<fhicl::ParameterSet>("SpacePointAlg"))
139  , fHIncChisq(0)
140  , fHPull(0)
141  , fNumEvent(0)
142 {
143  fHist = pset.get<bool>("Hist");
144  fUseClusterHits = pset.get<bool>("UseClusterHits");
145  fUsePFParticleHits = pset.get<bool>("UsePFParticleHits");
146  fUsePFParticleSeeds = pset.get<bool>("UsePFParticleSeeds");
147  fHitModuleLabel = pset.get<std::string>("HitModuleLabel");
148  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
149  fPFParticleModuleLabel = pset.get<std::string>("PFParticleModuleLabel");
150  if (fUseClusterHits && fUsePFParticleHits) {
151  throw cet::exception("Track3DKalmanHit")
152  << "Using input from both clustered and PFParticle hits.\n";
153  }
154 
155  produces<std::vector<recob::Track>>();
156  produces<std::vector<recob::SpacePoint>>();
158  recob::Hit>>(); // ****** REMEMBER to remove when FindMany improved ******
159  produces<art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta>>();
160  produces<art::Assns<recob::Track, recob::SpacePoint>>();
161  produces<art::Assns<recob::SpacePoint, recob::Hit>>();
162  produces<art::Assns<recob::PFParticle, recob::Track>>();
163 
164  // Report.
165 
166  mf::LogInfo("Track3DKalmanHit") << "Track3DKalmanHit configured with the following parameters:\n"
167  << " UseClusterHits = " << fUseClusterHits << "\n"
168  << " HitModuleLabel = " << fHitModuleLabel << "\n"
169  << " ClusterModuleLabel = " << fClusterModuleLabel;
170 }
171 
172 //----------------------------------------------------------------------------
173 /// Begin job method.
174 void
176 {
177  if (fHist) {
178  // Book histograms.
180  art::TFileDirectory dir = tfs->mkdir("hitkalman", "Track3DKalmanHit histograms");
181 
182  fHIncChisq = dir.make<TH1F>("IncChisq", "Incremental Chisquare", 100, 0., 20.);
183  fHPull = dir.make<TH1F>("Pull", "Hit Pull", 100, -10., 10.);
184  }
185 }
186 
187 //----------------------------------------------------------------------------
188 /// Produce method.
189 ///
190 void
192 {
193  ++fNumEvent;
194  auto tracks = std::make_unique<std::vector<recob::Track>>();
195  auto th_assn = std::make_unique<art::Assns<recob::Track, recob::Hit>>();
196  auto thm_assn = std::make_unique<art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta>>();
197  auto tsp_assn = std::make_unique<art::Assns<recob::Track, recob::SpacePoint>>();
198  auto pfPartTrack_assns = std::make_unique<art::Assns<recob::PFParticle, recob::Track>>();
199  auto spts = std::make_unique<std::vector<recob::SpacePoint>>();
200  auto sph_assn = std::make_unique<art::Assns<recob::SpacePoint, recob::Hit>>();
201 
202  prepareForInput();
203  auto inputs = getInput(evt);
204 
205  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
206  auto const detProp =
208  auto outputs = fTKHAlg.makeTracks(clockData, detProp, inputs);
209 
210  if (fHist) { fillHistograms(outputs); }
211 
212  createOutputs(evt,
213  detProp,
214  outputs,
215  inputs,
216  *tracks,
217  *spts,
218  *th_assn,
219  *thm_assn,
220  *tsp_assn,
221  *sph_assn,
222  *pfPartTrack_assns);
223 
225 
226  evt.put(std::move(tracks));
227  evt.put(std::move(spts));
228  evt.put(std::move(th_assn));
229  evt.put(std::move(thm_assn));
230  evt.put(std::move(tsp_assn));
231  evt.put(std::move(sph_assn));
232  evt.put(std::move(pfPartTrack_assns));
233 }
234 
235 //----------------------------------------------------------------------------
236 /// End job method.
237 void
239 {
240  mf::LogInfo("Track3DKalmanHit") << "Track3DKalmanHit statistics:\n"
241  << " Number of events = " << fNumEvent << "\n";
242 }
243 
244 //----------------------------------------------------------------------------
245 //A temporary method till we find an alternative
246 void
248 {
250 }
251 
252 //----------------------------------------------------------------------------
253 // There are three modes of operation:
254 // 1. Clustered hits (produces one hit collection).
255 // 2. PFParticle hits (products one hit collection for each PFParticle).
256 // 3. All hits (produces one hit collection).
257 
260 {
261  if (fUsePFParticleHits) return getPFParticleStuff(evt);
263 }
264 
265 //----------------------------------------------------------------------------
266 /// Fill a collection using clustered hits
269 {
270  Hits hits;
272  evt.getByLabel(fClusterModuleLabel, clusterh);
273 
274  if (!clusterh.isValid()) return hits;
275  // Get hits from all clusters.
276  art::FindManyP<recob::Hit> hitsbycluster(clusterh, evt, fClusterModuleLabel);
277 
278  for (size_t i = 0; i < clusterh->size(); ++i) {
279  std::vector<art::Ptr<recob::Hit>> clushits = hitsbycluster.at(i);
280  hits.insert(hits.end(), clushits.begin(), clushits.end());
281  }
282  return hits;
283 }
284 
285 //----------------------------------------------------------------------------
286 /// If both UseClusteredHits and UsePFParticles is false use this method to fill in hits
289 {
290  Hits hits;
292  evt.getByLabel(fHitModuleLabel, hith);
293  if (!hith.isValid()) return hits;
294  size_t nhits = hith->size();
295  hits.reserve(nhits);
296 
297  for (size_t i = 0; i < nhits; ++i) {
298  hits.push_back(art::Ptr<recob::Hit>(hith, i));
299  }
300  return hits;
301 }
302 
303 //----------------------------------------------------------------------------
304 /// If UsePFParticles is true use this method to fill in hits
307 {
309  // Our program is to drive the track creation/fitting off the PFParticles in the data store
310  // We'll use the hits associated to the PFParticles for each track - and only those hits.
311  // Without a valid collection of PFParticles there is nothing to do here
312  // We need a handle to the collection of clusters in the data store so we can
313  // handle associations to hits.
315  evt.getByLabel(fPFParticleModuleLabel, pfParticleHandle);
316  if (!pfParticleHandle.isValid()) return inputs;
317 
319  evt.getByLabel(fClusterModuleLabel, clusterHandle);
320 
321  // If there are no clusters then something is really wrong
322  if (!clusterHandle.isValid()) return inputs;
323 
324  // Recover the collection of associations between PFParticles and clusters, this will
325  // be the mechanism by which we actually deal with clusters
326  art::FindManyP<recob::Cluster> clusterAssns(pfParticleHandle, evt, fPFParticleModuleLabel);
327 
328  // Associations to seeds.
329  art::FindManyP<recob::Seed> seedAssns(pfParticleHandle, evt, fPFParticleModuleLabel);
330 
331  // Likewise, recover the collection of associations to hits
332  art::FindManyP<recob::Hit> clusterHitAssns(clusterHandle, evt, fClusterModuleLabel);
333 
334  inputs.reserve(pfParticleHandle->size());
335 
336  // While PFParticle describes a hierarchal structure, for now we simply loop over the collection
337  for (size_t partIdx = 0; partIdx < pfParticleHandle->size(); partIdx++) {
338 
339  // Add a new empty hit collection.
340  inputs.emplace_back();
341  KalmanInput& kalman_input = inputs.back();
342  kalman_input.pfPartPtr = art::Ptr<recob::PFParticle>(pfParticleHandle, partIdx);
343  Hits& hits = kalman_input.hits;
344 
345  // Fill this hit vector by looping over associated clusters and finding the
346  // hits associated to them
347  std::vector<art::Ptr<recob::Cluster>> clusterVec = clusterAssns.at(partIdx);
348 
349  for (auto const& cluster : clusterVec) {
350  std::vector<art::Ptr<recob::Hit>> hitVec = clusterHitAssns.at(cluster.key());
351  hits.insert(hits.end(), hitVec.begin(), hitVec.end());
352  }
353 
354  // If requested, fill associated seeds.
355  if (!fUsePFParticleSeeds) continue;
356  art::PtrVector<recob::Seed>& seeds = kalman_input.seeds;
357  std::vector<art::Ptr<recob::Seed>> seedVec = seedAssns.at(partIdx);
358  seeds.insert(seeds.end(), seedVec.begin(), seedVec.end());
359  art::FindManyP<recob::Hit> seedHitAssns(seedVec, evt, fPFParticleModuleLabel);
360 
361  for (size_t seedIdx = 0; seedIdx < seedVec.size(); ++seedIdx) {
362  std::vector<art::Ptr<recob::Hit>> seedHitVec;
363  //SS: why seedIdx can have an invalid value?
364  try {
365  seedHitVec = seedHitAssns.at(seedIdx);
366  }
367  catch (art::Exception const&) {
368  seedHitVec.clear();
369  }
370  kalman_input.seedhits.emplace_back();
371  Hits& seedhits = kalman_input.seedhits.back();
372  seedhits.insert(seedhits.end(), seedHitVec.begin(), seedHitVec.end());
373  }
374  }
375  return inputs;
376 }
377 
378 //-------------------------------------------------------------------------------------
379 void
381  const art::Event& evt,
382  detinfo::DetectorPropertiesData const& detProp,
383  std::vector<KalmanOutput> const& outputs,
384  KalmanInputs const& inputs,
385  std::vector<recob::Track>& tracks,
386  std::vector<recob::SpacePoint>& spts,
392 {
393  if (outputs.size() != inputs.size()) return;
394 
395  size_t tracksSize(0);
396  for (auto const& kalman_output : outputs) {
397  tracksSize += kalman_output.tracks.size();
398  }
399  tracks.reserve(tracksSize);
400 
401  auto const tid = evt.getProductID<std::vector<recob::Track>>();
402  auto const tidgetter = evt.productGetter(tid);
403 
404  auto const spacepointId = evt.getProductID<std::vector<recob::SpacePoint>>();
405  auto const getter = evt.productGetter(spacepointId);
406 
407  for (size_t i = 0; i < outputs.size(); ++i) {
408  // Recover the kalman tracks double ended queue
409  const std::deque<KGTrack>& kalman_tracks = outputs[i].tracks;
410 
411  for (auto const& kalman_track : kalman_tracks) {
412 
413  // Add Track object to collection.
414  recob::Track track;
415  kalman_track.fillTrack(detProp, track, tracks.size());
416  if (track.NumberTrajectoryPoints() < 2) { continue; }
417  unsigned int numtrajpts = track.NumberTrajectoryPoints();
418  tracks.emplace_back(std::move(track));
419  // SS: tracks->size() does not change after this point in each iteration
420 
421  //fill hits from this track
422  Hits trhits;
423  std::vector<unsigned int> hittpindex; //hit-trajectory point index
424  kalman_track.fillHits(trhits, hittpindex);
425  if (hittpindex.back() >= numtrajpts) { //something is wrong
426  throw cet::exception("Track3DKalmanHit")
427  << "Last hit corresponds to trajectory point index " << hittpindex.back()
428  << " while the number of trajectory points is " << numtrajpts << '\n';
429  }
430 
431  // Make space points from this track.
432  auto nspt = spts.size();
433  fSpacePointAlg.fillSpacePoints(detProp, spts, kalman_track.TrackMap());
434 
435  std::vector<art::Ptr<recob::SpacePoint>> sptvec;
436  for (auto ispt = nspt; ispt < spts.size(); ++ispt) {
437  sptvec.emplace_back(spacepointId, ispt, getter);
438  // Associate newly created space points with hits.
439  // Make space point to hit associations.
440  auto const& sphits = fSpacePointAlg.getAssociatedHits((spts)[ispt]);
441  for (auto const& sphit : sphits) {
442  sph_assn.addSingle(sptvec.back(), sphit);
443  }
444  }
445 
446  art::Ptr<recob::Track> aptr(tid, tracks.size() - 1, tidgetter);
447 
448  // Make Track to Hit associations.
449  for (size_t h = 0; h < trhits.size(); ++h) {
450  th_assn.addSingle(aptr, trhits[h]);
451  recob::TrackHitMeta metadata(hittpindex[h], -1);
452  thm_assn.addSingle(aptr, trhits[h], metadata);
453  }
454 
455  // Make track to space point associations
456  for (auto const& spt : sptvec) {
457  tsp_assn.addSingle(aptr, spt);
458  }
459 
460  // Optionally fill track-to-PFParticle associations.
461  if (fUsePFParticleHits) { pfPartTrack_assns.addSingle(inputs[i].pfPartPtr, aptr); }
462  } // end of loop over a given collection
463  }
464 }
465 
466 //----------------------------------------------------------------------------
467 /// Fill Histograms method
468 //fHPull and fHIncChisq are private data members of the class Track3DKalmanHit
469 
470 void
472 {
473  for (auto const& output : outputs) {
474  const std::deque<KGTrack>& kalman_tracks = output.tracks;
475  for (size_t i = 0; i < kalman_tracks.size(); ++i) {
476  const KGTrack& trg = kalman_tracks[i];
477  // Loop over measurements in this track.
478  const std::multimap<double, KHitTrack>& trackmap = trg.getTrackMap();
480  ih != trackmap.end();
481  ++ih) {
482  const KHitTrack& trh = (*ih).second;
483  const std::shared_ptr<const KHitBase>& hit = trh.getHit();
484  double chisq = hit->getChisq();
485  fHIncChisq->Fill(chisq);
486  const KHit<1>* ph1 = dynamic_cast<const KHit<1>*>(&*hit);
487  if (ph1 != 0) {
488  double pull = ph1->getResVector()(0) / std::sqrt(ph1->getResError()(0, 0));
489  fHPull->Fill(pull);
490  }
491  }
492  }
493  }
494 }
art::Ptr< recob::PFParticle > pfPartPtr
void reserve(size_type n)
Definition: PtrVector.h:337
bool fUseClusterHits
Use clustered hits as input.
const std::shared_ptr< const KHitBase > & getHit() const
Measurement.
Definition: KHitTrack.h:53
std::string fHitModuleLabel
Unclustered Hits.
Track3DKalmanHitAlg fTKHAlg
Track3DKalmanHit algorithm.
bool fUsePFParticleSeeds
Use PFParticle seeds.
ProductID getProductID(std::string const &instance_name="") const
Definition: DataViewImpl.h:338
KalmanInputs getInput(const art::Event &evt) const
std::string string
Definition: nybbler.cc:12
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
void beginJob() override
Begin job method.
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
std::string fPFParticleModuleLabel
PFParticle label.
void produces(std::string const &instanceName={}, Persistable const persistable=Persistable::Yes)
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:102
Class to keep data related to recob::Hit associated with recob::Track.
std::vector< trkf::KalmanOutput > makeTracks(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, KalmanInputs &kalman_inputs)
std::string fClusterModuleLabel
Clustered Hits.
intermediate_table::const_iterator const_iterator
Data related to recob::Hit associated with recob::Track.The purpose is to collect several variables t...
Definition: TrackHitMeta.h:43
string dir
Cluster finding and building.
Track3DKalmanHit(fhicl::ParameterSet const &pset)
const KVector< N >::type & getResVector() const
Residual vector.
Definition: KHit.h:138
Track3DKalmanHit Algorithm.
const art::PtrVector< recob::Hit > & getAssociatedHits(const recob::SpacePoint &spt) const
SpacePointAlg fSpacePointAlg
Space point algorithm.
const std::multimap< double, KHitTrack > & getTrackMap() const
KHitTrack collection, indexed by path distance.
Definition: KGTrack.h:59
bool isValid() const noexcept
Definition: Handle.h:191
const KSymMatrix< N >::type & getResError() const
Residual error matrix.
Definition: KHit.h:145
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
Hits getClusteredHits(const art::Event &evt) const
Fill a collection using clustered hits.
art::PtrVector< recob::Hit > hits
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
Kalman filter measurement class template.
def move(depos, offset)
Definition: depos.py:107
EDProductGetter const * productGetter(ProductID const pid) const
iterator end()
Definition: PtrVector.h:231
Hits getAllHits(const art::Event &evt) const
If both UseClusteredHits and UsePFParticles is false use this method to fill in hits.
std::vector< art::PtrVector< recob::Hit > > seedhits
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
reference at(size_type n)
Definition: PtrVector.h:359
void fillHistograms(std::vector< KalmanOutput > &outputs)
Fill Histograms method.
size_type size() const
Definition: PtrVector.h:302
void fillSpacePoints(detinfo::DetectorPropertiesData const &detProp, std::vector< recob::SpacePoint > &spts, std::multimap< double, KHitTrack > const &trackMap) const
Fill a collection of space points.
Detector simulation of raw signals on wires.
int fNumEvent
Number of events seen.
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
Definition: tracks.py:1
std::vector< TrajPoint > seeds
Definition: DataStructs.cxx:13
void clearHitMap() const
iterator insert(iterator position, Ptr< U > const &p)
Declaration of signal hit object.
bool fUsePFParticleHits
Use PFParticle hits as input.
void produce(art::Event &e) override
void addSingle(Ptr< left_t > const &left, Ptr< right_t > const &right, data_t const &data)
Definition: Assns.h:546
Provides recob::Track data product.
TrackCollectionProxyElement< TrackCollProxy > Track
Proxy to an element of a proxy collection of recob::Track objects.
Definition: Track.h:1036
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
TCEvent evt
Definition: DataStructs.cxx:7
Algorithm for generating space points from hits.
void createOutputs(const art::Event &evt, detinfo::DetectorPropertiesData const &detProp, const std::vector< KalmanOutput > &outputs, const KalmanInputs &inputs, std::vector< recob::Track > &tracks, std::vector< recob::SpacePoint > &spts, art::Assns< recob::Track, recob::Hit > &th_assn, art::Assns< recob::Track, recob::Hit, recob::TrackHitMeta > &thm_assn, art::Assns< recob::Track, recob::SpacePoint > &tsp_assn, art::Assns< recob::SpacePoint, recob::Hit > &sph_assn, art::Assns< recob::PFParticle, recob::Track > &pfPartTrack_assns)
std::vector< KalmanInput > KalmanInputs
void endJob() override
End job method.
TH1F * fHIncChisq
Incremental chisquare.
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:49
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
art::PtrVector< recob::Seed > seeds
KalmanInputs getPFParticleStuff(const art::Event &evt) const
If UsePFParticles is true use this method to fill in hits.