LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
trkf::Track3DKalmanHit Class Reference
Inheritance diagram for trkf::Track3DKalmanHit:
art::EDProducer art::ProducerBase art::Consumer art::EngineCreator art::ProductRegistryHelper

Public Types

using ModuleType = EDProducer
 
using WorkerType = WorkerT< EDProducer >
 
template<typename UserConfig , typename KeysToIgnore = void>
using Table = ProducerBase::Table< UserConfig, KeysToIgnore >
 

Public Member Functions

 Track3DKalmanHit (fhicl::ParameterSet const &pset)
 
virtual void reconfigure (fhicl::ParameterSet const &pset)
 
virtual void produce (art::Event &e) override
 
virtual void beginJob () override
 Begin job method. More...
 
virtual void endJob () override
 End job method. More...
 
template<typename PROD , BranchType B = InEvent>
ProductID getProductID (std::string const &instanceName={}) const
 
template<typename PROD , BranchType B>
ProductID getProductID (ModuleDescription const &moduleDescription, std::string const &instanceName) const
 
bool modifiesEvent () const
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ProductToken< T > consumes (InputTag const &it)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ViewToken< T > consumesView (InputTag const &it)
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ProductToken< T > mayConsume (InputTag const &it)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , art::BranchType BT>
art::ViewToken< T > mayConsumeView (InputTag const &it)
 
base_engine_tcreateEngine (seed_t seed)
 
base_engine_tcreateEngine (seed_t seed, std::string const &kind_of_engine_to_make)
 
base_engine_tcreateEngine (seed_t seed, std::string const &kind_of_engine_to_make, label_t const &engine_label)
 
seed_t get_seed_value (fhicl::ParameterSet const &pset, char const key[]="seed", seed_t const implicit_seed=-1)
 

Static Public Member Functions

static cet::exempt_ptr< Consumernon_module_context ()
 

Protected Member Functions

CurrentProcessingContext const * currentContext () const
 
void validateConsumedProduct (BranchType const bt, ProductInfo const &pi)
 
void prepareForJob (fhicl::ParameterSet const &pset)
 
void showMissingConsumes () const
 

Private Member Functions

void prepareForInput ()
 
KalmanInputs getInput (const art::Event &evt) const
 
Hits getClusteredHits (const art::Event &evt) const
 Fill a collection using clustered hits. More...
 
KalmanInputs getPFParticleStuff (const art::Event &evt) const
 If UsePFParticles is true use this method to fill in hits. More...
 
Hits getAllHits (const art::Event &evt) const
 If both UseClusteredHits and UsePFParticles is false use this method to fill in hits. More...
 
void createOutputs (const art::Event &evt, 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)
 
void fillHistograms (std::vector< KalmanOutput > &outputs)
 Fill Histograms method. More...
 

Private Attributes

bool fHist
 Make histograms. More...
 
bool fUseClusterHits
 Use clustered hits as input. More...
 
bool fUsePFParticleHits
 Use PFParticle hits as input. More...
 
bool fUsePFParticleSeeds
 Use PFParticle seeds. More...
 
std::string fHitModuleLabel
 Unclustered Hits. More...
 
std::string fClusterModuleLabel
 Clustered Hits. More...
 
std::string fPFParticleModuleLabel
 PFParticle label. More...
 
Track3DKalmanHitAlg fTKHAlg
 Track3DKalmanHit algorithm. More...
 
SpacePointAlg fSpacePointAlg
 Space point algorithm. More...
 
TH1F * fHIncChisq
 Incremental chisquare. More...
 
TH1F * fHPull
 Hit pull. More...
 
int fNumEvent
 Number of events seen. More...
 

Detailed Description

Definition at line 69 of file Track3DKalmanHit_module.cc.

Member Typedef Documentation

using art::EDProducer::ModuleType = EDProducer
inherited

Definition at line 34 of file EDProducer.h.

template<typename UserConfig , typename KeysToIgnore = void>
using art::EDProducer::Table = ProducerBase::Table<UserConfig, KeysToIgnore>
inherited

Definition at line 43 of file EDProducer.h.

using art::EDProducer::WorkerType = WorkerT<EDProducer>
inherited

Definition at line 35 of file EDProducer.h.

Constructor & Destructor Documentation

trkf::Track3DKalmanHit::Track3DKalmanHit ( fhicl::ParameterSet const &  pset)
explicit

Constructor.

Arguments:

p - Fcl parameters.

Definition at line 131 of file Track3DKalmanHit_module.cc.

References fClusterModuleLabel, fHitModuleLabel, fUseClusterHits, and reconfigure().

