22 #include "art_root_io/TFileService.h" 24 #include "canvas/Persistency/Common/FindManyP.h" 41 #include "nug4/ParticleNavigation/ParticleList.h" 61 std::map<std::pair<int, int>, std::pair<double, double>>& g4RecoBaseIDToPurityEfficiency);
64 art::Handle<std::vector<recob::Cluster>>
const& clscol,
72 art::Handle<std::vector<recob::Shower>>
const& scol,
76 art::Handle<std::vector<recob::Vertex>>
const& vtxcol,
80 art::Handle<std::vector<recob::Event>>
const& evtcol,
88 std::map<std::pair<int, int>, std::pair<double, double>>
const& g4RecoBaseIDToPurityEfficiency,
89 std::map<
int,
std::vector<std::pair<
int, std::pair<double, double>>>>&
90 g4IDToRecoBasePurityEfficiency,
93 TH1D* purityEfficiency,
94 TH2D* purityEfficiency2D);
181 mf::LogWarning(
"RecoVetter") <<
"attempting to run MC truth check on " 182 <<
"real data, bail";
189 std::vector<art::Ptr<recob::Hit>> allhits;
236 fEventPurity = tfs->make<TH1D>(
"eventPurity",
";Purity;Events", 100, 0., 1.1);
237 fEventEfficiency = tfs->make<TH1D>(
"eventEfficiency",
";Efficiency;Events", 100, 0., 1.1);
239 tfs->make<TH1D>(
"eventPurityEfficiency",
";purityEfficiency;Events", 110, 0., 1.1);
242 fVertexPurity = tfs->make<TH1D>(
"vertexPurity",
";Purity;Vertices", 100, 0., 1.1);
243 fVertexEfficiency = tfs->make<TH1D>(
"vertexEfficiency",
";Efficiency;Vertices", 100, 0., 1.1);
245 tfs->make<TH1D>(
"vertexPurityEfficiency",
";purityEfficiency;Vertex", 110, 0., 1.1);
248 fTrackPurity = tfs->make<TH1D>(
"trackPurity",
";Purity;Tracks", 100, 0., 1.1);
249 fTrackEfficiency = tfs->make<TH1D>(
"trackEfficiency",
";Efficiency;Tracks", 100, 0., 1.1);
251 tfs->make<TH1D>(
"trackPurityEfficiency",
";purityEfficiency;Tracks", 110, 0., 1.1);
253 tfs->make<TH2D>(
"trackPurityEfficiency2D",
";purity;efficiency", 110, 0., 1.1, 110, 0., 1.1);
256 fShowerPurity = tfs->make<TH1D>(
"showerPurity",
";Purity;Showers", 100, 0., 1.1);
257 fShowerEfficiency = tfs->make<TH1D>(
"showerEfficiency",
";Efficiency;Showers", 100, 0., 1.1);
259 tfs->make<TH1D>(
"showerPurityEfficiency",
";purityEfficiency;Showers", 110, 0., 1.1);
261 tfs->make<TH2D>(
"showerPurityEfficiency2D",
";purity;efficiency", 110, 0., 1.1, 110, 0., 1.1);
264 fClusterPurity = tfs->make<TH1D>(
"clusterPurity",
";Purity;Clusters", 110, 0., 1.1);
265 fClusterEfficiency = tfs->make<TH1D>(
"clusterEfficiency",
";Efficiency;Clusters", 110, 0., 1.1);
267 tfs->make<TH1D>(
"clusterPurityEfficiency",
";purityEfficiency;Clusters", 110, 0., 1.1);
269 "clusterPurityEfficiency2D",
";purity;efficiency", 110, 0., 1.1, 110, 0., 1.1);
272 fTree = tfs->make<TTree>(
"cheatertree",
"cheater tree");
304 std::map<std::pair<int, int>, std::pair<double, double>>& g4RecoBaseIDToPurityEfficiency)
313 while (itr != trackIDs.end()) {
324 std::pair<double, double> pe(purity, efficiency);
327 std::pair<int, int> g4reco(*itr, colID);
330 g4RecoBaseIDToPurityEfficiency[g4reco] = pe;
343 art::Handle<std::vector<recob::Cluster>>
const& clscol,
347 art::FindManyP<recob::Hit> fmh(clscol, evt, label);
349 for (
size_t c = 0;
c < clscol->size(); ++
c) {
352 std::vector<art::Ptr<recob::Hit>> hits = fmh.at(
c);
365 art::Handle<std::vector<recob::Track>>
const& tcol,
369 art::FindManyP<recob::Hit> fmh(tcol, evt, label);
371 for (
size_t p = 0;
p < tcol->size(); ++
p) {
374 std::vector<art::Ptr<recob::Hit>> hits = fmh.at(
p);
387 art::Handle<std::vector<recob::Shower>>
const& scol,
391 art::FindManyP<recob::Hit> fmh(scol, evt, label);
393 for (
size_t p = 0;
p < scol->size(); ++
p) {
396 std::vector<art::Ptr<recob::Hit>> hits = fmh.at(
p);
411 art::Handle<std::vector<recob::Vertex>>
const& vtxcol,
416 std::vector<std::set<int>> ids(1);
421 for (
const auto& PartPair : plist) {
422 auto trackID = PartPair.first;
423 if (!plist.IsPrimary(trackID))
continue;
425 ids[0].insert(trackID);
434 art::FindManyP<recob::Hit> fmh(vtxcol, evt, label);
438 for (
size_t v = 0; v < vtxcol->size(); ++v) {
441 std::vector<art::Ptr<recob::Hit>> hits = fmh.at(v);
443 double maxPurity = -1.;
444 double maxEfficiency = -1.;
446 for (
size_t tv = 0; tv < ids.size(); ++tv) {
453 if (purity > maxPurity) maxPurity = purity;
454 if (efficiency > maxEfficiency) maxEfficiency = efficiency;
473 art::Handle<std::vector<recob::Event>>
const& evtcol,
481 for (
const auto& PartPair : plist) {
482 auto trackID = PartPair.first;
483 if (!plist.IsPrimary(trackID))
continue;
490 art::FindManyP<recob::Hit> fmh(evtcol, evt, label);
493 for (
size_t ev = 0; ev < evtcol->size(); ++ev) {
496 std::vector<art::Ptr<recob::Hit>> hits = fmh.at(ev);
515 std::map<std::pair<int, int>, std::pair<double, double>>
const& g4RecoBaseIDToPurityEfficiency,
516 std::map<
int,
std::vector<std::pair<
int, std::pair<double, double>>>>&
517 g4IDToRecoBasePurityEfficiency,
520 TH1D* purityEfficiency,
521 TH2D* purityEfficiency2D)
524 std::map<std::pair<int, int>, std::pair<double, double>>
::const_iterator rbItr =
525 g4RecoBaseIDToPurityEfficiency.begin();
528 std::map<int, std::pair<double, double>> recoBIDToPurityEfficiency;
529 std::map<int, std::pair<double, double>>
::iterator rbpeItr;
531 while (rbItr != g4RecoBaseIDToPurityEfficiency.end()) {
534 std::pair<int, int> g4cl = rbItr->first;
536 std::pair<double, double> pe = rbItr->second;
541 std::pair<int, std::pair<double, double>> clpe(g4cl.second, pe);
543 g4IDToRecoBasePurityEfficiency[g4cl.first].push_back(clpe);
547 rbpeItr = recoBIDToPurityEfficiency.find(g4cl.second);
548 if (rbpeItr != recoBIDToPurityEfficiency.end()) {
549 std::pair<double, double> curpe = rbpeItr->second;
550 if (pe.first > curpe.first) recoBIDToPurityEfficiency[g4cl.second] = pe;
553 recoBIDToPurityEfficiency[g4cl.second] = pe;
558 rbpeItr = recoBIDToPurityEfficiency.begin();
561 while (rbpeItr != recoBIDToPurityEfficiency.end()) {
562 purity->Fill(rbpeItr->second.first);
563 efficiency->Fill(rbpeItr->second.second);
564 purityEfficiency->Fill(rbpeItr->second.first * rbpeItr->second.second);
565 purityEfficiency2D->Fill(rbpeItr->second.first, rbpeItr->second.second);
578 std::map<int, double> g4IDToHitEnergy;
579 for (
size_t h = 0;
h < allhits.size(); ++
h) {
580 const std::vector<sim::TrackIDE> hitTrackIDs =
fBT->
HitToTrackIDEs(clockData, allhits[
h]);
581 for (
size_t e = 0;
e < hitTrackIDs.size(); ++
e) {
582 g4IDToHitEnergy[hitTrackIDs[
e].trackID] += hitTrackIDs[
e].energy;
588 std::map<int, std::vector<std::pair<int, std::pair<double, double>>>>
589 g4IDToClusterPurityEfficiency;
590 std::map<int, std::vector<std::pair<int, std::pair<double, double>>>>
591 g4IDToShowerPurityEfficiency;
592 std::map<int, std::vector<std::pair<int, std::pair<double, double>>>> g4IDToTrackPurityEfficiency;
593 std::map<int, std::vector<std::pair<int, std::pair<double, double>>>>
::iterator g4peItr;
597 g4IDToClusterPurityEfficiency,
604 g4IDToShowerPurityEfficiency,
611 g4IDToTrackPurityEfficiency,
623 while (trackItr != trackIDs.end()) {
633 double totalDep = 0.;
634 for (
size_t i = 0; i < ides.size(); ++i)
635 totalDep += ides[i]->
energy;
637 if (totalDep > 0.)
fhiteff = g4IDToHitEnergy[*trackItr] / totalDep;
639 std::vector<std::pair<int, std::pair<double, double>>> clVec;
640 std::vector<std::pair<int, std::pair<double, double>>> shVec;
641 std::vector<std::pair<int, std::pair<double, double>>> trVec;
643 if (g4IDToClusterPurityEfficiency.find(*trackItr) != g4IDToClusterPurityEfficiency.end())
644 clVec = g4IDToClusterPurityEfficiency.find(*trackItr)->second;
646 if (g4IDToShowerPurityEfficiency.find(*trackItr) != g4IDToShowerPurityEfficiency.end())
647 shVec = g4IDToShowerPurityEfficiency.find(*trackItr)->second;
649 if (g4IDToTrackPurityEfficiency.find(*trackItr) != g4IDToTrackPurityEfficiency.end())
650 trVec = g4IDToTrackPurityEfficiency.find(*trackItr)->second;
652 fnclu = clVec.size();
653 fnshw = shVec.size();
654 fntrk = trVec.size();
656 for (
size_t c = 0;
c < clVec.size(); ++
c) {
657 fcluid.push_back(clVec[
c].first);
662 for (
size_t s = 0;
s < shVec.size(); ++
s) {
663 fshwid.push_back(shVec[
s].first);
668 for (
size_t t = 0;
t < trVec.size(); ++
t) {
669 ftrkid.push_back(trVec[
t].first);
art::ServiceHandle< cheat::BackTrackerService const > fBT
the back tracker service
int fpdg
particle pdg code
int fnshw
number of showers for this particle
void analyze(art::Event const &e) override
TH1D * fTrackPurityEfficiency
histogram of track efficiency times purity
bool fCheckVertices
should we check the reconstruction of vertices?
std::vector< sim::TrackIDE > HitToTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
const simb::MCParticle * TrackIdToParticle_P(int id) const
double fhiteff
hitfinder efficiency for this particle
TH1D * fClusterPurityEfficiency
histogram of cluster efficiency times purity
int fnclu
number of clusters for this particle
double HitCollectionEfficiency(detinfo::DetectorClocksData const &clockData, std::set< int > const &trackIds, std::vector< art::Ptr< recob::Hit >> const &hits, std::vector< art::Ptr< recob::Hit >> const &allhits, geo::View_t const &view) const
void CheckReco(detinfo::DetectorClocksData const &clockData, int const &colID, std::vector< art::Ptr< recob::Hit >> const &allhits, std::vector< art::Ptr< recob::Hit >> const &colHits, std::map< std::pair< int, int >, std::pair< double, double >> &g4RecoBaseIDToPurityEfficiency)
void CheckRecoEvents(art::Event const &evt, std::string const &label, art::Handle< std::vector< recob::Event >> const &evtcol, std::vector< art::Ptr< recob::Hit >> const &allhits)
art::ServiceHandle< cheat::ParticleInventoryService const > fPI
the back tracker service
int fntrk
number of tracks for this particle
std::vector< const sim::IDE * > TrackIdToSimIDEs_Ps(int const &id) const
RecoCheckAna(fhicl::ParameterSet const &p)
TH1D * fEventPurity
histogram of event purity
bool fCheckEvents
should we check the reconstruction of events?
std::map< std::pair< int, int >, std::pair< double, double > > fG4ClusterIDToPurityEfficiency
TH2D * fTrackPurityEfficiency2D
scatter histogram of cluster purity and efficiency
void beginRun(art::Run const &r) override
EDAnalyzer(fhicl::ParameterSet const &pset)
TH1D * fEventEfficiency
histogram of event efficiency
int NumberDaughters() const
TH1D * fVertexPurity
histogram of vertex purity
int Daughter(const int i) const
TH1D * fEventPurityEfficiency
histogram of event efficiency times purity
3-dimensional objects, potentially hits, clusters, prongs, etc.
TH1D * fVertexPurityEfficiency
histogram of vertex efficiency times purity
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
void CheckRecoTracks(art::Event const &evt, std::string const &label, art::Handle< std::vector< recob::Track >> const &tcol, std::vector< art::Ptr< recob::Hit >> const &allhits)
#define DEFINE_ART_MODULE(klass)
std::set< int > GetSetOfTrackIds() const
std::vector< double > fclupur
cluster purities
std::vector< double > fshweff
shower efficiencies
double HitCollectionPurity(detinfo::DetectorClocksData const &clockData, std::set< int > const &trackIds, std::vector< art::Ptr< recob::Hit >> const &hits) const
bool fCheckClusters
should we check the reconstruction of clusters?
double P(const int i=0) const
T get(std::string const &key) const
std::vector< double > fshwpur
shower purities
TTree * fTree
TTree to save efficiencies.
std::vector< int > fshwid
shower IDs
TH1D * fVertexEfficiency
histogram of vertex efficiency
std::vector< int > ftrkid
track IDs
std::vector< double > ftrkeff
track efficiencies
bool fCheckTracks
should we check the reconstruction of tracks?
std::vector< double > fclueff
cluster efficiencies
std::string fTrackModuleLabel
label for module making the tracks
double fpmom
particle momentum
void FlattenMap(std::map< std::pair< int, int >, std::pair< double, double >> const &g4RecoBaseIDToPurityEfficiency, std::map< int, std::vector< std::pair< int, std::pair< double, double >>>> &g4IDToRecoBasePurityEfficiency, TH1D *purity, TH1D *efficiency, TH1D *purityEfficiency, TH2D *purityEfficiency2D)
std::vector< int > fcluid
cluster IDs
std::vector< double > ftrkpur
track purities
Definition of data types for geometry description.
bool fCheckShowers
should we check the reconstruction of showers?
const sim::ParticleList & ParticleList() const
std::map< std::pair< int, int >, std::pair< double, double > > fG4TrackIDToPurityEfficiency
void CheckRecoShowers(art::Event const &evt, std::string const &label, art::Handle< std::vector< recob::Shower >> const &scol, std::vector< art::Ptr< recob::Hit >> const &allhits)
Declaration of signal hit object.
code to link reconstructed objects back to the MC truth information
int ftrackid
geant track ID
Contains all timing reference information for the detector.
TH1D * fShowerPurityEfficiency
histogram of shower efficiency times purity
TH1D * fClusterPurity
histogram of cluster purity
std::map< std::pair< int, int >, std::pair< double, double > > fG4ShowerIDToPurityEfficiency
void FillResults(detinfo::DetectorClocksData const &clockData, std::vector< art::Ptr< recob::Hit >> const &allhits)
void CheckRecoVertices(art::Event const &evt, std::string const &label, art::Handle< std::vector< recob::Vertex >> const &vtxcol, std::vector< art::Ptr< recob::Hit >> const &allhits)
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
Provides recob::Track data product.
std::string fClusterModuleLabel
label for module making the clusters
EventNumber_t event() const
TH1D * fClusterEfficiency
histogram of cluster efficiency
TH1D * fShowerPurity
histogram of shower purity
void CheckRecoClusters(art::Event const &evt, std::string const &label, art::Handle< std::vector< recob::Cluster >> const &clscol, std::vector< art::Ptr< recob::Hit >> const &allhits)
std::string fShowerModuleLabel
label for module making the showers
TH1D * fShowerEfficiency
histogram of shower efficiency
std::string fEventModuleLabel
label for module making the events
TH1D * fTrackPurity
histogram of track purity
TH2D * fClusterPurityEfficiency2D
scatter histogram of cluster purity and efficiency
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
std::string fVertexModuleLabel
label for module making the vertices
second_as<> second
Type of time stored in seconds, in double precision.
TH1D * fTrackEfficiency
histogram of track efficiency
TH2D * fShowerPurityEfficiency2D
scatter histogram of cluster purity and efficiency
std::string fHitModuleLabel
label for module making the hits