LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
ClusteringValidation::ClusterAnalyser Class Reference

Public Member Functions

 ClusterAnalyser (std::string &label)
 
void Analyse (std::vector< art::Ptr< recob::Hit > > &hits, std::vector< art::Ptr< recob::Cluster > > &clusters, const art::FindManyP< recob::Hit > &fmh, int numHits)
 
TrackID FindTrackID (art::Ptr< recob::Hit > &hit)
 
TrackID FindTrueTrack (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::Geometrygeometry
 
art::ServiceHandle< cheat::ParticleInventoryServicepi_serv
 
art::ServiceHandle< cheat::BackTrackerServicebt_serv
 

Detailed Description

Definition at line 173 of file ClusteringValidation_module.cc.

Constructor & Destructor Documentation

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

Definition at line 213 of file ClusteringValidation_module.cc.

213  {
214 
215  fClusterLabel = clusterLabel;
216 
217  // Make the histograms
218  hCompleteness = new TH1D("Completeness",";Completeness;",101,0,1.01);
219  hCompletenessEnergy = new TProfile("CompletenessEnergy",";True Energy (GeV);Completeness",100,0,2);
220  hCompletenessAngle = new TProfile("CompletenessAngle",";True Angle (deg);Completeness;",100,0,180);
221  hCompletenessConversionDistance = new TProfile("CompletenessConversionDistance",";True Distance from Vertex (cm);Completeness",100,0,200);
222  hCompletenessConversionSeparation = new TProfile("CompletenessConversionSeparation",";True Conversion Separation (cm);Completeness",100,0,200);
223  hCleanliness = new TH1D("Cleanliness",";Cleanliness;",101,0,1.01);
224  hCleanlinessEnergy = new TProfile("CleanlinessEnergy",";True Energy (GeV);Cleanliness",100,0,2);
225  hCleanlinessAngle = new TProfile("CleanlinessAngle",";True Angle (deg);Cleanliness;",100,0,180);
226  hCleanlinessConversionDistance = new TProfile("CleanlinessConversionDistance",";True Distance from Vertex (cm);Cleanliness",100,0,200);
227  hCleanlinessConversionSeparation = new TProfile("CleanlinessConversionSeparation",";True Conversion Separation (cm);Cleanliness",100,0,200);
228  hComplCleanl = new TH1D("CompletenessCleanliness",";Completeness * Cleanliness;",101,0,1.01);
229  hComplCleanlEnergy = new TProfile("CompletenessCleanlinessEnergy",";True Energy (GeV);Completeness * Cleanliness",100,0,2);
230  hComplCleanlAngle = new TProfile("CompletenessCleanlinessAngle",";True Angle (deg);Completeness * Cleanliness;",100,0,180);
231  hComplCleanlConversionDistance = new TProfile("CompletenessCleanlinessConversionDistance",";True Distance from Vertex (cm);Completeness * Cleanliness",100,0,200);
232  hComplCleanlConversionSeparation = new TProfile("CompletenessCleanlinessConversionSeparation",";True Conversion Separation (cm);Completeness * Cleanliness",100,0,200);
233  hPi0Energy = new TH1D("Pi0EnergyCut",";True Energy (GeV);",25,0,2); hPi0Energy ->Sumw2();
234  hPi0Angle = new TH1D("Pi0AngleCut",";True Angle (deg);",25,0,180); hPi0Angle ->Sumw2();
235  hPi0ConversionDistance = new TH1D("Pi0ConversionDistanceCut",";True Distance from Vertex (cm);",25,0,200); hPi0ConversionDistance ->Sumw2();
236  hPi0ConversionSeparation = new TH1D("Pi0ConversionSeparationCut",";True Separation from Vertex (cm);",25,0,200); hPi0ConversionSeparation ->Sumw2();
237  hPi0EnergyCut = new TH1D("Pi0EnergyCut",";True Energy (GeV);Efficiency",25,0,2); hPi0EnergyCut ->Sumw2();
238  hPi0AngleCut = new TH1D("Pi0AngleCut",";True Angle (deg);Efficiency",25,0,180); hPi0AngleCut ->Sumw2();
239  hPi0ConversionDistanceCut = new TH1D("Pi0ConversionDistanceCut",";True Distance from Vertex (cm);Efficiency",25,0,200); hPi0ConversionDistanceCut ->Sumw2();
240  hPi0ConversionSeparationCut = new TH1D("Pi0ConversionSeparationCut",";True Separation from Vertex (cm);Efficiency",25,0,200); hPi0ConversionSeparationCut->Sumw2();
241  hNumHitsCompleteness = new TH2D("NumHitsCompleteness",";Completeness;Size of Cluster",101,0,1.01,100,0,100);
242  hNumHitsEnergy = new TH2D("NumHitsEnergy",";True Energy (GeV);Size of Cluster",100,0,2,100,0,100);
243 
244 }

Member Function Documentation

void ClusteringValidation::ClusterAnalyser::Analyse ( 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 246 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().

246  {
247 
248  // Make a map of cluster counters in TPC/plane space
249  for (unsigned int tpc = 0; tpc < geometry->NTPC(0); ++tpc) {
250  for (unsigned int plane = 0; plane < geometry->Nplanes(tpc,0); ++plane) {
251  clusterMap[tpc][plane] = (std::unique_ptr<ClusterCounter>) new ClusterCounter(tpc, plane);
252  }
253  }
254 
255  // Save preclustered hits
256  for (size_t hitIt = 0; hitIt < hits.size(); ++hitIt) {
257  art::Ptr<recob::Hit> hit = hits.at(hitIt);
258  TrackID trackID = FindTrackID(hit);
259  clusterMap[hit->WireID().TPC%2][hit->WireID().Plane]->AddHitPreClustering(trackID);
260  }
261 
262  // Save true tracks
263  trueParticles.clear();
264  const sim::ParticleList& particles = pi_serv->ParticleList();
265  for (sim::ParticleList::const_iterator particleIt = particles.begin(); particleIt != particles.end(); ++particleIt) {
266  const simb::MCParticle *particle = particleIt->second;
267  trueParticles[(TrackID)particle->TrackId()] = particle;
268  }
269 
270  // Save the clustered hits
271  for (size_t clusIt = 0; clusIt < clusters.size(); ++clusIt) {
272 
273  // Get cluster information
274  unsigned int tpc = clusters.at(clusIt)->Plane().TPC%2;
275  unsigned int plane = clusters.at(clusIt)->Plane().Plane;
276  ClusterID id = (ClusterID)clusters.at(clusIt)->ID();
277 
278  // Only analyse planes with enough hits in to fairly represent the clustering
279  if (clusterMap[tpc][plane]->GetNumberHitsInPlane() < minHits) continue;
280 
281  // Get the hits from the cluster
282  std::vector<art::Ptr<recob::Hit> > clusterHits = fmh.at(clusIt);
283 
284  if (clusterHits.size() < 10) continue;
285 
286  // Find which track this cluster belongs to
287  TrackID trueTrackID = FindTrueTrack(clusterHits);
288 
289  // Save the info for this cluster
290  clusterMap[tpc][plane]->AssociateClusterAndTrack(id, trueTrackID);
291  for (std::vector<art::Ptr<recob::Hit> >::iterator clusHitIt = clusterHits.begin(); clusHitIt != clusterHits.end(); ++clusHitIt) {
292  art::Ptr<recob::Hit> hit = *clusHitIt;
293  TrackID trackID = FindTrackID(hit);
294  if (trackID == trueTrackID) clusterMap[tpc][plane]->AddSignalHitPostClustering(id);
295  else clusterMap[tpc][plane]->AddNoiseHitPostClustering (id);
296  }
297  }
298 
299  this->MakeHistograms();
300 
301 }
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
list_type::const_iterator const_iterator
Definition: ParticleList.h:132
std::map< TrackID, const simb::MCParticle * > trueParticles
int TrackId() const
Definition: MCParticle.h:214
art::ServiceHandle< geo::Geometry > geometry
std::map< unsigned int, std::map< unsigned int, std::unique_ptr< ClusterCounter > > > clusterMap
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
iterator begin()
Definition: ParticleList.h:305
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
Detector simulation of raw signals on wires.
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
TrackID FindTrackID(art::Ptr< recob::Hit > &hit)
TrackID FindTrueTrack(std::vector< art::Ptr< recob::Hit > > &clusterHits)
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
art::ServiceHandle< cheat::ParticleInventoryService > pi_serv
double ClusteringValidation::ClusterAnalyser::FindPhotonAngle ( )

Definition at line 334 of file ClusteringValidation_module.cc.

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

334  {
335  const simb::MCParticle* pi0 = GetPi0();
336  if (pi0->NumberDaughters() != 2) return -999;
337  double angle = (trueParticles.at((TrackID)pi0->Daughter(0))->Momentum().Angle(trueParticles.at((TrackID)pi0->Daughter(1))->Momentum().Vect()) * (180/TMath::Pi()));
338  return angle;
339 }
std::map< TrackID, const simb::MCParticle * > trueParticles
int NumberDaughters() const
Definition: MCParticle.h:221
int Daughter(const int i) const
Definition: MCParticle.cxx:112
TrackID ClusteringValidation::ClusterAnalyser::FindTrackID ( art::Ptr< recob::Hit > &  hit)

Definition at line 303 of file ClusteringValidation_module.cc.

303  {
304  double particleEnergy = 0;
305  TrackID likelyTrackID = (TrackID)0;
306  std::vector<sim::TrackIDE> trackIDs = bt_serv->HitToTrackIDEs(hit);
307  for (unsigned int idIt = 0; idIt < trackIDs.size(); ++idIt) {
308  if (trackIDs.at(idIt).energy > particleEnergy) {
309  particleEnergy = trackIDs.at(idIt).energy;
310  likelyTrackID = (TrackID)TMath::Abs(trackIDs.at(idIt).trackID);
311  }
312  }
313  return likelyTrackID;
314 }
art::ServiceHandle< cheat::BackTrackerService > bt_serv
const std::vector< sim::TrackIDE > HitToTrackIDEs(recob::Hit const &hit)
TrackID ClusteringValidation::ClusterAnalyser::FindTrueTrack ( std::vector< art::Ptr< recob::Hit > > &  clusterHits)

Definition at line 316 of file ClusteringValidation_module.cc.

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

316  {
317  std::map<TrackID,double> trackMap;
318  for (std::vector<art::Ptr<recob::Hit> >::iterator clusHitIt = clusterHits.begin(); clusHitIt != clusterHits.end(); ++clusHitIt) {
319  art::Ptr<recob::Hit> hit = *clusHitIt;
320  TrackID trackID = FindTrackID(hit);
321  trackMap[trackID] += hit->Integral();
322  }
323  //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;
324  double highestCharge = 0;
325  TrackID clusterTrack = (TrackID)0;
326  for (std::map<TrackID,double>::iterator trackIt = trackMap.begin(); trackIt != trackMap.end(); ++trackIt)
327  if (trackIt->second > highestCharge) {
328  highestCharge = trackIt->second;
329  clusterTrack = trackIt->first;
330  }
331  return clusterTrack;
332 }
intermediate_table::iterator iterator
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:225
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
Detector simulation of raw signals on wires.
TrackID FindTrackID(art::Ptr< recob::Hit > &hit)
double ClusteringValidation::ClusterAnalyser::GetEndTrackDistance ( TrackID  id1,
TrackID  id2 
)

Definition at line 349 of file ClusteringValidation_module.cc.

349  {
350  return TMath::Sqrt(TMath::Power(trueParticles.at(id1)->EndPosition().X() - trueParticles.at(id2)->EndPosition().X(),2)+
351  TMath::Power(trueParticles.at(id1)->EndPosition().Y() - trueParticles.at(id2)->EndPosition().Y(),2)+
352  TMath::Power(trueParticles.at(id1)->EndPosition().Z() - trueParticles.at(id2)->EndPosition().Z(),2));
353 }
std::map< TrackID, const simb::MCParticle * > trueParticles
TObjArray ClusteringValidation::ClusterAnalyser::GetHistograms ( )

Definition at line 355 of file ClusteringValidation_module.cc.

355  {
356 
357  // Make efficiency histograms
358  hEfficiencyEnergy = new TEfficiency(*hPi0EnergyCut,*hPi0Energy);
359  hEfficiencyAngle = new TEfficiency(*hPi0AngleCut,*hPi0Angle);
362 
363  hEfficiencyEnergy ->SetName("EfficiencyEnergy");
364  hEfficiencyAngle ->SetName("EnergyAngle");
365  hEfficiencyConversionDistance ->SetName("EfficiencyConversionDistance");
366  hEfficiencyConversionSeparation->SetName("EfficiencyConversionSeparation");
367 
368  // Add all the histograms to the object array
374 
375  return fHistArray;
376 }
const simb::MCParticle * ClusteringValidation::ClusterAnalyser::GetPi0 ( )

Definition at line 341 of file ClusteringValidation_module.cc.

341  {
342  const simb::MCParticle* pi0 = nullptr;
343  for (std::map<TrackID,const simb::MCParticle*>::iterator particleIt = trueParticles.begin(); particleIt != trueParticles.end(); ++particleIt)
344  if (particleIt->second->PdgCode() == 111)
345  pi0 = particleIt->second;
346  return pi0;
347 }
intermediate_table::iterator iterator
std::map< TrackID, const simb::MCParticle * > trueParticles
void ClusteringValidation::ClusterAnalyser::MakeHistograms ( )

Definition at line 378 of file ClusteringValidation_module.cc.

References E, and min.

378  {
379 
380  // Loop over the tpcs and planes in the geometry
381  for (unsigned int tpc = 0; tpc < geometry->NTPC(0); ++tpc) {
382  for (unsigned int plane = 0; plane < geometry->Nplanes(tpc,0); ++plane) {
383 
384  ClusterIDs clusterIDs = clusterMap[tpc][plane]->GetListOfClusterIDs();
385 
386  // Fill histograms for the efficiency
387  if (clusterMap[tpc][plane]->GetPhotons().size() == 2) {
388  hPi0Angle ->Fill(FindPhotonAngle());
389  hPi0Energy ->Fill(GetPi0()->Momentum().E());
390  hPi0ConversionDistance ->Fill(std::min(GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first, (TrackID)GetPi0()->TrackId()),
391  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(1).first, (TrackID)GetPi0()->TrackId())));
392  hPi0ConversionSeparation->Fill(GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first, clusterMap[tpc][plane]->GetPhotons().at(1).first));
393  if (clusterMap[tpc][plane]->PassesCut()) {
394  hPi0AngleCut ->Fill(FindPhotonAngle());
395  hPi0EnergyCut ->Fill(GetPi0()->Momentum().E());
396  hPi0ConversionDistanceCut ->Fill(std::min(GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first, (TrackID)GetPi0()->TrackId()),
397  GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(1).first, (TrackID)GetPi0()->TrackId())));
398  hPi0ConversionSeparationCut->Fill(GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first, clusterMap[tpc][plane]->GetPhotons().at(1).first));
399  }
400  else
401  std::cout << "TPC " << tpc << ", Plane " << plane << " fails the cut" << std::endl;
402  }
403 
404  // Look at all the clusters
405  for (unsigned int cluster = 0; cluster < clusterIDs.size(); ++cluster) {
406 
407  ClusterID clusID = clusterIDs.at(cluster);
408  double completeness = clusterMap[tpc][plane]->GetCompleteness(clusID);
409  double cleanliness = clusterMap[tpc][plane]->GetCleanliness(clusID);
410  int numClusterHits = clusterMap[tpc][plane]->GetNumberHitsInCluster(clusID);
411 
412  // Fill histograms for this cluster
413  hCompleteness ->Fill(completeness, numClusterHits);
414  hCleanliness ->Fill(cleanliness, numClusterHits);
415  hComplCleanl ->Fill(completeness*cleanliness, numClusterHits);
416  hNumHitsCompleteness ->Fill(completeness, numClusterHits);
417 
418  // Is this cluster doesn't correspond to a true particle continue
419  if (clusterMap[tpc][plane]->IsNoise(clusID)) continue;
420 
421  double pi0Energy = GetPi0()->Momentum().E();
422  double pi0DecayAngle = FindPhotonAngle();
423  double conversionDistance = GetEndTrackDistance(clusterMap[tpc][plane]->GetTrack(clusID), (TrackID)GetPi0()->TrackId());
424 
425  hCompletenessEnergy ->Fill(pi0Energy, completeness, numClusterHits);
426  hCompletenessAngle ->Fill(pi0DecayAngle, completeness, numClusterHits);
427  hCompletenessConversionDistance->Fill(conversionDistance, completeness, numClusterHits);
428  hCleanlinessEnergy ->Fill(pi0Energy, cleanliness, numClusterHits);
429  hCleanlinessAngle ->Fill(pi0DecayAngle, cleanliness, numClusterHits);
430  hCleanlinessConversionDistance ->Fill(conversionDistance, cleanliness, numClusterHits);
431  hComplCleanlEnergy ->Fill(pi0Energy, cleanliness * completeness, numClusterHits);
432  hComplCleanlAngle ->Fill(pi0DecayAngle, cleanliness * completeness, numClusterHits);
433  hComplCleanlConversionDistance ->Fill(conversionDistance, cleanliness * completeness, numClusterHits);
434  hNumHitsEnergy ->Fill(pi0Energy, numClusterHits);
435 
436  // Continue if there are not two photons in the view
437  if (clusterMap[tpc][plane]->GetPhotons().size() != 2) continue;
438 
439  double conversionSeparation = GetEndTrackDistance(clusterMap[tpc][plane]->GetPhotons().at(0).first, clusterMap[tpc][plane]->GetPhotons().at(1).first);
440 
441  hCompletenessConversionSeparation->Fill(conversionSeparation, completeness, numClusterHits);
442  hCleanlinessConversionSeparation ->Fill(conversionSeparation, cleanliness, numClusterHits);
443  }
444 
445  }
446  }
447 
448 }
Float_t E
Definition: plot.C:23
double GetEndTrackDistance(TrackID id1, TrackID id2)
Cluster finding and building.
art::ServiceHandle< geo::Geometry > geometry
std::map< unsigned int, std::map< unsigned int, std::unique_ptr< ClusterCounter > > > clusterMap
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
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
Int_t min
Definition: plot.C:26
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:224
void ClusteringValidation::ClusterAnalyser::WriteFile ( )