131  :
132 fHist(false),
133 fUseClusterHits(false),
134 fUsePFParticleHits(false),
135 fUsePFParticleSeeds(false),
136 fTKHAlg(pset.get<fhicl::ParameterSet>("Track3DKalmanHitAlg")),
137 fSpacePointAlg(pset.get<fhicl::ParameterSet>("SpacePointAlg")),
138 fHIncChisq(0),
139 fHPull(0),
140 fNumEvent(0)
141 {
142  reconfigure(pset);
143  produces<std::vector<recob::Track> >();
144  produces<std::vector<recob::SpacePoint> >();
145  produces<art::Assns<recob::Track, recob::Hit> >(); // ****** REMEMBER to remove when FindMany improved ******
146  produces<art::Assns<recob::Track, recob::Hit, recob::TrackHitMeta> >();
147  produces<art::Assns<recob::Track, recob::SpacePoint> >();
148  produces<art::Assns<recob::SpacePoint, recob::Hit> >();
149  produces<art::Assns<recob::PFParticle, recob::Track> >();
150 
151  // Report.
152 
153  mf::LogInfo("Track3DKalmanHit")
154  << "Track3DKalmanHit configured with the following parameters:\n"
155  << " UseClusterHits = " << fUseClusterHits << "\n"
156  << " HitModuleLabel = " << fHitModuleLabel << "\n"
157  << " ClusterModuleLabel = " << fClusterModuleLabel;
158 }
bool fUseClusterHits
Use clustered hits as input.
std::string fHitModuleLabel
Unclustered Hits.
Track3DKalmanHitAlg fTKHAlg
Track3DKalmanHit algorithm.
bool fUsePFParticleSeeds
Use PFParticle seeds.
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::string fClusterModuleLabel
Clustered Hits.
SpacePointAlg fSpacePointAlg
Space point algorithm.
int fNumEvent
Number of events seen.
bool fUsePFParticleHits
Use PFParticle hits as input.
virtual void reconfigure(fhicl::ParameterSet const &pset)
TH1F * fHIncChisq
Incremental chisquare.

Member Function Documentation

void trkf::Track3DKalmanHit::beginJob ( )
overridevirtual

Begin job method.

Reimplemented from art::EDProducer.

Definition at line 187 of file Track3DKalmanHit_module.cc.

References dir, fHIncChisq, fHist, fHPull, art::TFileDirectory::make(), and art::TFileDirectory::mkdir().

188 {
189  if(fHist) {
190  // Book histograms.
192  art::TFileDirectory dir = tfs->mkdir("hitkalman", "Track3DKalmanHit histograms");
193 
194  fHIncChisq = dir.make<TH1F>("IncChisq", "Incremental Chisquare", 100, 0., 20.);
195  fHPull = dir.make<TH1F>("Pull", "Hit Pull", 100, -10., 10.);
196  }
197 }
TFileDirectory mkdir(std::string const &dir, std::string const &descr="")
T * make(ARGS...args) const
TDirectory * dir
Definition: macro.C:5
TH1F * fHIncChisq
Incremental chisquare.
template<typename T , BranchType = InEvent>
ProductToken<T> art::Consumer::consumes ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ProductToken<T> art::Consumer::consumes ( InputTag const &  it)
inherited

Definition at line 147 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

148 {
149  if (!moduleContext_)
150  return ProductToken<T>::invalid();
151 
152  consumables_[BT].emplace_back(ConsumableType::Product,
153  TypeID{typeid(T)},
154  it.label(),
155  it.instance(),
156  it.process());
157  return ProductToken<T>{it};
158 }
static ProductToken< T > invalid()
Definition: ProductToken.h:47
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename T , art::BranchType BT>
void art::Consumer::consumesMany ( )
inherited

Definition at line 162 of file Consumer.h.

163 {
164  if (!moduleContext_)
165  return;
166 
167  consumables_[BT].emplace_back(ConsumableType::Many, TypeID{typeid(T)});
168 }
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename Element , BranchType = InEvent>
ViewToken<Element> art::Consumer::consumesView ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ViewToken<T> art::Consumer::consumesView ( InputTag const &  it)
inherited

Definition at line 172 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

173 {
174  if (!moduleContext_)
175  return ViewToken<T>::invalid();
176 
177  consumables_[BT].emplace_back(ConsumableType::ViewElement,
178  TypeID{typeid(T)},
179  it.label(),
180  it.instance(),
181  it.process());
182  return ViewToken<T>{it};
183 }
static ViewToken< Element > invalid()
Definition: ProductToken.h:75
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
EngineCreator::base_engine_t & EngineCreator::createEngine ( seed_t  seed,
std::string const &  kind_of_engine_to_make 
)
inherited

Definition at line 32 of file EngineCreator.cc.

References art::EngineCreator::rng().

34 {
35  return rng()->createEngine(
36  placeholder_schedule_id(), seed, kind_of_engine_to_make);
37 }
long seed
Definition: chem4.cc:68
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
EngineCreator::base_engine_t & EngineCreator::createEngine ( seed_t  seed,
std::string const &  kind_of_engine_to_make,
label_t const &  engine_label 
)
inherited

Definition at line 40 of file EngineCreator.cc.

References art::EngineCreator::rng().

