Pandora
Pandora source code navigator
Loading...
Searching...
No Matches
PFParticleMonitoring_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
14
15#include <string>
16
17//------------------------------------------------------------------------------------------------------------------------------------------
18
19namespace lar_pandora {
20
24 class PFParticleMonitoring : public art::EDAnalyzer {
25 public:
31 PFParticleMonitoring(fhicl::ParameterSet const& pset);
32
36 virtual ~PFParticleMonitoring();
37
38 void beginJob();
39 void endJob();
40 void analyze(const art::Event& evt);
41 void reconfigure(fhicl::ParameterSet const& pset);
42
43 private:
44 typedef std::set<art::Ptr<recob::PFParticle>> PFParticleSet;
45 typedef std::set<art::Ptr<simb::MCParticle>> MCParticleSet;
46 typedef std::set<art::Ptr<simb::MCTruth>> MCTruthSet;
47
56 void BuildTrueNeutrinoHitMaps(const MCTruthToMCParticles& truthToParticles,
57 const MCParticlesToHits& trueParticlesToHits,
58 MCTruthToHits& trueNeutrinosToHits,
59 HitsToMCTruth& trueHitsToNeutrinos) const;
60
69 void BuildRecoNeutrinoHitMaps(const PFParticleMap& recoParticleMap,
70 const PFParticlesToHits& recoParticlesToHits,
71 PFParticlesToHits& recoNeutrinosToHits,
72 HitsToPFParticles& recoHitsToNeutrinos) const;
73
82 void GetRecoToTrueMatches(const PFParticlesToHits& recoNeutrinosToHits,
83 const HitsToMCTruth& trueHitsToNeutrinos,
84 MCTruthToPFParticles& matchedNeutrinos,
85 MCTruthToHits& matchedNeutrinoHits) const;
86
97 void GetRecoToTrueMatches(const PFParticlesToHits& recoNeutrinosToHits,
98 const HitsToMCTruth& trueHitsToNeutrinos,
99 MCTruthToPFParticles& matchedNeutrinos,
100 MCTruthToHits& matchedNeutrinoHits,
101 PFParticleSet& recoVeto,
102 MCTruthSet& trueVeto) const;
103
112 void GetRecoToTrueMatches(const PFParticlesToHits& recoParticlesToHits,
113 const HitsToMCParticles& trueHitsToParticles,
114 MCParticlesToPFParticles& matchedParticles,
115 MCParticlesToHits& matchedHits) const;
116
127 void GetRecoToTrueMatches(const PFParticlesToHits& recoParticlesToHits,
128 const HitsToMCParticles& trueHitsToParticles,
129 MCParticlesToPFParticles& matchedParticles,
130 MCParticlesToHits& matchedHits,
131 PFParticleSet& recoVeto,
132 MCParticleSet& trueVeto) const;
133
140 int CountHitsByType(const int view, const HitVector& hitVector) const;
141
149 void GetStartAndEndPoints(const art::Ptr<simb::MCParticle> trueParticle,
150 int& startT,
151 int& endT) const;
152
160 double GetLength(const art::Ptr<simb::MCParticle> trueParticle,
161 const int startT,
162 const int endT) const;
163
164 TTree* m_pRecoTree;
165
166 int m_run;
169
174
183
191
194 double m_pfoVtxX;
195 double m_pfoVtxY;
196 double m_pfoVtxZ;
197 double m_pfoEndX;
198 double m_pfoEndY;
199 double m_pfoEndZ;
200 double m_pfoDirX;
201 double m_pfoDirY;
202 double m_pfoDirZ;
203 double m_pfoLength;
205
207 double m_mcVtxX;
208 double m_mcVtxY;
209 double m_mcVtxZ;
210 double m_mcEndX;
211 double m_mcEndY;
212 double m_mcEndZ;
213 double m_mcDirX;
214 double m_mcDirY;
215 double m_mcDirZ;
216 double m_mcEnergy;
217 double m_mcLength;
219
221 double m_purity;
222
226
230
234
238
239 int
241 int
243
246
247 std::string m_hitfinderLabel;
248 std::string m_trackLabel;
249 std::string m_particleLabel;
250 std::string m_backtrackerLabel;
251 std::string m_geantModuleLabel;
252
257
260 bool
262 };
263
264 DEFINE_ART_MODULE(PFParticleMonitoring)
265
266} // namespace lar_pandora
267
268//------------------------------------------------------------------------------------------------------------------------------------------
269// implementation follows
270
271#include "art/Framework/Principal/Event.h"
272#include "art/Framework/Principal/Handle.h"
273#include "art/Framework/Services/Registry/ServiceHandle.h"
274#include "art_root_io/TFileDirectory.h"
275#include "art_root_io/TFileService.h"
276#include "fhiclcpp/ParameterSet.h"
277#include "messagefacility/MessageLogger/MessageLogger.h"
278
279#include "larcore/Geometry/Geometry.h"
280#include "larcorealg/Geometry/CryostatGeo.h"
281#include "larcorealg/Geometry/PlaneGeo.h"
282#include "larcorealg/Geometry/TPCGeo.h"
283#include "lardata/DetectorInfoServices/DetectorPropertiesService.h"
284#include "lardata/DetectorInfoServices/LArPropertiesService.h"
285#include "lardata/Utilities/AssociationUtil.h"
286#include "lardataobj/AnalysisBase/T0.h"
287#include "lardataobj/RecoBase/Cluster.h"
288#include "lardataobj/RecoBase/Hit.h"
289#include "lardataobj/RecoBase/PFParticle.h"
290#include "lardataobj/RecoBase/SpacePoint.h"
291#include "lardataobj/RecoBase/Track.h"
292#include "lardataobj/RecoBase/Vertex.h"
293#include "nusimdata/SimulationBase/MCParticle.h"
294#include "nusimdata/SimulationBase/MCTruth.h"
295
296#include <iostream>
297
298namespace lar_pandora {
299
300 PFParticleMonitoring::PFParticleMonitoring(fhicl::ParameterSet const& pset)
301 : art::EDAnalyzer(pset)
302 {
303 this->reconfigure(pset);
304 }
305
306 //------------------------------------------------------------------------------------------------------------------------------------------
307
309
310 //------------------------------------------------------------------------------------------------------------------------------------------
311
312 void PFParticleMonitoring::reconfigure(fhicl::ParameterSet const& pset)
313 {
314 m_trackLabel = pset.get<std::string>("TrackModule", "pandoraTracks");
315 m_particleLabel = pset.get<std::string>("PFParticleModule", "pandora");
316 m_hitfinderLabel = pset.get<std::string>("HitFinderModule", "gaushit");
317 m_backtrackerLabel = pset.get<std::string>("BackTrackerModule", "gaushitTruthMatch");
318 m_geantModuleLabel = pset.get<std::string>("GeantModule", "largeant");
319
320 m_useDaughterPFParticles = pset.get<bool>("UseDaughterPFParticles", false);
321 m_useDaughterMCParticles = pset.get<bool>("UseDaughterMCParticles", true);
322 m_addDaughterPFParticles = pset.get<bool>("AddDaughterPFParticles", true);
323 m_addDaughterMCParticles = pset.get<bool>("AddDaughterMCParticles", true);
324
325 m_recursiveMatching = pset.get<bool>("RecursiveMatching", false);
326 m_printDebug = pset.get<bool>("PrintDebug", false);
327 m_disableRealDataCheck = pset.get<bool>("DisableRealDataCheck", false);
328 }
329
330 //------------------------------------------------------------------------------------------------------------------------------------------
331
333 {
334 mf::LogDebug("LArPandora") << " *** PFParticleMonitoring::beginJob() *** " << std::endl;
335
336 //
337 art::ServiceHandle<art::TFileService const> tfs;
338
339 m_pRecoTree = tfs->make<TTree>("pandora", "LAr Reco vs True");
340 m_pRecoTree->Branch("run", &m_run, "run/I");
341 m_pRecoTree->Branch("event", &m_event, "event/I");
342 m_pRecoTree->Branch("index", &m_index, "index/I");
343 m_pRecoTree->Branch("nMCParticles", &m_nMCParticles, "nMCParticles/I");
344 m_pRecoTree->Branch("nNeutrinoPfos", &m_nNeutrinoPfos, "nNeutrinoPfos/I");
345 m_pRecoTree->Branch("nPrimaryPfos", &m_nPrimaryPfos, "nPrimaryPfos/I");
346 m_pRecoTree->Branch("nDaughterPfos", &m_nDaughterPfos, "nDaughterPfos/I");
347 m_pRecoTree->Branch("mcPdg", &m_mcPdg, "mcPdg/I");
348 m_pRecoTree->Branch("mcNuPdg", &m_mcNuPdg, "mcNuPdg/I");
349 m_pRecoTree->Branch("mcParentPdg", &m_mcParentPdg, "mcParentPdg/I");
350 m_pRecoTree->Branch("mcPrimaryPdg", &m_mcPrimaryPdg, "mcPrimaryPdg/I");
351 m_pRecoTree->Branch("mcIsNeutrino", &m_mcIsNeutrino, "mcIsNeutrino/I");
352 m_pRecoTree->Branch("mcIsPrimary", &m_mcIsPrimary, "mcIsPrimary/I");
353 m_pRecoTree->Branch("mcIsDecay", &m_mcIsDecay, "mcIsDecay/I");
354 m_pRecoTree->Branch("mcIsCC", &m_mcIsCC, "mcIsCC/I");
355 m_pRecoTree->Branch("pfoPdg", &m_pfoPdg, "pfoPdg/I");
356 m_pRecoTree->Branch("pfoNuPdg", &m_pfoNuPdg, "pfoNuPdg/I");
357 m_pRecoTree->Branch("pfoParentPdg", &m_pfoParentPdg, "pfoParentPdg/I");
358 m_pRecoTree->Branch("pfoPrimaryPdg", &m_pfoPrimaryPdg, "pfoPrimaryPdg/I");
359 m_pRecoTree->Branch("pfoIsNeutrino", &m_pfoIsNeutrino, "pfoIsNeutrino/I");
360 m_pRecoTree->Branch("pfoIsPrimary", &m_pfoIsPrimary, "pfoIsPrimary/I");
361 m_pRecoTree->Branch("pfoIsStitched", &m_pfoIsStitched, "pfoIsStitched/I");
362 m_pRecoTree->Branch("pfoTrack", &m_pfoTrack, "pfoTrack/I");
363 m_pRecoTree->Branch("pfoVertex", &m_pfoVertex, "pfoVertex/I");
364 m_pRecoTree->Branch("pfoVtxX", &m_pfoVtxX, "pfoVtxX/D");
365 m_pRecoTree->Branch("pfoVtxY", &m_pfoVtxY, "pfoVtxY/D");
366 m_pRecoTree->Branch("pfoVtxZ", &m_pfoVtxZ, "pfoVtxZ/D");
367 m_pRecoTree->Branch("pfoEndX", &m_pfoEndX, "pfoEndX/D");
368 m_pRecoTree->Branch("pfoEndY", &m_pfoEndY, "pfoEndY/D");
369 m_pRecoTree->Branch("pfoEndZ", &m_pfoEndZ, "pfoEndZ/D");
370 m_pRecoTree->Branch("pfoDirX", &m_pfoDirX, "pfoDirX/D");
371 m_pRecoTree->Branch("pfoDirY", &m_pfoDirY, "pfoDirY/D");
372 m_pRecoTree->Branch("pfoDirZ", &m_pfoDirZ, "pfoDirZ/D");
373 m_pRecoTree->Branch("pfoLength", &m_pfoLength, "pfoLength/D");
374 m_pRecoTree->Branch("pfoStraightLength", &m_pfoStraightLength, "pfoStraightLength/D");
375 m_pRecoTree->Branch("mcVertex", &m_mcVertex, "mcVertex/I");
376 m_pRecoTree->Branch("mcVtxX", &m_mcVtxX, "mcVtxX/D");
377 m_pRecoTree->Branch("mcVtxY", &m_mcVtxY, "mcVtxY/D");
378 m_pRecoTree->Branch("mcVtxZ", &m_mcVtxZ, "mcVtxZ/D");
379 m_pRecoTree->Branch("mcEndX", &m_mcEndX, "mcEndX/D");
380 m_pRecoTree->Branch("mcEndY", &m_mcEndY, "mcEndY/D");
381 m_pRecoTree->Branch("mcEndZ", &m_mcEndZ, "mcEndZ/D");
382 m_pRecoTree->Branch("mcDirX", &m_mcDirX, "mcDirX/D");
383 m_pRecoTree->Branch("mcDirY", &m_mcDirY, "mcDirY/D");
384 m_pRecoTree->Branch("mcDirZ", &m_mcDirZ, "mcDirZ/D");
385 m_pRecoTree->Branch("mcEnergy", &m_mcEnergy, "mcEnergy/D");
386 m_pRecoTree->Branch("mcLength", &m_mcLength, "mcLength/D");
387 m_pRecoTree->Branch("mcStraightLength", &m_mcStraightLength, "mcStraightLength/D");
388 m_pRecoTree->Branch("completeness", &m_completeness, "completeness/D");
389 m_pRecoTree->Branch("purity", &m_purity, "purity/D");
390 m_pRecoTree->Branch("nMCHits", &m_nMCHits, "nMCHits/I");
391 m_pRecoTree->Branch("nPfoHits", &m_nPfoHits, "nPfoHits/I");
392 m_pRecoTree->Branch("nMatchedHits", &m_nMatchedHits, "nMatchedHits/I");
393 m_pRecoTree->Branch("nMCHitsU", &m_nMCHitsU, "nMCHitsU/I");
394 m_pRecoTree->Branch("nMCHitsV", &m_nMCHitsV, "nMCHitsV/I");
395 m_pRecoTree->Branch("nMCHitsW", &m_nMCHitsW, "nMCHitsW/I");
396 m_pRecoTree->Branch("nPfoHitsU", &m_nPfoHitsU, "nPfoHitsU/I");
397 m_pRecoTree->Branch("nPfoHitsV", &m_nPfoHitsV, "nPfoHitsV/I");
398 m_pRecoTree->Branch("nPfoHitsW", &m_nPfoHitsW, "nPfoHitsW/I");
399 m_pRecoTree->Branch("nMatchedHitsU", &m_nMatchedHitsU, "nMatchedHitsU/I");
400 m_pRecoTree->Branch("nMatchedHitsV", &m_nMatchedHitsV, "nMatchedHitsV/I");
401 m_pRecoTree->Branch("nMatchedHitsW", &m_nMatchedHitsW, "nMatchedHitsW/I");
402 m_pRecoTree->Branch("nTrueWithoutRecoHits", &m_nTrueWithoutRecoHits, "nTrueWithoutRecoHits/I");
403 m_pRecoTree->Branch("nRecoWithoutTrueHits", &m_nRecoWithoutTrueHits, "nRecoWithoutTrueHits/I");
404 m_pRecoTree->Branch("spacepointsMinX", &m_spacepointsMinX, "spacepointsMinX/D");
405 m_pRecoTree->Branch("spacepointsMaxX", &m_spacepointsMaxX, "spacepointsMaxX/D");
406 }
407
408 //------------------------------------------------------------------------------------------------------------------------------------------
409
411
412 //------------------------------------------------------------------------------------------------------------------------------------------
413
414 void PFParticleMonitoring::analyze(const art::Event& evt)
415 {
416 if (m_printDebug) std::cout << " *** PFParticleMonitoring::analyze(...) *** " << std::endl;
417
418 m_run = evt.run();
419 m_event = evt.id().event();
420 m_index = 0;
421
422 m_nMCParticles = 0;
423 m_nNeutrinoPfos = 0;
424 m_nPrimaryPfos = 0;
425 m_nDaughterPfos = 0;
426
427 m_mcPdg = 0;
428 m_mcNuPdg = 0;
429 m_mcParentPdg = 0;
430 m_mcPrimaryPdg = 0;
431 m_mcIsNeutrino = 0;
432 m_mcIsPrimary = 0;
433 m_mcIsDecay = 0;
434 m_mcIsCC = 0;
435
436 m_pfoPdg = 0;
437 m_pfoNuPdg = 0;
438 m_pfoParentPdg = 0;
439 m_pfoPrimaryPdg = 0;
440 m_pfoIsNeutrino = 0;
441 m_pfoIsPrimary = 0;
442 m_pfoIsStitched = 0;
443 m_pfoTrack = 0;
444 m_pfoVertex = 0;
445 m_pfoVtxX = 0.0;
446 m_pfoVtxY = 0.0;
447 m_pfoVtxZ = 0.0;
448 m_pfoEndX = 0.0;
449 m_pfoEndY = 0.0;
450 m_pfoEndZ = 0.0;
451 m_pfoDirX = 0.0;
452 m_pfoDirY = 0.0;
453 m_pfoDirZ = 0.0;
454 m_pfoLength = 0.0;
456
457 m_mcVertex = 0;
458 m_mcVtxX = 0.0;
459 m_mcVtxY = 0.0;
460 m_mcVtxZ = 0.0;
461 m_mcEndX = 0.0;
462 m_mcEndY = 0.0;
463 m_mcEndZ = 0.0;
464 m_mcDirX = 0.0;
465 m_mcDirY = 0.0;
466 m_mcDirZ = 0.0;
467 m_mcEnergy = 0.0;
468 m_mcLength = 0.0;
469 m_mcStraightLength = 0.0;
470
471 m_completeness = 0.0;
472 m_purity = 0.0;
473
474 m_nMCHits = 0;
475 m_nPfoHits = 0;
476 m_nMatchedHits = 0;
477 m_nMCHitsU = 0;
478 m_nMCHitsV = 0;
479 m_nMCHitsW = 0;
480 m_nPfoHitsU = 0;
481 m_nPfoHitsV = 0;
482 m_nPfoHitsW = 0;
483 m_nMatchedHitsU = 0;
484 m_nMatchedHitsV = 0;
485 m_nMatchedHitsW = 0;
486
489
490 m_spacepointsMinX = 0.0;
491 m_spacepointsMaxX = 0.0;
492
493 if (m_printDebug) {
494 std::cout << " Run: " << m_run << std::endl;
495 std::cout << " Event: " << m_event << std::endl;
496 }
497
498 // Collect Hits
499 // ============
500 HitVector hitVector;
502
503 if (m_printDebug) std::cout << " Hits: " << hitVector.size() << std::endl;
504
505 // Collect SpacePoints and SpacePoint <-> Hit Associations
506 // =======================================================
507 SpacePointVector spacePointVector;
508 SpacePointsToHits spacePointsToHits;
509 HitsToSpacePoints hitsToSpacePoints;
511 evt, m_particleLabel, spacePointVector, spacePointsToHits, hitsToSpacePoints);
512
513 if (m_printDebug) std::cout << " SpacePoints: " << spacePointVector.size() << std::endl;
514
515 // Collect Tracks and PFParticle <-> Track Associations
516 // ====================================================
517 TrackVector recoTrackVector;
518 PFParticlesToTracks recoParticlesToTracks;
519 LArPandoraHelper::CollectTracks(evt, m_trackLabel, recoTrackVector, recoParticlesToTracks);
520
521 if (m_printDebug) std::cout << " Tracks: " << recoTrackVector.size() << std::endl;
522
523 // Collect TOs and PFParticle <-> T0 Associations
524 // ==============================================
525 T0Vector t0Vector;
526 PFParticlesToT0s particlesToT0s;
527 LArPandoraHelper::CollectT0s(evt, m_particleLabel, t0Vector, particlesToT0s);
528
529 // Collect Vertices and PFParticle <-> Vertex Associations
530 // =======================================================
531 VertexVector recoVertexVector;
532 PFParticlesToVertices recoParticlesToVertices;
534 evt, m_particleLabel, recoVertexVector, recoParticlesToVertices);
535
536 if (m_printDebug) std::cout << " Vertices: " << recoVertexVector.size() << std::endl;
537
538 // Collect PFParticles and match Reco Particles to Hits
539 // ====================================================
540 PFParticleVector recoParticleVector;
541 PFParticleVector recoNeutrinoVector;
542 PFParticlesToHits recoParticlesToHits;
543 HitsToPFParticles recoHitsToParticles;
544
546 LArPandoraHelper::SelectNeutrinoPFParticles(recoParticleVector, recoNeutrinoVector);
548 evt,
550 recoParticlesToHits,
551 recoHitsToParticles,
555
556 if (m_printDebug) {
557 std::cout << " RecoNeutrinos: " << recoNeutrinoVector.size() << std::endl;
558 std::cout << " RecoParticles: " << recoParticleVector.size() << std::endl;
559 }
560
561 // Collect MCParticles and match True Particles to Hits
562 // ====================================================
563 MCParticleVector trueParticleVector;
564 MCTruthToMCParticles truthToParticles;
565 MCParticlesToMCTruth particlesToTruth;
566 MCParticlesToHits trueParticlesToHits;
567 HitsToMCParticles trueHitsToParticles;
568
569 if (m_disableRealDataCheck || !evt.isRealData()) {
572 evt, m_geantModuleLabel, truthToParticles, particlesToTruth);
573
575 evt,
577 hitVector,
578 trueParticlesToHits,
579 trueHitsToParticles,
583
584 if (trueHitsToParticles.empty()) {
585 if (m_backtrackerLabel.empty())
586 throw cet::exception("LArPandora") << " PFParticleMonitoring::analyze - no sim channels "
587 "found, backtracker module must be set in FHiCL "
588 << std::endl;
589
591 evt,
595 trueParticlesToHits,
596 trueHitsToParticles,
600 }
601 }
602
603 if (m_printDebug) {
604 std::cout << " TrueParticles: " << particlesToTruth.size() << std::endl;
605 std::cout << " TrueEvents: " << truthToParticles.size() << std::endl;
606 std::cout << " MatchedParticles: " << trueParticlesToHits.size() << std::endl;
607 }
608
609 if (trueParticlesToHits.empty()) {
610 m_pRecoTree->Fill();
611 return;
612 }
613
614 // Build Reco and True Particle Maps (for Parent/Daughter Navigation)
615 // =================================================================
616 MCParticleMap trueParticleMap;
617 PFParticleMap recoParticleMap;
618
619 LArPandoraHelper::BuildMCParticleMap(trueParticleVector, trueParticleMap);
620 LArPandoraHelper::BuildPFParticleMap(recoParticleVector, recoParticleMap);
621
622 m_nMCParticles = trueParticlesToHits.size();
623 m_nNeutrinoPfos = 0;
624 m_nPrimaryPfos = 0;
625 m_nDaughterPfos = 0;
626
627 // Count reconstructed particles
628 for (PFParticleVector::const_iterator iter = recoParticleVector.begin(),
629 iterEnd = recoParticleVector.end();
630 iter != iterEnd;
631 ++iter) {
632 const art::Ptr<recob::PFParticle> recoParticle = *iter;
633
634 if (LArPandoraHelper::IsNeutrino(recoParticle)) { m_nNeutrinoPfos++; }
635 else if (LArPandoraHelper::IsFinalState(recoParticleMap, recoParticle)) {
637 }
638 else {
640 }
641 }
642
643 // Match Reco Neutrinos to True Neutrinos
644 // ======================================
645 PFParticlesToHits recoNeutrinosToHits;
646 HitsToPFParticles recoHitsToNeutrinos;
647 HitsToMCTruth trueHitsToNeutrinos;
648 MCTruthToHits trueNeutrinosToHits;
650 recoParticleMap, recoParticlesToHits, recoNeutrinosToHits, recoHitsToNeutrinos);
652 truthToParticles, trueParticlesToHits, trueNeutrinosToHits, trueHitsToNeutrinos);
653
654 MCTruthToPFParticles matchedNeutrinos;
655 MCTruthToHits matchedNeutrinoHits;
657 recoNeutrinosToHits, trueHitsToNeutrinos, matchedNeutrinos, matchedNeutrinoHits);
658
659 for (MCTruthToHits::const_iterator iter = trueNeutrinosToHits.begin(),
660 iterEnd = trueNeutrinosToHits.end();
661 iter != iterEnd;
662 ++iter) {
663 const art::Ptr<simb::MCTruth> trueEvent = iter->first;
664 const HitVector& trueHitVector = iter->second;
665
666 if (trueHitVector.empty()) continue;
667
668 if (!trueEvent->NeutrinoSet()) continue;
669
670 const simb::MCNeutrino trueNeutrino(trueEvent->GetNeutrino());
671 const simb::MCParticle trueParticle(trueNeutrino.Nu());
672
673 m_mcIsCC = ((simb::kCC == trueNeutrino.CCNC()) ? 1 : 0);
674 m_mcPdg = trueParticle.PdgCode();
676 m_mcParentPdg = 0;
677 m_mcPrimaryPdg = 0;
678 m_mcIsNeutrino = 1;
679 m_mcIsPrimary = 0;
680 m_mcIsDecay = 0;
681
682 m_mcVertex = 1;
683 m_mcVtxX = trueParticle.Vx();
684 m_mcVtxY = trueParticle.Vy();
685 m_mcVtxZ = trueParticle.Vz();
689 m_mcDirX = trueParticle.Px() / trueParticle.P();
690 m_mcDirY = trueParticle.Py() / trueParticle.P();
691 m_mcDirZ = trueParticle.Pz() / trueParticle.P();
692 m_mcEnergy = trueParticle.E();
693 m_mcLength = 0.0;
694 m_mcStraightLength = 0.0;
695
696 m_nMCHits = trueHitVector.size();
697 m_nMCHitsU = this->CountHitsByType(geo::kU, trueHitVector);
698 m_nMCHitsV = this->CountHitsByType(geo::kV, trueHitVector);
699 m_nMCHitsW = this->CountHitsByType(geo::kW, trueHitVector);
700
701 m_pfoPdg = 0;
702 m_pfoNuPdg = 0;
703 m_pfoParentPdg = 0;
704 m_pfoPrimaryPdg = 0;
705 m_pfoIsNeutrino = 0;
706 m_pfoIsPrimary = 0;
707 m_pfoIsStitched = 0;
708 m_pfoTrack = 0;
709 m_pfoVertex = 0;
710 m_pfoVtxX = 0.0;
711 m_pfoVtxY = 0.0;
712 m_pfoVtxZ = 0.0;
713 m_pfoEndX = 0.0;
714 m_pfoEndY = 0.0;
715 m_pfoEndZ = 0.0;
716 m_pfoDirX = 0.0;
717 m_pfoDirY = 0.0;
718 m_pfoDirZ = 0.0;
719 m_pfoLength = 0.0;
721
722 m_nPfoHits = 0;
723 m_nPfoHitsU = 0;
724 m_nPfoHitsV = 0;
725 m_nPfoHitsW = 0;
726
727 m_nMatchedHits = 0;
728 m_nMatchedHitsU = 0;
729 m_nMatchedHitsV = 0;
730 m_nMatchedHitsW = 0;
731
734
735 m_spacepointsMinX = 0.0;
736 m_spacepointsMaxX = 0.0;
737
738 m_completeness = 0.0;
739 m_purity = 0.0;
740
741 for (HitVector::const_iterator hIter1 = trueHitVector.begin(),
742 hIterEnd1 = trueHitVector.end();
743 hIter1 != hIterEnd1;
744 ++hIter1) {
745 if (recoHitsToNeutrinos.find(*hIter1) == recoHitsToNeutrinos.end())
747 }
748
749 MCTruthToPFParticles::const_iterator pIter1 = matchedNeutrinos.find(trueEvent);
750 if (matchedNeutrinos.end() != pIter1) {
751 const art::Ptr<recob::PFParticle> recoParticle = pIter1->second;
752 m_pfoPdg = recoParticle->PdgCode();
755 m_pfoPrimaryPdg = 0;
756 m_pfoIsNeutrino = 1;
757 m_pfoIsPrimary = 0;
758
759 if (!LArPandoraHelper::IsNeutrino(recoParticle))
760 std::cout << " Warning: Found neutrino with an invalid PDG code " << std::endl;
761
762 PFParticlesToHits::const_iterator pIter2 = recoNeutrinosToHits.find(recoParticle);
763 if (recoParticlesToHits.end() != pIter2) {
764 const HitVector& recoHitVector = pIter2->second;
765
766 for (HitVector::const_iterator hIter2 = recoHitVector.begin(),
767 hIterEnd2 = recoHitVector.end();
768 hIter2 != hIterEnd2;
769 ++hIter2) {
770 if (trueHitsToNeutrinos.find(*hIter2) == trueHitsToNeutrinos.end())
772 }
773
774 MCTruthToHits::const_iterator pIter3 = matchedNeutrinoHits.find(trueEvent);
775 if (matchedNeutrinoHits.end() != pIter3) {
776 const HitVector& matchedHitVector = pIter3->second;
777
778 m_nPfoHits = recoHitVector.size();
779 m_nPfoHitsU = this->CountHitsByType(geo::kU, recoHitVector);
780 m_nPfoHitsV = this->CountHitsByType(geo::kV, recoHitVector);
781 m_nPfoHitsW = this->CountHitsByType(geo::kW, recoHitVector);
782
783 m_nMatchedHits = matchedHitVector.size();
784 m_nMatchedHitsU = this->CountHitsByType(geo::kU, matchedHitVector);
785 m_nMatchedHitsV = this->CountHitsByType(geo::kV, matchedHitVector);
786 m_nMatchedHitsW = this->CountHitsByType(geo::kW, matchedHitVector);
787
788 PFParticlesToVertices::const_iterator pIter4 =
789 recoParticlesToVertices.find(recoParticle);
790 if (recoParticlesToVertices.end() != pIter4) {
791 const VertexVector& vertexVector = pIter4->second;
792 if (!vertexVector.empty()) {
793 if (vertexVector.size() != 1 && m_printDebug)
794 std::cout << " Warning: Found particle with more than one associated vertex "
795 << std::endl;
796
797 const art::Ptr<recob::Vertex> recoVertex = *(vertexVector.begin());
798 double xyz[3] = {0.0, 0.0, 0.0};
799 recoVertex->XYZ(xyz);
800
801 m_pfoVertex = 1;
802 m_pfoVtxX = xyz[0];
803 m_pfoVtxY = xyz[1];
804 m_pfoVtxZ = xyz[2];
805 }
806 }
807 }
808 }
809 }
810
811 m_purity =
812 ((m_nPfoHits == 0) ? 0.0 :
813 static_cast<double>(m_nMatchedHits) / static_cast<double>(m_nPfoHits));
815 ((m_nPfoHits == 0) ? 0.0 :
816 static_cast<double>(m_nMatchedHits) / static_cast<double>(m_nMCHits));
817
818 if (m_printDebug)
819 std::cout << " MCNeutrino [" << m_index << "]"
820 << " trueNu=" << m_mcNuPdg << ", truePdg=" << m_mcPdg
821 << ", recoNu=" << m_pfoNuPdg << ", recoPdg=" << m_pfoPdg
822 << ", mcHits=" << m_nMCHits << ", pfoHits=" << m_nPfoHits
823 << ", matchedHits=" << m_nMatchedHits
824 << ", availableHits=" << m_nTrueWithoutRecoHits << std::endl;
825
826 m_pRecoTree->Fill();
827 ++m_index; // Increment index number
828 }
829
830 // Match Reco Particles to True Particles
831 // ======================================
832 MCParticlesToPFParticles matchedParticles;
833 MCParticlesToHits matchedParticleHits;
835 recoParticlesToHits, trueHitsToParticles, matchedParticles, matchedParticleHits);
836
837 // Compare true and reconstructed particles
838 for (MCParticlesToHits::const_iterator iter = trueParticlesToHits.begin(),
839 iterEnd = trueParticlesToHits.end();
840 iter != iterEnd;
841 ++iter) {
842 const art::Ptr<simb::MCParticle> trueParticle = iter->first;
843 const HitVector& trueHitVector = iter->second;
844
845 if (trueHitVector.empty()) continue;
846
847 m_mcPdg = trueParticle->PdgCode();
848 m_mcNuPdg = 0;
849 m_mcParentPdg = 0;
850 m_mcPrimaryPdg = 0;
851 m_mcIsNeutrino = 0;
852 m_mcIsPrimary = 0;
853 m_mcIsDecay = 0;
854 m_mcIsCC = 0;
855
856 m_pfoPdg = 0;
857 m_pfoNuPdg = 0;
858 m_pfoParentPdg = 0;
859 m_pfoPrimaryPdg = 0;
860 m_pfoIsNeutrino = 0;
861 m_pfoIsPrimary = 0;
862 m_pfoIsStitched = 0;
863 m_pfoTrack = 0;
864 m_pfoVertex = 0;
865 m_pfoVtxX = 0.0;
866 m_pfoVtxY = 0.0;
867 m_pfoVtxZ = 0.0;
868 m_pfoEndX = 0.0;
869 m_pfoEndY = 0.0;
870 m_pfoEndZ = 0.0;
871 m_pfoDirX = 0.0;
872 m_pfoDirY = 0.0;
873 m_pfoDirZ = 0.0;
874 m_pfoLength = 0.0;
876
877 m_mcVertex = 0;
878 m_mcVtxX = 0.0;
879 m_mcVtxY = 0.0;
880 m_mcVtxZ = 0.0;
881 m_mcEndX = 0.0;
882 m_mcEndY = 0.0;
883 m_mcEndZ = 0.0;
884 m_mcDirX = 0.0;
885 m_mcDirY = 0.0;
886 m_mcDirZ = 0.0;
887 m_mcEnergy = 0.0;
888 m_mcLength = 0.0;
889 m_mcStraightLength = 0.0;
890
891 m_completeness = 0.0;
892 m_purity = 0.0;
893
894 m_nMCHits = 0;
895 m_nMCHitsU = 0;
896 m_nMCHitsV = 0;
897 m_nMCHitsW = 0;
898
899 m_nPfoHits = 0;
900 m_nPfoHitsU = 0;
901 m_nPfoHitsV = 0;
902 m_nPfoHitsW = 0;
903
904 m_nMatchedHits = 0;
905 m_nMatchedHitsU = 0;
906 m_nMatchedHitsV = 0;
907 m_nMatchedHitsW = 0;
908
911
912 m_spacepointsMinX = 0.0;
913 m_spacepointsMaxX = 0.0;
914
915 // Set true properties
916 try {
917 int startT(-1);
918 int endT(-1);
919 this->GetStartAndEndPoints(trueParticle, startT, endT);
920
921 // vertex and end positions
922 m_mcVertex = 1;
923 m_mcVtxX = trueParticle->Vx(startT);
924 m_mcVtxY = trueParticle->Vy(startT);
925 m_mcVtxZ = trueParticle->Vz(startT);
926 m_mcEndX = trueParticle->Vx(endT);
927 m_mcEndY = trueParticle->Vy(endT);
928 m_mcEndZ = trueParticle->Vz(endT);
929
930 const double dx(m_mcEndX - m_mcVtxX);
931 const double dy(m_mcEndY - m_mcVtxY);
932 const double dz(m_mcEndZ - m_mcVtxZ);
933
934 m_mcStraightLength = std::sqrt(dx * dx + dy * dy + dz * dz);
935 m_mcLength = this->GetLength(trueParticle, startT, endT);
936
937 // energy and momentum
938 const double Ptot(trueParticle->P(startT));
939
940 if (Ptot > 0.0) {
941 m_mcDirX = trueParticle->Px(startT) / Ptot;
942 m_mcDirY = trueParticle->Py(startT) / Ptot;
943 m_mcDirZ = trueParticle->Pz(startT) / Ptot;
944 m_mcEnergy = trueParticle->E(startT);
945 }
946 }
947 catch (cet::exception& e) {
948 }
949
950 // Get the true parent neutrino
951 MCParticlesToMCTruth::const_iterator nuIter = particlesToTruth.find(trueParticle);
952 if (particlesToTruth.end() == nuIter)
953 throw cet::exception("LArPandora") << " PFParticleMonitoring::analyze --- Found a true "
954 "particle without any ancestry information ";
955
956 const art::Ptr<simb::MCTruth> trueEvent = nuIter->second;
957
958 if (trueEvent->NeutrinoSet()) {
959 const simb::MCNeutrino neutrino(trueEvent->GetNeutrino());
960 m_mcNuPdg = neutrino.Nu().PdgCode();
961 m_mcIsCC = ((simb::kCC == neutrino.CCNC()) ? 1 : 0);
962 }
963
964 // Get the true 'parent' and 'primary' particles
965 try {
966 const art::Ptr<simb::MCParticle> parentParticle(
967 LArPandoraHelper::GetParentMCParticle(trueParticleMap, trueParticle));
968 const art::Ptr<simb::MCParticle> primaryParticle(
969 LArPandoraHelper::GetFinalStateMCParticle(trueParticleMap, trueParticle));
970 m_mcParentPdg = ((parentParticle != trueParticle) ? parentParticle->PdgCode() : 0);
971 m_mcPrimaryPdg = primaryParticle->PdgCode();
972 m_mcIsPrimary = (primaryParticle == trueParticle);
973 m_mcIsDecay = ("Decay" == trueParticle->Process());
974 }
975 catch (cet::exception& e) {
976 }
977
978 // Find min and max X positions of space points
979 bool foundSpacePoints(false);
980
981 for (HitVector::const_iterator hIter1 = trueHitVector.begin(),
982 hIterEnd1 = trueHitVector.end();
983 hIter1 != hIterEnd1;
984 ++hIter1) {
985 const art::Ptr<recob::Hit> hit = *hIter1;
986
987 HitsToSpacePoints::const_iterator hIter2 = hitsToSpacePoints.find(hit);
988 if (hitsToSpacePoints.end() == hIter2) continue;
989
990 const art::Ptr<recob::SpacePoint> spacepoint = hIter2->second;
991 const double X(spacepoint->XYZ()[0]);
992
993 if (!foundSpacePoints) {
996 foundSpacePoints = true;
997 }
998 else {
1000 m_spacepointsMaxX = std::max(m_spacepointsMaxX, X);
1001 }
1002 }
1003
1004 // Count number of available hits
1005 for (HitVector::const_iterator hIter1 = trueHitVector.begin(),
1006 hIterEnd1 = trueHitVector.end();
1007 hIter1 != hIterEnd1;
1008 ++hIter1) {
1009 if (recoHitsToParticles.find(*hIter1) == recoHitsToParticles.end())
1011 }
1012
1013 // Match true and reconstructed hits
1014 m_nMCHits = trueHitVector.size();
1015 m_nMCHitsU = this->CountHitsByType(geo::kU, trueHitVector);
1016 m_nMCHitsV = this->CountHitsByType(geo::kV, trueHitVector);
1017 m_nMCHitsW = this->CountHitsByType(geo::kW, trueHitVector);
1018
1019 MCParticlesToPFParticles::const_iterator pIter1 = matchedParticles.find(trueParticle);
1020 if (matchedParticles.end() != pIter1) {
1021 const art::Ptr<recob::PFParticle> recoParticle = pIter1->second;
1022 m_pfoPdg = recoParticle->PdgCode();
1023 m_pfoNuPdg = LArPandoraHelper::GetParentNeutrino(recoParticleMap, recoParticle);
1024 m_pfoIsPrimary = LArPandoraHelper::IsFinalState(recoParticleMap, recoParticle);
1025
1026 const art::Ptr<recob::PFParticle> parentParticle =
1027 LArPandoraHelper::GetParentPFParticle(recoParticleMap, recoParticle);
1028 m_pfoParentPdg = parentParticle->PdgCode();
1029
1030 const art::Ptr<recob::PFParticle> primaryParticle =
1031 LArPandoraHelper::GetFinalStatePFParticle(recoParticleMap, recoParticle);
1032 m_pfoPrimaryPdg = primaryParticle->PdgCode();
1033
1034 PFParticlesToHits::const_iterator pIter2 = recoParticlesToHits.find(recoParticle);
1035 if (recoParticlesToHits.end() == pIter2)
1036 throw cet::exception("LArPandora")
1037 << " PFParticleMonitoring::analyze --- Found a reco particle without any hits ";
1038
1039 const HitVector& recoHitVector = pIter2->second;
1040
1041 for (HitVector::const_iterator hIter2 = recoHitVector.begin(),
1042 hIterEnd2 = recoHitVector.end();
1043 hIter2 != hIterEnd2;
1044 ++hIter2) {
1045 if (trueHitsToParticles.find(*hIter2) == trueHitsToParticles.end())
1047 }
1048
1049 MCParticlesToHits::const_iterator pIter3 = matchedParticleHits.find(trueParticle);
1050 if (matchedParticleHits.end() == pIter3)
1051 throw cet::exception("LArPandora") << " PFParticleMonitoring::analyze --- Found a "
1052 "matched true particle without matched hits ";
1053
1054 const HitVector& matchedHitVector = pIter3->second;
1055
1056 m_nPfoHits = recoHitVector.size();
1057 m_nPfoHitsU = this->CountHitsByType(geo::kU, recoHitVector);
1058 m_nPfoHitsV = this->CountHitsByType(geo::kV, recoHitVector);
1059 m_nPfoHitsW = this->CountHitsByType(geo::kW, recoHitVector);
1060
1061 m_nMatchedHits = matchedHitVector.size();
1062 m_nMatchedHitsU = this->CountHitsByType(geo::kU, matchedHitVector);
1063 m_nMatchedHitsV = this->CountHitsByType(geo::kV, matchedHitVector);
1064 m_nMatchedHitsW = this->CountHitsByType(geo::kW, matchedHitVector);
1065
1066 PFParticlesToVertices::const_iterator pIter4 = recoParticlesToVertices.find(recoParticle);
1067 if (recoParticlesToVertices.end() != pIter4) {
1068 const VertexVector& vertexVector = pIter4->second;
1069 if (!vertexVector.empty()) {
1070 if (vertexVector.size() != 1 && m_printDebug)
1071 std::cout << " Warning: Found particle with more than one associated vertex "
1072 << std::endl;
1073
1074 const art::Ptr<recob::Vertex> recoVertex = *(vertexVector.begin());
1075 double xyz[3] = {0.0, 0.0, 0.0};
1076 recoVertex->XYZ(xyz);
1077
1078 m_pfoVertex = 1;
1079 m_pfoVtxX = xyz[0];
1080 m_pfoVtxY = xyz[1];
1081 m_pfoVtxZ = xyz[2];
1082 }
1083 }
1084
1085 PFParticlesToTracks::const_iterator pIter5 = recoParticlesToTracks.find(recoParticle);
1086 if (recoParticlesToTracks.end() != pIter5) {
1087 const TrackVector& trackVector = pIter5->second;
1088 if (!trackVector.empty()) {
1089 if (trackVector.size() != 1 && m_printDebug)
1090 std::cout << " Warning: Found particle with more than one associated track "
1091 << std::endl;
1092
1093 const art::Ptr<recob::Track> recoTrack = *(trackVector.begin());
1094 const auto& vtxPosition = recoTrack->Vertex();
1095 const auto& endPosition = recoTrack->End();
1096 const auto& vtxDirection = recoTrack->VertexDirection();
1097
1098 m_pfoTrack = 1;
1099 m_pfoVtxX = vtxPosition.x();
1100 m_pfoVtxY = vtxPosition.y();
1101 m_pfoVtxZ = vtxPosition.z();
1102 m_pfoEndX = endPosition.x();
1103 m_pfoEndY = endPosition.y();
1104 m_pfoEndZ = endPosition.z();
1105 m_pfoDirX = vtxDirection.x();
1106 m_pfoDirY = vtxDirection.y();
1107 m_pfoDirZ = vtxDirection.z();
1108 m_pfoStraightLength = (endPosition - vtxPosition).R();
1109 m_pfoLength = recoTrack->Length();
1110 }
1111 }
1112
1113 m_pfoIsStitched = (particlesToT0s.end() != particlesToT0s.find(recoParticle));
1114 }
1115
1116 m_purity =
1117 ((m_nPfoHits == 0) ? 0.0 :
1118 static_cast<double>(m_nMatchedHits) / static_cast<double>(m_nPfoHits));
1120 ((m_nPfoHits == 0) ? 0.0 :
1121 static_cast<double>(m_nMatchedHits) / static_cast<double>(m_nMCHits));
1122
1123 if (m_printDebug)
1124 std::cout << " MCParticle [" << m_index << "]"
1125 << " trueNu=" << m_mcNuPdg << ", truePdg=" << m_mcPdg
1126 << ", recoNu=" << m_pfoNuPdg << ", recoPdg=" << m_pfoPdg
1127 << ", mcHits=" << m_nMCHits << ", pfoHits=" << m_nPfoHits
1128 << ", matchedHits=" << m_nMatchedHits
1129 << ", availableHits=" << m_nTrueWithoutRecoHits << std::endl;
1130
1131 m_pRecoTree->Fill();
1132 ++m_index; // Increment index number
1133 }
1134 }
1135
1136 //------------------------------------------------------------------------------------------------------------------------------------------
1137
1139 const MCParticlesToHits& trueParticlesToHits,
1140 MCTruthToHits& trueNeutrinosToHits,
1141 HitsToMCTruth& trueHitsToNeutrinos) const
1142 {
1143 for (MCTruthToMCParticles::const_iterator iter1 = truthToParticles.begin(),
1144 iterEnd1 = truthToParticles.end();
1145 iter1 != iterEnd1;
1146 ++iter1) {
1147 const art::Ptr<simb::MCTruth> trueNeutrino = iter1->first;
1148 const MCParticleVector& trueParticleVector = iter1->second;
1149
1150 for (MCParticleVector::const_iterator iter2 = trueParticleVector.begin(),
1151 iterEnd2 = trueParticleVector.end();
1152 iter2 != iterEnd2;
1153 ++iter2) {
1154 const MCParticlesToHits::const_iterator iter3 = trueParticlesToHits.find(*iter2);
1155 if (trueParticlesToHits.end() == iter3) continue;
1156
1157 const HitVector& hitVector = iter3->second;
1158
1159 for (HitVector::const_iterator iter4 = hitVector.begin(), iterEnd4 = hitVector.end();
1160 iter4 != iterEnd4;
1161 ++iter4) {
1162 const art::Ptr<recob::Hit> hit = *iter4;
1163 trueHitsToNeutrinos[hit] = trueNeutrino;
1164 trueNeutrinosToHits[trueNeutrino].push_back(hit);
1165 }
1166 }
1167 }
1168 }
1169
1170 //------------------------------------------------------------------------------------------------------------------------------------------
1171
1173 const PFParticlesToHits& recoParticlesToHits,
1174 PFParticlesToHits& recoNeutrinosToHits,
1175 HitsToPFParticles& recoHitsToNeutrinos) const
1176 {
1177 for (PFParticleMap::const_iterator iter1 = recoParticleMap.begin(),
1178 iterEnd1 = recoParticleMap.end();
1179 iter1 != iterEnd1;
1180 ++iter1) {
1181 const art::Ptr<recob::PFParticle> recoParticle = iter1->second;
1182 const art::Ptr<recob::PFParticle> recoNeutrino =
1183 LArPandoraHelper::GetParentPFParticle(recoParticleMap, recoParticle);
1184
1185 if (!LArPandoraHelper::IsNeutrino(recoNeutrino)) continue;
1186
1187 const PFParticlesToHits::const_iterator iter2 = recoParticlesToHits.find(recoParticle);
1188 if (recoParticlesToHits.end() == iter2) continue;
1189
1190 const HitVector& hitVector = iter2->second;
1191
1192 for (HitVector::const_iterator iter3 = hitVector.begin(), iterEnd3 = hitVector.end();
1193 iter3 != iterEnd3;
1194 ++iter3) {
1195 const art::Ptr<recob::Hit> hit = *iter3;
1196 recoHitsToNeutrinos[hit] = recoNeutrino;
1197 recoNeutrinosToHits[recoNeutrino].push_back(hit);
1198 }
1199 }
1200 }
1201
1202 //------------------------------------------------------------------------------------------------------------------------------------------
1203
1205 const HitsToMCTruth& trueHitsToNeutrinos,
1206 MCTruthToPFParticles& matchedNeutrinos,
1207 MCTruthToHits& matchedNeutrinoHits) const
1208 {
1209 PFParticleSet recoVeto;
1210 MCTruthSet trueVeto;
1211
1212 this->GetRecoToTrueMatches(recoNeutrinosToHits,
1213 trueHitsToNeutrinos,
1214 matchedNeutrinos,
1215 matchedNeutrinoHits,
1216 recoVeto,
1217 trueVeto);
1218 }
1219
1220 //------------------------------------------------------------------------------------------------------------------------------------------
1221
1223 const HitsToMCTruth& trueHitsToNeutrinos,
1224 MCTruthToPFParticles& matchedNeutrinos,
1225 MCTruthToHits& matchedNeutrinoHits,
1226 PFParticleSet& vetoReco,
1227 MCTruthSet& vetoTrue) const
1228 {
1229 bool foundMatches(false);
1230
1231 for (PFParticlesToHits::const_iterator iter1 = recoNeutrinosToHits.begin(),
1232 iterEnd1 = recoNeutrinosToHits.end();
1233 iter1 != iterEnd1;
1234 ++iter1) {
1235 const art::Ptr<recob::PFParticle> recoNeutrino = iter1->first;
1236 if (vetoReco.count(recoNeutrino) > 0) continue;
1237
1238 const HitVector& hitVector = iter1->second;
1239
1240 MCTruthToHits truthContributionMap;
1241
1242 for (HitVector::const_iterator iter2 = hitVector.begin(), iterEnd2 = hitVector.end();
1243 iter2 != iterEnd2;
1244 ++iter2) {
1245 const art::Ptr<recob::Hit> hit = *iter2;
1246
1247 HitsToMCTruth::const_iterator iter3 = trueHitsToNeutrinos.find(hit);
1248 if (trueHitsToNeutrinos.end() == iter3) continue;
1249
1250 const art::Ptr<simb::MCTruth> trueNeutrino = iter3->second;
1251 if (vetoTrue.count(trueNeutrino) > 0) continue;
1252
1253 truthContributionMap[trueNeutrino].push_back(hit);
1254 }
1255
1256 MCTruthToHits::const_iterator mIter = truthContributionMap.end();
1257
1258 for (MCTruthToHits::const_iterator iter4 = truthContributionMap.begin(),
1259 iterEnd4 = truthContributionMap.end();
1260 iter4 != iterEnd4;
1261 ++iter4) {
1262 if ((truthContributionMap.end() == mIter) ||
1263 (iter4->second.size() > mIter->second.size())) {
1264 mIter = iter4;
1265 }
1266 }
1267
1268 if (truthContributionMap.end() != mIter) {
1269 const art::Ptr<simb::MCTruth> trueNeutrino = mIter->first;
1270
1271 MCTruthToHits::const_iterator iter5 = matchedNeutrinoHits.find(trueNeutrino);
1272
1273 if ((matchedNeutrinoHits.end() == iter5) || (mIter->second.size() > iter5->second.size())) {
1274 matchedNeutrinos[trueNeutrino] = recoNeutrino;
1275 matchedNeutrinoHits[trueNeutrino] = mIter->second;
1276 foundMatches = true;
1277 }
1278 }
1279 }
1280
1281 if (!foundMatches) return;
1282
1283 for (MCTruthToPFParticles::const_iterator pIter = matchedNeutrinos.begin(),
1284 pIterEnd = matchedNeutrinos.end();
1285 pIter != pIterEnd;
1286 ++pIter) {
1287 vetoTrue.insert(pIter->first);
1288 vetoReco.insert(pIter->second);
1289 }
1290
1292 this->GetRecoToTrueMatches(recoNeutrinosToHits,
1293 trueHitsToNeutrinos,
1294 matchedNeutrinos,
1295 matchedNeutrinoHits,
1296 vetoReco,
1297 vetoTrue);
1298 }
1299
1300 //------------------------------------------------------------------------------------------------------------------------------------------
1301
1303 const HitsToMCParticles& trueHitsToParticles,
1304 MCParticlesToPFParticles& matchedParticles,
1305 MCParticlesToHits& matchedHits) const
1306 {
1307 PFParticleSet recoVeto;
1308 MCParticleSet trueVeto;
1309
1311 recoParticlesToHits, trueHitsToParticles, matchedParticles, matchedHits, recoVeto, trueVeto);
1312 }
1313
1314 //------------------------------------------------------------------------------------------------------------------------------------------
1315
1317 const HitsToMCParticles& trueHitsToParticles,
1318 MCParticlesToPFParticles& matchedParticles,
1319 MCParticlesToHits& matchedHits,
1320 PFParticleSet& vetoReco,
1321 MCParticleSet& vetoTrue) const
1322 {
1323 bool foundMatches(false);
1324
1325 for (PFParticlesToHits::const_iterator iter1 = recoParticlesToHits.begin(),
1326 iterEnd1 = recoParticlesToHits.end();
1327 iter1 != iterEnd1;
1328 ++iter1) {
1329 const art::Ptr<recob::PFParticle> recoParticle = iter1->first;
1330 if (vetoReco.count(recoParticle) > 0) continue;
1331
1332 const HitVector& hitVector = iter1->second;
1333
1334 MCParticlesToHits truthContributionMap;
1335
1336 for (HitVector::const_iterator iter2 = hitVector.begin(), iterEnd2 = hitVector.end();
1337 iter2 != iterEnd2;
1338 ++iter2) {
1339 const art::Ptr<recob::Hit> hit = *iter2;
1340
1341 HitsToMCParticles::const_iterator iter3 = trueHitsToParticles.find(hit);
1342 if (trueHitsToParticles.end() == iter3) continue;
1343
1344 const art::Ptr<simb::MCParticle> trueParticle = iter3->second;
1345 if (vetoTrue.count(trueParticle) > 0) continue;
1346
1347 truthContributionMap[trueParticle].push_back(hit);
1348 }
1349
1350 MCParticlesToHits::const_iterator mIter = truthContributionMap.end();
1351
1352 for (MCParticlesToHits::const_iterator iter4 = truthContributionMap.begin(),
1353 iterEnd4 = truthContributionMap.end();
1354 iter4 != iterEnd4;
1355 ++iter4) {
1356 if ((truthContributionMap.end() == mIter) ||
1357 (iter4->second.size() > mIter->second.size())) {
1358 mIter = iter4;
1359 }
1360 }
1361
1362 if (truthContributionMap.end() != mIter) {
1363 const art::Ptr<simb::MCParticle> trueParticle = mIter->first;
1364
1365 MCParticlesToHits::const_iterator iter5 = matchedHits.find(trueParticle);
1366
1367 if ((matchedHits.end() == iter5) || (mIter->second.size() > iter5->second.size())) {
1368 matchedParticles[trueParticle] = recoParticle;
1369 matchedHits[trueParticle] = mIter->second;
1370 foundMatches = true;
1371 }
1372 }
1373 }
1374
1375 if (!foundMatches) return;
1376
1377 for (MCParticlesToPFParticles::const_iterator pIter = matchedParticles.begin(),
1378 pIterEnd = matchedParticles.end();
1379 pIter != pIterEnd;
1380 ++pIter) {
1381 vetoTrue.insert(pIter->first);
1382 vetoReco.insert(pIter->second);
1383 }
1384
1386 this->GetRecoToTrueMatches(recoParticlesToHits,
1387 trueHitsToParticles,
1388 matchedParticles,
1389 matchedHits,
1390 vetoReco,
1391 vetoTrue);
1392 }
1393
1394 //------------------------------------------------------------------------------------------------------------------------------------------
1395
1396 int PFParticleMonitoring::CountHitsByType(const int view, const HitVector& hitVector) const
1397 {
1398 int nHits(0);
1399
1400 for (HitVector::const_iterator iter = hitVector.begin(), iterEnd = hitVector.end();
1401 iter != iterEnd;
1402 ++iter) {
1403 const art::Ptr<recob::Hit> hit = *iter;
1404 if (hit->View() == view) ++nHits;
1405 }
1406
1407 return nHits;
1408 }
1409
1410 //------------------------------------------------------------------------------------------------------------------------------------------
1411
1412 void PFParticleMonitoring::GetStartAndEndPoints(const art::Ptr<simb::MCParticle> particle,
1413 int& startT,
1414 int& endT) const
1415 {
1416 art::ServiceHandle<geo::Geometry const> theGeometry;
1417
1418 bool foundStartPosition(false);
1419
1420 const int numTrajectoryPoints(static_cast<int>(particle->NumberTrajectoryPoints()));
1421
1422 for (int nt = 0; nt < numTrajectoryPoints; ++nt) {
1423 geo::Point_t const pos{particle->Vx(nt), particle->Vy(nt), particle->Vz(nt)};
1424 if (theGeometry->PositionToTPCptr(pos) == nullptr) { continue; }
1425
1426 // TODO: Apply fiducial cut due to readout window
1427
1428 endT = nt;
1429 if (!foundStartPosition) {
1430 startT = endT;
1431 foundStartPosition = true;
1432 }
1433 }
1434
1435 if (!foundStartPosition) throw cet::exception("LArPandora");
1436 }
1437
1438 //------------------------------------------------------------------------------------------------------------------------------------------
1439
1440 double PFParticleMonitoring::GetLength(const art::Ptr<simb::MCParticle> particle,
1441 const int startT,
1442 const int endT) const
1443 {
1444 if (endT <= startT) return 0.0;
1445
1446 double length(0.0);
1447
1448 for (int nt = startT; nt < endT; ++nt) {
1449 const double dx(particle->Vx(nt + 1) - particle->Vx(nt));
1450 const double dy(particle->Vy(nt + 1) - particle->Vy(nt));
1451 const double dz(particle->Vz(nt + 1) - particle->Vz(nt));
1452 length += sqrt(dx * dx + dy * dy + dz * dz);
1453 }
1454
1455 return length;
1456 }
1457
1458} //namespace lar_pandora
helper function for LArPandoraInterface producer module
static void SelectNeutrinoPFParticles(const PFParticleVector &inputParticles, PFParticleVector &outputParticles)
Select reconstructed neutrino particles from a list of all reconstructed particles.
static void CollectHits(const art::Event &evt, const std::string &label, HitVector &hitVector)
Collect the reconstructed Hits from the ART event record.
static void BuildMCParticleMap(const MCParticleVector &particleVector, MCParticleMap &particleMap)
Build particle maps for true particles.
static void CollectVertices(const art::Event &evt, const std::string &label, VertexVector &vertexVector, PFParticlesToVertices &particlesToVertices)
Collect the reconstructed PFParticles and associated Vertices from the ART event record.
static art::Ptr< simb::MCParticle > GetFinalStateMCParticle(const MCParticleMap &particleMap, const art::Ptr< simb::MCParticle > daughterParticle)
Return the final-state parent particle by navigating up the chain of parent/daughter associations.
static void BuildPFParticleMap(const PFParticleVector &particleVector, PFParticleMap &particleMap)
Build particle maps for reconstructed particles.
static art::Ptr< recob::PFParticle > GetFinalStatePFParticle(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Return the final-state parent particle by navigating up the chain of parent/daughter associations.
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 art::Ptr< simb::MCParticle > GetParentMCParticle(const MCParticleMap &particleMap, const art::Ptr< simb::MCParticle > daughterParticle)
Return the top-level parent particle by navigating up the chain of parent/daughter associations.
static void CollectT0s(const art::Event &evt, const std::string &label, T0Vector &t0Vector, PFParticlesToT0s &particlesToT0s)
Collect a vector of T0s from the ART event record.
static void BuildMCParticleHitMaps(const art::Event &evt, const HitVector &hitVector, const SimChannelVector &simChannelVector, HitsToTrackIDEs &hitsToTrackIDEs)
Collect the links from reconstructed hits to their true energy deposits.
static art::Ptr< recob::PFParticle > GetParentPFParticle(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Return the top-level parent particle by navigating up the chain of parent/daughter associations.
static int GetParentNeutrino(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Return the parent neutrino PDG code (or zero for cosmics) for a given reconstructed particle.
static void CollectPFParticles(const art::Event &evt, const std::string &label, PFParticleVector &particleVector)
Collect the reconstructed PFParticles from the ART event record.
static bool IsFinalState(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Determine whether a particle has been reconstructed as a final-state particle.
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.
static void CollectMCParticles(const art::Event &evt, const std::string &label, MCParticleVector &particleVector)
Collect a vector of MCParticle objects from the ART event record.
std::set< art::Ptr< simb::MCParticle > > MCParticleSet
std::set< art::Ptr< recob::PFParticle > > PFParticleSet
int m_nRecoWithoutTrueHits
Reconstructed hits which don't belong to any true particle - "missing".
void GetRecoToTrueMatches(const PFParticlesToHits &recoNeutrinosToHits, const HitsToMCTruth &trueHitsToNeutrinos, MCTruthToPFParticles &matchedNeutrinos, MCTruthToHits &matchedNeutrinoHits) const
Perform matching between true and reconstructed neutrino events.
void reconfigure(fhicl::ParameterSet const &pset)
void BuildRecoNeutrinoHitMaps(const PFParticleMap &recoParticleMap, const PFParticlesToHits &recoParticlesToHits, PFParticlesToHits &recoNeutrinosToHits, HitsToPFParticles &recoHitsToNeutrinos) const
Build mapping from reconstructed neutrinos to hits.
void GetStartAndEndPoints(const art::Ptr< simb::MCParticle > trueParticle, int &startT, int &endT) const
Find the start and end points of the true particle in the active region of detector.
std::set< art::Ptr< simb::MCTruth > > MCTruthSet
int CountHitsByType(const int view, const HitVector &hitVector) const
Count the number of reconstructed hits in a given wire plane.
bool m_printDebug
switch for print statements (TODO: use message service!)
PFParticleMonitoring(fhicl::ParameterSet const &pset)
Constructor.
double GetLength(const art::Ptr< simb::MCParticle > trueParticle, const int startT, const int endT) const
Find the length of the true particle trajectory through the active region of the detector.
bool m_disableRealDataCheck
Whether to check if the input file contains real data before accessing MC information.
void BuildTrueNeutrinoHitMaps(const MCTruthToMCParticles &truthToParticles, const MCParticlesToHits &trueParticlesToHits, MCTruthToHits &trueNeutrinosToHits, HitsToMCTruth &trueHitsToNeutrinos) const
Build mapping from true neutrinos to hits.
int m_nTrueWithoutRecoHits
True hits which don't belong to any reconstructed particle - "available".
std::map< art::Ptr< recob::PFParticle >, T0Vector > PFParticlesToT0s
std::vector< art::Ptr< anab::T0 > > T0Vector
std::map< int, art::Ptr< recob::PFParticle > > PFParticleMap
std::vector< art::Ptr< recob::Track > > TrackVector
std::map< int, art::Ptr< simb::MCParticle > > MCParticleMap
std::map< art::Ptr< simb::MCParticle >, art::Ptr< simb::MCTruth > > MCParticlesToMCTruth
std::vector< art::Ptr< recob::Vertex > > VertexVector
std::map< art::Ptr< simb::MCTruth >, MCParticleVector > MCTruthToMCParticles
std::map< art::Ptr< simb::MCTruth >, art::Ptr< recob::PFParticle > > MCTruthToPFParticles
std::map< art::Ptr< recob::PFParticle >, HitVector > PFParticlesToHits
std::vector< art::Ptr< recob::Hit > > HitVector
std::vector< art::Ptr< recob::PFParticle > > PFParticleVector
std::vector< art::Ptr< simb::MCParticle > > MCParticleVector
std::map< art::Ptr< simb::MCParticle >, art::Ptr< recob::PFParticle > > MCParticlesToPFParticles
std::map< art::Ptr< recob::PFParticle >, TrackVector > PFParticlesToTracks
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > SpacePointsToHits
std::map< art::Ptr< simb::MCParticle >, HitVector > MCParticlesToHits
std::map< art::Ptr< simb::MCTruth >, HitVector > MCTruthToHits
std::vector< art::Ptr< recob::SpacePoint > > SpacePointVector
std::map< art::Ptr< recob::Hit >, art::Ptr< simb::MCTruth > > HitsToMCTruth
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::PFParticle > > HitsToPFParticles
std::map< art::Ptr< recob::Hit >, art::Ptr< simb::MCParticle > > HitsToMCParticles
std::map< art::Ptr< recob::PFParticle >, VertexVector > PFParticlesToVertices
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > HitsToSpacePoints