LArSoft  v10_04_05
Liquid Argon Software toolkit - https://larsoft.org/
ClusteringValidation::ClusterAnalyser Class Reference

Public Member Functions

 ClusterAnalyser (std::string &label)
 
void Analyse (detinfo::DetectorClocksData const &clockData, std::vector< art::Ptr< recob::Hit >> &hits, std::vector< art::Ptr< recob::Cluster >> &clusters, const art::FindManyP< recob::Hit > &fmh, int numHits)
 
TrackID FindTrackID (detinfo::DetectorClocksData const &clockData, art::Ptr< recob::Hit > &hit)
 
TrackID FindTrueTrack (detinfo::DetectorClocksData const &clockData, std::vector< art::Ptr< recob::Hit >> &clusterHits)
 
double FindPhotonAngle ()
 
double GetEndTrackDistance (TrackID id1, TrackID id2)
 
const simb::MCParticleGetPi0 ()
 
TObjArray GetHistograms ()
 
void MakeHistograms ()
 
void WriteFile ()
 

Private Attributes

std::string fClusterLabel
 
TH1 * hCompleteness
 
TH1 * hCleanliness
 
TH1 * hComplCleanl
 
TH1 * hPi0Angle
 
TH1 * hPi0Energy
 
TH1 * hPi0ConversionDistance
 
TH1 * hPi0ConversionSeparation
 
TH1 * hPi0AngleCut
 
TH1 * hPi0EnergyCut
 
TH1 * hPi0ConversionDistanceCut
 
TH1 * hPi0ConversionSeparationCut
 
TH2 * hNumHitsCompleteness
 
TH2 * hNumHitsEnergy
 
TProfile * hCompletenessEnergy
 
TProfile * hCompletenessAngle
 
TProfile * hCompletenessConversionDistance
 
TProfile * hCompletenessConversionSeparation
 
TProfile * hCleanlinessEnergy
 
TProfile * hCleanlinessAngle
 
TProfile * hCleanlinessConversionDistance
 
TProfile * hCleanlinessConversionSeparation
 
TProfile * hComplCleanlEnergy
 
TProfile * hComplCleanlAngle
 
TProfile * hComplCleanlConversionDistance
 
TProfile * hComplCleanlConversionSeparation
 
TEfficiency * hEfficiencyAngle
 
TEfficiency * hEfficiencyEnergy
 
TEfficiency * hEfficiencyConversionDistance
 
TEfficiency * hEfficiencyConversionSeparation
 
TObjArray fHistArray
 
std::map< unsigned int, std::map< unsigned int, std::unique_ptr< ClusterCounter > > > clusterMap
 
std::map< TrackID, const simb::MCParticle * > trueParticles
 
geo::WireReadoutGeom const & wireReadoutGeom = art::ServiceHandle<geo::WireReadout>()->Get()
 
art::ServiceHandle< cheat::ParticleInventoryService const > pi_serv
 
art::ServiceHandle< cheat::BackTrackerService const > bt_serv
 

Detailed Description

Definition at line 214 of file ClusteringValidation_module.cc.

Constructor & Destructor Documentation

ClusteringValidation::ClusterAnalyser::ClusterAnalyser ( std::string &  label)
explicit

Definition at line 278 of file ClusteringValidation_module.cc.