43 {
44  return rng()->createEngine(
45  placeholder_schedule_id(), seed, kind_of_engine_to_make, engine_label);
46 }
long seed
Definition: chem4.cc:68
static art::ServiceHandle< art::RandomNumberGenerator > & rng()
void trkf::Track3DKalmanHit::createOutputs ( const art::Event evt,
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 
)
private

Definition at line 372 of file Track3DKalmanHit_module.cc.

References art::Assns< L, R, D >::addSingle(), trkf::SpacePointAlg::fillSpacePoints(), fSpacePointAlg, fUsePFParticleHits, trkf::SpacePointAlg::getAssociatedHits(), recob::Track::NumberTrajectoryPoints(), art::Event::productGetter(), art::PtrVector< T >::size(), and track.

Referenced by produce().

382 {
383  if(outputs.size()!= inputs.size()) return;
384 
385  size_t tracksSize(0);
386  for(auto const &kalman_output : outputs) {
387  tracksSize += kalman_output.tracks.size();
388  }
389  tracks.reserve(tracksSize);
390 
391  auto const tid = getProductID<std::vector<recob::Track> >();
392  auto const tidgetter = evt.productGetter(tid);
393 
394  auto const spacepointId = getProductID<std::vector<recob::SpacePoint> >();
395  auto const getter = evt.productGetter(spacepointId);
396  for (size_t i = 0; i<outputs.size();++i) {
397  // Recover the kalman tracks double ended queue
398  const std::deque<KGTrack>& kalman_tracks = outputs[i].tracks;
399 
400  for(auto const& kalman_track:kalman_tracks) {
401 
402  // Add Track object to collection.
404  kalman_track.fillTrack(track, tracks.size());
405  if(track.NumberTrajectoryPoints() < 2) {
406  continue;
407  }
408  unsigned int numtrajpts = track.NumberTrajectoryPoints();
409  tracks.emplace_back(std::move(track));
410  // SS: tracks->size() does not change after this point in each iteration
411 
412  //fill hits from this track
413  Hits trhits;
414  std::vector<unsigned int> hittpindex; //hit-trajectory point index
415  kalman_track.fillHits(trhits, hittpindex);
416  if (hittpindex.back()>=numtrajpts){ //something is wrong
417  throw cet::exception("Track3DKalmanHit")
418  << "Last hit corresponds to trajectory point index "<<hittpindex.back()<<" while the number of trajectory points is "<<numtrajpts<<'\n';
419  }
420 
421  // Make space points from this track.
422  auto nspt = spts.size();
423  fSpacePointAlg.fillSpacePoints(spts, kalman_track.TrackMap());
424 
425  std::vector<art::Ptr<recob::SpacePoint>> sptvec;
426  for(auto ispt = nspt; ispt < spts.size(); ++ispt) {
427  sptvec.emplace_back(spacepointId, ispt, getter);
428  // Associate newly created space points with hits.
429  // Make space point to hit associations.
430  auto const &sphits = fSpacePointAlg.getAssociatedHits((spts)[ispt]);
431  for(auto const& sphit: sphits) {
432  sph_assn.addSingle(sptvec.back(), sphit);
433  }
434  }
435 
436  art::Ptr<recob::Track> aptr(tid, tracks.size()-1, tidgetter);
437 
438  // Make Track to Hit associations.
439  for (size_t h = 0; h< trhits.size(); ++h){
440  th_assn.addSingle(aptr, trhits[h]);
441  recob::TrackHitMeta metadata(hittpindex[h], -1);
442  thm_assn.addSingle(aptr, trhits[h], metadata);
443  }
444 
445  // Make track to space point associations
446  for (auto const& spt: sptvec) {
447  tsp_assn.addSingle(aptr, spt);
448  }
449 
450  // Optionally fill track-to-PFParticle associations.
451  if (fUsePFParticleHits) {
452  pfPartTrack_assns.addSingle(inputs[i].pfPartPtr, aptr);
453  }
454  } // end of loop over a given collection
455  }
456 }
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:109
Data related to recob::Hit associated with recob::Track.The purpose is to collect several variables t...
Definition: TrackHitMeta.h:43
EDProductGetter const * productGetter(ProductID const) const
Definition: Event.cc:64
const art::PtrVector< recob::Hit > & getAssociatedHits(const recob::SpacePoint &spt) const
SpacePointAlg fSpacePointAlg
Space point algorithm.
void fillSpacePoints(std::vector< recob::SpacePoint > &spts, std::multimap< double, KHitTrack > const &trackMap) const
Fill a collection of space points.
art::PtrVector< recob::Hit > Hits
bool fUsePFParticleHits
Use PFParticle hits as input.
void addSingle(Ptr< left_t > const &left, Ptr< right_t > const &right, data_t const &data)
Definition: Assns.h:489
Float_t track
Definition: plot.C:34
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:51
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
CurrentProcessingContext const * art::EDProducer::currentContext ( ) const
protectedinherited

Definition at line 120 of file EDProducer.cc.

References art::EDProducer::current_context_.

