Public Member Functions | Private Attributes | List of all members
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 240 of file ClusteringValidation_module.cc.

Constructor & Destructor Documentation

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

Definition at line 287 of file ClusteringValidation_module.cc.

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

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 358 of file ClusteringValidation_module.cc.

363 {
364 
365  // Make a map of cluster counters in TPC/plane space
366  for (unsigned int tpc = 0; tpc < geometry->NTPC(0); ++tpc) {
367  for (unsigned int plane = 0; plane < geometry->Nplanes(tpc, 0); ++plane) {
368  clusterMap[tpc][plane] = (std::unique_ptr<ClusterCounter>)new ClusterCounter(tpc, plane);
369  }
370  }
371 
372  // Save preclustered hits
373  for (size_t hitIt = 0; hitIt < hits.size(); ++hitIt) {
374  art::Ptr<recob::Hit> hit = hits.at(hitIt);
375  TrackID trackID = FindTrackID(clockData, hit);
376  clusterMap[hit->WireID().TPC % 2][hit->WireID().Plane]->AddHitPreClustering(trackID);
377  }
378 
379  // Save true tracks
380  trueParticles.clear();
381  const sim::ParticleList& particles = pi_serv->ParticleList();
382  for (sim::ParticleList::const_iterator particleIt = particles.begin();
383  particleIt != particles.end();
384  ++particleIt) {
385  const simb::MCParticle* particle = particleIt->second;
386  trueParticles[(TrackID)particle->TrackId()] = particle;
387  }
388 
389  // Save the clustered hits
390  for (size_t clusIt = 0; clusIt < clusters.size(); ++clusIt) {
391 
392  // Get cluster information
393  unsigned int tpc = clusters.at(clusIt)->Plane().TPC % 2;
394  unsigned int plane = clusters.at(clusIt)->Plane().Plane;
395  ClusterID id = (ClusterID)clusters.at(clusIt)->ID();
396 
397  // Only analyse planes with enough hits in to fairly represent the clustering
398  if (clusterMap[tpc][plane]->GetNumberHitsInPlane() < minHits) continue;
399 
400  // Get the hits from the cluster
401  std::vector<art::Ptr<recob::Hit>> clusterHits = fmh.at(clusIt);
402 
403  if (clusterHits.size() < 10) continue;
404 
405  // Find which track this cluster belongs to
406  TrackID trueTrackID = FindTrueTrack(clockData, clusterHits);
407 
408  // Save the info for this cluster
409  clusterMap[tpc][plane]->AssociateClusterAndTrack(id, trueTrackID);
410  for (std::vector<art::Ptr<recob::Hit>>::iterator clusHitIt = clusterHits.begin();
411  clusHitIt != clusterHits.end();
412  ++clusHitIt) {
413  art::Ptr<recob::Hit> hit = *clusHitIt;
414  TrackID trackID = FindTrackID(clockData, hit);
415  if (trackID == trueTrackID)
416  clusterMap[tpc][plane]->AddSignalHitPostClustering(id);
417  else
418  clusterMap[tpc][plane]->AddNoiseHitPostClustering(id);
419  }
420  }
421 
422  this->MakeHistograms();
423 }
geo::WireID WireID() const
Definition: Hit.h:233
struct vector vector
intermediate_table::const_iterator const_iterator
int TrackId() const
Definition: MCParticle.h:210
TrackID FindTrackID(detinfo::DetectorClocksData const &clockData, art::Ptr< recob::Hit > &hit)
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
art::ServiceHandle< cheat::ParticleInventoryService const > pi_serv
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:493
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
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
std::map< TrackID, const simb::MCParticle * > trueParticles
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
art::ServiceHandle< geo::Geometry const > geometry
double ClusteringValidation::ClusterAnalyser::FindPhotonAngle ( )

Definition at line 466 of file ClusteringValidation_module.cc.