279 {
280 
281  fClusterLabel = clusterLabel;
282 
283  // Make the histograms
284  hCompleteness = new TH1D("Completeness", ";Completeness;", 101, 0, 1.01);
286  new TProfile("CompletenessEnergy", ";True Energy (GeV);Completeness", 100, 0, 2);
288  new TProfile("CompletenessAngle", ";True Angle (deg);Completeness;", 100, 0, 180);
289  hCompletenessConversionDistance = new TProfile(
290  "CompletenessConversionDistance", ";True Distance from Vertex (cm);Completeness", 100, 0, 200);
291  hCompletenessConversionSeparation = new TProfile("CompletenessConversionSeparation",
292  ";True Conversion Separation (cm);Completeness",
293  100,
294  0,
295  200);
296  hCleanliness = new TH1D("Cleanliness", ";Cleanliness;", 101, 0, 1.01);
298  new TProfile("CleanlinessEnergy", ";True Energy (GeV);Cleanliness", 100, 0, 2);
300  new TProfile("CleanlinessAngle", ";True Angle (deg);Cleanliness;", 100, 0, 180);
301  hCleanlinessConversionDistance = new TProfile(
302  "CleanlinessConversionDistance", ";True Distance from Vertex (cm);Cleanliness", 100, 0, 200);
303  hCleanlinessConversionSeparation = new TProfile(
304  "CleanlinessConversionSeparation", ";True Conversion Separation (cm);Cleanliness", 100, 0, 200);
305  hComplCleanl = new TH1D("CompletenessCleanliness", ";Completeness * Cleanliness;", 101, 0, 1.01);
306  hComplCleanlEnergy = new TProfile(
307  "CompletenessCleanlinessEnergy", ";True Energy (GeV);Completeness * Cleanliness", 100, 0, 2);
308  hComplCleanlAngle = new TProfile(
309  "CompletenessCleanlinessAngle", ";True Angle (deg);Completeness * Cleanliness;", 100, 0, 180);
311  new TProfile("CompletenessCleanlinessConversionDistance",
312  ";True Distance from Vertex (cm);Completeness * Cleanliness",
313  100,
314  0,
315  200);
317  new TProfile("CompletenessCleanlinessConversionSeparation",
318  ";True Conversion Separation (cm);Completeness * Cleanliness",
319  100,
320  0,
321  200);
322  hPi0Energy = new TH1D("Pi0EnergyCut", ";True Energy (GeV);", 25, 0, 2);
323  hPi0Energy->Sumw2();
324  hPi0Angle = new TH1D("Pi0AngleCut", ";True Angle (deg);", 25, 0, 180);
325  hPi0Angle->Sumw2();
327  new TH1D("Pi0ConversionDistanceCut", ";True Distance from Vertex (cm);", 25, 0, 200);
328  hPi0ConversionDistance->Sumw2();
330  new TH1D("Pi0ConversionSeparationCut", ";True Separation from Vertex (cm);", 25, 0, 200);
331  hPi0ConversionSeparation->Sumw2();
332  hPi0EnergyCut = new TH1D("Pi0EnergyCut", ";True Energy (GeV);Efficiency", 25, 0, 2);
333  hPi0EnergyCut->Sumw2();
334  hPi0AngleCut = new TH1D("Pi0AngleCut", ";True Angle (deg);Efficiency", 25, 0, 180);
335  hPi0AngleCut->Sumw2();
337  new TH1D("Pi0ConversionDistanceCut", ";True Distance from Vertex (cm);Efficiency", 25, 0, 200);
338  hPi0ConversionDistanceCut->Sumw2();
339  hPi0ConversionSeparationCut = new TH1D(
340  "Pi0ConversionSeparationCut", ";True Separation from Vertex (cm);Efficiency", 25, 0, 200);
343  new TH2D("NumHitsCompleteness", ";Completeness;Size of Cluster", 101, 0, 1.01, 100, 0, 100);
345  new TH2D("NumHitsEnergy", ";True Energy (GeV);Size of Cluster", 100, 0, 2, 100, 0, 100);
346 }

Member Function Documentation

void ClusteringValidation::ClusterAnalyser::Analyse ( detinfo::DetectorClocksData const &  clockData,
std::vector< art::Ptr< recob::Hit >> &  hits,
std::vector< art::Ptr< recob::Cluster >> &  clusters,
const art::FindManyP< recob::Hit > &  fmh,
int  numHits 
)

Definition at line 348 of file ClusteringValidation_module.cc.

References sim::ParticleList::begin(), sim::ParticleList::end(), hits(), sim::ParticleList::ParticleList(), geo::PlaneID::Plane, geo::TPCID::TPC, simb::MCParticle::TrackId(), lar::dump::vector(), and recob::Hit::WireID().