121  {
122  return current_context_.get();
123  }
CPC_exempt_ptr current_context_
Definition: EDProducer.h:116
void trkf::Track3DKalmanHit::endJob ( )
overridevirtual

End job method.

Reimplemented from art::EDProducer.

Definition at line 237 of file Track3DKalmanHit_module.cc.

References fNumEvent.

238 {
239  mf::LogInfo("Track3DKalmanHit")
240  << "Track3DKalmanHit statistics:\n"
241  << " Number of events = " << fNumEvent << "\n";
242 }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
int fNumEvent
Number of events seen.
void trkf::Track3DKalmanHit::fillHistograms ( std::vector< KalmanOutput > &  outputs)
private

Fill Histograms method.

Definition at line 462 of file Track3DKalmanHit_module.cc.

References fHIncChisq, fHPull, trkf::KHitTrack::getHit(), trkf::KHit< N >::getResError(), trkf::KHit< N >::getResVector(), and trkf::KGTrack::getTrackMap().

Referenced by produce().

463 {
464  for(auto const &output : outputs) {
465  const std::deque<KGTrack>& kalman_tracks = output.tracks;
466  for (size_t i = 0; i < kalman_tracks.size(); ++i) {
467  const KGTrack& trg = kalman_tracks[i];
468  // Loop over measurements in this track.
469  const std::multimap<double, KHitTrack>& trackmap = trg.getTrackMap();
470  for(std::multimap<double, KHitTrack>::const_iterator ih = trackmap.begin();
471  ih != trackmap.end(); ++ih) {
472  const KHitTrack& trh = (*ih).second;
473  const std::shared_ptr<const KHitBase>& hit = trh.getHit();
474  double chisq = hit->getChisq();
475  fHIncChisq->Fill(chisq);
476  const KHit<1>* ph1 = dynamic_cast<const KHit<1>*>(&*hit);
477  if(ph1 != 0) {
478  double pull = ph1->getResVector()(0) / std::sqrt(ph1->getResError()(0, 0));
479  fHPull->Fill(pull);
480  }
481  }
482  }
483  }
484 }
intermediate_table::const_iterator const_iterator
Detector simulation of raw signals on wires.
TH1F * fHIncChisq
Incremental chisquare.
EngineCreator::seed_t EngineCreator::get_seed_value ( fhicl::ParameterSet const &  pset,
char const  key[] = "seed",
seed_t const  implicit_seed = -1 
)
inherited

Definition at line 49 of file EngineCreator.cc.

References fhicl::ParameterSet::get().

Referenced by art::MixFilter< T >::initEngine_().

52 {
53  auto const& explicit_seeds = pset.get<std::vector<int>>(key, {});
54  return explicit_seeds.empty() ? implicit_seed : explicit_seeds.front();
55 }
trkf::Hits trkf::Track3DKalmanHit::getAllHits ( const art::Event evt) const
private

If both UseClusteredHits and UsePFParticles is false use this method to fill in hits.

Definition at line 284 of file Track3DKalmanHit_module.cc.

References fHitModuleLabel, art::DataViewImpl::getByLabel(), hits(), art::Handle< T >::isValid(), art::PtrVector< T >::push_back(), and art::PtrVector< T >::reserve().

Referenced by getInput().

284  {
285  Hits hits;
287  evt.getByLabel(fHitModuleLabel, hith);
288  if(!hith.isValid()) return hits;
289  size_t nhits = hith->size();
290  hits.reserve(nhits);
291 
292  for(size_t i = 0; i < nhits; ++i) {
293  hits.push_back(art::Ptr<recob::Hit>(hith, i));
294  }
295  return hits;
296 }
std::string fHitModuleLabel
Unclustered Hits.
bool isValid() const
Definition: Handle.h:190
void hits()
Definition: readHits.C:15
art::PtrVector< recob::Hit > Hits
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
trkf::Hits trkf::Track3DKalmanHit::getClusteredHits ( const art::Event evt) const
private

Fill a collection using clustered hits.

Definition at line 266 of file Track3DKalmanHit_module.cc.

References art::PtrVector< T >::end(), fClusterModuleLabel, art::DataViewImpl::getByLabel(), hits(), art::PtrVector< T >::insert(), and art::Handle< T >::isValid().

Referenced by getInput().

266  {
267  Hits hits;
269  evt.getByLabel(fClusterModuleLabel, clusterh);
270 
271  if (!clusterh.isValid()) return hits;
272  // Get hits from all clusters.
273  art::FindManyP<recob::Hit> hitsbycluster(clusterh, evt, fClusterModuleLabel);
274 
275  for(size_t i = 0; i < clusterh->size(); ++i) {
276  std::vector< art::Ptr<recob::Hit> > clushits = hitsbycluster.at(i);
277  hits.insert(hits.end(), clushits.begin(), clushits.end());
278  }
279  return hits;
280 }
std::string fClusterModuleLabel
Clustered Hits.
bool isValid() const
Definition: Handle.h:190
void hits()
Definition: readHits.C:15
art::PtrVector< recob::Hit > Hits
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
trkf::KalmanInputs trkf::Track3DKalmanHit::getInput ( const art::Event evt) const
private

