Pandora
Pandora source code navigator
Loading...
Searching...
No Matches
PFParticleHitDumper_module.cc
Go to the documentation of this file.
1
8#include "art/Framework/Core/EDAnalyzer.h"
9#include "art/Framework/Core/ModuleMacros.h"
10
11#include "TTree.h"
12
13#include "larcoreobj/SimpleTypesAndConstants/geo_types.h"
14#include "lardata/DetectorInfoServices/DetectorClocksService.h"
16
17#include <string>
18
19//------------------------------------------------------------------------------------------------------------------------------------------
20
21namespace lar_pandora {
22
26 class PFParticleHitDumper : public art::EDAnalyzer {
27 public:
33 PFParticleHitDumper(fhicl::ParameterSet const& pset);
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
133 TTree* m_pReco3D;
134 TTree* m_pReco2D;
136 TTree* m_pRecoWire;
137
138 int m_run;
143
145 int m_tpc;
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
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
188#include "larcore/Geometry/Geometry.h"
189#include "larcorealg/Geometry/CryostatGeo.h"
190#include "larcorealg/Geometry/PlaneGeo.h"
191#include "larcorealg/Geometry/TPCGeo.h"
192#include "lardata/DetectorInfoServices/DetectorPropertiesService.h"
193#include "lardata/DetectorInfoServices/LArPropertiesService.h"
194#include "lardata/Utilities/AssociationUtil.h"
195#include "lardataobj/RecoBase/Cluster.h"
196#include "lardataobj/RecoBase/Hit.h"
197#include "lardataobj/RecoBase/PFParticle.h"
198#include "lardataobj/RecoBase/Shower.h"
199#include "lardataobj/RecoBase/SpacePoint.h"
200#include "lardataobj/RecoBase/Track.h"
201#include "lardataobj/RecoBase/Wire.h"
202#include "nusimdata/SimulationBase/MCParticle.h"
203#include "nusimdata/SimulationBase/MCTruth.h"
204
205#include <iostream>
206
207namespace 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 PFParticleHitDumper::reconfigure(fhicl::ParameterSet const& pset)
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 //
240 art::ServiceHandle<art::TFileService const> tfs;
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
307 void PFParticleHitDumper::analyze(const art::Event& evt)
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
337 art::ServiceHandle<geo::Geometry const> theGeometry;
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,
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;
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
608 void PFParticleHitDumper::FillReco2D(const art::Event& e,
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
659 void PFParticleHitDumper::FillRecoWires(const art::Event& e, const WireVector& wireVector)
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
666 art::ServiceHandle<geo::Geometry const> theGeometry;
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
721 art::ServiceHandle<geo::Geometry const> theGeometry;
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
750 art::ServiceHandle<geo::Geometry const> theGeometry;
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
763 art::ServiceHandle<geo::Geometry const> theGeometry;
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
helper function for LArPandoraInterface producer module
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.
static void CollectHits(const art::Event &evt, const std::string &label, HitVector &hitVector)
Collect the reconstructed Hits from the ART event record.
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.
static void CollectWires(const art::Event &evt, const std::string &label, WireVector &wireVector)
Collect the reconstructed wires 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.
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 bool IsNeutrino(const art::Ptr< recob::PFParticle > particle)
Determine whether a particle has been reconstructed as a neutrino.
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 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.
void FillRecoTracks(const PFParticlesToTracks &particlesToTracks)
Store 3D track hits.
void FillReco2D(const art::Event &event, const HitVector &hitVector, const HitsToPFParticles &hitsToParticles)
Store 2D hits.
void FillRecoWires(const art::Event &event, const WireVector &wireVector)
Store raw data.
void FillReco3D(const PFParticleVector &particleVector, const PFParticlesToSpacePoints &particlesToSpacePoints, const SpacePointsToHits &spacePointsToHits)
Store 3D hits.
PFParticleHitDumper(fhicl::ParameterSet const &pset)
Constructor.
bool m_printDebug
switch for print statements (TODO: use message service!)
double YZtoU(const unsigned int cstat, const unsigned int tpc, const double y, const double z) const
Convert from (Y,Z) to U coordinate.
void reconfigure(fhicl::ParameterSet const &pset)
double YZtoV(const unsigned int cstat, const unsigned int tpc, const double y, const double z) const
Convert from (Y,Z) to V coordinate.
double GetUVW(const geo::WireID &wireID) const
Conversion from wire ID to U/V/W coordinate.
std::map< int, art::Ptr< recob::PFParticle > > PFParticleMap
std::vector< art::Ptr< recob::Track > > TrackVector
std::vector< art::Ptr< recob::Shower > > ShowerVector
std::map< art::Ptr< recob::PFParticle >, ShowerVector > PFParticlesToShowers
std::map< art::Ptr< recob::PFParticle >, SpacePointVector > PFParticlesToSpacePoints
std::map< art::Ptr< recob::PFParticle >, HitVector > PFParticlesToHits
std::vector< art::Ptr< recob::Hit > > HitVector
std::vector< art::Ptr< recob::PFParticle > > PFParticleVector
std::map< art::Ptr< recob::Shower >, HitVector > ShowersToHits
std::map< art::Ptr< recob::PFParticle >, TrackVector > PFParticlesToTracks
std::map< art::Ptr< recob::Track >, HitVector > TracksToHits
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > SpacePointsToHits
std::vector< art::Ptr< recob::SpacePoint > > SpacePointVector
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::PFParticle > > HitsToPFParticles
std::vector< art::Ptr< recob::Wire > > WireVector