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