Definition at line 257 of file Track3DKalmanHit_module.cc.

References fUseClusterHits, fUsePFParticleHits, getAllHits(), getClusteredHits(), and getPFParticleStuff().

Referenced by produce().

257  {
258  if (fUsePFParticleHits) return getPFParticleStuff(evt);
259  return KalmanInputs(1, KalmanInput(fUseClusterHits ?
260  getClusteredHits(evt):
261  getAllHits(evt)));
262 }
bool fUseClusterHits
Use clustered hits as input.
Hits getClusteredHits(const art::Event &evt) const
Fill a collection using clustered hits.
Hits getAllHits(const art::Event &evt) const
If both UseClusteredHits and UsePFParticles is false use this method to fill in hits.
bool fUsePFParticleHits
Use PFParticle hits as input.
std::vector< KalmanInput > KalmanInputs
KalmanInputs getPFParticleStuff(const art::Event &evt) const
If UsePFParticles is true use this method to fill in hits.
trkf::KalmanInputs trkf::Track3DKalmanHit::getPFParticleStuff ( const art::Event evt) const
private

If UsePFParticles is true use this method to fill in hits.

Definition at line 300 of file Track3DKalmanHit_module.cc.

References art::PtrVector< T >::end(), fClusterModuleLabel, fPFParticleModuleLabel, fUsePFParticleSeeds, art::DataViewImpl::getByLabel(), hits(), trkf::KalmanInput::hits, art::PtrVector< T >::insert(), art::Handle< T >::isValid(), trkf::KalmanInput::pfPartPtr, trkf::KalmanInput::seedhits, trkf::KalmanInput::seeds, and x.

Referenced by getInput().

300  {
301  KalmanInputs inputs;
302  // Our program is to drive the track creation/fitting off the PFParticles in the data store
303  // We'll use the hits associated to the PFParticles for each track - and only those hits.
304  // Without a valid collection of PFParticles there is nothing to do here
305  // We need a handle to the collection of clusters in the data store so we can
306  // handle associations to hits.
308  evt.getByLabel(fPFParticleModuleLabel, pfParticleHandle);
309  if (!pfParticleHandle.isValid()) return inputs;
310 
312  evt.getByLabel(fClusterModuleLabel, clusterHandle);
313 
314  // If there are no clusters then something is really wrong
315  if (!clusterHandle.isValid()) return inputs;
316 
317  // Recover the collection of associations between PFParticles and clusters, this will
318  // be the mechanism by which we actually deal with clusters
319  art::FindManyP<recob::Cluster> clusterAssns(pfParticleHandle, evt, fPFParticleModuleLabel);
320 
321  // Associations to seeds.
322  art::FindManyP<recob::Seed> seedAssns(pfParticleHandle, evt, fPFParticleModuleLabel);
323 
324  // Likewise, recover the collection of associations to hits
325  art::FindManyP<recob::Hit> clusterHitAssns(clusterHandle, evt, fClusterModuleLabel);
326 
327  inputs.reserve(pfParticleHandle->size());
328 
329  // While PFParticle describes a hierarchal structure, for now we simply loop over the collection
330  for(size_t partIdx = 0; partIdx < pfParticleHandle->size(); partIdx++) {
331 
332  // Add a new empty hit collection.
333  inputs.emplace_back();
334  KalmanInput& kalman_input = inputs.back();
335  kalman_input.pfPartPtr = art::Ptr<recob::PFParticle>(pfParticleHandle, partIdx);
336  Hits& hits = kalman_input.hits;
337 
338  // Fill this hit vector by looping over associated clusters and finding the
339  // hits associated to them
340  std::vector<art::Ptr<recob::Cluster> > clusterVec = clusterAssns.at(partIdx);
341 
342  for(auto const &cluster : clusterVec) {
343  std::vector<art::Ptr<recob::Hit> > hitVec = clusterHitAssns.at(cluster.key());
344  hits.insert(hits.end(), hitVec.begin(), hitVec.end());
345  }
346 
347  // If requested, fill associated seeds.
348  if(!fUsePFParticleSeeds) continue;
349  art::PtrVector<recob::Seed>& seeds = kalman_input.seeds;
350  std::vector<art::Ptr<recob::Seed> > seedVec = seedAssns.at(partIdx);
351  seeds.insert(seeds.end(), seedVec.begin(), seedVec.end());
352  art::FindManyP<recob::Hit> seedHitAssns(seedVec, evt, fPFParticleModuleLabel);
353 
354  for(size_t seedIdx = 0; seedIdx < seedVec.size(); ++seedIdx) {
355  std::vector<art::Ptr<recob::Hit> > seedHitVec;
356  //SS: why seedIdx can have an invalid value?
357  try {
358  seedHitVec = seedHitAssns.at(seedIdx);
359  }
360  catch(art::Exception x) {
361  seedHitVec.clear();
362  }
363  kalman_input.seedhits.emplace_back();
364  Hits& seedhits = kalman_input.seedhits.back();
365  seedhits.insert(seedhits.end(), seedHitVec.begin(), seedHitVec.end());
366  }
367  }
368  return inputs;
369 }
Float_t x
Definition: compare.C:6
bool fUsePFParticleSeeds
Use PFParticle seeds.
std::string fPFParticleModuleLabel
PFParticle label.
std::string fClusterModuleLabel
Clustered Hits.
Cluster finding and building.
bool isValid() const
Definition: Handle.h:190
reference back()
Definition: PtrVector.h:393
void hits()
Definition: readHits.C:15
iterator end()
Definition: PtrVector.h:237
reference at(size_type n)
Definition: PtrVector.h:365
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
art::PtrVector< recob::Hit > Hits
iterator insert(iterator position, Ptr< U > const &p)
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
std::vector< KalmanInput > KalmanInputs
template<typename PROD , BranchType B>
ProductID art::EDProducer::getProductID ( std::string const &  instanceName = {}) const
inlineinherited