Definition at line 450 of file ClusteringValidation_module.cc.

450  {
451 
452  // Average completeness/cleanliness
453  double avCompleteness = hCompleteness->GetMean();
454  double avCleanliness = hCleanliness ->GetMean();
455 
456  // Write file
457  std::ofstream outFile("effpur");
458  outFile << avCompleteness << " " << avCleanliness;
459  outFile.close();
460 
461 }

Member Data Documentation

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

Definition at line 209 of file ClusteringValidation_module.cc.

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

Definition at line 203 of file ClusteringValidation_module.cc.

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

Definition at line 191 of file ClusteringValidation_module.cc.

TObjArray ClusteringValidation::ClusterAnalyser::fHistArray
private

Definition at line 201 of file ClusteringValidation_module.cc.

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

Definition at line 207 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hCleanliness
private

Definition at line 194 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCleanlinessAngle
private

Definition at line 198 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCleanlinessConversionDistance
private

Definition at line 198 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCleanlinessConversionSeparation
private

Definition at line 198 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCleanlinessEnergy
private

Definition at line 198 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hComplCleanl
private

Definition at line 194 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hComplCleanlAngle
private

Definition at line 199 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hComplCleanlConversionDistance
private

Definition at line 199 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hComplCleanlConversionSeparation
private