467 {
468  const simb::MCParticle* pi0 = GetPi0();
469  if (pi0->NumberDaughters() != 2) return -999;
470  double angle = (trueParticles.at((TrackID)pi0->Daughter(0))
471  ->Momentum()
472  .Angle(trueParticles.at((TrackID)pi0->Daughter(1))->Momentum().Vect()) *
473  (180 / TMath::Pi()));
474  return angle;
475 }
int NumberDaughters() const
Definition: MCParticle.h:217
int Daughter(const int i) const
Definition: MCParticle.cxx:112
std::map< TrackID, const simb::MCParticle * > trueParticles
TrackID ClusteringValidation::ClusterAnalyser::FindTrackID ( detinfo::DetectorClocksData const &  clockData,
art::Ptr< recob::Hit > &  hit 
)

Definition at line 426 of file ClusteringValidation_module.cc.

428 {
429  double particleEnergy = 0;
430  TrackID likelyTrackID = (TrackID)0;
431  std::vector<sim::TrackIDE> trackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
432  for (unsigned int idIt = 0; idIt < trackIDs.size(); ++idIt) {
433  if (trackIDs.at(idIt).energy > particleEnergy) {
434  particleEnergy = trackIDs.at(idIt).energy;
435  likelyTrackID = (TrackID)TMath::Abs(trackIDs.at(idIt).trackID);
436  }
437  }
438  return likelyTrackID;
439 }
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 442 of file ClusteringValidation_module.cc.

444 {
445  std::map<TrackID, double> trackMap;
446  for (std::vector<art::Ptr<recob::Hit>>::iterator clusHitIt = clusterHits.begin();
447  clusHitIt != clusterHits.end();
448  ++clusHitIt) {
449  art::Ptr<recob::Hit> hit = *clusHitIt;
450  TrackID trackID = FindTrackID(clockData, hit);
451  trackMap[trackID] += hit->Integral();
452  }
453  //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;
454  double highestCharge = 0;
455  TrackID clusterTrack = (TrackID)0;
456  for (std::map<TrackID, double>::iterator trackIt = trackMap.begin(); trackIt != trackMap.end();
457  ++trackIt)
458  if (trackIt->second > highestCharge) {
459  highestCharge = trackIt->second;
460  clusterTrack = trackIt->first;
461  }
462  return clusterTrack;
463 }
intermediate_table::iterator iterator
struct vector vector
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:224
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 489 of file ClusteringValidation_module.cc.

490 {
491  return TMath::Sqrt(
492  TMath::Power(
493  trueParticles.at(id1)->EndPosition().X() - trueParticles.at(id2)->EndPosition().X(), 2) +
494  TMath::Power(
495  trueParticles.at(id1)->EndPosition().Y() - trueParticles.at(id2)->EndPosition().Y(), 2) +
496  TMath::Power(
497  trueParticles.at(id1)->EndPosition().Z() - trueParticles.at(id2)->EndPosition().Z(), 2));
498 }
std::map< TrackID, const simb::MCParticle * > trueParticles
TObjArray ClusteringValidation::ClusterAnalyser::GetHistograms ( )

Definition at line 501 of file ClusteringValidation_module.cc.

502 {
503 
504  // Make efficiency histograms
505  hEfficiencyEnergy = new TEfficiency(*hPi0EnergyCut, *hPi0Energy);
506  hEfficiencyAngle = new TEfficiency(*hPi0AngleCut, *hPi0Angle);
511 
512  hEfficiencyEnergy->SetName("EfficiencyEnergy");
513  hEfficiencyAngle->SetName("EnergyAngle");
514  hEfficiencyConversionDistance->SetName("EfficiencyConversionDistance");
515  hEfficiencyConversionSeparation->SetName("EfficiencyConversionSeparation");
516 
517  // Add all the histograms to the object array
539 
540  return fHistArray;
541 }
const simb::MCParticle * ClusteringValidation::ClusterAnalyser::GetPi0 ( )

Definition at line 478 of file ClusteringValidation_module.cc.