Definition at line 123 of file EDProducer.h.

References art::EDProducer::moduleDescription_.

124  {
125  return ProducerBase::getProductID<PROD, B>(moduleDescription_,
126  instanceName);
127  }
ModuleDescription moduleDescription_
Definition: EDProducer.h:115
template<typename PROD , BranchType B>
ProductID art::ProducerBase::getProductID ( ModuleDescription const &  moduleDescription,
std::string const &  instanceName 
) const
inherited

Definition at line 56 of file ProducerBase.h.

References B, and art::ModuleDescription::moduleLabel().

Referenced by art::ProducerBase::modifiesEvent().

58  {
59  auto const& pd =
60  get_ProductDescription<PROD>(B, md.moduleLabel(), instanceName);
61  return pd.productID();
62  }
Int_t B
Definition: plot.C:25
template<typename T , BranchType = InEvent>
ProductToken<T> art::Consumer::mayConsume ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ProductToken<T> art::Consumer::mayConsume ( InputTag const &  it)
inherited

Definition at line 190 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

191 {
192  if (!moduleContext_)
193  return ProductToken<T>::invalid();
194 
195  consumables_[BT].emplace_back(ConsumableType::Product,
196  TypeID{typeid(T)},
197  it.label(),
198  it.instance(),
199  it.process());
200  return ProductToken<T>{it};
201 }
static ProductToken< T > invalid()
Definition: ProductToken.h:47
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename T , art::BranchType BT>
void art::Consumer::mayConsumeMany ( )
inherited

Definition at line 205 of file Consumer.h.

206 {
207  if (!moduleContext_)
208  return;
209 
210  consumables_[BT].emplace_back(ConsumableType::Many, TypeID{typeid(T)});
211 }
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
template<typename Element , BranchType = InEvent>
ViewToken<Element> art::Consumer::mayConsumeView ( InputTag const &  )
inherited
template<typename T , art::BranchType BT>
art::ViewToken<T> art::Consumer::mayConsumeView ( InputTag const &  it)
inherited

Definition at line 215 of file Consumer.h.

References art::InputTag::instance(), art::InputTag::label(), and art::InputTag::process().

216 {
217  if (!moduleContext_)
218  return ViewToken<T>::invalid();
219 
220  consumables_[BT].emplace_back(ConsumableType::ViewElement,
221  TypeID{typeid(T)},
222  it.label(),
223  it.instance(),
224  it.process());
225  return ViewToken<T>{it};
226 }
static ViewToken< Element > invalid()
Definition: ProductToken.h:75
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
bool art::ProducerBase::modifiesEvent ( ) const
inlineinherited

Definition at line 40 of file ProducerBase.h.

References art::ProducerBase::getProductID().

41  {
42  return true;
43  }
void trkf::Track3DKalmanHit::prepareForInput ( )
private

Definition at line 247 of file Track3DKalmanHit_module.cc.

References trkf::SpacePointAlg::clearHitMap(), and fSpacePointAlg.

Referenced by produce().

247  {
249 }
SpacePointAlg fSpacePointAlg
Space point algorithm.
void clearHitMap() const
void art::Consumer::prepareForJob ( fhicl::ParameterSet const &  pset)
protectedinherited

Definition at line 89 of file Consumer.cc.

References fhicl::ParameterSet::get_if_present().

Referenced by art::EDProducer::doBeginJob(), art::EDFilter::doBeginJob(), and art::EDAnalyzer::doBeginJob().

90 {
91  if (!moduleContext_)
92  return;
93 
94  pset.get_if_present("errorOnMissingConsumes", requireConsumes_);
95  for (auto& consumablesPerBranch : consumables_) {
96  cet::sort_all(consumablesPerBranch);
97  }
98 }
bool requireConsumes_
Definition: Consumer.h:137
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
void trkf::Track3DKalmanHit::produce ( art::Event evt)
overridevirtual

