LArSoft  v09_90_00
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
 
art::ServiceHandle< geo::Geometry const > geometry
 
art::ServiceHandle< cheat::ParticleInventoryService const > pi_serv
 
art::ServiceHandle< cheat::BackTrackerService const > bt_serv
 

Detailed Description

Definition at line 215 of file ClusteringValidation_module.cc.

Constructor & Destructor Documentation

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

Definition at line 279 of file ClusteringValidation_module.cc.

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

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 349 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().

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

Definition at line 455 of file ClusteringValidation_module.cc.

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

456 {
457  const simb::MCParticle* pi0 = GetPi0();
458  if (pi0->NumberDaughters() != 2) return -999;
459  double angle = (trueParticles.at((TrackID)pi0->Daughter(0))
460  ->Momentum()
461  .Angle(trueParticles.at((TrackID)pi0->Daughter(1))->Momentum().Vect()) *
462  (180 / TMath::Pi()));
463  return angle;
464 }
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 415 of file ClusteringValidation_module.cc.

418 {
419  double particleEnergy = 0;
420  TrackID likelyTrackID = (TrackID)0;
421  std::vector<sim::TrackIDE> trackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
422  for (unsigned int idIt = 0; idIt < trackIDs.size(); ++idIt) {
423  if (trackIDs.at(idIt).energy > particleEnergy) {
424  particleEnergy = trackIDs.at(idIt).energy;
425  likelyTrackID = (TrackID)TMath::Abs(trackIDs.at(idIt).trackID);
426  }
427  }
428  return likelyTrackID;
429 }
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 431 of file ClusteringValidation_module.cc.

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

434 {
435  std::map<TrackID, double> trackMap;
436  for (std::vector<art::Ptr<recob::Hit>>::iterator clusHitIt = clusterHits.begin();
437  clusHitIt != clusterHits.end();
438  ++clusHitIt) {
439  art::Ptr<recob::Hit> hit = *clusHitIt;
440  TrackID trackID = FindTrackID(clockData, hit);
441  trackMap[trackID] += hit->Integral();
442  }
443  //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;
444  double highestCharge = 0;
445  TrackID clusterTrack = (TrackID)0;
446  for (std::map<TrackID, double>::iterator trackIt = trackMap.begin(); trackIt != trackMap.end();
447  ++trackIt)
448  if (trackIt->second > highestCharge) {
449  highestCharge = trackIt->second;
450  clusterTrack = trackIt->first;
451  }
452  return clusterTrack;
453 }
intermediate_table::iterator iterator
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:244
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 476 of file ClusteringValidation_module.cc.

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

Definition at line 487 of file ClusteringValidation_module.cc.

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

Definition at line 466 of file ClusteringValidation_module.cc.

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

Definition at line 529 of file ClusteringValidation_module.cc.

References util::counter(), and E.

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

Definition at line 606 of file ClusteringValidation_module.cc.

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

Member Data Documentation

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

Definition at line 276 of file ClusteringValidation_module.cc.

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

Definition at line 270 of file ClusteringValidation_module.cc.

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

Definition at line 236 of file ClusteringValidation_module.cc.

TObjArray ClusteringValidation::ClusterAnalyser::fHistArray
private

Definition at line 268 of file ClusteringValidation_module.cc.

art::ServiceHandle<geo::Geometry const> ClusteringValidation::ClusterAnalyser::geometry
private

Definition at line 274 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hCleanliness
private

Definition at line 240 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessAngle
private

Definition at line 257 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessConversionDistance
private

Definition at line 258 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessConversionSeparation
private

Definition at line 259 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessEnergy
private

Definition at line 256 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hComplCleanl
private

Definition at line 241 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlAngle
private

Definition at line 261 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlConversionDistance
private

Definition at line 262 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlConversionSeparation
private

Definition at line 263 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlEnergy
private

Definition at line 260 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hCompleteness
private

Definition at line 239 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessAngle
private

Definition at line 253 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessConversionDistance
private

Definition at line 254 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessConversionSeparation
private

Definition at line 255 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessEnergy
private

Definition at line 252 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyAngle
private

Definition at line 264 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyConversionDistance
private

Definition at line 266 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyConversionSeparation
private

Definition at line 267 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyEnergy
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH2* ClusteringValidation::ClusterAnalyser::hNumHitsCompleteness
private

Definition at line 250 of file ClusteringValidation_module.cc.

TH2* ClusteringValidation::ClusterAnalyser::hNumHitsEnergy
private

Definition at line 251 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0Angle
private

Definition at line 242 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0AngleCut
private

Definition at line 246 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionDistance
private

Definition at line 244 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionDistanceCut
private

Definition at line 248 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparation
private

Definition at line 245 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparationCut
private

Definition at line 249 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0Energy
private

Definition at line 243 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0EnergyCut
private

Definition at line 247 of file ClusteringValidation_module.cc.

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

Definition at line 275 of file ClusteringValidation_module.cc.

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

Definition at line 271 of file ClusteringValidation_module.cc.


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