479 {
480  const simb::MCParticle* pi0 = nullptr;
482  particleIt != trueParticles.end();
483  ++particleIt)
484  if (particleIt->second->PdgCode() == 111) pi0 = particleIt->second;
485  return pi0;
486 }
intermediate_table::iterator iterator
std::map< TrackID, const simb::MCParticle * > trueParticles
void ClusteringValidation::ClusterAnalyser::MakeHistograms ( )

Definition at line 544 of file ClusteringValidation_module.cc.

545 {
546 
547  // Loop over the tpcs and planes in the geometry
548  for (unsigned int tpc = 0; tpc < geometry->NTPC(0); ++tpc) {
549  for (unsigned int plane = 0; plane < geometry->Nplanes(tpc, 0); ++plane) {
550 
551  ClusterIDs clusterIDs = clusterMap[tpc][plane]->GetListOfClusterIDs();
552 
553  // Fill histograms for the efficiency
554  if (clusterMap[tpc][plane]->GetPhotons().size() == 2) {
555  hPi0Angle->Fill(FindPhotonAngle());
556  hPi0Energy->Fill(GetPi0()->Momentum().E());
558  std::min(GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first,
559  (TrackID)GetPi0()->TrackId()),
560  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(1).first,
561  (TrackID)GetPi0()->TrackId())));
563  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first,
564  clusterMap[tpc][plane]->GetPhotons().at(1).first));
565  if (clusterMap[tpc][plane]->PassesCut()) {
566  hPi0AngleCut->Fill(FindPhotonAngle());
567  hPi0EnergyCut->Fill(GetPi0()->Momentum().E());
569  std::min(GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first,
570  (TrackID)GetPi0()->TrackId()),
571  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(1).first,
572  (TrackID)GetPi0()->TrackId())));
574  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first,
575  clusterMap[tpc][plane]->GetPhotons().at(1).first));
576  }
577  else
578  std::cout << "TPC " << tpc << ", Plane " << plane << " fails the cut" << std::endl;
579  }
580 
581  // Look at all the clusters
582  for (unsigned int cluster = 0; cluster < clusterIDs.size(); ++cluster) {
583 
584  ClusterID clusID = clusterIDs.at(cluster);
585  double completeness = clusterMap[tpc][plane]->GetCompleteness(clusID);
586  double cleanliness = clusterMap[tpc][plane]->GetCleanliness(clusID);
587  int numClusterHits = clusterMap[tpc][plane]->GetNumberHitsInCluster(clusID);
588 
589  // Fill histograms for this cluster
590  hCompleteness->Fill(completeness, numClusterHits);
591  hCleanliness->Fill(cleanliness, numClusterHits);
592  hComplCleanl->Fill(completeness * cleanliness, numClusterHits);
593  hNumHitsCompleteness->Fill(completeness, numClusterHits);
594 
595  // Is this cluster doesn't correspond to a true particle continue
596  if (clusterMap[tpc][plane]->IsNoise(clusID)) continue;
597 
598  double pi0Energy = GetPi0()->Momentum().E();
599  double pi0DecayAngle = FindPhotonAngle();
600  double conversionDistance = GetEndTrackDistance(clusterMap[tpc][plane]->GetTrack(clusID),
601  (TrackID)GetPi0()->TrackId());
602 
603  hCompletenessEnergy->Fill(pi0Energy, completeness, numClusterHits);
604  hCompletenessAngle->Fill(pi0DecayAngle, completeness, numClusterHits);
605  hCompletenessConversionDistance->Fill(conversionDistance, completeness, numClusterHits);
606  hCleanlinessEnergy->Fill(pi0Energy, cleanliness, numClusterHits);
607  hCleanlinessAngle->Fill(pi0DecayAngle, cleanliness, numClusterHits);
608  hCleanlinessConversionDistance->Fill(conversionDistance, cleanliness, numClusterHits);
609  hComplCleanlEnergy->Fill(pi0Energy, cleanliness * completeness, numClusterHits);
610  hComplCleanlAngle->Fill(pi0DecayAngle, cleanliness * completeness, numClusterHits);
612  conversionDistance, cleanliness * completeness, numClusterHits);
613  hNumHitsEnergy->Fill(pi0Energy, numClusterHits);
614 
615  // Continue if there are not two photons in the view
616  if (clusterMap[tpc][plane]->GetPhotons().size() != 2) continue;
617 
618  double conversionSeparation =
619  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first,
620  clusterMap[tpc][plane]->GetPhotons().at(1).first);
621 
622  hCompletenessConversionSeparation->Fill(conversionSeparation, completeness, numClusterHits);
623  hCleanlinessConversionSeparation->Fill(conversionSeparation, cleanliness, numClusterHits);
624  }
625  }
626  }
627 }
void cluster(In first, In last, Out result, Pred *pred)
Definition: NNClusters.h:41
double GetEndTrackDistance(TrackID id1, TrackID id2)
Cluster finding and building.
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
std::vector< ClusterID > ClusterIDs
std::map< unsigned int, std::map< unsigned int, std::unique_ptr< ClusterCounter > > > clusterMap
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
E
Definition: 018_def.c:13
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:220
QTextStream & endl(QTextStream &s)
art::ServiceHandle< geo::Geometry const > geometry
void ClusteringValidation::ClusterAnalyser::WriteFile ( )