Produce method.

Implements art::EDProducer.

Definition at line 202 of file Track3DKalmanHit_module.cc.

References trkf::SpacePointAlg::clearHitMap(), createOutputs(), fHist, fillHistograms(), fNumEvent, fSpacePointAlg, fTKHAlg, getInput(), trkf::Track3DKalmanHitAlg::makeTracks(), prepareForInput(), and art::Event::put().

203 {
204  ++fNumEvent;
205  auto tracks = std::make_unique<std::vector<recob::Track>>();
206  auto th_assn = std::make_unique<art::Assns<recob::Track, recob::Hit>>();
207  auto thm_assn = std::make_unique< art::Assns< recob::Track, recob::Hit, recob::TrackHitMeta > >();
208  auto tsp_assn = std::make_unique<art::Assns<recob::Track, recob::SpacePoint>>();
209  auto pfPartTrack_assns = std::make_unique<art::Assns<recob::PFParticle, recob::Track>>();
210  auto spts = std::make_unique<std::vector<recob::SpacePoint>>();
211  auto sph_assn = std::make_unique<art::Assns<recob::SpacePoint, recob::Hit>>();
212 
213  prepareForInput();
214  auto inputs = getInput(evt);
215 
216  auto outputs = fTKHAlg.makeTracks(inputs);
217 
218  if(fHist) {
219  fillHistograms(outputs);
220  }
221 
222  createOutputs(evt, outputs, inputs, *tracks, *spts, *th_assn, *thm_assn, *tsp_assn, *sph_assn, *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 }
Track3DKalmanHitAlg fTKHAlg
Track3DKalmanHit algorithm.
KalmanInputs getInput(const art::Event &evt) const
std::vector< trkf::KalmanOutput > makeTracks(KalmanInputs &kalman_inputs)
SpacePointAlg fSpacePointAlg
Space point algorithm.
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
void fillHistograms(std::vector< KalmanOutput > &outputs)
Fill Histograms method.
void createOutputs(const art::Event &evt, 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)
int fNumEvent
Number of events seen.
void clearHitMap() const
void trkf::Track3DKalmanHit::reconfigure ( fhicl::ParameterSet const &  pset)
virtual

Reconfigure method.

Arguments:

p - Fcl parameter set.

Definition at line 168 of file Track3DKalmanHit_module.cc.

References fClusterModuleLabel, fHist, fHitModuleLabel, fPFParticleModuleLabel, fSpacePointAlg, fTKHAlg, fUseClusterHits, fUsePFParticleHits, fUsePFParticleSeeds, fhicl::ParameterSet::get(), trkf::Track3DKalmanHitAlg::reconfigure(), and trkf::SpacePointAlg::reconfigure().

Referenced by Track3DKalmanHit().

169 {
170  fHist = pset.get<bool>("Hist");
171  fTKHAlg.reconfigure(pset.get<fhicl::ParameterSet>("Track3DKalmanHitAlg"));
172  fSpacePointAlg.reconfigure(pset.get<fhicl::ParameterSet>("SpacePointAlg"));
173  fUseClusterHits = pset.get<bool>("UseClusterHits");
174  fUsePFParticleHits = pset.get<bool>("UsePFParticleHits");
175  fUsePFParticleSeeds = pset.get<bool>("UsePFParticleSeeds");
176  fHitModuleLabel = pset.get<std::string>("HitModuleLabel");
177  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
178  fPFParticleModuleLabel = pset.get<std::string>("PFParticleModuleLabel");
179  if(fUseClusterHits && fUsePFParticleHits) {
180  throw cet::exception("Track3DKalmanHit")
181  << "Using input from both clustered and PFParticle hits.\n";
182  }
183 }
bool fUseClusterHits
Use clustered hits as input.
std::string fHitModuleLabel
Unclustered Hits.
Track3DKalmanHitAlg fTKHAlg
Track3DKalmanHit algorithm.
bool fUsePFParticleSeeds
Use PFParticle seeds.
std::string fPFParticleModuleLabel
PFParticle label.
std::string fClusterModuleLabel
Clustered Hits.
SpacePointAlg fSpacePointAlg
Space point algorithm.
void reconfigure(const fhicl::ParameterSet &pset)
Reconfigure method.
void reconfigure(const fhicl::ParameterSet &pset)
bool fUsePFParticleHits
Use PFParticle hits as input.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void art::Consumer::showMissingConsumes ( ) const
protectedinherited

Definition at line 125 of file Consumer.cc.

Referenced by art::EDProducer::doEndJob(), art::EDFilter::doEndJob(), art::EDAnalyzer::doEndJob(), and art::RootOutput::endJob().

126 {
127  if (!moduleContext_)
128  return;
129 
130  // If none of the branches have missing consumes statements, exit early.
131  if (std::all_of(cbegin(missingConsumes_),
132  cend(missingConsumes_),
133  [](auto const& perBranch) { return perBranch.empty(); }))
134  return;
135 
136  constexpr cet::HorizontalRule rule{60};
137  mf::LogPrint log{"MTdiagnostics"};
138  log << '\n'
139  << rule('=') << '\n'
140  << "The following consumes (or mayConsume) statements are missing from\n"
141  << module_context(moduleDescription_) << '\n'
142  << rule('-') << '\n';
143 
144  cet::for_all_with_index(
145  missingConsumes_, [&log](std::size_t const i, auto const& perBranch) {
146  for (auto const& pi : perBranch) {
147  log << " "
148  << assemble_consumes_statement(static_cast<BranchType>(i), pi)
149  << '\n';
150  }
151  });
152  log << rule('=');
153 }
cet::exempt_ptr< ModuleDescription const > moduleDescription_
Definition: Consumer.h:140
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
bool moduleContext_
Definition: Consumer.h:136
ConsumableProductSets missingConsumes_
Definition: Consumer.h:139
void art::Consumer::validateConsumedProduct ( BranchType const  bt,
ProductInfo const &  pi 
)
protectedinherited

Definition at line 101 of file Consumer.cc.

References art::errors::ProductRegistrationFailure.

103 {
104  // Early exits if consumes tracking has been disabled or if the
105  // consumed product is an allowed consumable.
106  if (!moduleContext_)
107  return;
108 
109  if (cet::binary_search_all(consumables_[bt], pi))
110  return;
111 
112  if (requireConsumes_) {
114  "Consumer: an error occurred during validation of a "
115  "retrieved product\n\n")
116  << "The following consumes (or mayConsume) statement is missing from\n"
117  << module_context(moduleDescription_) << ":\n\n"
118  << " " << assemble_consumes_statement(bt, pi) << "\n\n";
119  }
120 
121  missingConsumes_[bt].insert(pi);
122 }
cet::exempt_ptr< ModuleDescription const > moduleDescription_
Definition: Consumer.h:140
bool requireConsumes_
Definition: Consumer.h:137
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
ConsumableProducts consumables_
Definition: Consumer.h:138
bool moduleContext_
Definition: Consumer.h:136
ConsumableProductSets missingConsumes_
Definition: Consumer.h:139

Member Data Documentation

std::string trkf::Track3DKalmanHit::fClusterModuleLabel
private

Clustered Hits.

Definition at line 109 of file Track3DKalmanHit_module.cc.

Referenced by getClusteredHits(), getPFParticleStuff(), reconfigure(), and Track3DKalmanHit().

TH1F* trkf::Track3DKalmanHit::fHIncChisq
private

Incremental chisquare.

Definition at line 115 of file Track3DKalmanHit_module.cc.

Referenced by beginJob(), and fillHistograms().

bool trkf::Track3DKalmanHit::fHist
private

Make histograms.

Definition at line 104 of file Track3DKalmanHit_module.cc.

Referenced by beginJob(), produce(), and reconfigure().

std::string trkf::Track3DKalmanHit::fHitModuleLabel
private

Unclustered Hits.

Definition at line 108 of file Track3DKalmanHit_module.cc.

Referenced by getAllHits(), reconfigure(), and Track3DKalmanHit().

TH1F* trkf::Track3DKalmanHit::fHPull
private

Hit pull.

Definition at line 116 of file Track3DKalmanHit_module.cc.

Referenced by beginJob(), and fillHistograms().

int trkf::Track3DKalmanHit::fNumEvent
private

Number of events seen.

Definition at line 119 of file Track3DKalmanHit_module.cc.

Referenced by endJob(), and produce().

std::string trkf::Track3DKalmanHit::fPFParticleModuleLabel
private

PFParticle label.

Definition at line 110 of file Track3DKalmanHit_module.cc.

Referenced by getPFParticleStuff(), and reconfigure().

SpacePointAlg trkf::Track3DKalmanHit::fSpacePointAlg
private

Space point algorithm.

Definition at line 112 of file Track3DKalmanHit_module.cc.

Referenced by createOutputs(), prepareForInput(), produce(), and reconfigure().

Track3DKalmanHitAlg trkf::Track3DKalmanHit::fTKHAlg
private

Track3DKalmanHit algorithm.

Definition at line 111 of file Track3DKalmanHit_module.cc.

Referenced by produce(), and reconfigure().

bool trkf::Track3DKalmanHit::fUseClusterHits
private

Use clustered hits as input.

Definition at line 105 of file Track3DKalmanHit_module.cc.

Referenced by getInput(), reconfigure(), and Track3DKalmanHit().

bool trkf::Track3DKalmanHit::fUsePFParticleHits
private

Use PFParticle hits as input.

Definition at line 106 of file Track3DKalmanHit_module.cc.

Referenced by createOutputs(), getInput(), and reconfigure().

bool trkf::Track3DKalmanHit::fUsePFParticleSeeds
private

Use PFParticle seeds.

Definition at line 107 of file Track3DKalmanHit_module.cc.

Referenced by getPFParticleStuff(), and reconfigure().


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