353 {
354 
355  // Make a map of cluster counters in TPC/plane space
356  for (auto const& id : wireReadoutGeom.Iterate<geo::PlaneID>(geo::CryostatID{0})) {
357  auto const [tpc, plane] = std::make_tuple(id.TPC, id.Plane);
358  clusterMap[tpc][plane] = std::make_unique<ClusterCounter>();
359  }
360 
361  // Save preclustered hits
362  for (size_t hitIt = 0; hitIt < hits.size(); ++hitIt) {
363  art::Ptr<recob::Hit> hit = hits.at(hitIt);
364  TrackID trackID = FindTrackID(clockData, hit);
365  clusterMap[hit->WireID().TPC % 2][hit->WireID().Plane]->AddHitPreClustering(trackID);
366  }
367 
368  // Save true tracks
369  trueParticles.clear();
370  const sim::ParticleList& particles = pi_serv->ParticleList();
371  for (sim::ParticleList::const_iterator particleIt = particles.begin();
372  particleIt != particles.end();
373  ++particleIt) {
374  const simb::MCParticle* particle = particleIt->second;
375  trueParticles[(TrackID)particle->TrackId()] = particle;
376  }
377 
378  // Save the clustered hits
379  for (size_t clusIt = 0; clusIt < clusters.size(); ++clusIt) {
380 
381  // Get cluster information
382  unsigned int tpc = clusters.at(clusIt)->Plane().TPC % 2;
383  unsigned int plane = clusters.at(clusIt)->Plane().Plane;
384  ClusterID id = (ClusterID)clusters.at(clusIt)->ID();
385 
386  // Only analyse planes with enough hits in to fairly represent the clustering
387  if (clusterMap[tpc][plane]->GetNumberHitsInPlane() < minHits) continue;
388 
389  // Get the hits from the cluster
390  std::vector<art::Ptr<recob::Hit>> clusterHits = fmh.at(clusIt);
391 
392  if (clusterHits.size() < 10) continue;
393 
394  // Find which track this cluster belongs to
395  TrackID trueTrackID = FindTrueTrack(clockData, clusterHits);
396 
397  // Save the info for this cluster
398  clusterMap[tpc][plane]->AssociateClusterAndTrack(id, trueTrackID);
399  for (std::vector<art::Ptr<recob::Hit>>::iterator clusHitIt = clusterHits.begin();
400  clusHitIt != clusterHits.end();
401  ++clusHitIt) {
402  art::Ptr<recob::Hit> hit = *clusHitIt;
403  TrackID trackID = FindTrackID(clockData, hit);
404  if (trackID == trueTrackID)
405  clusterMap[tpc][plane]->AddSignalHitPostClustering(id);
406  else
407  clusterMap[tpc][plane]->AddNoiseHitPostClustering(id);
408  }
409  }
410 
411  this->MakeHistograms();
412 }
list_type::const_iterator const_iterator
Definition: ParticleList.h:132
The data type to uniquely identify a Plane.
Definition: geo_types.h:364
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:290
int TrackId() const
Definition: MCParticle.h:211
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
TrackID FindTrackID(detinfo::DetectorClocksData const &clockData, art::Ptr< recob::Hit > &hit)
art::ServiceHandle< cheat::ParticleInventoryService const > pi_serv
iterator begin()
Definition: ParticleList.h:305
TrackID FindTrueTrack(detinfo::DetectorClocksData const &clockData, std::vector< art::Ptr< recob::Hit >> &clusterHits)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:373
std::map< unsigned int, std::map< unsigned int, std::unique_ptr< ClusterCounter > > > clusterMap
Detector simulation of raw signals on wires.
const sim::ParticleList & ParticleList() const
std::map< TrackID, const simb::MCParticle * > trueParticles
range_type< T > Iterate() const
Definition: Iterable.h:121
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:315
recob::tracking::Plane Plane
Definition: TrackState.h:17
The data type to uniquely identify a cryostat.
Definition: geo_types.h:187
double ClusteringValidation::ClusterAnalyser::FindPhotonAngle ( )

Definition at line 454 of file ClusteringValidation_module.cc.

References simb::MCParticle::Daughter(), and simb::MCParticle::NumberDaughters().

455 {
456  const simb::MCParticle* pi0 = GetPi0();
457  if (pi0->NumberDaughters() != 2) return -999;
458  double angle = (trueParticles.at((TrackID)pi0->Daughter(0))
459  ->Momentum()
460  .Angle(trueParticles.at((TrackID)pi0->Daughter(1))->Momentum().Vect()) *
461  (180 / TMath::Pi()));
462  return angle;
463 }
int NumberDaughters() const
Definition: MCParticle.h:218
int Daughter(const int i) const
Definition: MCParticle.cxx:118
std::map< TrackID, const simb::MCParticle * > trueParticles
TrackID ClusteringValidation::ClusterAnalyser::FindTrackID ( detinfo::DetectorClocksData const &  clockData,
art::Ptr< recob::Hit > &  hit 
)

Definition at line 414 of file ClusteringValidation_module.cc.