Definition at line 630 of file ClusteringValidation_module.cc.

631 {
632 
633  // Average completeness/cleanliness
634  double avCompleteness = hCompleteness->GetMean();
635  double avCleanliness = hCleanliness->GetMean();
636 
637  // Write file
638  std::ofstream outFile("effpur");
639  outFile << avCompleteness << " " << avCleanliness;
640  outFile.close();
641 }
TFile * outFile
Definition: makeDST.cxx:36

Member Data Documentation

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

Definition at line 284 of file ClusteringValidation_module.cc.

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

Definition at line 278 of file ClusteringValidation_module.cc.

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

Definition at line 261 of file ClusteringValidation_module.cc.

TObjArray ClusteringValidation::ClusterAnalyser::fHistArray
private

Definition at line 276 of file ClusteringValidation_module.cc.

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

Definition at line 282 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hCleanliness
private

Definition at line 264 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCleanlinessAngle
private

Definition at line 270 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCleanlinessConversionDistance
private

Definition at line 270 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCleanlinessConversionSeparation
private

Definition at line 270 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessEnergy
private

Definition at line 270 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hComplCleanl
private

Definition at line 264 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hComplCleanlAngle
private

Definition at line 272 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hComplCleanlConversionDistance
private

Definition at line 272 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hComplCleanlConversionSeparation
private

Definition at line 272 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlEnergy
private

Definition at line 272 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hCompleteness
private

Definition at line 264 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCompletenessAngle
private

Definition at line 268 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCompletenessConversionDistance
private

Definition at line 268 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCompletenessConversionSeparation
private

Definition at line 268 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessEnergy
private

Definition at line 268 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyAngle
private

Definition at line 274 of file ClusteringValidation_module.cc.

TEfficiency * ClusteringValidation::ClusterAnalyser::hEfficiencyConversionDistance
private

Definition at line 274 of file ClusteringValidation_module.cc.

TEfficiency * ClusteringValidation::ClusterAnalyser::hEfficiencyConversionSeparation
private

Definition at line 274 of file ClusteringValidation_module.cc.

TEfficiency * ClusteringValidation::ClusterAnalyser::hEfficiencyEnergy
private

Definition at line 274 of file ClusteringValidation_module.cc.

TH2* ClusteringValidation::ClusterAnalyser::hNumHitsCompleteness
private

Definition at line 267 of file ClusteringValidation_module.cc.

TH2 * ClusteringValidation::ClusterAnalyser::hNumHitsEnergy
private

Definition at line 267 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0Angle
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0AngleCut
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionDistance
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionDistanceCut
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparation
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparationCut
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0Energy
private

Definition at line 265 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0EnergyCut
private

Definition at line 265 of file ClusteringValidation_module.cc.

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

Definition at line 283 of file ClusteringValidation_module.cc.

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

Definition at line 279 of file ClusteringValidation_module.cc.


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