8#include "art/Framework/Core/EDAnalyzer.h"
9#include "art/Framework/Core/ModuleMacros.h"
11#include "larcoreobj/SimpleTypesAndConstants/geo_types.h"
39 void analyze(
const art::Event& evt);
117 typedef std::map<SimpleMCPrimary, SimpleMatchedPfoList>
339#include "art/Framework/Principal/Event.h"
341#include "fhiclcpp/ParameterSet.h"
343#include "nusimdata/SimulationBase/MCParticle.h"
344#include "nusimdata/SimulationBase/MCTruth.h"
346#include "lardataobj/RecoBase/Hit.h"
347#include "lardataobj/RecoBase/PFParticle.h"
355 :
art::EDAnalyzer(pset)
368 m_particleLabel = pset.get<std::string>(
"PFParticleModule",
"pandoraNu");
414 if (hitsToMCParticles.empty()) {
416 throw cet::exception(
"LArPandora") <<
" PFParticleValidation::analyze - no sim channels "
417 "found, backtracker module must be set in FHiCL "
431 pfParticlesToHits, mcParticlesToHits, hitsToMCParticles, mcParticleMatchingMap);
435 evt, mcParticlesToHits, hitsToMCParticles, mcParticleMatchingMap, simpleMCPrimaryList);
439 simpleMCPrimaryList, mcParticleMatchingMap, pfParticlesToHits, mcPrimaryMatchingMap);
448 this->
PrintAllOutput(mcTruthVector, recoNeutrinoVector, mcPrimaryMatchingMap);
466 for (
const MCParticlesToHits::value_type& mcParticleToHitsEntry : mcParticlesToHits) {
467 if (!mcParticleToHitsEntry.second.empty())
468 (
void)mcParticleMatchingMap.insert(MCParticleMatchingMap::value_type(
473 for (
const PFParticlesToHits::value_type& recoParticleToHits : pfParticlesToHits) {
474 const art::Ptr<recob::PFParticle> pRecoParticle(recoParticleToHits.first);
475 const HitVector& hitVector(recoParticleToHits.second);
477 for (
const art::Ptr<recob::Hit> pHit : hitVector) {
478 HitsToMCParticles::const_iterator mcParticleIter = hitsToMCParticles.find(pHit);
480 if (hitsToMCParticles.end() == mcParticleIter)
continue;
482 const art::Ptr<simb::MCParticle> pTrueParticle = mcParticleIter->second;
483 mcParticleMatchingMap[pTrueParticle][pRecoParticle].push_back(pHit);
491 const art::Event& evt,
502 for (
const MCParticlesToHits::value_type& mapEntry : mcParticlesToHits) {
503 const art::Ptr<simb::MCParticle> pMCPrimary(mapEntry.first);
510 simpleMCPrimary.
m_pAddress = pMCPrimary.get();
511 simpleMCPrimary.
m_pdgCode = pMCPrimary->PdgCode();
512 simpleMCPrimary.
m_energy = pMCPrimary->E();
514 MCParticlesToHits::const_iterator trueHitsIter = mcParticlesToHits.find(pMCPrimary);
516 if (mcParticlesToHits.end() != trueHitsIter) {
517 const HitVector& hitVector(trueHitsIter->second);
524 MCParticleMatchingMap::const_iterator matchedPfoIter = mcParticleMatchingMap.find(pMCPrimary);
526 if (mcParticleMatchingMap.end() != matchedPfoIter)
529 simpleMCPrimaryList.push_back(simpleMCPrimary);
532 std::sort(simpleMCPrimaryList.begin(),
533 simpleMCPrimaryList.end(),
538 simpleMCPrimary.m_id = mcPrimaryId++;
544 const art::Ptr<simb::MCParticle> pMCParticle,
547 MCParticlesToMCTruth::const_iterator iter = artMCParticlesToMCTruth.find(pMCParticle);
549 if (artMCParticlesToMCTruth.end() == iter)
return false;
551 const art::Ptr<simb::MCTruth> pMCTruth = iter->second;
552 return (simb::kBeamNeutrino == pMCTruth->Origin());
565 MCParticleMatchingMap::const_iterator matchedPfoIter = mcParticleMatchingMap.end();
568 for (MCParticleMatchingMap::const_iterator iter = mcParticleMatchingMap.begin(),
569 iterEnd = mcParticleMatchingMap.end();
572 if (simpleMCPrimary.m_pAddress == iter->first.get()) {
573 matchedPfoIter = iter;
578 if (mcParticleMatchingMap.end() != matchedPfoIter) {
579 for (
const PFParticleToMatchedHits::value_type& contribution : matchedPfoIter->second) {
580 const art::Ptr<recob::PFParticle> pMatchedPfo(contribution.first);
581 const HitVector& matchedHitVector(contribution.second);
584 simpleMatchedPfo.
m_pAddress = pMatchedPfo.get();
585 simpleMatchedPfo.
m_id = pMatchedPfo->Self();
588 PFParticlesToHits::const_iterator parentPfoIter = pfParticlesToHits.end();
591 for (PFParticlesToHits::const_iterator iter = pfParticlesToHits.begin(),
592 iterEnd = pfParticlesToHits.end();
595 if (pMatchedPfo->Parent() == iter->first->Self()) {
596 parentPfoIter = iter;
601 if ((pfParticlesToHits.end() != parentPfoIter) &&
603 simpleMatchedPfo.
m_parentId = parentPfoIter->first->Self();
605 simpleMatchedPfo.
m_pdgCode = pMatchedPfo->PdgCode();
611 PFParticlesToHits::const_iterator pfoHitsIter = pfParticlesToHits.find(pMatchedPfo);
613 if (pfParticlesToHits.end() == pfoHitsIter)
614 throw cet::exception(
"LArPandora")
615 <<
" PFParticleValidation::analyze --- Presence of PFParticle in map mandatory.";
617 const HitVector& pfoHitVector(pfoHitsIter->second);
624 simpleMatchedPfoList.push_back(simpleMatchedPfo);
629 std::sort(simpleMatchedPfoList.begin(),
630 simpleMatchedPfoList.end(),
633 if (!mcPrimaryMatchingMap
634 .insert(MCPrimaryMatchingMap::value_type(simpleMCPrimary, simpleMatchedPfoList))
636 throw cet::exception(
"LArPandora")
637 <<
" PFParticleValidation::analyze --- Double-counting MC primaries.";
650 for (
const auto& mapEntry : artMCTruthToMCParticles) {
651 const art::Ptr<simb::MCTruth> truth = mapEntry.first;
653 if (!truth->NeutrinoSet())
continue;
655 if (mcTruthVector.end() == std::find(mcTruthVector.begin(), mcTruthVector.end(), truth))
656 mcTruthVector.push_back(truth);
676 std::cout <<
"---RAW-MATCHING-OUTPUT-----------------------------------------------------------"
680 for (
const art::Ptr<simb::MCTruth> pMCTruth : mcTruthVector) {
681 std::cout <<
"MCNeutrino, PDG " << pMCTruth->GetNeutrino().Nu().PdgCode()
682 <<
", InteractionType " << pMCTruth->GetNeutrino().InteractionType() << std::endl;
685 for (
const art::Ptr<recob::PFParticle> pPfo : recoNeutrinoVector) {
686 std::cout <<
"RecoNeutrino, PDG " << pPfo->PdgCode() <<
", nDaughters "
687 << pPfo->NumDaughters() << std::endl;
690 for (
const MCPrimaryMatchingMap::value_type& mapValue : mcPrimaryMatchingMap) {
693 std::cout << std::endl
694 <<
"Primary " << simpleMCPrimary.
m_id <<
", PDG " << simpleMCPrimary.
m_pdgCode
697 << simpleMCPrimary.
m_nMCHitsW <<
")" << std::endl;
700 std::cout <<
"-MatchedPfo " << simpleMatchedPfo.m_id;
702 if (simpleMatchedPfo.m_parentId >= 0)
703 std::cout <<
", ParentPfo " << simpleMatchedPfo.m_parentId;
705 std::cout <<
", PDG " << simpleMatchedPfo.m_pdgCode <<
", nMatchedHits "
706 << simpleMatchedPfo.m_nMatchedHitsTotal <<
" ("
707 << simpleMatchedPfo.m_nMatchedHitsU <<
", " << simpleMatchedPfo.m_nMatchedHitsV
708 <<
", " << simpleMatchedPfo.m_nMatchedHitsW <<
")"
709 <<
", nPfoHits " << simpleMatchedPfo.m_nPfoHitsTotal <<
" ("
710 << simpleMatchedPfo.m_nPfoHitsU <<
", " << simpleMatchedPfo.m_nPfoHitsV <<
", "
711 << simpleMatchedPfo.m_nPfoHitsW <<
")" << std::endl;
715 std::cout <<
"---------------------------------------------------------------------------------"
726 IntSet usedMCIds, usedPfoIds;
727 while (
GetStrongestPfoMatch(mcPrimaryMatchingMap, usedMCIds, usedPfoIds, matchingDetailsMap)) {}
740 int bestPfoMatchId(-1);
743 for (
const MCPrimaryMatchingMap::value_type& mapValue : mcPrimaryMatchingMap) {
748 if (usedMCIds.count(simpleMCPrimary.
m_id))
continue;
751 if (usedPfoIds.count(simpleMatchedPfo.m_id))
continue;
753 if (!this->
IsGoodMatch(simpleMCPrimary, simpleMatchedPfo))
continue;
755 if (simpleMatchedPfo.m_nMatchedHitsTotal > bestMatchingDetails.
m_nMatchedHits) {
756 bestPfoMatchId = simpleMatchedPfo.m_id;
758 bestMatchingDetails.
m_nMatchedHits = simpleMatchedPfo.m_nMatchedHitsTotal;
760 static_cast<float>(simpleMatchedPfo.m_nMatchedHitsTotal) /
766 if (bestPfoMatchId > -1) {
767 matchingDetailsMap[bestPfoMatchId] = bestMatchingDetails;
769 usedPfoIds.insert(bestPfoMatchId);
783 for (
const MCPrimaryMatchingMap::value_type& mapValue : mcPrimaryMatchingMap) {
789 if (usedPfoIds.count(simpleMatchedPfo.m_id))
continue;
791 MatchingDetails& matchingDetails(matchingDetailsMap[simpleMatchedPfo.m_id]);
793 if (simpleMatchedPfo.m_nMatchedHitsTotal > matchingDetails.
m_nMatchedHits) {
795 matchingDetails.
m_nMatchedHits = simpleMatchedPfo.m_nMatchedHitsTotal;
797 static_cast<float>(simpleMatchedPfo.m_nMatchedHitsTotal) /
809 std::cout <<
"---PROCESSED-MATCHING-OUTPUT-----------------------------------------------------"
819 bool isCorrect(
true), isCalculable(
false);
821 for (
const MCPrimaryMatchingMap::value_type& mapValue : mcPrimaryMatchingMap) {
823 const bool hasMatch(this->
HasMatch(simpleMCPrimary, mapValue.second, matchingDetailsMap));
827 if (!hasMatch && !isTargetPrimary)
continue;
829 std::cout << std::endl
830 << (!isTargetPrimary ?
"(Non target) " :
"") <<
"Primary " << simpleMCPrimary.
m_id
831 <<
", PDG " << simpleMCPrimary.
m_pdgCode <<
", nMCHits "
836 if (2112 != simpleMCPrimary.
m_pdgCode) isCalculable =
true;
838 unsigned int nMatches(0);
841 if (matchingDetailsMap.count(simpleMatchedPfo.m_id) &&
842 (simpleMCPrimary.
m_id ==
843 matchingDetailsMap.at(simpleMatchedPfo.m_id).m_matchedPrimaryId)) {
844 const bool isGoodMatch(this->
IsGoodMatch(simpleMCPrimary, simpleMatchedPfo));
846 if (isGoodMatch) ++nMatches;
847 std::cout <<
"-" << (!isGoodMatch ?
"(Below threshold) " :
"") <<
"MatchedPfo "
848 << simpleMatchedPfo.m_id;
850 if (simpleMatchedPfo.m_parentId >= 0)
851 std::cout <<
", ParentPfo " << simpleMatchedPfo.m_parentId;
853 std::cout <<
", PDG " << simpleMatchedPfo.m_pdgCode <<
", nMatchedHits "
854 << simpleMatchedPfo.m_nMatchedHitsTotal <<
" ("
855 << simpleMatchedPfo.m_nMatchedHitsU <<
", " << simpleMatchedPfo.m_nMatchedHitsV
856 <<
", " << simpleMatchedPfo.m_nMatchedHitsW <<
")"
857 <<
", nPfoHits " << simpleMatchedPfo.m_nPfoHitsTotal <<
" ("
858 << simpleMatchedPfo.m_nPfoHitsU <<
", " << simpleMatchedPfo.m_nPfoHitsV <<
", "
859 << simpleMatchedPfo.m_nPfoHitsW <<
")" << std::endl;
863 if (isTargetPrimary && (1 != nMatches)) isCorrect =
false;
866 std::cout << std::endl <<
"Is correct? " << (isCorrect && isCalculable) << std::endl;
867 std::cout <<
"---------------------------------------------------------------------------------"
895 if (matchingDetailsMap.count(simpleMatchedPfo.m_id) &&
896 (simpleMCPrimary.
m_id == matchingDetailsMap.at(simpleMatchedPfo.m_id).m_matchedPrimaryId))
926 unsigned int nHitsOfSpecifiedType(0);
928 for (
const art::Ptr<recob::Hit> pHit : hitVector) {
929 if (view == pHit->View()) ++nHitsOfSpecifiedType;
932 return nHitsOfSpecifiedType;
971 , m_pAddress(nullptr)
978 if (
this == &rhs)
return false;
980 return (m_id < rhs.
m_id);
994 , m_nMatchedHitsTotal(0)
998 , m_pAddress(nullptr)
1005 : m_matchedPrimaryId(-1), m_nMatchedHits(0), m_completeness(0.f)
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 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 void CollectPFParticles(const art::Event &evt, const std::string &label, PFParticleVector &particleVector)
Collect the reconstructed PFParticles 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.
float m_completeness
The completeness of the match.
int m_matchedPrimaryId
The total number of occurences.
int m_nMatchedHits
The number of times the primary has 0 pfo matches.
MatchingDetails()
Default constructor.
int m_pdgCode
The pdg code.
int m_nMCHitsU
The number of u mc hits.
int m_nMCHitsV
The number of v mc hits.
bool operator<(const SimpleMCPrimary &rhs) const
operator <
int m_nMCHitsTotal
The total number of mc hits.
int m_nMCHitsW
The number of w mc hits.
const simb::MCParticle * m_pAddress
The address of the mc primary.
int m_id
The unique identifier.
float m_energy
The energy.
SimpleMCPrimary()
Constructor.
int m_nMatchedPfos
The number of matched pfos.
int m_nMatchedHitsV
The number of v matched hits.
int m_nMatchedHitsU
The number of u matched hits.
int m_nPfoHitsV
The number of v pfo hits.
int m_nPfoHitsW
The number of w pfo hits.
int m_nPfoHitsU
The number of u pfo hits.
int m_nMatchedHitsTotal
The total number of matched hits.
int m_parentId
The unique identifier of the parent pfo (-1 if no parent set)
int m_id
The unique identifier.
const recob::PFParticle * m_pAddress
The address of the pf primary.
int m_nPfoHitsTotal
The total number of pfo hits.
int m_nMatchedHitsW
The number of w matched hits.
SimpleMatchedPfo()
Constructor.
int m_pdgCode
The pdg code.
PFParticleValidation class.
std::string m_hitfinderLabel
The name/label of the hit producer module.
int m_matchingMinHitsForGoodView
The minimum number of good mc primary hits in given view to declare view to be good.
std::map< SimpleMCPrimary, SimpleMatchedPfoList > MCPrimaryMatchingMap
unsigned int CountHitsByType(const geo::View_t view, const HitVector &hitVector) const
Count the number of hits, in a provided vector, of a specified view.
bool m_useSmallPrimaries
Whether to consider matches to mc primaries with fewer than m_matchingMinPrimaryHits.
void GetMCTruth(const art::Event &evt, MCTruthVector &mcTruthVector) const
Obtain a vector of mc truth.
int m_matchingMinPrimaryHits
The minimum number of good mc primary hits used in matching scheme.
void PrintMatchingOutput(const MCPrimaryMatchingMap &mcPrimaryMatchingMap, const MatchingDetailsMap &matchingDetailsMap) const
Print the results of the matching procedure.
std::map< art::Ptr< simb::MCParticle >, PFParticleToMatchedHits > MCParticleMatchingMap
int m_matchingMinPrimaryGoodViews
The minimum number of good views for a mc primary.
void GetRecoNeutrinos(const art::Event &evt, PFParticleVector &recoNeutrinoVector) const
Obtain a vector of reco neutrinos.
std::vector< SimpleMCPrimary > SimpleMCPrimaryList
void GetMCPrimaryMatchingMap(const SimpleMCPrimaryList &simpleMCPrimaryList, const MCParticleMatchingMap &mcParticleMatchingMap, const PFParticlesToHits &pfParticlesToHits, MCPrimaryMatchingMap &mcPrimaryMatchingMap) const
Obtain a sorted list of matched pfos for each mc primary.
bool HasMatch(const SimpleMCPrimary &simpleMCPrimary, const SimpleMatchedPfoList &simpleMatchedPfoList, const MatchingDetailsMap &matchingDetailsMap) const
Whether a provided mc primary has a match, of any quality (use simple matched pfo list and informatio...
void PerformMatching(const MCPrimaryMatchingMap &mcPrimaryMatchingMap, MatchingDetailsMap &matchingDetailsMap) const
Apply a well-defined matching procedure to the comprehensive matches in the provided mc primary match...
float m_matchingMinCompleteness
The minimum particle completeness to declare a match.
float m_matchingMinPurity
The minimum particle purity to declare a match.
bool m_neutrinoInducedOnly
Whether to consider only mc particles that were neutrino induced.
bool GetStrongestPfoMatch(const MCPrimaryMatchingMap &mcPrimaryMatchingMap, IntSet &usedMCIds, IntSet &usedPfoIds, MatchingDetailsMap &matchingDetailsMap) const
Get the strongest pfo match (most matched hits) between an available mc primary and an available pfo.
std::vector< SimpleMatchedPfo > SimpleMatchedPfoList
std::map< art::Ptr< recob::PFParticle >, HitVector > PFParticleToMatchedHits
void analyze(const art::Event &evt)
static bool SortSimpleMCPrimaries(const SimpleMCPrimary &lhs, const SimpleMCPrimary &rhs)
Sort simple mc primaries by number of mc hits.
void GetMCParticleMatchingMap(const PFParticlesToHits &recoParticlesToHits, const MCParticlesToHits &trueParticlesToHits, const HitsToMCParticles &hitsToTrueParticles, MCParticleMatchingMap &mcParticleMatchingMap) const
Performing matching between true and reconstructed particles.
bool m_printAllToScreen
Whether to print all/raw matching details to screen.
std::string m_backtrackerLabel
The name/label of the back-tracker module.
static bool SortSimpleMatchedPfos(const SimpleMatchedPfo &lhs, const SimpleMatchedPfo &rhs)
Sort simple matched pfos by number of matched hits.
std::string m_particleLabel
The name/label of the particle producer module.
virtual ~PFParticleValidation()
Destructor.
void GetSimpleMCPrimaryList(const art::Event &evt, const MCParticlesToHits &mcParticlesToHits, const HitsToMCParticles &hitsToMCParticles, const MCParticleMatchingMap &mcParticleMatchingMap, SimpleMCPrimaryList &simpleMCPrimaryList) const
Extract details of each mc primary (ordered by number of true hits)
bool IsGoodMCPrimary(const SimpleMCPrimary &simpleMCPrimary) const
Whether a provided mc primary passes selection, based on number of "good" hits.
void PrintAllOutput(const MCTruthVector &mcTruthVector, const PFParticleVector &recoNeutrinoVector, const MCPrimaryMatchingMap &mcPrimaryMatchingMap) const
Print all the raw matching output to screen.
void reconfigure(fhicl::ParameterSet const &pset)
std::string m_geantModuleLabel
The name/label of the geant module.
bool IsGoodMatch(const SimpleMCPrimary &simpleMCPrimary, const SimpleMatchedPfo &simpleMatchedPfo) const
Whether a provided mc primary and pfo are deemed to be a good match.
void GetRemainingPfoMatches(const MCPrimaryMatchingMap &mcPrimaryMatchingMap, const IntSet &usedPfoIds, MatchingDetailsMap &matchingDetailsMap) const
Get the best matches for any pfos left-over after the strong matching procedure.
std::map< int, MatchingDetails > MatchingDetailsMap
PFParticleValidation(fhicl::ParameterSet const &pset)
Constructor.
bool m_printMatchingToScreen
Whether to print matching output to screen.
int m_matchingMinSharedHits
The minimum number of shared hits used in matching scheme.
bool IsNeutrinoInduced(const art::Ptr< simb::MCParticle > pMCParticle, const MCParticlesToMCTruth &artMCParticlesToMCTruth) const
Whether a mc particle is neutrino induced.
std::map< art::Ptr< simb::MCParticle >, art::Ptr< simb::MCTruth > > MCParticlesToMCTruth
std::map< art::Ptr< simb::MCTruth >, MCParticleVector > MCTruthToMCParticles
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::MCTruth > > MCTruthVector
std::map< art::Ptr< simb::MCParticle >, HitVector > MCParticlesToHits
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::PFParticle > > HitsToPFParticles
std::map< art::Ptr< recob::Hit >, art::Ptr< simb::MCParticle > > HitsToMCParticles