417 {
418  double particleEnergy = 0;
419  TrackID likelyTrackID = (TrackID)0;
420  std::vector<sim::TrackIDE> trackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
421  for (unsigned int idIt = 0; idIt < trackIDs.size(); ++idIt) {
422  if (trackIDs.at(idIt).energy > particleEnergy) {
423  particleEnergy = trackIDs.at(idIt).energy;
424  likelyTrackID = (TrackID)TMath::Abs(trackIDs.at(idIt).trackID);
425  }
426  }
427  return likelyTrackID;
428 }
art::ServiceHandle< cheat::BackTrackerService const > bt_serv
std::vector< sim::TrackIDE > HitToTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
TrackID ClusteringValidation::ClusterAnalyser::FindTrueTrack ( detinfo::DetectorClocksData const &  clockData,
std::vector< art::Ptr< recob::Hit >> &  clusterHits 
)

Definition at line 430 of file ClusteringValidation_module.cc.

References recob::Hit::Integral(), and lar::dump::vector().

433 {
434  std::map<TrackID, double> trackMap;
435  for (std::vector<art::Ptr<recob::Hit>>::iterator clusHitIt = clusterHits.begin();
436  clusHitIt != clusterHits.end();
437  ++clusHitIt) {
438  art::Ptr<recob::Hit> hit = *clusHitIt;
439  TrackID trackID = FindTrackID(clockData, hit);
440  trackMap[trackID] += hit->Integral();
441  }
442  //return std::max_element(trackMap.begin(), trackMap.end(), [](const std::pair<int,double>& p1, const std::pair<int,double>& p2) {return p1.second < p2.second;} )->first;
443  double highestCharge = 0;
444  TrackID clusterTrack = (TrackID)0;
445  for (std::map<TrackID, double>::iterator trackIt = trackMap.begin(); trackIt != trackMap.end();
446  ++trackIt)
447  if (trackIt->second > highestCharge) {
448  highestCharge = trackIt->second;
449  clusterTrack = trackIt->first;
450  }
451  return clusterTrack;
452 }
intermediate_table::iterator iterator
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:254
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
TrackID FindTrackID(detinfo::DetectorClocksData const &clockData, art::Ptr< recob::Hit > &hit)
Detector simulation of raw signals on wires.
double ClusteringValidation::ClusterAnalyser::GetEndTrackDistance ( TrackID  id1,
TrackID  id2 
)

Definition at line 475 of file ClusteringValidation_module.cc.

476 {
477  return TMath::Sqrt(
478  TMath::Power(
479  trueParticles.at(id1)->EndPosition().X() - trueParticles.at(id2)->EndPosition().X(), 2) +
480  TMath::Power(
481  trueParticles.at(id1)->EndPosition().Y() - trueParticles.at(id2)->EndPosition().Y(), 2) +
482  TMath::Power(
483  trueParticles.at(id1)->EndPosition().Z() - trueParticles.at(id2)->EndPosition().Z(), 2));
484 }
std::map< TrackID, const simb::MCParticle * > trueParticles
TObjArray ClusteringValidation::ClusterAnalyser::GetHistograms ( )

Definition at line 486 of file ClusteringValidation_module.cc.

487 {
488 
489  // Make efficiency histograms
490  hEfficiencyEnergy = new TEfficiency(*hPi0EnergyCut, *hPi0Energy);
491  hEfficiencyAngle = new TEfficiency(*hPi0AngleCut, *hPi0Angle);
496 
497  hEfficiencyEnergy->SetName("EfficiencyEnergy");
498  hEfficiencyAngle->SetName("EnergyAngle");
499  hEfficiencyConversionDistance->SetName("EfficiencyConversionDistance");
500  hEfficiencyConversionSeparation->SetName("EfficiencyConversionSeparation");
501 
502  // Add all the histograms to the object array
524 
525  return fHistArray;
526 }
const simb::MCParticle * ClusteringValidation::ClusterAnalyser::GetPi0 ( )

Definition at line 465 of file ClusteringValidation_module.cc.

466 {
467  const simb::MCParticle* pi0 = nullptr;
469  particleIt != trueParticles.end();
470  ++particleIt)
471  if (particleIt->second->PdgCode() == 111) pi0 = particleIt->second;
472  return pi0;
473 }
intermediate_table::iterator iterator
std::map< TrackID, const simb::MCParticle * > trueParticles
void ClusteringValidation::ClusterAnalyser::MakeHistograms ( )

