All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PFParticleHitDumper_module.cc
Go to the documentation of this file.
1 /**
2  * @file larpandora/LArPandoraAnalysis/PFParticleHitDumper_module.cc
3  *
4  * @brief Analysis module for created particles
5  *
6  */
7 
8 #include "art/Framework/Core/EDAnalyzer.h"
9 #include "art/Framework/Core/ModuleMacros.h"
10 
11 #include "TTree.h"
12 
16 
17 #include <string>
18 
19 //------------------------------------------------------------------------------------------------------------------------------------------
20 
21 namespace lar_pandora {
22 
23  /**
24  * @brief PFParticleHitDumper class
25  */
26  class PFParticleHitDumper : public art::EDAnalyzer {
27  public:
28  /**
29  * @brief Constructor
30  *
31  * @param pset
32  */
33  PFParticleHitDumper(fhicl::ParameterSet const& pset);
34 
35  /**
36  * @brief Destructor
37  */
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  /**
47  * @brief Store 3D track hits
48  *
49  * @param particlesToTracks mapping between 3D track hits and PFParticles
50  */
51  void FillRecoTracks(const PFParticlesToTracks& particlesToTracks);
52 
53  /**
54  * @brief Store 3D hits
55  *
56  * @param particleVector the input vector of PFParticles
57  * @param particlesToSpacePoints mapping between 3D hits and PFParticles
58  * @param spacePointsToHits mapping between 3D hits and 2D hits
59  */
60  void FillReco3D(const PFParticleVector& particleVector,
61  const PFParticlesToSpacePoints& particlesToSpacePoints,
62  const SpacePointsToHits& spacePointsToHits);
63 
64  /**
65  * @brief Store 2D hits
66  *
67  * @param hitVector the input vector of 2D hits
68  * @param hitsToParticles mapping between 2D hits and PFParticles
69  */
70  void FillReco2D(const art::Event& event,
71  const HitVector& hitVector,
72  const HitsToPFParticles& hitsToParticles);
73 
74  /**
75  * @brief Store number of 2D hits associated to PFParticle in different ways
76  *
77  * @param particleVector the input vector of PFParticles
78  * @param particlesToHits mapping between PFParticles and 2D hits through space points
79  * @param particlesToHitsClusters mapping between PFParticles and 2D hits through clusters
80  * @param particlesToTracks mapping between PFParticles and tracks
81  * @param particlesToShowers mapping between PFParticles and showers
82  */
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 
92  /**
93  * @brief Store raw data
94  *
95  * @param wireVector the input vector of reconstructed wires
96  */
97  void FillRecoWires(const art::Event& event, const WireVector& wireVector);
98 
99  /**
100  * @brief Conversion from wire ID to U/V/W coordinate
101  *
102  * @param wireID the input wire ID
103  */
104  double GetUVW(const geo::WireID& wireID) const;
105 
106  /**
107  * @brief Convert from (Y,Z) to U coordinate
108  *
109  * @param cstat the cryostat
110  * @param tpc the tpc
111  * @param y the y coordinate
112  * @param z the z coordinate
113  */
114  double YZtoU(const unsigned int cstat,
115  const unsigned int tpc,
116  const double y,
117  const double z) const;
118 
119  /**
120  * @brief Convert from (Y,Z) to V coordinate
121  *
122  * @param cstat the crystat
123  * @param tpc the tpc
124  * @param y the y coordinate
125  * @param z the z coordinate
126  */
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; ///<
135  TTree* m_pRecoComparison; ///<
136  TTree* m_pRecoWire; ///<
137 
138  int m_run; ///<
139  int m_event; ///<
140  int m_particle; ///<
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 
169  bool m_storeWires; ///<
170  bool m_printDebug; ///< switch for print statements (TODO: use message service!)
171  };
172 
173  DEFINE_ART_MODULE(PFParticleHitDumper)
174 
175 } // namespace lar_pandora
176 
177 //------------------------------------------------------------------------------------------------------------------------------------------
178 // implementation follows
179 
180 #include "art/Framework/Principal/Event.h"
181 #include "art/Framework/Principal/Handle.h"
182 #include "art/Framework/Services/Registry/ServiceHandle.h"
183 #include "art_root_io/TFileDirectory.h"
184 #include "art_root_io/TFileService.h"
185 #include "fhiclcpp/ParameterSet.h"
186 #include "messagefacility/MessageLogger/MessageLogger.h"
187 
196 #include "lardataobj/RecoBase/Hit.h"
202 #include "nusimdata/SimulationBase/MCParticle.h"
203 #include "nusimdata/SimulationBase/MCTruth.h"
204 
205 #include <iostream>
206 
207 namespace lar_pandora {
208 
209  PFParticleHitDumper::PFParticleHitDumper(fhicl::ParameterSet const& pset) : art::EDAnalyzer(pset)
210  {
211  this->reconfigure(pset);
212  }
213 
214  //------------------------------------------------------------------------------------------------------------------------------------------
215 
217 
218  //------------------------------------------------------------------------------------------------------------------------------------------
219 
220  void
221  PFParticleHitDumper::reconfigure(fhicl::ParameterSet const& pset)
222  {
223  m_storeWires = pset.get<bool>("StoreWires", false);
224  m_trackLabel = pset.get<std::string>("TrackModule", "pandoraTrack");
225  m_showerLabel = pset.get<std::string>("ShowerModule", "pandoraShower");
226  m_particleLabel = pset.get<std::string>("PFParticleModule", "pandora");
227  m_spacepointLabel = pset.get<std::string>("SpacePointModule", "pandora");
228  m_clusterLabel = pset.get<std::string>("ClusterModule", "pandora");
229  m_hitfinderLabel = pset.get<std::string>("HitFinderModule", "gaushit");
230  m_calwireLabel = pset.get<std::string>("CalWireModule", "caldata");
231  m_printDebug = pset.get<bool>("PrintDebug", false);
232  }
233 
234  //------------------------------------------------------------------------------------------------------------------------------------------
235 
236  void
238  {
239  mf::LogDebug("LArPandora") << " *** PFParticleHitDumper::beginJob() *** " << std::endl;
240 
241  //
242  art::ServiceHandle<art::TFileService const> tfs;
243 
244  m_pRecoTracks = tfs->make<TTree>("pandoraTracks", "LAr Reco Tracks");
245  m_pRecoTracks->Branch("run", &m_run, "run/I");
246  m_pRecoTracks->Branch("event", &m_event, "event/I");
247  m_pRecoTracks->Branch("particle", &m_particle, "particle/I");
248  m_pRecoTracks->Branch("x", &m_x, "x/D");
249  m_pRecoTracks->Branch("y", &m_y, "y/D");
250  m_pRecoTracks->Branch("z", &m_z, "z/D");
251 
252  m_pReco3D = tfs->make<TTree>("pandora3D", "LAr Reco 3D");
253  m_pReco3D->Branch("run", &m_run, "run/I");
254  m_pReco3D->Branch("event", &m_event, "event/I");
255  m_pReco3D->Branch("particle", &m_particle, "particle/I");
256  m_pReco3D->Branch("primary", &m_primary, "primary/I");
257  m_pReco3D->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
258  m_pReco3D->Branch("cstat", &m_cstat, "cstat/I");
259  m_pReco3D->Branch("tpc", &m_tpc, "tpc/I");
260  m_pReco3D->Branch("plane", &m_plane, "plane/I");
261  m_pReco3D->Branch("x", &m_x, "x/D");
262  m_pReco3D->Branch("y", &m_y, "y/D");
263  m_pReco3D->Branch("u", &m_u, "u/D");
264  m_pReco3D->Branch("v", &m_v, "v/D");
265  m_pReco3D->Branch("z", &m_z, "z/D");
266 
267  m_pReco2D = tfs->make<TTree>("pandora2D", "LAr Reco 2D");
268  m_pReco2D->Branch("run", &m_run, "run/I");
269  m_pReco2D->Branch("event", &m_event, "event/I");
270  m_pReco2D->Branch("particle", &m_particle, "particle/I");
271  m_pReco2D->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
272  m_pReco2D->Branch("cstat", &m_cstat, "cstat/I");
273  m_pReco2D->Branch("tpc", &m_tpc, "tpc/I");
274  m_pReco2D->Branch("plane", &m_plane, "plane/I");
275  m_pReco2D->Branch("wire", &m_wire, "wire/I");
276  m_pReco2D->Branch("x", &m_x, "x/D");
277  m_pReco2D->Branch("w", &m_w, "w/D");
278  m_pReco2D->Branch("q", &m_q, "q/D");
279 
280  m_pRecoComparison = tfs->make<TTree>("pandora2Dcomparison", "LAr Reco 2D (comparison)");
281  m_pRecoComparison->Branch("run", &m_run, "run/I");
282  m_pRecoComparison->Branch("event", &m_event, "event/I");
283  m_pRecoComparison->Branch("particle", &m_particle, "particle/I");
284  m_pRecoComparison->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
285  m_pRecoComparison->Branch(
286  "hitsFromSpacePoints", &m_hitsFromSpacePoints, "hitsFromSpacePoints/I");
287  m_pRecoComparison->Branch("hitsFromClusters", &m_hitsFromClusters, "hitsFromClusters/I");
288  m_pRecoComparison->Branch(
289  "hitsFromTrackOrShower", &m_hitsFromTrackOrShower, "hitsFromTrackOrShower/I");
290 
291  m_pRecoWire = tfs->make<TTree>("rawdata", "LAr Reco Wires");
292  m_pRecoWire->Branch("run", &m_run, "run/I");
293  m_pRecoWire->Branch("event", &m_event, "event/I");
294  m_pRecoWire->Branch("cstat", &m_cstat, "cstat/I");
295  m_pRecoWire->Branch("tpc", &m_tpc, "tpc/I");
296  m_pRecoWire->Branch("plane", &m_plane, "plane/I");
297  m_pRecoWire->Branch("wire", &m_wire, "wire/I");
298  m_pRecoWire->Branch("x", &m_x, "x/D");
299  m_pRecoWire->Branch("w", &m_w, "w/D");
300  m_pRecoWire->Branch("q", &m_q, "q/D");
301  }
302 
303  //------------------------------------------------------------------------------------------------------------------------------------------
304 
305  void
307  {}
308 
309  //------------------------------------------------------------------------------------------------------------------------------------------
310 
311  void
313  {
314  if (m_printDebug) std::cout << " *** PFParticleHitDumper::analyze(...) *** " << std::endl;
315 
316  m_run = evt.run();
317  m_event = evt.id().event();
318 
319  m_particle = -1;
320  m_primary = 0;
321  m_pdgcode = 0;
322 
323  m_cstat = 0;
324  m_tpc = 0;
325  m_plane = 0;
326  m_wire = 0;
327 
328  m_x = 0.0;
329  m_y = 0.0;
330  m_u = 0.0;
331  m_v = 0.0;
332  m_z = 0.0;
333  m_w = 0.0;
334  m_q = 0.0;
335 
336  if (m_printDebug) {
337  std::cout << " Run: " << m_run << std::endl;
338  std::cout << " Event: " << m_event << std::endl;
339  }
340 
341  // Need geometry service to convert channel to wire ID
342  art::ServiceHandle<geo::Geometry const> theGeometry;
343 
344  // Get particles, tracks, space points, hits (and wires)
345  // ====================================================
346  TrackVector trackVector, trackVectorExtra;
347  ShowerVector showerVector, showerVectorExtra;
348  PFParticleVector particleVector;
349  SpacePointVector spacePointVector;
350  HitVector hitVector;
351  WireVector wireVector;
352 
353  PFParticlesToTracks particlesToTracks;
354  PFParticlesToShowers particlesToShowers;
355  PFParticlesToSpacePoints particlesToSpacePoints;
356  PFParticlesToHits particlesToHits, particlesToHitsClusters;
357  TracksToHits tracksToHits;
358  ShowersToHits showersToHits;
359  HitsToPFParticles hitsToParticles, hitsToParticlesClusters;
360  SpacePointsToHits spacePointsToHits;
361 
364  evt, m_spacepointLabel, spacePointVector, spacePointsToHits);
365  LArPandoraHelper::CollectTracks(evt, m_trackLabel, trackVector, particlesToTracks);
366  LArPandoraHelper::CollectTracks(evt, m_trackLabel, trackVectorExtra, tracksToHits);
367  LArPandoraHelper::CollectShowers(evt, m_showerLabel, showerVector, particlesToShowers);
368  LArPandoraHelper::CollectShowers(evt, m_showerLabel, showerVectorExtra, showersToHits);
370  evt, m_particleLabel, particleVector, particlesToSpacePoints);
374  particlesToHits,
375  hitsToParticles,
376  LArPandoraHelper::DaughterMode::kUseDaughters,
377  false);
379  evt, m_particleLabel, m_clusterLabel, particlesToHitsClusters, hitsToParticlesClusters);
380 
382 
383  if (m_printDebug) std::cout << " PFParticles: " << particleVector.size() << std::endl;
384 
385  // Loop over Tracks (Fill 3D Track Tree)
386  // =====================================
387  if (m_printDebug) std::cout << " PFParticleHitDumper::FillRecoTracks(...) " << std::endl;
388  this->FillRecoTracks(particlesToTracks);
389 
390  // Loop over PFParticles (Fill 3D Reco Tree)
391  // =========================================
392  if (m_printDebug) std::cout << " PFParticleHitDumper::FillReco3D(...) " << std::endl;
393  this->FillReco3D(particleVector, particlesToSpacePoints, spacePointsToHits);
394 
395  // Loop over Hits (Fill 2D Reco Tree)
396  // ==================================
397  if (m_printDebug) std::cout << " PFParticleHitDumper::FillReco2D(...) " << std::endl;
398  this->FillReco2D(evt, hitVector, hitsToParticles);
399 
400  // Loop over Hits (Fill Associated 2D Hits Tree)
401  // =============================================
402  if (m_printDebug)
403  std::cout << " PFParticleHitDumper::FillAssociated2DHits(...) " << std::endl;
404  this->FillAssociated2DHits(evt,
405  particleVector,
406  particlesToHits,
407  particlesToHitsClusters,
408  particlesToTracks,
409  tracksToHits,
410  particlesToShowers,
411  showersToHits);
412 
413  // Loop over Wires (Fill Reco Wire Tree)
414  // =====================================
415  if (m_printDebug) std::cout << " PFParticleHitDumper::FillRecoWires(...) " << std::endl;
416  this->FillRecoWires(evt, wireVector);
417  }
418 
419  //------------------------------------------------------------------------------------------------------------------------------------------
420 
421  void
423  {
424  // Initialise variables
425  m_particle = -1;
426  m_x = 0.0;
427  m_y = 0.0;
428  m_z = 0.0;
429 
430  // Create dummy entry if there are no particles
431  if (particlesToTracks.empty()) { m_pRecoTracks->Fill(); }
432 
433  // Loop over tracks
434  for (PFParticlesToTracks::const_iterator iter = particlesToTracks.begin(),
435  iterEnd = particlesToTracks.end();
436  iter != iterEnd;
437  ++iter) {
438  const art::Ptr<recob::PFParticle> particle = iter->first;
439  const TrackVector& trackVector = iter->second;
440 
441  m_particle = particle->Self();
442 
443  if (!trackVector.empty()) {
444  if (trackVector.size() != 1 && m_printDebug)
445  std::cout << " Warning: Found particle with more than one associated track " << std::endl;
446 
447  const art::Ptr<recob::Track> track = *(trackVector.begin());
448 
449  if (m_printDebug)
450  std::cout << " PFPARTICLE [" << m_particle << "] (" << track->NumberTrajectoryPoints()
451  << " Trajectory Points)" << std::endl;
452 
453  for (unsigned int p = 0; p < track->NumberTrajectoryPoints(); ++p) {
454  const auto position(track->LocationAtPoint(p));
455  m_x = position.x();
456  m_y = position.y();
457  m_z = position.z();
458 
459  m_pRecoTracks->Fill();
460  }
461  }
462  }
463  }
464 
465  //------------------------------------------------------------------------------------------------------------------------------------------
466 
467  void
469  const PFParticlesToSpacePoints& particlesToSpacePoints,
470  const SpacePointsToHits& spacePointsToHits)
471  {
472  // Initialise variables
473  m_particle = -1;
474  m_primary = 0;
475  m_pdgcode = 0;
476  m_cstat = 0;
477  m_tpc = 0;
478  m_plane = 0;
479  m_x = 0.0;
480  m_u = 0.0;
481  m_v = 0.0;
482  m_y = 0.0;
483  m_z = 0.0;
484 
485  // Create dummy entry if there are no particles
486  if (particleVector.empty()) { m_pReco3D->Fill(); }
487 
488  // Store associations between particle and particle ID
489  PFParticleMap theParticleMap;
490 
491  for (unsigned int i = 0; i < particleVector.size(); ++i) {
492  const art::Ptr<recob::PFParticle> particle = particleVector.at(i);
493  theParticleMap[particle->Self()] = particle;
494  }
495 
496  // Loop over particles
497  for (PFParticlesToSpacePoints::const_iterator iter1 = particlesToSpacePoints.begin(),
498  iterEnd1 = particlesToSpacePoints.end();
499  iter1 != iterEnd1;
500  ++iter1) {
501  const art::Ptr<recob::PFParticle> particle = iter1->first;
502  const SpacePointVector& spacepoints = iter1->second;
503 
504  m_particle = particle->Self();
505  m_pdgcode = particle->PdgCode();
506  m_primary = 0;
507 
508  if (particle->IsPrimary()) { m_primary = 1; }
509  else {
510  const size_t parentID(particle->Parent());
511  PFParticleMap::const_iterator pIter = theParticleMap.find(parentID);
512  if (theParticleMap.end() == pIter)
513  throw cet::exception("LArPandora")
514  << " PFParticleHitDumper::analyze --- Found particle with ID code";
515 
516  const art::Ptr<recob::PFParticle> particleParent = pIter->second;
517  if (LArPandoraHelper::IsNeutrino(particleParent)) m_primary = 1;
518  }
519 
520  if (m_printDebug)
521  std::cout << " PFPARTICLE [" << m_particle << "] [Primary=" << m_primary << "] ("
522  << spacepoints.size() << " Space Points)" << std::endl;
523 
524  for (unsigned int j = 0; j < spacepoints.size(); ++j) {
525  const art::Ptr<recob::SpacePoint> spacepoint = spacepoints.at(j);
526 
527  m_x = spacepoint->XYZ()[0];
528  m_y = spacepoint->XYZ()[1];
529  m_z = spacepoint->XYZ()[2];
530 
531  SpacePointsToHits::const_iterator iter2 = spacePointsToHits.find(spacepoint);
532  if (spacePointsToHits.end() == iter2)
533  throw cet::exception("LArPandora")
534  << " PFParticleHitDumper::analyze --- Found space point without associated hit";
535 
536  const art::Ptr<recob::Hit> hit = iter2->second;
537  const geo::WireID& wireID(hit->WireID());
538 
539  m_cstat = wireID.Cryostat;
540  m_tpc = wireID.TPC;
541  m_plane = wireID.Plane;
542 
543  m_u = this->YZtoU(m_cstat, m_tpc, m_y, m_z);
544  m_v = this->YZtoV(m_cstat, m_tpc, m_y, m_z);
545 
546  m_pReco3D->Fill();
547  }
548  }
549  }
550 
551  //------------------------------------------------------------------------------------------------------------------------------------------
552 
553  void
555  const PFParticleVector& particleVector,
556  const PFParticlesToHits& particlesToHits,
557  const PFParticlesToHits& particlesToHitsClusters,
558  const PFParticlesToTracks& particlesToTracks,
559  const TracksToHits& tracksToHits,
560  const PFParticlesToShowers& particlesToShowers,
561  const ShowersToHits& showersToHits)
562  {
563  // Create dummy entry if there are no 2D hits
564  if (particleVector.empty()) { m_pRecoComparison->Fill(); }
565 
566  for (unsigned int i = 0; i < particleVector.size(); ++i) {
567  //initialise variables
568  m_particle = -1;
569  m_pdgcode = 0;
571  m_hitsFromClusters = 0;
573 
574  const art::Ptr<recob::PFParticle> particle = particleVector.at(i);
575  m_particle = particle->Self();
576  m_pdgcode = particle->PdgCode();
577 
578  PFParticlesToHits::const_iterator pIter = particlesToHits.find(particle);
579  PFParticlesToHits::const_iterator pIter2 = particlesToHitsClusters.find(particle);
580  if (particlesToHits.end() != pIter) m_hitsFromSpacePoints = pIter->second.size();
581  if (particlesToHitsClusters.end() != pIter2) m_hitsFromClusters = pIter2->second.size();
582 
583  if (m_pdgcode == 13) {
584  PFParticlesToTracks::const_iterator iter = particlesToTracks.find(particle);
585  const art::Ptr<recob::Track> track = *(iter->second.begin());
586 
587  TracksToHits::const_iterator iter2 = tracksToHits.find(track);
588  if (tracksToHits.end() != iter2) {
589  const HitVector& hitVector = iter2->second;
590  m_hitsFromTrackOrShower = hitVector.size();
591  }
592  }
593  else if (m_pdgcode == 11) {
594  PFParticlesToShowers::const_iterator iter = particlesToShowers.find(particle);
595  const art::Ptr<recob::Shower> shower = *(iter->second.begin());
596 
597  ShowersToHits::const_iterator iter2 = showersToHits.find(shower);
598  if (showersToHits.end() != iter2) {
599  const HitVector& hitVector = iter2->second;
600  m_hitsFromTrackOrShower = hitVector.size();
601  }
602  }
603 
604  if (m_printDebug)
605  std::cout << " PFParticle " << m_particle << " (PDG " << m_pdgcode << ") has "
606  << m_hitsFromSpacePoints << " hits from space points and " << m_hitsFromClusters
607  << " hits from clusters, and its recob::Track/Shower has "
608  << m_hitsFromTrackOrShower << " associated hits " << std::endl;
609 
610  m_pRecoComparison->Fill();
611  }
612  }
613 
614  //------------------------------------------------------------------------------------------------------------------------------------------
615 
616  void
618  const HitVector& hitVector,
619  const HitsToPFParticles& hitsToParticles)
620  {
621  // Initialise variables
622  m_particle = -1;
623  m_pdgcode = 0;
624  m_cstat = 0;
625  m_tpc = 0;
626  m_plane = 0;
627  m_wire = 0;
628  m_x = 0.0;
629  m_w = 0.0;
630  m_q = 0.0;
631 
632  // Create dummy entry if there are no 2D hits
633  if (hitVector.empty()) { m_pReco2D->Fill(); }
634 
635  // Need DetectorProperties service to convert from ticks to X
636  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(e);
637 
638  // Loop over 2D hits
639  for (unsigned int i = 0; i < hitVector.size(); ++i) {
640  const art::Ptr<recob::Hit> hit = hitVector.at(i);
641 
642  m_particle = -1;
643  m_pdgcode = 0;
644 
645  HitsToPFParticles::const_iterator pIter = hitsToParticles.find(hit);
646  if (hitsToParticles.end() != pIter) {
647  const art::Ptr<recob::PFParticle> particle = pIter->second;
648  m_particle = particle->Self();
649  m_pdgcode = particle->PdgCode();
650  }
651 
652  const geo::WireID& wireID(hit->WireID());
653  m_cstat = wireID.Cryostat;
654  m_tpc = wireID.TPC;
655  m_plane = wireID.Plane;
656  m_wire = wireID.Wire;
657 
658  m_q = hit->Integral();
659  m_x = detProp.ConvertTicksToX(hit->PeakTime(), wireID.Plane, wireID.TPC, wireID.Cryostat);
660  m_w = this->GetUVW(wireID);
661 
662  m_pReco2D->Fill();
663  }
664  }
665 
666  //------------------------------------------------------------------------------------------------------------------------------------------
667 
668  void
669  PFParticleHitDumper::FillRecoWires(const art::Event& e, const WireVector& wireVector)
670  {
671 
672  // Create dummy entry if there are no wires
673  if (wireVector.empty()) { m_pRecoWire->Fill(); }
674 
675  // Need geometry service to convert channel to wire ID
676  art::ServiceHandle<geo::Geometry const> theGeometry;
677 
678  // Need DetectorProperties service to convert from ticks to X
679  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(e);
680 
681  // Loop over wires
682  int signalCounter(0);
683 
684  for (unsigned int i = 0; i < wireVector.size(); ++i) {
685  const art::Ptr<recob::Wire> wire = wireVector.at(i);
686 
687  const std::vector<float>& signals(wire->Signal());
688  const std::vector<geo::WireID> wireIds = theGeometry->ChannelToWire(wire->Channel());
689 
690  if ((signalCounter++) < 10 && m_printDebug)
691  std::cout << " numWires=" << wireVector.size() << " numSignals=" << signals.size()
692  << std::endl;
693 
694  double time(0.0);
695 
696  m_q = 0.0;
697 
698  for (std::vector<float>::const_iterator tIter = signals.begin(), tIterEnd = signals.end();
699  tIter != tIterEnd;
700  ++tIter) {
701  time += 1.0;
702  m_q = *tIter;
703 
704  if (m_q < 2.0) // seems to remove most noise
705  continue;
706 
707  for (std::vector<geo::WireID>::const_iterator wIter = wireIds.begin(),
708  wIterEnd = wireIds.end();
709  wIter != wIterEnd;
710  ++wIter) {
711  const geo::WireID& wireID = *wIter;
712  m_cstat = wireID.Cryostat;
713  m_tpc = wireID.TPC;
714  m_plane = wireID.Plane;
715  m_wire = wireID.Wire;
716 
717  m_x = detProp.ConvertTicksToX(time, wireID.Plane, wireID.TPC, wireID.Cryostat);
718  m_w = this->GetUVW(wireID);
719 
720  m_pRecoWire->Fill();
721  }
722  }
723  }
724  }
725 
726  //------------------------------------------------------------------------------------------------------------------------------------------
727 
728  double
730  {
731  // define UVW as closest distance from (0,0) to wire axis
732  art::ServiceHandle<geo::Geometry const> theGeometry;
733 
734  double xyzStart[3];
735  theGeometry->Cryostat(wireID.Cryostat)
736  .TPC(wireID.TPC)
737  .Plane(wireID.Plane)
738  .Wire(wireID.Wire)
739  .GetStart(xyzStart);
740  const double ay(xyzStart[1]);
741  const double az(xyzStart[2]);
742 
743  double xyzEnd[3];
744  theGeometry->Cryostat(wireID.Cryostat)
745  .TPC(wireID.TPC)
746  .Plane(wireID.Plane)
747  .Wire(wireID.Wire)
748  .GetEnd(xyzEnd);
749  const double by(xyzEnd[1]);
750  const double bz(xyzEnd[2]);
751 
752  const double ny(by - ay);
753  const double nz(bz - az);
754  const double N2(ny * ny + nz * nz);
755 
756  const double ry(ay - (ay * ny + az * nz) * ny / N2);
757  const double rz(az - (ay * ny + az * nz) * nz / N2);
758  const double sign((rz > 0.0) ? +1.0 : -1.0);
759 
760  return sign * std::sqrt(ry * ry + rz * rz);
761  }
762 
763  //------------------------------------------------------------------------------------------------------------------------------------------
764 
765  double
766  PFParticleHitDumper::YZtoU(const unsigned int cstat,
767  const unsigned int tpc,
768  const double y,
769  const double z) const
770  {
771  // TODO: Check that this stills works in DUNE
772  art::ServiceHandle<geo::Geometry const> theGeometry;
773  const double m_theta(theGeometry->WireAngleToVertical(geo::kU, tpc, cstat));
774  return z * std::sin(m_theta) - y * std::cos(m_theta);
775  }
776 
777  //------------------------------------------------------------------------------------------------------------------------------------------
778 
779  double
780  PFParticleHitDumper::YZtoV(const unsigned int cstat,
781  const unsigned int tpc,
782  const double y,
783  const double z) const
784  {
785  // TODO; Check that this still works in DUNE
786  art::ServiceHandle<geo::Geometry const> theGeometry;
787  const double m_theta(theGeometry->WireAngleToVertical(geo::kV, tpc, cstat));
788  return z * std::sin(m_theta) - y * std::cos(m_theta);
789  }
790 
791  //------------------------------------------------------------------------------------------------------------------------------------------
792 
793 } //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.
process_name opflash particleana ie ie ie z
double YZtoU(const unsigned int cstat, const unsigned int tpc, const double y, const double z) const
Convert from (Y,Z) to U coordinate.
Encapsulate the construction of a single cyostat.
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:130
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.
pdgs p
Definition: selectors.fcl:22
static void CollectWires(const art::Event &evt, const std::string &label, WireVector &wireVector)
Collect the reconstructed wires from the ART event record.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:212
std::vector< art::Ptr< recob::Shower > > ShowerVector
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:580
process_name use argoneut_mc_hitfinder track
void FillRecoWires(const art::Event &event, const WireVector &wireVector)
Store raw data.
process_name hit
Definition: cheaterreco.fcl:51
process_name shower
Definition: cheaterreco.fcl:51
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:129
bool m_printDebug
switch for print statements (TODO: use message service!)
std::map< art::Ptr< recob::PFParticle >, ShowerVector > PFParticlesToShowers
process_name opflash particleana ie ie y
void reconfigure(fhicl::ParameterSet const &pset)
std::vector< art::Ptr< recob::PFParticle > > PFParticleVector
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.
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.
std::vector< art::Ptr< recob::Track > > TrackVector
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
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.
Provides recob::Track data product.
static void CollectHits(const art::Event &evt, const std::string &label, HitVector &hitVector)
Collect the reconstructed Hits from the ART event record.
void FillRecoTracks(const PFParticlesToTracks &particlesToTracks)
Store 3D track hits.
int sign(double val)
Definition: UtilFunc.cxx:104
std::vector< art::Ptr< recob::Hit > > HitVector
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
do i e
Declaration of basic channel signal object.
art::ServiceHandle< art::TFileService > tfs
TCEvent evt
Definition: DataStructs.cxx:8
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
helper function for LArPandoraInterface producer module
art framework interface to geometry description
BEGIN_PROLOG could also be cout
auto const detProp
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > SpacePointsToHits
Encapsulate the construction of a single detector plane.