LArSoft  v10_04_05
Liquid Argon Software toolkit - https://larsoft.org/
PFParticleCosmicAna_module.cc
Go to the documentation of this file.
1 
9 
11 
12 #include "TTree.h"
13 
14 #include <string>
15 
16 //------------------------------------------------------------------------------------------------------------------------------------------
17 
18 namespace lar_pandora {
19 
24  public:
31 
32  void beginJob();
33  void endJob();
34  void analyze(const art::Event& evt);
35  void reconfigure(fhicl::ParameterSet const& pset);
36 
37  private:
45  void FillRecoTree(const PFParticlesToHits& recoParticlesToHits,
46  const PFParticlesToTracks& recoParticlesToTracks,
47  const TracksToCosmicTags& recoTracksToCosmicTags);
48 
59  void FillTrueTree(const HitVector& hitVector,
60  const HitsToMCParticles& trueHitsToParticles,
61  const HitsToPFParticles& recoHitsToParticles,
62  const MCParticlesToMCTruth& particlesToTruth,
63  const PFParticlesToTracks& particlesToTracks,
64  const TracksToCosmicTags& tracksToCosmicTags);
65 
73  float GetCosmicScore(const art::Ptr<recob::PFParticle> particle,
74  const PFParticlesToTracks& recoParticlesToTracks,
75  const TracksToCosmicTags& recoTracksToCosmicTags) const;
76 
77  TTree* m_pRecoTree;
78  TTree* m_pTrueTree;
79 
80  int m_run;
81  int m_event;
82  int m_index;
83 
84  int m_self;
85  int m_pdgCode;
88  float m_cosmicScore;
89  int m_nTracks;
90  int m_nHits;
91 
92  float m_trackVtxX;
93  float m_trackVtxY;
94  float m_trackVtxZ;
95  float m_trackEndX;
96  float m_trackEndY;
97  float m_trackEndZ;
110 
113 
120 
127 
128  std::string m_hitfinderLabel;
129  std::string m_trackfitLabel;
130  std::string m_particleLabel;
131  std::string m_cosmicLabel;
132  std::string m_geantModuleLabel;
133 
136 
138  };
139 
141 
142 } // namespace lar_pandora
143 
144 //------------------------------------------------------------------------------------------------------------------------------------------
145 // implementation follows
146 
150 #include "art_root_io/TFileDirectory.h"
151 #include "art_root_io/TFileService.h"
154 #include "fhiclcpp/ParameterSet.h"
156 
162 
163 #include <iostream>
164 
165 namespace lar_pandora {
166 
168  {
169  m_cosmicLabel = pset.get<std::string>("CosmicTagModule", "cosmictagger");
170  m_particleLabel = pset.get<std::string>("PFParticleModule", "pandora");
171  m_trackfitLabel = pset.get<std::string>("TrackFitModule", "trackfit");
172  m_hitfinderLabel = pset.get<std::string>("HitFinderModule", "gaushit");
173  m_geantModuleLabel = pset.get<std::string>("GeantModule", "largeant");
174 
175  m_useDaughterPFParticles = pset.get<bool>("UseDaughterPFParticles", true);
176  m_useDaughterMCParticles = pset.get<bool>("UseDaughterMCParticles", true);
177 
178  m_cosmicContainmentCut = pset.get<double>("CosmicContainmentCut", 5.0);
179  }
180 
181  //------------------------------------------------------------------------------------------------------------------------------------------
182 
184  {
185  mf::LogDebug("LArPandora") << " *** PFParticleCosmicAna::beginJob() *** " << std::endl;
186 
187  //
189 
190  m_pRecoTree = tfs->make<TTree>("recoTree", "LAr Cosmic Reco Tree");
191  m_pRecoTree->Branch("run", &m_run, "run/I");
192  m_pRecoTree->Branch("event", &m_event, "event/I");
193  m_pRecoTree->Branch("index", &m_index, "index/I");
194  m_pRecoTree->Branch("self", &m_self, "self/I");
195  m_pRecoTree->Branch("pdgCode", &m_pdgCode, "pdgCode/I");
196  m_pRecoTree->Branch("isTrackLike", &m_isTrackLike, "isTrackLike/I");
197  m_pRecoTree->Branch("isPrimary", &m_isPrimary, "isPrimary/I");
198  m_pRecoTree->Branch("cosmicScore", &m_cosmicScore, "cosmicScore/F");
199  m_pRecoTree->Branch("trackVtxX", &m_trackVtxX, "trackVtxX/F");
200  m_pRecoTree->Branch("trackVtxY", &m_trackVtxY, "trackVtxY/F");
201  m_pRecoTree->Branch("trackVtxZ", &m_trackVtxZ, "trackVtxZ/F");
202  m_pRecoTree->Branch("trackEndX", &m_trackEndX, "trackEndX/F");
203  m_pRecoTree->Branch("trackEndY", &m_trackEndY, "trackEndY/F");
204  m_pRecoTree->Branch("trackEndZ", &m_trackEndZ, "trackEndZ/F");
205  m_pRecoTree->Branch("trackVtxDirX", &m_trackVtxDirX, "trackVtxDirX/F");
206  m_pRecoTree->Branch("trackVtxDirY", &m_trackVtxDirY, "trackVtxDirY/F");
207  m_pRecoTree->Branch("trackVtxDirZ", &m_trackVtxDirZ, "trackVtxDirZ/F");
208  m_pRecoTree->Branch("trackEndDirX", &m_trackEndDirX, "trackEndDirX/F");
209  m_pRecoTree->Branch("trackEndDirY", &m_trackEndDirY, "trackEndDirY/F");
210  m_pRecoTree->Branch("trackEndDirZ", &m_trackEndDirZ, "trackEndDirZ/F");
211  m_pRecoTree->Branch("trackLength", &m_trackLength, "trackLength/F");
212  m_pRecoTree->Branch("trackWidthX", &m_trackWidthX, "trackWidthX/F");
213  m_pRecoTree->Branch("trackWidthY", &m_trackWidthY, "trackWidthY/F");
214  m_pRecoTree->Branch("trackWidthZ", &m_trackWidthZ, "trackWidthZ/F");
215  m_pRecoTree->Branch("trackVtxDeltaYZ", &m_trackVtxDeltaYZ, "trackVtxDeltaYZ/F");
216  m_pRecoTree->Branch("trackEndDeltaYZ", &m_trackEndDeltaYZ, "trackEndDeltaYZ/F");
217  m_pRecoTree->Branch("trackVtxContained", &m_trackVtxContained, "trackVtxContained/I");
218  m_pRecoTree->Branch("trackEndContained", &m_trackEndContained, "trackEndContained/I");
219  m_pRecoTree->Branch("nTracks", &m_nTracks, "nTracks/I");
220  m_pRecoTree->Branch("nHits", &m_nHits, "nHits/I");
221 
222  m_pTrueTree = tfs->make<TTree>("trueTree", "LAr Cosmic True Tree");
223  m_pTrueTree->Branch("run", &m_run, "run/I");
224  m_pTrueTree->Branch("event", &m_event, "event/I");
225  m_pTrueTree->Branch("nHits", &m_nHits, "nHits/I");
226  m_pTrueTree->Branch("nNeutrinoHits", &m_nNeutrinoHits, "nNeutrinoHits/I");
227  m_pTrueTree->Branch(
228  "nNeutrinoHitsFullyTagged", &m_nNeutrinoHitsFullyTagged, "nNeutrinoHitsFullyTagged/I");
229  m_pTrueTree->Branch(
230  "nNeutrinoHitsSemiTagged", &m_nNeutrinoHitsSemiTagged, "nNeutrinoHitsSemiTagged/I");
231  m_pTrueTree->Branch(
232  "nNeutrinoHitsNotTagged", &m_nNeutrinoHitsNotTagged, "nNeutrinoHitsNotTagged/I");
233  m_pTrueTree->Branch("nNeutrinoHitsNotReconstructed",
235  "nNeutrinoHitsNotReconstructed/I");
236  m_pTrueTree->Branch(
237  "nNeutrinoHitsReconstructed", &m_nNeutrinoHitsReconstructed, "nNeutrinoHitsReconstructed/I");
238  m_pTrueTree->Branch("nCosmicHits", &m_nCosmicHits, "nCosmicHits/I");
239  m_pTrueTree->Branch(
240  "nCosmicHitsFullyTagged", &m_nCosmicHitsFullyTagged, "nCosmicHitsFullyTagged/I");
241  m_pTrueTree->Branch(
242  "nCosmicHitsSemiTagged", &m_nCosmicHitsSemiTagged, "nCosmicHitsSemiTagged/I");
243  m_pTrueTree->Branch("nCosmicHitsNotTagged", &m_nCosmicHitsNotTagged, "nCosmicHitsNotTagged/I");
244  m_pTrueTree->Branch("nCosmicHitsNotReconstructed",
246  "nCosmicHitsNotReconstructed/I");
247  m_pTrueTree->Branch(
248  "nCosmicHitsReconstructed", &m_nCosmicHitsReconstructed, "nCosmicHitsReconstructed/I");
249  }
250 
251  //------------------------------------------------------------------------------------------------------------------------------------------
252 
254 
255  //------------------------------------------------------------------------------------------------------------------------------------------
256 
258  {
259  std::cout << " *** PFParticleCosmicAna::analyze(...) *** " << std::endl;
260 
261  //
262  // Note: I've made this is MicroBooNE-only module
263  //
264 
265  m_run = evt.run();
266  m_event = evt.id().event();
267 
268  std::cout << " Run: " << m_run << std::endl;
269  std::cout << " Event: " << m_event << std::endl;
270 
271  // Collect True Particles
272  // ======================
273  HitVector hitVector;
274  MCTruthToMCParticles truthToParticles;
275  MCParticlesToMCTruth particlesToTruth;
276  MCParticlesToHits trueParticlesToHits;
277  HitsToMCParticles trueHitsToParticles;
278 
281  evt, m_geantModuleLabel, truthToParticles, particlesToTruth);
284  hitVector,
285  trueParticlesToHits,
286  trueHitsToParticles,
290 
291  // Collect Reco Particles
292  // ======================
293  PFParticleVector recoParticleVector;
294  PFParticlesToHits recoParticlesToHits;
295  HitsToPFParticles recoHitsToParticles;
296 
297  LArPandoraHelper::CollectPFParticles(evt, m_particleLabel, recoParticleVector);
301  recoParticlesToHits,
302  recoHitsToParticles,
305  LArPandoraHelper::kIgnoreDaughters));
306 
307  std::cout << " PFParticles: " << recoParticleVector.size() << std::endl;
308 
309  // Collect Reco Tracks
310  // ===================
311  TrackVector recoTrackVector;
312  PFParticlesToTracks recoParticlesToTracks;
313  LArPandoraHelper::CollectTracks(evt, m_trackfitLabel, recoTrackVector, recoParticlesToTracks);
314 
315  // Collect Cosmic Tags
316  // =====================
317  CosmicTagVector recoCosmicTagVector;
318  TracksToCosmicTags recoTracksToCosmicTags;
320  evt, m_cosmicLabel, recoCosmicTagVector, recoTracksToCosmicTags);
321 
322  // Analyse Reconstructed Particles
323  // ===============================
324  this->FillRecoTree(recoParticlesToHits, recoParticlesToTracks, recoTracksToCosmicTags);
325 
326  // Analyse True Hits
327  // =================
328  this->FillTrueTree(hitVector,
329  trueHitsToParticles,
330  recoHitsToParticles,
331  particlesToTruth,
332  recoParticlesToTracks,
333  recoTracksToCosmicTags);
334  }
335 
336  //------------------------------------------------------------------------------------------------------------------------------------------
337 
338  void PFParticleCosmicAna::FillRecoTree(const PFParticlesToHits& recoParticlesToHits,
339  const PFParticlesToTracks& recoParticlesToTracks,
340  const TracksToCosmicTags& recoTracksToCosmicTags)
341  {
342  // Set up Geometry Service
343  // =======================
344  auto const& tpc = art::ServiceHandle<geo::Geometry const>()->TPC({0, 0});
345 
346  const double xmin(0.0);
347  const double xmax(2.0 * tpc.HalfWidth());
348  const double ymin(-tpc.HalfHeight());
349  const double ymax(+tpc.HalfHeight());
350  const double zmin(0.0);
351  const double zmax(tpc.Length());
352  const double xyzCut(m_cosmicContainmentCut);
353 
354  m_index = 0;
355 
356  m_self = 0;
357  m_pdgCode = 0;
358  m_isTrackLike = 0;
359  m_isPrimary = 0;
360  m_cosmicScore = 0.f;
361 
362  m_trackVtxX = 0.f;
363  m_trackVtxY = 0.f;
364  m_trackVtxZ = 0.f;
365  m_trackEndX = 0.f;
366  m_trackEndY = 0.f;
367  m_trackEndZ = 0.f;
368  m_trackVtxDirX = 0.f;
369  m_trackVtxDirY = 0.f;
370  m_trackVtxDirZ = 0.f;
371  m_trackEndDirX = 0.f;
372  m_trackEndDirY = 0.f;
373  m_trackEndDirZ = 0.f;
374  m_trackLength = 0.f;
375  m_trackWidthX = 0.f;
376  m_trackWidthY = 0.f;
377  m_trackWidthZ = 0.f;
378  m_trackVtxDeltaYZ = 0.f;
379  m_trackEndDeltaYZ = 0.f;
380 
383 
384  m_nTracks = 0;
385  m_nHits = 0;
386 
387  // Loop over Reco Particles
388  // ========================
389  for (PFParticlesToHits::const_iterator iter1 = recoParticlesToHits.begin(),
390  iterEnd1 = recoParticlesToHits.end();
391  iter1 != iterEnd1;
392  ++iter1) {
393  const art::Ptr<recob::PFParticle> recoParticle = iter1->first;
394 
395  const HitVector& hitVector = iter1->second;
396  if (hitVector.empty()) continue;
397 
398  PFParticlesToTracks::const_iterator iter2 = recoParticlesToTracks.find(recoParticle);
399  if (recoParticlesToTracks.end() == iter2) continue;
400 
401  const TrackVector& trackVector = iter2->second;
402  if (trackVector.empty()) continue;
403 
404  m_nHits = hitVector.size();
405  m_nTracks = trackVector.size();
406 
407  m_self = recoParticle->Self();
408  m_pdgCode = recoParticle->PdgCode();
409  m_isPrimary = recoParticle->IsPrimary();
411  m_cosmicScore =
412  this->GetCosmicScore(recoParticle, recoParticlesToTracks, recoTracksToCosmicTags);
413 
414  m_trackVtxX = 0.f;
415  m_trackVtxY = 0.f;
416  m_trackVtxZ = 0.f;
417  m_trackEndX = 0.f;
418  m_trackEndY = 0.f;
419  m_trackEndZ = 0.f;
420  m_trackVtxDirX = 0.f;
421  m_trackVtxDirY = 0.f;
422  m_trackVtxDirZ = 0.f;
423  m_trackEndDirX = 0.f;
424  m_trackEndDirY = 0.f;
425  m_trackEndDirZ = 0.f;
426  m_trackLength = 0.f;
427  m_trackWidthX = 0.f;
428  m_trackWidthY = 0.f;
429  m_trackWidthZ = 0.f;
430  m_trackVtxDeltaYZ = 0.f;
431  m_trackEndDeltaYZ = 0.f;
432 
435 
436  for (TrackVector::const_iterator iter3 = trackVector.begin(), iterEnd3 = trackVector.end();
437  iter3 != iterEnd3;
438  ++iter3) {
439  const art::Ptr<recob::Track> track = *iter3;
440  const float trackLength(track->Length());
441 
442  if (trackLength < m_trackLength) continue;
443 
444  m_trackLength = trackLength;
445 
446  const auto& trackVtxPosition = track->Vertex();
447  const auto& trackVtxDirection = track->VertexDirection();
448  const auto& trackEndPosition = track->End();
449  const auto& trackEndDirection = track->EndDirection();
450 
451  m_trackVtxX = trackVtxPosition.x();
452  m_trackVtxY = trackVtxPosition.y();
453  m_trackVtxZ = trackVtxPosition.z();
454  m_trackVtxDirX = trackVtxDirection.x();
455  m_trackVtxDirY = trackVtxDirection.y();
456  m_trackVtxDirZ = trackVtxDirection.z();
457  m_trackEndX = trackEndPosition.x();
458  m_trackEndY = trackEndPosition.y();
459  m_trackEndZ = trackEndPosition.z();
460  m_trackEndDirX = trackEndDirection.x();
461  m_trackEndDirY = trackEndDirection.y();
462  m_trackEndDirZ = trackEndDirection.z();
463 
464  m_trackWidthX = std::fabs(m_trackEndX - m_trackVtxX);
465  m_trackWidthY = std::fabs(m_trackEndY - m_trackVtxY);
466  m_trackWidthZ = std::fabs(m_trackEndZ - m_trackVtxZ);
467 
469  std::min((ymax - m_trackVtxY), std::min((m_trackVtxZ - zmin), (zmax - m_trackVtxZ)));
471  std::min((m_trackEndY - ymin), std::min((m_trackEndZ - zmin), (zmax - m_trackEndZ)));
472 
473  m_trackVtxContained = ((m_trackVtxX > xmin + xyzCut && m_trackVtxX < xmax - xyzCut) &&
474  (m_trackVtxY > ymin + xyzCut && m_trackVtxY < ymax - xyzCut) &&
475  (m_trackVtxZ > zmin + xyzCut && m_trackVtxZ < zmax - xyzCut));
476  m_trackEndContained = ((m_trackEndX > xmin + xyzCut && m_trackEndX < xmax - xyzCut) &&
477  (m_trackEndY > ymin + xyzCut && m_trackEndY < ymax - xyzCut) &&
478  (m_trackEndZ > zmin + xyzCut && m_trackEndZ < zmax - xyzCut));
479  }
480 
481  std::cout << " PFParticle: [" << m_index << "] nHits=" << m_nHits
482  << ", nTracks=" << m_nTracks << ", cosmicScore=" << m_cosmicScore << std::endl;
483 
484  m_pRecoTree->Fill();
485  ++m_index;
486  }
487  }
488 
489  //------------------------------------------------------------------------------------------------------------------------------------------
490 
492  const HitsToMCParticles& trueHitsToParticles,
493  const HitsToPFParticles& recoHitsToParticles,
494  const MCParticlesToMCTruth& particlesToTruth,
495  const PFParticlesToTracks& particlesToTracks,
496  const TracksToCosmicTags& tracksToCosmicTags)
497  {
498  m_nHits = 0;
499 
500  m_nNeutrinoHits = 0;
506 
507  m_nCosmicHits = 0;
513 
514  for (HitVector::const_iterator iter2 = hitVector.begin(), iterEnd2 = hitVector.end();
515  iter2 != iterEnd2;
516  ++iter2) {
517  const art::Ptr<recob::Hit> hit = *iter2;
518 
519  HitsToMCParticles::const_iterator iter3 = trueHitsToParticles.find(hit);
520  if (trueHitsToParticles.end() == iter3) continue;
521 
522  const art::Ptr<simb::MCParticle> trueParticle = iter3->second;
523 
524  MCParticlesToMCTruth::const_iterator iter4 = particlesToTruth.find(trueParticle);
525  if (particlesToTruth.end() == iter4)
526  throw cet::exception("LArPandora") << " PFParticleCosmicAna::analyze --- Found a true "
527  "particle without any ancestry information ";
528 
529  const art::Ptr<simb::MCTruth> truth = iter4->second;
530 
531  float cosmicScore(-0.2);
532 
533  HitsToPFParticles::const_iterator iter5 = recoHitsToParticles.find(hit);
534  if (recoHitsToParticles.end() != iter5) {
535  const art::Ptr<recob::PFParticle> particle = iter5->second;
536  cosmicScore = this->GetCosmicScore(particle, particlesToTracks, tracksToCosmicTags);
537  }
538 
539  ++m_nHits;
540 
541  if (truth->NeutrinoSet()) {
542  ++m_nNeutrinoHits;
543 
544  if (cosmicScore >= 0)
546  else
548 
549  if (cosmicScore > 0.51)
551  else if (cosmicScore > 0.39)
553  else
555  }
556  else {
557  ++m_nCosmicHits;
558 
559  if (cosmicScore >= 0)
561  else
563 
564  if (cosmicScore > 0.51)
566  else if (cosmicScore > 0.39)
568  else
570  }
571  }
572 
573  m_pTrueTree->Fill();
574  }
575 
576  //------------------------------------------------------------------------------------------------------------------------------------------
577 
579  const PFParticlesToTracks& recoParticlesToTracks,
580  const TracksToCosmicTags& recoTracksToCosmicTags) const
581  {
582  float cosmicScore(0.f);
583 
584  // Get cosmic tags associated with this particle
585  PFParticlesToTracks::const_iterator iter2 = recoParticlesToTracks.find(particle);
586  if (recoParticlesToTracks.end() != iter2) {
587  for (TrackVector::const_iterator iter3 = iter2->second.begin(),
588  iterEnd3 = iter2->second.end();
589  iter3 != iterEnd3;
590  ++iter3) {
591  const art::Ptr<recob::Track> track = *iter3;
592 
593  TracksToCosmicTags::const_iterator iter4 = recoTracksToCosmicTags.find(track);
594  if (recoTracksToCosmicTags.end() != iter4) {
595  for (CosmicTagVector::const_iterator iter5 = iter4->second.begin(),
596  iterEnd5 = iter4->second.end();
597  iter5 != iterEnd5;
598  ++iter5) {
599  const art::Ptr<anab::CosmicTag> cosmicTag = *iter5;
600  if (cosmicTag->CosmicScore() > cosmicScore) cosmicScore = cosmicTag->CosmicScore();
601  }
602  }
603  }
604  }
605 
606  return cosmicScore;
607  }
608 
609 } //namespace lar_pandora
static void BuildPFParticleHitMaps(const PFParticleVector &particleVector, const PFParticlesToSpacePoints &particlesToSpacePoints, const SpacePointsToHits &spacePointsToHits, PFParticlesToHits &particlesToHits, HitsToPFParticles &hitsToParticles, const DaughterMode daughterMode=kUseDaughters)
Build mapping between PFParticles and Hits using PFParticle/SpacePoint/Hit maps.
std::map< art::Ptr< recob::Hit >, art::Ptr< simb::MCParticle > > HitsToMCParticles
size_t Self() const
Returns the index of this particle.
Definition: PFParticle.h:88
void FillTrueTree(const HitVector &hitVector, const HitsToMCParticles &trueHitsToParticles, const HitsToPFParticles &recoHitsToParticles, const MCParticlesToMCTruth &particlesToTruth, const PFParticlesToTracks &particlesToTracks, const TracksToCosmicTags &tracksToCosmicTags)
Fill track-level variables using input maps between reconstructed objects.
std::map< art::Ptr< simb::MCParticle >, art::Ptr< simb::MCTruth > > MCParticlesToMCTruth
Vector_t VertexDirection() const
Access to track direction at different points.
Definition: Track.h:166
int PdgCode() const
Return the type of particle as a PDG ID.
Definition: PFParticle.h:79
intermediate_table::const_iterator const_iterator
std::vector< art::Ptr< anab::CosmicTag > > CosmicTagVector
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
std::map< art::Ptr< recob::PFParticle >, TrackVector > PFParticlesToTracks
float & CosmicScore()
Definition: CosmicTag.h:54
std::map< art::Ptr< recob::Track >, CosmicTagVector > TracksToCosmicTags
TFile f
Definition: plotHisto.C:6
void reconfigure(fhicl::ParameterSet const &pset)
double Length(size_t p=0) const
Access to various track properties.
Definition: Track.h:207
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
std::map< art::Ptr< simb::MCParticle >, HitVector > MCParticlesToHits
std::vector< art::Ptr< recob::PFParticle > > PFParticleVector
std::map< art::Ptr< simb::MCTruth >, MCParticleVector > MCTruthToMCParticles
T get(std::string const &key) const
Definition: ParameterSet.h:314
std::map< art::Ptr< recob::PFParticle >, HitVector > PFParticlesToHits
Point_t const & Vertex() const
Access to track position at different points.
Definition: Track.h:158
bool IsPrimary() const
Returns whether the particle is the root of the flow.
Definition: PFParticle.h:82
std::vector< art::Ptr< recob::Track > > TrackVector
static void CollectTracks(const art::Event &evt, const std::string &label, TrackVector &trackVector, PFParticlesToTracks &particlesToTracks)
Collect the reconstructed PFParticles and associated Tracks from the ART event record.
static void CollectPFParticles(const art::Event &evt, const std::string &label, PFParticleVector &particleVector)
Collect the reconstructed PFParticles from the ART event record.
void FillRecoTree(const PFParticlesToHits &recoParticlesToHits, const PFParticlesToTracks &recoParticlesToTracks, const TracksToCosmicTags &recoTracksToCosmicTags)
Fill event-level variables using input maps between reconstructed objects.
static void CollectHits(const art::Event &evt, const std::string &label, HitVector &hitVector)
Collect the reconstructed Hits from the ART event record.
Detector simulation of raw signals on wires.
std::vector< art::Ptr< recob::Hit > > HitVector
float GetCosmicScore(const art::Ptr< recob::PFParticle > particle, const PFParticlesToTracks &recoParticlesToTracks, const TracksToCosmicTags &recoTracksToCosmicTags) const
Get cosmic score for a PFParticle using track-level information.
Provides recob::Track data product.
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::PFParticle > > HitsToPFParticles
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
Vector_t EndDirection() const
Access to track direction at different points.
Definition: Track.h:167
static bool IsTrack(const art::Ptr< recob::PFParticle > particle)
Determine whether a particle has been reconstructed as track-like.
Definition: MVAAlg.h:12
EventNumber_t event() const
Definition: EventID.h:116
Point_t const & End() const
Access to track position at different points.
Definition: Track.h:159
bool NeutrinoSet() const
Definition: MCTruth.h:78
TCEvent evt
Definition: DataStructs.cxx:8
RunNumber_t run() const
Definition: Event.cc:29
static void BuildMCParticleHitMaps(const art::Event &evt, const HitVector &hitVector, const SimChannelVector &simChannelVector, HitsToTrackIDEs &hitsToTrackIDEs)
Collect the links from reconstructed hits to their true energy deposits.
static void CollectMCParticles(const art::Event &evt, const std::string &label, MCParticleVector &particleVector)
Collect a vector of MCParticle objects from the ART event record.
Float_t track
Definition: plot.C:35
helper function for LArPandoraInterface producer module
PFParticleCosmicAna(fhicl::ParameterSet const &pset)
Constructor.
EventID id() const
Definition: Event.cc:23
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
static void CollectCosmicTags(const art::Event &evt, const std::string &label, CosmicTagVector &cosmicTagVector, TracksToCosmicTags &tracksToCosmicTags)
Collect a vector of cosmic tags from the ART event record.