Definition at line 528 of file ClusteringValidation_module.cc.

References util::counter(), and E.

529 {
530 
531  // Loop over the tpcs and planes in the geometry
532  for (auto const& id : wireReadoutGeom.Iterate<geo::PlaneID>(geo::CryostatID{0})) {
533  auto const [tpc, plane] = std::make_tuple(id.TPC, id.Plane);
534  auto& counter = clusterMap[tpc][plane];
535  ClusterIDs clusterIDs = counter->GetListOfClusterIDs();
536 
537  // Fill histograms for the efficiency
538  if (counter->GetPhotons().size() == 2) {
539  hPi0Angle->Fill(FindPhotonAngle());
540  hPi0Energy->Fill(GetPi0()->Momentum().E());
541  hPi0ConversionDistance->Fill(std::min(
542  GetEndTrackDistance(counter->GetPhotons().at(0).first, (TrackID)GetPi0()->TrackId()),
543  GetEndTrackDistance(counter->GetPhotons().at(1).first, (TrackID)GetPi0()->TrackId())));
545  GetEndTrackDistance(counter->GetPhotons().at(0).first, counter->GetPhotons().at(1).first));
546  if (counter->PassesCut()) {
547  hPi0AngleCut->Fill(FindPhotonAngle());
548  hPi0EnergyCut->Fill(GetPi0()->Momentum().E());
549  hPi0ConversionDistanceCut->Fill(std::min(
550  GetEndTrackDistance(counter->GetPhotons().at(0).first, (TrackID)GetPi0()->TrackId()),
551  GetEndTrackDistance(counter->GetPhotons().at(1).first, (TrackID)GetPi0()->TrackId())));
552  hPi0ConversionSeparationCut->Fill(GetEndTrackDistance(counter->GetPhotons().at(0).first,
553  counter->GetPhotons().at(1).first));
554  }
555  else
556  std::cout << "TPC " << tpc << ", Plane " << plane << " fails the cut" << std::endl;
557  }
558 
559  // Look at all the clusters
560  for (unsigned int cluster = 0; cluster < clusterIDs.size(); ++cluster) {
561 
562  ClusterID clusID = clusterIDs.at(cluster);
563  double completeness = counter->GetCompleteness(clusID);
564  double cleanliness = counter->GetCleanliness(clusID);
565  int numClusterHits = counter->GetNumberHitsInCluster(clusID);
566 
567  // Fill histograms for this cluster
568  hCompleteness->Fill(completeness, numClusterHits);
569  hCleanliness->Fill(cleanliness, numClusterHits);
570  hComplCleanl->Fill(completeness * cleanliness, numClusterHits);
571  hNumHitsCompleteness->Fill(completeness, numClusterHits);
572 
573  // Is this cluster doesn't correspond to a true particle continue
574  if (counter->IsNoise(clusID)) continue;
575 
576  double pi0Energy = GetPi0()->Momentum().E();
577  double pi0DecayAngle = FindPhotonAngle();
578  double conversionDistance =
579  GetEndTrackDistance(counter->GetTrack(clusID), (TrackID)GetPi0()->TrackId());
580 
581  hCompletenessEnergy->Fill(pi0Energy, completeness, numClusterHits);
582  hCompletenessAngle->Fill(pi0DecayAngle, completeness, numClusterHits);
583  hCompletenessConversionDistance->Fill(conversionDistance, completeness, numClusterHits);
584  hCleanlinessEnergy->Fill(pi0Energy, cleanliness, numClusterHits);
585  hCleanlinessAngle->Fill(pi0DecayAngle, cleanliness, numClusterHits);
586  hCleanlinessConversionDistance->Fill(conversionDistance, cleanliness, numClusterHits);
587  hComplCleanlEnergy->Fill(pi0Energy, cleanliness * completeness, numClusterHits);
588  hComplCleanlAngle->Fill(pi0DecayAngle, cleanliness * completeness, numClusterHits);
590  conversionDistance, cleanliness * completeness, numClusterHits);
591  hNumHitsEnergy->Fill(pi0Energy, numClusterHits);
592 
593  // Continue if there are not two photons in the view
594  if (counter->GetPhotons().size() != 2) continue;
595 
596  double conversionSeparation =
597  GetEndTrackDistance(counter->GetPhotons()[0].first, counter->GetPhotons()[1].first);
598 
599  hCompletenessConversionSeparation->Fill(conversionSeparation, completeness, numClusterHits);
600  hCleanlinessConversionSeparation->Fill(conversionSeparation, cleanliness, numClusterHits);
601  }
602  }
603 }
double GetEndTrackDistance(TrackID id1, TrackID id2)
The data type to uniquely identify a Plane.
Definition: geo_types.h:364
Cluster finding and building.
int TrackId() const
Definition: MCParticle.h:211
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:292
Float_t E
Definition: plot.C:20
std::vector< ClusterID > ClusterIDs
std::map< unsigned int, std::map< unsigned int, std::unique_ptr< ClusterCounter > > > clusterMap
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:221
range_type< T > Iterate() const
Definition: Iterable.h:121
recob::tracking::Plane Plane
Definition: TrackState.h:17
The data type to uniquely identify a cryostat.
Definition: geo_types.h:187
void ClusteringValidation::ClusterAnalyser::WriteFile ( )