Definition at line 199 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hComplCleanlEnergy
private

Definition at line 199 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hCompleteness
private

Definition at line 194 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCompletenessAngle
private

Definition at line 197 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCompletenessConversionDistance
private

Definition at line 197 of file ClusteringValidation_module.cc.

TProfile * ClusteringValidation::ClusterAnalyser::hCompletenessConversionSeparation
private

Definition at line 197 of file ClusteringValidation_module.cc.

TProfile* ClusteringValidation::ClusterAnalyser::hCompletenessEnergy
private

Definition at line 197 of file ClusteringValidation_module.cc.

TEfficiency* ClusteringValidation::ClusterAnalyser::hEfficiencyAngle
private

Definition at line 200 of file ClusteringValidation_module.cc.

TEfficiency * ClusteringValidation::ClusterAnalyser::hEfficiencyConversionDistance
private

Definition at line 200 of file ClusteringValidation_module.cc.

TEfficiency * ClusteringValidation::ClusterAnalyser::hEfficiencyConversionSeparation
private

Definition at line 200 of file ClusteringValidation_module.cc.

TEfficiency * ClusteringValidation::ClusterAnalyser::hEfficiencyEnergy
private

Definition at line 200 of file ClusteringValidation_module.cc.

TH2* ClusteringValidation::ClusterAnalyser::hNumHitsCompleteness
private

Definition at line 196 of file ClusteringValidation_module.cc.

TH2 * ClusteringValidation::ClusterAnalyser::hNumHitsEnergy
private

Definition at line 196 of file ClusteringValidation_module.cc.

TH1* ClusteringValidation::ClusterAnalyser::hPi0Angle
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0AngleCut
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionDistance
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionDistanceCut
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparation
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0ConversionSeparationCut
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0Energy
private

Definition at line 195 of file ClusteringValidation_module.cc.

TH1 * ClusteringValidation::ClusterAnalyser::hPi0EnergyCut
private

Definition at line 195 of file ClusteringValidation_module.cc.

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

Definition at line 208 of file ClusteringValidation_module.cc.

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

Definition at line 204 of file ClusteringValidation_module.cc.


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