LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
PFParticleHitDumper_module.cc
Go to the documentation of this file.
1 
10 
11 #include "TTree.h"
12 
16 
17 #include <string>
18 
19 //------------------------------------------------------------------------------------------------------------------------------------------
20 
21 namespace lar_pandora {
22 
27  public:
34 
38  virtual ~PFParticleHitDumper();
39 
40  void beginJob();
41  void endJob();
42  void analyze(const art::Event& evt);
43  void reconfigure(fhicl::ParameterSet const& pset);
44 
45  private:
51  void FillRecoTracks(const PFParticlesToTracks& particlesToTracks);
52 
60  void FillReco3D(const PFParticleVector& particleVector,
61  const PFParticlesToSpacePoints& particlesToSpacePoints,
62  const SpacePointsToHits& spacePointsToHits);
63 
70  void FillReco2D(const art::Event& event,
71  const HitVector& hitVector,
72  const HitsToPFParticles& hitsToParticles);
73 
83  void FillAssociated2DHits(const art::Event& evt,
84  const PFParticleVector& particleVector,
85  const PFParticlesToHits& particlesToHits,
86  const PFParticlesToHits& particlesToHitsClusters,
87  const PFParticlesToTracks& particlesToTracks,
88  const TracksToHits& tracksToHits,
89  const PFParticlesToShowers& particlesToShowers,
90  const ShowersToHits& showersToHits);
91 
97  void FillRecoWires(const art::Event& event, const WireVector& wireVector);
98 
104  double GetUVW(const geo::WireID& wireID) const;
105 
114  double YZtoU(const unsigned int cstat,
115  const unsigned int tpc,
116  const double y,
117  const double z) const;
118 
127  double YZtoV(const unsigned int cstat,
128  const unsigned int tpc,
129  const double y,
130  const double z) const;
131 
132  TTree* m_pRecoTracks;
133  TTree* m_pReco3D;
134  TTree* m_pReco2D;
136  TTree* m_pRecoWire;
137 
138  int m_run;
139  int m_event;
141  int m_primary;
142  int m_pdgcode;
143 
144  int m_cstat;
145  int m_tpc;
146  int m_plane;
147  int m_wire;
148 
149  double m_u;
150  double m_v;
151  double m_w;
152  double m_x;
153  double m_y;
154  double m_z;
155  double m_q;
156 
160 
161  std::string m_calwireLabel;
162  std::string m_hitfinderLabel;
163  std::string m_clusterLabel;
164  std::string m_spacepointLabel;
165  std::string m_particleLabel;
166  std::string m_trackLabel;
167  std::string m_showerLabel;
168 
171  };
172 
174 
175 } // namespace lar_pandora
176 
177 //------------------------------------------------------------------------------------------------------------------------------------------
178 // implementation follows
179 
183 #include "art_root_io/TFileDirectory.h"
184 #include "art_root_io/TFileService.h"
185 #include "fhiclcpp/ParameterSet.h"
187 
196 #include "lardataobj/RecoBase/Hit.h"
204 
205 #include <iostream>
206 
207 namespace lar_pandora {
208 
210  {
211  this->reconfigure(pset);
212  }
213 
214  //------------------------------------------------------------------------------------------------------------------------------------------
215 
217 
218  //------------------------------------------------------------------------------------------------------------------------------------------
219 
221  {
222  m_storeWires = pset.get<bool>("StoreWires", false);
223  m_trackLabel = pset.get<std::string>("TrackModule", "pandoraTrack");
224  m_showerLabel = pset.get<std::string>("ShowerModule", "pandoraShower");
225  m_particleLabel = pset.get<std::string>("PFParticleModule", "pandora");
226  m_spacepointLabel = pset.get<std::string>("SpacePointModule", "pandora");
227  m_clusterLabel = pset.get<std::string>("ClusterModule", "pandora");
228  m_hitfinderLabel = pset.get<std::string>("HitFinderModule", "gaushit");
229  m_calwireLabel = pset.get<std::string>("CalWireModule", "caldata");
230  m_printDebug = pset.get<bool>("PrintDebug", false);
231  }
232 
233  //------------------------------------------------------------------------------------------------------------------------------------------
234 
236  {
237  mf::LogDebug("LArPandora") << " *** PFParticleHitDumper::beginJob() *** " << std::endl;
238 
239  //
241 
242  m_pRecoTracks = tfs->make<TTree>("pandoraTracks", "LAr Reco Tracks");
243  m_pRecoTracks->Branch("run", &m_run, "run/I");
244  m_pRecoTracks->Branch("event", &m_event, "event/I");
245  m_pRecoTracks->Branch("particle", &m_particle, "particle/I");
246  m_pRecoTracks->Branch("x", &m_x, "x/D");
247  m_pRecoTracks->Branch("y", &m_y, "y/D");
248  m_pRecoTracks->Branch("z", &m_z, "z/D");
249 
250  m_pReco3D = tfs->make<TTree>("pandora3D", "LAr Reco 3D");
251  m_pReco3D->Branch("run", &m_run, "run/I");
252  m_pReco3D->Branch("event", &m_event, "event/I");
253  m_pReco3D->Branch("particle", &m_particle, "particle/I");
254  m_pReco3D->Branch("primary", &m_primary, "primary/I");
255  m_pReco3D->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
256  m_pReco3D->Branch("cstat", &m_cstat, "cstat/I");
257  m_pReco3D->Branch("tpc", &m_tpc, "tpc/I");
258  m_pReco3D->Branch("plane", &m_plane, "plane/I");
259  m_pReco3D->Branch("x", &m_x, "x/D");
260  m_pReco3D->Branch("y", &m_y, "y/D");
261  m_pReco3D->Branch("u", &m_u, "u/D");
262  m_pReco3D->Branch("v", &m_v, "v/D");
263  m_pReco3D->Branch("z", &m_z, "z/D");
264 
265  m_pReco2D = tfs->make<TTree>("pandora2D", "LAr Reco 2D");
266  m_pReco2D->Branch("run", &m_run, "run/I");
267  m_pReco2D->Branch("event", &m_event, "event/I");
268  m_pReco2D->Branch("particle", &m_particle, "particle/I");
269  m_pReco2D->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
270  m_pReco2D->Branch("cstat", &m_cstat, "cstat/I");
271  m_pReco2D->Branch("tpc", &m_tpc, "tpc/I");
272  m_pReco2D->Branch("plane", &m_plane, "plane/I");
273  m_pReco2D->Branch("wire", &m_wire, "wire/I");
274  m_pReco2D->Branch("x", &m_x, "x/D");
275  m_pReco2D->Branch("w", &m_w, "w/D");
276  m_pReco2D->Branch("q", &m_q, "q/D");
277 
278  m_pRecoComparison = tfs->make<TTree>("pandora2Dcomparison", "LAr Reco 2D (comparison)");
279  m_pRecoComparison->Branch("run", &m_run, "run/I");
280  m_pRecoComparison->Branch("event", &m_event, "event/I");
281  m_pRecoComparison->Branch("particle", &m_particle, "particle/I");
282  m_pRecoComparison->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
283  m_pRecoComparison->Branch(
284  "hitsFromSpacePoints", &m_hitsFromSpacePoints, "hitsFromSpacePoints/I");
285  m_pRecoComparison->Branch("hitsFromClusters", &m_hitsFromClusters, "hitsFromClusters/I");
286  m_pRecoComparison->Branch(
287  "hitsFromTrackOrShower", &m_hitsFromTrackOrShower, "hitsFromTrackOrShower/I");
288 
289  m_pRecoWire = tfs->make<TTree>("rawdata", "LAr Reco Wires");
290  m_pRecoWire->Branch("run", &m_run, "run/I");
291  m_pRecoWire->Branch("event", &m_event, "event/I");
292  m_pRecoWire->Branch("cstat", &m_cstat, "cstat/I");
293  m_pRecoWire->Branch("tpc", &m_tpc, "tpc/I");
294  m_pRecoWire->Branch("plane", &m_plane, "plane/I");
295  m_pRecoWire->Branch("wire", &m_wire, "wire/I");
296  m_pRecoWire->Branch("x", &m_x, "x/D");
297  m_pRecoWire->Branch("w", &m_w, "w/D");
298  m_pRecoWire->Branch("q", &m_q, "q/D");
299  }
300 
301  //------------------------------------------------------------------------------------------------------------------------------------------
302 
304 
305  //------------------------------------------------------------------------------------------------------------------------------------------
306 
308  {
309  if (m_printDebug) std::cout << " *** PFParticleHitDumper::analyze(...) *** " << std::endl;
310 
311  m_run = evt.run();
312  m_event = evt.id().event();
313 
314  m_particle = -1;
315  m_primary = 0;
316  m_pdgcode = 0;
317 
318  m_cstat = 0;
319  m_tpc = 0;
320  m_plane = 0;
321  m_wire = 0;
322 
323  m_x = 0.0;
324  m_y = 0.0;
325  m_u = 0.0;
326  m_v = 0.0;
327  m_z = 0.0;
328  m_w = 0.0;
329  m_q = 0.0;
330 
331  if (m_printDebug) {
332  std::cout << " Run: " << m_run << std::endl;
333  std::cout << " Event: " << m_event << std::endl;
334  }
335 
336  // Need geometry service to convert channel to wire ID
338 
339  // Get particles, tracks, space points, hits (and wires)
340  // ====================================================
341  TrackVector trackVector, trackVectorExtra;
342  ShowerVector showerVector, showerVectorExtra;
343  PFParticleVector particleVector;
344  SpacePointVector spacePointVector;
345  HitVector hitVector;
346  WireVector wireVector;
347 
348  PFParticlesToTracks particlesToTracks;
349  PFParticlesToShowers particlesToShowers;
350  PFParticlesToSpacePoints particlesToSpacePoints;
351  PFParticlesToHits particlesToHits, particlesToHitsClusters;
352  TracksToHits tracksToHits;
353  ShowersToHits showersToHits;
354  HitsToPFParticles hitsToParticles, hitsToParticlesClusters;
355  SpacePointsToHits spacePointsToHits;
356 
359  evt, m_spacepointLabel, spacePointVector, spacePointsToHits);
360  LArPandoraHelper::CollectTracks(evt, m_trackLabel, trackVector, particlesToTracks);
361  LArPandoraHelper::CollectTracks(evt, m_trackLabel, trackVectorExtra, tracksToHits);
362  LArPandoraHelper::CollectShowers(evt, m_showerLabel, showerVector, particlesToShowers);
363  LArPandoraHelper::CollectShowers(evt, m_showerLabel, showerVectorExtra, showersToHits);
365  evt, m_particleLabel, particleVector, particlesToSpacePoints);
369  particlesToHits,
370  hitsToParticles,
371  LArPandoraHelper::DaughterMode::kUseDaughters,
372  false);
374  evt, m_particleLabel, m_clusterLabel, particlesToHitsClusters, hitsToParticlesClusters);
375 
377 
378  if (m_printDebug) std::cout << " PFParticles: " << particleVector.size() << std::endl;
379 
380  // Loop over Tracks (Fill 3D Track Tree)
381  // =====================================
382  if (m_printDebug) std::cout << " PFParticleHitDumper::FillRecoTracks(...) " << std::endl;
383  this->FillRecoTracks(particlesToTracks);
384 
385  // Loop over PFParticles (Fill 3D Reco Tree)
386  // =========================================
387  if (m_printDebug) std::cout << " PFParticleHitDumper::FillReco3D(...) " << std::endl;
388  this->FillReco3D(particleVector, particlesToSpacePoints, spacePointsToHits);
389 
390  // Loop over Hits (Fill 2D Reco Tree)
391  // ==================================
392  if (m_printDebug) std::cout << " PFParticleHitDumper::FillReco2D(...) " << std::endl;
393  this->FillReco2D(evt, hitVector, hitsToParticles);
394 
395  // Loop over Hits (Fill Associated 2D Hits Tree)
396  // =============================================
397  if (m_printDebug)
398  std::cout << " PFParticleHitDumper::FillAssociated2DHits(...) " << std::endl;
399  this->FillAssociated2DHits(evt,
400  particleVector,
401  particlesToHits,
402  particlesToHitsClusters,
403  particlesToTracks,
404  tracksToHits,
405  particlesToShowers,
406  showersToHits);
407 
408  // Loop over Wires (Fill Reco Wire Tree)
409  // =====================================
410  if (m_printDebug) std::cout << " PFParticleHitDumper::FillRecoWires(...) " << std::endl;
411  this->FillRecoWires(evt, wireVector);
412  }
413 
414  //------------------------------------------------------------------------------------------------------------------------------------------
415 
417  {
418  // Initialise variables
419  m_particle = -1;
420  m_x = 0.0;
421  m_y = 0.0;
422  m_z = 0.0;
423 
424  // Create dummy entry if there are no particles
425  if (particlesToTracks.empty()) { m_pRecoTracks->Fill(); }
426 
427  // Loop over tracks
428  for (PFParticlesToTracks::const_iterator iter = particlesToTracks.begin(),
429  iterEnd = particlesToTracks.end();
430  iter != iterEnd;
431  ++iter) {
432  const art::Ptr<recob::PFParticle> particle = iter->first;
433  const TrackVector& trackVector = iter->second;
434 
435  m_particle = particle->Self();
436 
437  if (!trackVector.empty()) {
438  if (trackVector.size() != 1 && m_printDebug)
439  std::cout << " Warning: Found particle with more than one associated track " << std::endl;
440 
441  const art::Ptr<recob::Track> track = *(trackVector.begin());
442 
443  if (m_printDebug)
444  std::cout << " PFPARTICLE [" << m_particle << "] (" << track->NumberTrajectoryPoints()
445  << " Trajectory Points)" << std::endl;
446 
447  for (unsigned int p = 0; p < track->NumberTrajectoryPoints(); ++p) {
448  const auto position(track->LocationAtPoint(p));
449  m_x = position.x();
450  m_y = position.y();
451  m_z = position.z();
452 
453  m_pRecoTracks->Fill();
454  }
455  }
456  }
457  }
458 
459  //------------------------------------------------------------------------------------------------------------------------------------------
460 
462  const PFParticlesToSpacePoints& particlesToSpacePoints,
463  const SpacePointsToHits& spacePointsToHits)
464  {
465  // Initialise variables
466  m_particle = -1;
467  m_primary = 0;
468  m_pdgcode = 0;
469  m_cstat = 0;
470  m_tpc = 0;
471  m_plane = 0;
472  m_x = 0.0;
473  m_u = 0.0;
474  m_v = 0.0;
475  m_y = 0.0;
476  m_z = 0.0;
477 
478  // Create dummy entry if there are no particles
479  if (particleVector.empty()) { m_pReco3D->Fill(); }
480 
481  // Store associations between particle and particle ID
482  PFParticleMap theParticleMap;
483 
484  for (unsigned int i = 0; i < particleVector.size(); ++i) {
485  const art::Ptr<recob::PFParticle> particle = particleVector.at(i);
486  theParticleMap[particle->Self()] = particle;
487  }
488 
489  // Loop over particles
490  for (PFParticlesToSpacePoints::const_iterator iter1 = particlesToSpacePoints.begin(),
491  iterEnd1 = particlesToSpacePoints.end();
492  iter1 != iterEnd1;
493  ++iter1) {
494  const art::Ptr<recob::PFParticle> particle = iter1->first;
495  const SpacePointVector& spacepoints = iter1->second;
496 
497  m_particle = particle->Self();
498  m_pdgcode = particle->PdgCode();
499  m_primary = 0;
500 
501  if (particle->IsPrimary()) { m_primary = 1; }
502  else {
503  const size_t parentID(particle->Parent());
504  PFParticleMap::const_iterator pIter = theParticleMap.find(parentID);
505  if (theParticleMap.end() == pIter)
506  throw cet::exception("LArPandora")
507  << " PFParticleHitDumper::analyze --- Found particle with ID code";
508 
509  const art::Ptr<recob::PFParticle> particleParent = pIter->second;
510  if (LArPandoraHelper::IsNeutrino(particleParent)) m_primary = 1;
511  }
512 
513  if (m_printDebug)
514  std::cout << " PFPARTICLE [" << m_particle << "] [Primary=" << m_primary << "] ("
515  << spacepoints.size() << " Space Points)" << std::endl;
516 
517  for (unsigned int j = 0; j < spacepoints.size(); ++j) {
518  const art::Ptr<recob::SpacePoint> spacepoint = spacepoints.at(j);
519 
520  m_x = spacepoint->XYZ()[0];
521  m_y = spacepoint->XYZ()[1];
522  m_z = spacepoint->XYZ()[2];
523 
524  SpacePointsToHits::const_iterator iter2 = spacePointsToHits.find(spacepoint);
525  if (spacePointsToHits.end() == iter2)
526  throw cet::exception("LArPandora")
527  << " PFParticleHitDumper::analyze --- Found space point without associated hit";
528 
529  const art::Ptr<recob::Hit> hit = iter2->second;
530  const geo::WireID& wireID(hit->WireID());
531 
532  m_cstat = wireID.Cryostat;
533  m_tpc = wireID.TPC;
534  m_plane = wireID.Plane;
535 
536  m_u = this->YZtoU(m_cstat, m_tpc, m_y, m_z);
537  m_v = this->YZtoV(m_cstat, m_tpc, m_y, m_z);
538 
539  m_pReco3D->Fill();
540  }
541  }
542  }
543 
544  //------------------------------------------------------------------------------------------------------------------------------------------
545 
547  const PFParticleVector& particleVector,
548  const PFParticlesToHits& particlesToHits,
549  const PFParticlesToHits& particlesToHitsClusters,
550  const PFParticlesToTracks& particlesToTracks,
551  const TracksToHits& tracksToHits,
552  const PFParticlesToShowers& particlesToShowers,
553  const ShowersToHits& showersToHits)
554  {
555  // Create dummy entry if there are no 2D hits
556  if (particleVector.empty()) { m_pRecoComparison->Fill(); }
557 
558  for (unsigned int i = 0; i < particleVector.size(); ++i) {
559  //initialise variables
560  m_particle = -1;
561  m_pdgcode = 0;
563  m_hitsFromClusters = 0;
565 
566  const art::Ptr<recob::PFParticle> particle = particleVector.at(i);
567  m_particle = particle->Self();
568  m_pdgcode = particle->PdgCode();
569 
570  PFParticlesToHits::const_iterator pIter = particlesToHits.find(particle);
571  PFParticlesToHits::const_iterator pIter2 = particlesToHitsClusters.find(particle);
572  if (particlesToHits.end() != pIter) m_hitsFromSpacePoints = pIter->second.size();
573  if (particlesToHitsClusters.end() != pIter2) m_hitsFromClusters = pIter2->second.size();
574 
575  if (m_pdgcode == 13) {
576  PFParticlesToTracks::const_iterator iter = particlesToTracks.find(particle);
577  const art::Ptr<recob::Track> track = *(iter->second.begin());
578 
579  TracksToHits::const_iterator iter2 = tracksToHits.find(track);
580  if (tracksToHits.end() != iter2) {
581  const HitVector& hitVector = iter2->second;
582  m_hitsFromTrackOrShower = hitVector.size();
583  }
584  }
585  else if (m_pdgcode == 11) {
586  PFParticlesToShowers::const_iterator iter = particlesToShowers.find(particle);
587  const art::Ptr<recob::Shower> shower = *(iter->second.begin());
588 
589  ShowersToHits::const_iterator iter2 = showersToHits.find(shower);
590  if (showersToHits.end() != iter2) {
591  const HitVector& hitVector = iter2->second;
592  m_hitsFromTrackOrShower = hitVector.size();
593  }
594  }
595 
596  if (m_printDebug)
597  std::cout << " PFParticle " << m_particle << " (PDG " << m_pdgcode << ") has "
598  << m_hitsFromSpacePoints << " hits from space points and " << m_hitsFromClusters
599  << " hits from clusters, and its recob::Track/Shower has "
600  << m_hitsFromTrackOrShower << " associated hits " << std::endl;
601 
602  m_pRecoComparison->Fill();
603  }
604  }
605 
606  //------------------------------------------------------------------------------------------------------------------------------------------
607 
609  const HitVector& hitVector,
610  const HitsToPFParticles& hitsToParticles)
611  {
612  // Initialise variables
613  m_particle = -1;
614  m_pdgcode = 0;
615  m_cstat = 0;
616  m_tpc = 0;
617  m_plane = 0;
618  m_wire = 0;
619  m_x = 0.0;
620  m_w = 0.0;
621  m_q = 0.0;
622 
623  // Create dummy entry if there are no 2D hits
624  if (hitVector.empty()) { m_pReco2D->Fill(); }
625 
626  // Need DetectorProperties service to convert from ticks to X
627  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(e);
628 
629  // Loop over 2D hits
630  for (unsigned int i = 0; i < hitVector.size(); ++i) {
631  const art::Ptr<recob::Hit> hit = hitVector.at(i);
632 
633  m_particle = -1;
634  m_pdgcode = 0;
635 
636  HitsToPFParticles::const_iterator pIter = hitsToParticles.find(hit);
637  if (hitsToParticles.end() != pIter) {
638  const art::Ptr<recob::PFParticle> particle = pIter->second;
639  m_particle = particle->Self();
640  m_pdgcode = particle->PdgCode();
641  }
642 
643  const geo::WireID& wireID(hit->WireID());
644  m_cstat = wireID.Cryostat;
645  m_tpc = wireID.TPC;
646  m_plane = wireID.Plane;
647  m_wire = wireID.Wire;
648 
649  m_q = hit->Integral();
650  m_x = detProp.ConvertTicksToX(hit->PeakTime(), wireID.Plane, wireID.TPC, wireID.Cryostat);
651  m_w = this->GetUVW(wireID);
652 
653  m_pReco2D->Fill();
654  }
655  }
656 
657  //------------------------------------------------------------------------------------------------------------------------------------------
658 
660  {
661 
662  // Create dummy entry if there are no wires
663  if (wireVector.empty()) { m_pRecoWire->Fill(); }
664 
665  // Need geometry service to convert channel to wire ID
667 
668  // Need DetectorProperties service to convert from ticks to X
669  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(e);
670 
671  // Loop over wires
672  int signalCounter(0);
673 
674  for (unsigned int i = 0; i < wireVector.size(); ++i) {
675  const art::Ptr<recob::Wire> wire = wireVector.at(i);
676 
677  const std::vector<float>& signals(wire->Signal());
678  const std::vector<geo::WireID> wireIds = theGeometry->ChannelToWire(wire->Channel());
679 
680  if ((signalCounter++) < 10 && m_printDebug)
681  std::cout << " numWires=" << wireVector.size() << " numSignals=" << signals.size()
682  << std::endl;
683 
684  double time(0.0);
685 
686  m_q = 0.0;
687 
688  for (std::vector<float>::const_iterator tIter = signals.begin(), tIterEnd = signals.end();
689  tIter != tIterEnd;
690  ++tIter) {
691  time += 1.0;
692  m_q = *tIter;
693 
694  if (m_q < 2.0) // seems to remove most noise
695  continue;
696 
697  for (std::vector<geo::WireID>::const_iterator wIter = wireIds.begin(),
698  wIterEnd = wireIds.end();
699  wIter != wIterEnd;
700  ++wIter) {
701  const geo::WireID& wireID = *wIter;
702  m_cstat = wireID.Cryostat;
703  m_tpc = wireID.TPC;
704  m_plane = wireID.Plane;
705  m_wire = wireID.Wire;
706 
707  m_x = detProp.ConvertTicksToX(time, wireID.Plane, wireID.TPC, wireID.Cryostat);
708  m_w = this->GetUVW(wireID);
709 
710  m_pRecoWire->Fill();
711  }
712  }
713  }
714  }
715 
716  //------------------------------------------------------------------------------------------------------------------------------------------
717 
718  double PFParticleHitDumper::GetUVW(const geo::WireID& wireID) const
719  {
720  // define UVW as closest distance from (0,0) to wire axis
722 
723  auto const xyzStart = theGeometry->Wire(wireID).GetStart();
724  const double ay(xyzStart.Y());
725  const double az(xyzStart.Z());
726 
727  auto const xyzEnd = theGeometry->Wire(wireID).GetEnd();
728  const double by(xyzEnd.Y());
729  const double bz(xyzEnd.Z());
730 
731  const double ny(by - ay);
732  const double nz(bz - az);
733  const double N2(ny * ny + nz * nz);
734 
735  const double ry(ay - (ay * ny + az * nz) * ny / N2);
736  const double rz(az - (ay * ny + az * nz) * nz / N2);
737  const double sign((rz > 0.0) ? +1.0 : -1.0);
738 
739  return sign * std::sqrt(ry * ry + rz * rz);
740  }
741 
742  //------------------------------------------------------------------------------------------------------------------------------------------
743 
744  double PFParticleHitDumper::YZtoU(const unsigned int cstat,
745  const unsigned int tpc,
746  const double y,
747  const double z) const
748  {
749  // TODO: Check that this stills works in DUNE
751  const double m_theta(theGeometry->WireAngleToVertical(geo::kU, geo::TPCID{cstat, tpc}));
752  return z * std::sin(m_theta) - y * std::cos(m_theta);
753  }
754 
755  //------------------------------------------------------------------------------------------------------------------------------------------
756 
757  double PFParticleHitDumper::YZtoV(const unsigned int cstat,
758  const unsigned int tpc,
759  const double y,
760  const double z) const
761  {
762  // TODO; Check that this still works in DUNE
764  const double m_theta(theGeometry->WireAngleToVertical(geo::kV, geo::TPCID{cstat, tpc}));
765  return z * std::sin(m_theta) - y * std::cos(m_theta);
766  }
767 
768  //------------------------------------------------------------------------------------------------------------------------------------------
769 
770 } //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.
void FillReco2D(const art::Event &event, const HitVector &hitVector, const HitsToPFParticles &hitsToParticles)
Store 2D hits.
double YZtoU(const unsigned int cstat, const unsigned int tpc, const double y, const double z) const
Convert from (Y,Z) to U coordinate.
std::vector< WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
Encapsulate the construction of a single cyostat.
size_t Self() const
Returns the index of this particle.
Definition: PFParticle.h:88
static bool IsNeutrino(const art::Ptr< recob::PFParticle > particle)
Determine whether a particle has been reconstructed as a neutrino.
Planes which measure V.
Definition: geo_types.h:136
Declaration of signal hit object.
double YZtoV(const unsigned int cstat, const unsigned int tpc, const double y, const double z) const
Convert from (Y,Z) to V coordinate.
Float_t y
Definition: compare.C:6
Point_t const & LocationAtPoint(size_t i) const
Access to track position at different points.
Definition: Track.h:160
Double_t z
Definition: plot.C:276
static void CollectWires(const art::Event &evt, const std::string &label, WireVector &wireVector)
Collect the reconstructed wires from the ART event record.
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:136
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
int PdgCode() const
Return the type of particle as a PDG ID.
Definition: PFParticle.h:79
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:244
std::vector< art::Ptr< recob::Shower > > ShowerVector
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
intermediate_table::const_iterator const_iterator
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
void FillRecoWires(const art::Event &event, const WireVector &wireVector)
Store raw data.
WireGeo const & Wire(WireID const &wireid) const
Returns the specified wire.
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
Point_t GetStart() const
Returns the world coordinate of one end of the wire [cm].
Definition: WireGeo.h:226
Particle class.
Point_t GetEnd() const
Returns the world coordinate of one end of the wire [cm].
Definition: WireGeo.h:231
std::map< art::Ptr< recob::PFParticle >, TrackVector > PFParticlesToTracks
double GetUVW(const geo::WireID &wireID) const
Conversion from wire ID to U/V/W coordinate.
std::map< int, art::Ptr< recob::PFParticle > > PFParticleMap
static void CollectSpacePoints(const art::Event &evt, const std::string &label, SpacePointVector &spacePointVector, SpacePointsToHits &spacePointsToHits)
Collect the reconstructed SpacePoints and associated hits from the ART event record.
void FillAssociated2DHits(const art::Event &evt, const PFParticleVector &particleVector, const PFParticlesToHits &particlesToHits, const PFParticlesToHits &particlesToHitsClusters, const PFParticlesToTracks &particlesToTracks, const TracksToHits &tracksToHits, const PFParticlesToShowers &particlesToShowers, const ShowersToHits &showersToHits)
Store number of 2D hits associated to PFParticle in different ways.
Planes which measure U.
Definition: geo_types.h:135
double WireAngleToVertical(View_t view, TPCID const &tpcid) const
Returns the angle of the wires in the specified view from vertical.
bool m_printDebug
switch for print statements (TODO: use message service!)
std::map< art::Ptr< recob::PFParticle >, ShowerVector > PFParticlesToShowers
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
raw::ChannelID_t Channel() const
Returns the ID of the channel (or InvalidChannelID)
Definition: Wire.h:223
void reconfigure(fhicl::ParameterSet const &pset)
std::vector< art::Ptr< recob::PFParticle > > PFParticleVector
size_t Parent() const
Definition: PFParticle.h:92
PFParticleHitDumper(fhicl::ParameterSet const &pset)
Constructor.
std::map< art::Ptr< recob::Shower >, HitVector > ShowersToHits
static void CollectShowers(const art::Event &evt, const std::string &label, ShowerVector &showerVector, PFParticlesToShowers &particlesToShowers)
Collect the reconstructed PFParticles and associated Showers from the ART event record.
T get(std::string const &key) const
Definition: ParameterSet.h:314
std::vector< art::Ptr< recob::SpacePoint > > SpacePointVector
std::map< art::Ptr< recob::PFParticle >, HitVector > PFParticlesToHits
void FillReco3D(const PFParticleVector &particleVector, const PFParticlesToSpacePoints &particlesToSpacePoints, const SpacePointsToHits &spacePointsToHits)
Store 3D hits.
Provides recob::Track data product.
bool IsPrimary() const
Returns whether the particle is the root of the flow.
Definition: PFParticle.h:82
const Double32_t * XYZ() const
Definition: SpacePoint.h:78
std::vector< art::Ptr< recob::Track > > TrackVector
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:481
std::map< art::Ptr< recob::PFParticle >, SpacePointVector > PFParticlesToSpacePoints
Declaration of cluster object.
std::map< art::Ptr< recob::Track >, HitVector > TracksToHits
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.
Definition of data types for geometry description.
static void CollectPFParticles(const art::Event &evt, const std::string &label, PFParticleVector &particleVector)
Collect the reconstructed PFParticles from the ART event record.
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< float > Signal() const
Return a zero-padded full length vector filled with RoI signal.
Definition: Wire.cxx:30
void FillRecoTracks(const PFParticlesToTracks &particlesToTracks)
Store 3D track hits.
int sign(double val)
Definition: UtilFunc.cxx:104
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
std::vector< art::Ptr< recob::Hit > > HitVector
Utility object to perform functions of association.
Encapsulate the construction of a single detector plane.
std::vector< art::Ptr< recob::Wire > > WireVector
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::PFParticle > > HitsToPFParticles
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
Definition: MVAAlg.h:12
EventNumber_t event() const
Definition: EventID.h:116
Declaration of basic channel signal object.
TCEvent evt
Definition: DataStructs.cxx:8
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
Float_t e
Definition: plot.C:35
RunNumber_t run() const
Definition: Event.cc:29
Float_t track
Definition: plot.C:35
helper function for LArPandoraInterface producer module
EventID id() const
Definition: Event.cc:23
art framework interface to geometry description
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > SpacePointsToHits
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Event finding and building.
Encapsulate the construction of a single detector plane.