Definition at line 605 of file ClusteringValidation_module.cc.

606 {
607 
608  // Average completeness/cleanliness
609  double avCompleteness = hCompleteness->GetMean();
610  double avCleanliness = hCleanliness->GetMean();
611 
612  // Write file
613  std::ofstream outFile("effpur");
614  outFile << avCompleteness << " " << avCleanliness;
615  outFile.close();
616 }

Member Data Documentation

art::ServiceHandle<cheat::BackTrackerService const> ClusteringValidation::ClusterAnalyser::bt_serv
private

Definition at line 275 of file ClusteringValidation_module.cc.

std::map<unsigned int, std::map<unsigned int, std::unique_ptr<ClusterCounter> > > ClusteringValidation::ClusterAnalyser::clusterMap
private

Definition at line 269 of file ClusteringValidation_module.cc.

std::string ClusteringValidation::ClusterAnalyser::fClusterLabel
private

Definition at line 235 of file ClusteringValidation_module.cc.

TObjArray ClusteringValidation::ClusterAnalyser::fHistArray
private

Definition at line 267 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hCleanliness
private

Definition at line 239 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessAngle
private

Definition at line 256 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessConversionDistance
private

Definition at line 257 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessConversionSeparation
private

Definition at line 258 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessEnergy
private

Definition at line 255 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hComplCleanl
private

Definition at line 240 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlAngle
private

Definition at line 260 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlConversionDistance
private

Definition at line 261 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlConversionSeparation
private

Definition at line 262 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlEnergy
private

Definition at line 259 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hCompleteness
private

Definition at line 238 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessAngle
private

Definition at line 252 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessConversionDistance
private

Definition at line 253 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessConversionSeparation
private

Definition at line 254 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessEnergy
private

Definition at line 251 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyAngle
private

Definition at line 263 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyConversionDistance
private

Definition at line 265 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyConversionSeparation
private

Definition at line 266 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyEnergy
private

Definition at line 264 of file ClusteringValidation_module.cc.

TH2* ClusteringValidation::ClusterAnalyser::hNumHitsCompleteness
private

Definition at line 249 of file ClusteringValidation_module.cc.

TH2* ClusteringValidation::ClusterAnalyser::hNumHitsEnergy
private

Definition at line 250 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0Angle
private

Definition at line 241 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0AngleCut
private

Definition at line 245 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionDistance
private

Definition at line 243 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionDistanceCut
private

Definition at line 247 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparation
private

Definition at line 244 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparationCut
private

Definition at line 248 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0Energy
private

Definition at line 242 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0EnergyCut
private

Definition at line 246 of file ClusteringValidation_module.cc.

art::ServiceHandle<cheat::ParticleInventoryService const> ClusteringValidation::ClusterAnalyser::pi_serv
private

Definition at line 274 of file ClusteringValidation_module.cc.

std::map<TrackID, const simb::MCParticle*> ClusteringValidation::ClusterAnalyser::trueParticles
private

Definition at line 270 of file ClusteringValidation_module.cc.

geo::WireReadoutGeom const& ClusteringValidation::ClusterAnalyser::wireReadoutGeom = art::ServiceHandle<geo::WireReadout>()->Get()
private

Definition at line 273 of file ClusteringValidation_module.cc.


The documentation for this class was generated from the following file: