LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
lar_content::ParticleRecoveryAlgorithm Class Reference

ParticleRecoveryAlgorithm class. More...

#include "ParticleRecoveryAlgorithm.h"

Inheritance diagram for lar_content::ParticleRecoveryAlgorithm:

Classes

class  SimpleOverlapTensor
 SimpleOverlapTensor class. More...
 

Public Member Functions

 ParticleRecoveryAlgorithm ()
 Default constructor. More...
 

Private Member Functions

pandora::StatusCode Run ()
 
void GetInputClusters (pandora::ClusterList &inputClusterListU, pandora::ClusterList &inputClusterListV, pandora::ClusterList &inputClusterListW) const
 Get the input cluster lists for processing in this algorithm. More...
 
void SelectInputClusters (const pandora::ClusterList &inputClusterList, pandora::ClusterList &selectedClusterList) const
 Select a subset of input clusters for processing in this algorithm. More...
 
void StandardClusterSelection (const pandora::ClusterList &inputClusterList, pandora::ClusterList &selectedClusterList) const
 Select a subset of input clusters for processing in this algorithm. More...
 
void VertexClusterSelection (const pandora::ClusterList &inputClusterList, pandora::ClusterList &selectedClusterList) const
 Select a subset of input clusters nodally associated with the vertices of existing particles. More...
 
void FindOverlaps (const pandora::ClusterList &clusterList1, const pandora::ClusterList &clusterList2, SimpleOverlapTensor &overlapTensor) const
 Find cluster overlaps and record these in the overlap tensor. More...
 
bool IsOverlap (const pandora::Cluster *const pCluster1, const pandora::Cluster *const pCluster2) const
 Whether two clusters overlap convincingly in x. More...
 
void CalculateEffectiveOverlapFractions (const pandora::Cluster *const pCluster1, const float xMin1, const float xMax1, const pandora::Cluster *const pCluster2, const float xMin2, const float xMax2, float &xOverlapFraction1, float &xOverlapFraction2) const
 Calculate effective overlap fractions taking into account gaps. More...
 
void CalculateEffectiveSpan (const pandora::Cluster *const pCluster, const float xMin, const float xMax, float &xMinEff, float &xMaxEff) const
 Calculate effective span for a given clsuter taking gaps into account. More...
 
void ExamineTensor (const SimpleOverlapTensor &overlapTensor) const
 Identify unambiguous cluster overlaps and resolve ambiguous overlaps, creating new track particles. More...
 
bool CheckConsistency (const pandora::Cluster *const pClusterU, const pandora::Cluster *const pClusterV, const pandora::Cluster *const pClusterW) const
 Whether a trio of clusters are consistent with representing projections of the same 3d trajectory. More...
 
void CreateTrackParticle (const pandora::ClusterList &clusterList) const
 Create and save a track particle containing the provided clusters. More...
 
pandora::StatusCode ReadSettings (const pandora::TiXmlHandle xmlHandle)
 

Private Attributes

pandora::StringVector m_inputClusterListNames
 The list of cluster list names. More...
 
std::string m_outputPfoListName
 The output pfo list name. More...
 
bool m_checkGaps
 Whether to check for gaps in the calculation of the overlap. More...
 
unsigned int m_minClusterCaloHits
 The min number of hits in base cluster selection method. More...
 
float m_minClusterLengthSquared
 The min length (squared) in base cluster selection method. More...
 
float m_minClusterXSpan
 The min x span required in order to consider a cluster. More...
 
bool m_vertexClusterMode
 Whether to demand clusters are associated with vertices of existing particles. More...
 
float m_minVertexLongitudinalDistance
 Vertex association check: min longitudinal distance cut. More...
 
float m_maxVertexTransverseDistance
 Vertex association check: max transverse distance cut. More...
 
float m_minXOverlapFraction
 The min x overlap fraction required in order to id overlapping clusters. More...
 
float m_minXOverlapFractionGaps
 The min x overlap fraction when there are gaps involved. More...
 
float m_sampleStepSize
 The sampling step size used in association checks, units cm. More...
 
unsigned int m_slidingFitHalfWindow
 The half window for the fit sliding result constructor. More...
 
float m_pseudoChi2Cut
 The selection cut on the matched chi2. More...
 

Detailed Description

ParticleRecoveryAlgorithm class.

Definition at line 21 of file ParticleRecoveryAlgorithm.h.

Constructor & Destructor Documentation

lar_content::ParticleRecoveryAlgorithm::ParticleRecoveryAlgorithm ( )

Default constructor.

Definition at line 26 of file ParticleRecoveryAlgorithm.cc.

26  :
27  m_checkGaps(true),
30  m_minClusterXSpan(0.25f),
31  m_vertexClusterMode(false),
36  m_sampleStepSize(0.5f),
39 {
40 }
unsigned int m_minClusterCaloHits
The min number of hits in base cluster selection method.
float m_minClusterLengthSquared
The min length (squared) in base cluster selection method.
float m_maxVertexTransverseDistance
Vertex association check: max transverse distance cut.
float m_minXOverlapFractionGaps
The min x overlap fraction when there are gaps involved.
TFile f
Definition: plotHisto.C:6
float m_minVertexLongitudinalDistance
Vertex association check: min longitudinal distance cut.
float m_minClusterXSpan
The min x span required in order to consider a cluster.
float m_pseudoChi2Cut
The selection cut on the matched chi2.
unsigned int m_slidingFitHalfWindow
The half window for the fit sliding result constructor.
float m_minXOverlapFraction
The min x overlap fraction required in order to id overlapping clusters.
bool m_checkGaps
Whether to check for gaps in the calculation of the overlap.
bool m_vertexClusterMode
Whether to demand clusters are associated with vertices of existing particles.
float m_sampleStepSize
The sampling step size used in association checks, units cm.

Member Function Documentation

void lar_content::ParticleRecoveryAlgorithm::CalculateEffectiveOverlapFractions ( const pandora::Cluster *const  pCluster1,
const float  xMin1,
const float  xMax1,
const pandora::Cluster *const  pCluster2,
const float  xMin2,
const float  xMax2,
float &  xOverlapFraction1,
float &  xOverlapFraction2 
) const
private

Calculate effective overlap fractions taking into account gaps.

Parameters
pCluster1address of the first cluster
xMin1min x value of the first cluster
xMax1max x value of the first cluster
pCluster2address of the second cluster
xMin2min x value of the second cluster
xMax2max x value of the second cluster
xOverlapFraction1to receive the effective overlap fraction for the first cluster
xOverlapFraction2to receive the effective overlap fraction for the second cluster

Definition at line 230 of file ParticleRecoveryAlgorithm.cc.

References CalculateEffectiveSpan(), f, max, and min.

Referenced by IsOverlap().

232 {
233  if (PandoraContentApi::GetGeometry(*this)->GetDetectorGapList().empty())
234  return;
235 
236  const float xMin(std::min(xMin1, xMin2));
237  const float xMax(std::max(xMax1, xMax2));
238  float xMinEff1(xMin1), xMaxEff1(xMax1), xMinEff2(xMin2), xMaxEff2(xMax2);
239 
240  this->CalculateEffectiveSpan(pCluster1, xMin, xMax, xMinEff1, xMaxEff1);
241  this->CalculateEffectiveSpan(pCluster2, xMin, xMax, xMinEff2, xMaxEff2);
242 
243  const float effectiveXSpan1(xMaxEff1 - xMinEff1), effectiveXSpan2(xMaxEff2 - xMinEff2);
244  const float effectiveXOverlapSpan(std::min(xMaxEff1, xMaxEff2) - std::max(xMinEff1, xMinEff2));
245 
246  if ((effectiveXSpan1 > std::numeric_limits<float>::epsilon()) && (effectiveXSpan2 > std::numeric_limits<float>::epsilon()))
247  {
248  xOverlapFraction1 = std::min(1.f, (effectiveXOverlapSpan / effectiveXSpan1));
249  xOverlapFraction2 = std::min(1.f, (effectiveXOverlapSpan / effectiveXSpan2));
250  }
251 }
TFile f
Definition: plotHisto.C:6
void CalculateEffectiveSpan(const pandora::Cluster *const pCluster, const float xMin, const float xMax, float &xMinEff, float &xMaxEff) const
Calculate effective span for a given clsuter taking gaps into account.
Int_t max
Definition: plot.C:27
Int_t min
Definition: plot.C:26
void lar_content::ParticleRecoveryAlgorithm::CalculateEffectiveSpan ( const pandora::Cluster *const  pCluster,
const float  xMin,
const float  xMax,
float &  xMinEff,
float &  xMaxEff 
) const
private

Calculate effective span for a given clsuter taking gaps into account.

Parameters
pClusteraddress of the cluster
xMinthe min x value above which checks for gaps will be performed
xMaxthe max x value below which checks for gaps will be performed
xMinEffto receive the effective min x value for the cluster, including adjacent gaps
xMaxEffto receive the effective max x value for the cluster, including adjacent gaps

Definition at line 255 of file ParticleRecoveryAlgorithm.cc.

References f, lar_content::LArGeometryHelper::GetWireZPitch(), lar_content::LArGeometryHelper::IsXSamplingPointInGap(), m_sampleStepSize, m_slidingFitHalfWindow, max, and min.

Referenced by CalculateEffectiveOverlapFractions().

256 {
257  // TODO cache sliding linear fit results and optimise protection against exceptions from TwoDSlidingFitResult and IsXSamplingPointInGap
258  try
259  {
260  const float slidingFitPitch(LArGeometryHelper::GetWireZPitch(this->GetPandora()));
261 
262  const TwoDSlidingFitResult slidingFitResult(pCluster, m_slidingFitHalfWindow, slidingFitPitch);
263 
264  const int nSamplingPointsLeft(1 + static_cast<int>((xMinEff - xMin) / m_sampleStepSize));
265  const int nSamplingPointsRight(1 + static_cast<int>((xMax - xMaxEff) / m_sampleStepSize));
266  float dxMin(0.f), dxMax(0.f);
267 
268  for (int iSample = 1; iSample <= nSamplingPointsLeft; ++iSample)
269  {
270  const float xSample(std::max(xMin, xMinEff - static_cast<float>(iSample) * m_sampleStepSize));
271 
272  if (!LArGeometryHelper::IsXSamplingPointInGap(this->GetPandora(), xSample, slidingFitResult, m_sampleStepSize))
273  break;
274 
275  dxMin = xMinEff - xSample;
276  }
277 
278  for (int iSample = 1; iSample <= nSamplingPointsRight; ++iSample)
279  {
280  const float xSample(std::min(xMax, xMaxEff + static_cast<float>(iSample) * m_sampleStepSize));
281 
282  if (!LArGeometryHelper::IsXSamplingPointInGap(this->GetPandora(), xSample, slidingFitResult, m_sampleStepSize))
283  break;
284 
285  dxMax = xSample - xMaxEff;
286  }
287 
288  xMinEff = xMinEff - dxMin;
289  xMaxEff = xMaxEff + dxMax;
290  }
291  catch (StatusCodeException &) {}
292 }
static float GetWireZPitch(const pandora::Pandora &pandora, const float maxWirePitchDiscrepancy=0.01)
Return the wire pitch.
TFile f
Definition: plotHisto.C:6
Int_t max
Definition: plot.C:27
static bool IsXSamplingPointInGap(const pandora::Pandora &pandora, const float xSample, const TwoDSlidingFitResult &slidingFitResult, const float gapTolerance=0.f)
Whether there is a gap in a cluster (described via its sliding fit result) at a specified x sampling ...
unsigned int m_slidingFitHalfWindow
The half window for the fit sliding result constructor.
Int_t min
Definition: plot.C:26
float m_sampleStepSize
The sampling step size used in association checks, units cm.
bool lar_content::ParticleRecoveryAlgorithm::CheckConsistency ( const pandora::Cluster *const  pClusterU,
const pandora::Cluster *const  pClusterV,
const pandora::Cluster *const  pClusterW 
) const
private

Whether a trio of clusters are consistent with representing projections of the same 3d trajectory.

Parameters
pClusterUthe address of cluster u
pClusterVthe address of cluster v
pClusterWthe address of cluster w
Returns
boolean

Definition at line 334 of file ParticleRecoveryAlgorithm.cc.

References f, lar_content::LArClusterHelper::GetAverageZ(), lar_content::LArClusterHelper::GetClusterSpanX(), m_pseudoChi2Cut, max, lar_content::LArGeometryHelper::MergeTwoPositions(), min, and w.

Referenced by ExamineTensor().

335 {
336  // Requirements on X matching
337  float xMinU(0.f), xMinV(0.f), xMinW(0.f), xMaxU(0.f), xMaxV(0.f), xMaxW(0.f);
338  LArClusterHelper::GetClusterSpanX(pClusterU, xMinU, xMaxU);
339  LArClusterHelper::GetClusterSpanX(pClusterV, xMinV, xMaxV);
340  LArClusterHelper::GetClusterSpanX(pClusterW, xMinW, xMaxW);
341 
342  const float xMin(std::max(xMinU, std::max(xMinV, xMinW)));
343  const float xMax(std::min(xMaxU, std::min(xMaxV, xMaxW)));
344  const float xOverlap(xMax - xMin);
345 
346  if (xOverlap < std::numeric_limits<float>::epsilon())
347  return false;
348 
349  // Requirements on 3D matching
351 
352  if ((STATUS_CODE_SUCCESS != LArClusterHelper::GetAverageZ(pClusterU, xMin, xMax, u)) ||
353  (STATUS_CODE_SUCCESS != LArClusterHelper::GetAverageZ(pClusterV, xMin, xMax, v)) ||
354  (STATUS_CODE_SUCCESS != LArClusterHelper::GetAverageZ(pClusterW, xMin, xMax, w)))
355  {
356  return false;
357  }
358 
359  const float uv2w(LArGeometryHelper::MergeTwoPositions(this->GetPandora(), TPC_VIEW_U, TPC_VIEW_V, u, v));
360  const float vw2u(LArGeometryHelper::MergeTwoPositions(this->GetPandora(), TPC_VIEW_V, TPC_VIEW_W, v, w));
361  const float wu2v(LArGeometryHelper::MergeTwoPositions(this->GetPandora(), TPC_VIEW_W, TPC_VIEW_U, w, u));
362 
363  const float pseudoChi2(((u - vw2u) * (u - vw2u) + (v - wu2v) * (v - wu2v) + (w - uv2w) * (w - uv2w)) / 3.f);
364 
365  if (pseudoChi2 > m_pseudoChi2Cut)
366  return false;
367 
368  return true;
369 }
static pandora::StatusCode GetAverageZ(const pandora::Cluster *const pCluster, const float xmin, const float xmax, float &averageZ)
Get average Z positions of the calo hits in a cluster in range xmin to xmax.
TFile f
Definition: plotHisto.C:6
Int_t max
Definition: plot.C:27
static float MergeTwoPositions(const pandora::Pandora &pandora, const pandora::HitType view1, const pandora::HitType view2, const float position1, const float position2)
Merge two views (U,V) to give a third view (Z).
float m_pseudoChi2Cut
The selection cut on the matched chi2.
Int_t min
Definition: plot.C:26
static void GetClusterSpanX(const pandora::Cluster *const pCluster, float &xmin, float &xmax)
Get minimum and maximum X positions of the calo hits in a cluster.
Float_t w
Definition: plot.C:23
void lar_content::ParticleRecoveryAlgorithm::CreateTrackParticle ( const pandora::ClusterList &  clusterList) const
private

Create and save a track particle containing the provided clusters.

Parameters
clusterListthe cluster list

Definition at line 373 of file ParticleRecoveryAlgorithm.cc.

References f, and m_outputPfoListName.

Referenced by ExamineTensor().

374 {
375  if (clusterList.size() < 2)
376  throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
377 
378  const PfoList *pPfoList = NULL; std::string pfoListName;
379  PANDORA_THROW_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::CreateTemporaryListAndSetCurrent(*this, pPfoList, pfoListName));
380 
381  // TODO Correct these placeholder parameters
382  PandoraContentApi::ParticleFlowObject::Parameters pfoParameters;
383  pfoParameters.m_particleId = MU_MINUS;
384  pfoParameters.m_charge = PdgTable::GetParticleCharge(pfoParameters.m_particleId.Get());
385  pfoParameters.m_mass = PdgTable::GetParticleMass(pfoParameters.m_particleId.Get());
386  pfoParameters.m_energy = 0.f;
387  pfoParameters.m_momentum = CartesianVector(0.f, 0.f, 0.f);
388  pfoParameters.m_clusterList = clusterList;
389 
390  const ParticleFlowObject *pPfo(NULL);
391  PANDORA_THROW_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::ParticleFlowObject::Create(*this, pfoParameters, pPfo));
392 
393  if (!pPfoList->empty())
394  {
395  PANDORA_THROW_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::SaveList<Pfo>(*this, m_outputPfoListName));
396  PANDORA_THROW_RESULT_IF(STATUS_CODE_SUCCESS, !=, PandoraContentApi::ReplaceCurrentList<Pfo>(*this, m_outputPfoListName));
397  }
398 }
TFile f
Definition: plotHisto.C:6
std::string m_outputPfoListName
The output pfo list name.
void lar_content::ParticleRecoveryAlgorithm::ExamineTensor ( const SimpleOverlapTensor overlapTensor) const
private

Identify unambiguous cluster overlaps and resolve ambiguous overlaps, creating new track particles.

Parameters
overlapTensorthe overlap tensor

Definition at line 296 of file ParticleRecoveryAlgorithm.cc.

References CheckConsistency(), CreateTrackParticle(), lar_content::ParticleRecoveryAlgorithm::SimpleOverlapTensor::GetConnectedElements(), lar_content::ParticleRecoveryAlgorithm::SimpleOverlapTensor::GetKeyClusters(), and lar_content::LArClusterHelper::SortByNHits().

Referenced by Run().

297 {
298  ClusterVector sortedKeyClusters(overlapTensor.GetKeyClusters().begin(), overlapTensor.GetKeyClusters().end());
299  std::sort(sortedKeyClusters.begin(), sortedKeyClusters.end(), LArClusterHelper::SortByNHits);
300 
301  for (const Cluster *const pKeyCluster : sortedKeyClusters)
302  {
303  ClusterList clusterListU, clusterListV, clusterListW;
304 
305  overlapTensor.GetConnectedElements(pKeyCluster, true, clusterListU, clusterListV, clusterListW);
306  const unsigned int nU(clusterListU.size()), nV(clusterListV.size()), nW(clusterListW.size());
307 
308  if ((0 == nU * nV) && (0 == nV * nW) && (0 == nW * nU))
309  continue;
310 
311  ClusterList clusterListAll;
312  clusterListAll.insert(clusterListAll.end(), clusterListU.begin(), clusterListU.end());
313  clusterListAll.insert(clusterListAll.end(), clusterListV.begin(), clusterListV.end());
314  clusterListAll.insert(clusterListAll.end(), clusterListW.begin(), clusterListW.end());
315 
316  if ((1 == nU * nV * nW) && this->CheckConsistency(*(clusterListU.begin()), *(clusterListV.begin()), *(clusterListW.begin())))
317  {
318  this->CreateTrackParticle(clusterListAll);
319  }
320  else if ((0 == nU * nV * nW) && ((1 == nU && 1 == nV) || (1 == nV && 1 == nW) || (1 == nW && 1 == nU)))
321  {
322  this->CreateTrackParticle(clusterListAll);
323  }
324  else
325  {
326  // TODO - check here whether there is a gap in the 2 in one view when 1:2:0
327  // May later choose to resolve simple ambiguities, e.g. of form nU:nV:nW == 1:2:0
328  }
329  }
330 }
static bool SortByNHits(const pandora::Cluster *const pLhs, const pandora::Cluster *const pRhs)
Sort clusters by number of hits, then layer span, then inner layer, then position, then pulse-height.
bool CheckConsistency(const pandora::Cluster *const pClusterU, const pandora::Cluster *const pClusterV, const pandora::Cluster *const pClusterW) const
Whether a trio of clusters are consistent with representing projections of the same 3d trajectory...
void CreateTrackParticle(const pandora::ClusterList &clusterList) const
Create and save a track particle containing the provided clusters.
std::vector< art::Ptr< recob::Cluster > > ClusterVector
void lar_content::ParticleRecoveryAlgorithm::FindOverlaps ( const pandora::ClusterList &  clusterList1,
const pandora::ClusterList &  clusterList2,
SimpleOverlapTensor overlapTensor 
) const
private

Find cluster overlaps and record these in the overlap tensor.

Parameters
clusterList1the first cluster list
clusterList2the second cluster list
overlapTensorthe overlap tensor

Definition at line 184 of file ParticleRecoveryAlgorithm.cc.

References lar_content::ParticleRecoveryAlgorithm::SimpleOverlapTensor::AddAssociation(), and IsOverlap().

Referenced by Run().

185 {
186  for (ClusterList::const_iterator iter1 = clusterList1.begin(), iter1End = clusterList1.end(); iter1 != iter1End; ++iter1)
187  {
188  for (ClusterList::const_iterator iter2 = clusterList2.begin(), iter2End = clusterList2.end(); iter2 != iter2End; ++iter2)
189  {
190  if (this->IsOverlap(*iter1, *iter2))
191  overlapTensor.AddAssociation(*iter1, *iter2);
192  }
193  }
194 }
intermediate_table::const_iterator const_iterator
bool IsOverlap(const pandora::Cluster *const pCluster1, const pandora::Cluster *const pCluster2) const
Whether two clusters overlap convincingly in x.
void lar_content::ParticleRecoveryAlgorithm::GetInputClusters ( pandora::ClusterList &  inputClusterListU,
pandora::ClusterList &  inputClusterListV,
pandora::ClusterList &  inputClusterListW 
) const
private

Get the input cluster lists for processing in this algorithm.

Parameters
inputClusterListUto receive the list of clusters in the u view
inputClusterListUto receive the list of clusters in the v view
inputClusterListUto receive the list of clusters in the w view

Definition at line 65 of file ParticleRecoveryAlgorithm.cc.

References lar_content::LArClusterHelper::GetClusterHitType(), and m_inputClusterListNames.

Referenced by Run().

66 {
67  for (StringVector::const_iterator iter = m_inputClusterListNames.begin(), iterEnd = m_inputClusterListNames.end(); iter != iterEnd; ++iter)
68  {
69  const ClusterList *pClusterList(NULL);
70  PANDORA_THROW_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_INITIALIZED, !=, PandoraContentApi::GetList(*this, *iter, pClusterList));
71 
72  if (!pClusterList || pClusterList->empty())
73  {
74  if (PandoraContentApi::GetSettings(*this)->ShouldDisplayAlgorithmInfo())
75  std::cout << "ParticleRecoveryAlgorithm: unable to find cluster list " << *iter << std::endl;
76 
77  continue;
78  }
79 
80  for (ClusterList::const_iterator cIter = pClusterList->begin(), cIterEnd = pClusterList->end(); cIter != cIterEnd; ++cIter)
81  {
82  const Cluster *const pCluster(*cIter);
83  const HitType hitType(LArClusterHelper::GetClusterHitType(pCluster));
84 
85  if ((TPC_VIEW_U != hitType) && (TPC_VIEW_V != hitType) && (TPC_VIEW_W != hitType))
86  throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
87 
88  ClusterList &clusterList((TPC_VIEW_U == hitType) ? inputClusterListU : (TPC_VIEW_V == hitType) ? inputClusterListV : inputClusterListW);
89  clusterList.push_back(pCluster);
90  }
91  }
92 }
pandora::StringVector m_inputClusterListNames
The list of cluster list names.
static pandora::HitType GetClusterHitType(const pandora::Cluster *const pCluster)
Get the hit type associated with a two dimensional cluster.
intermediate_table::const_iterator const_iterator
bool lar_content::ParticleRecoveryAlgorithm::IsOverlap ( const pandora::Cluster *const  pCluster1,
const pandora::Cluster *const  pCluster2 
) const
private

Whether two clusters overlap convincingly in x.

Parameters
pCluster1address of the first cluster
pCluster2address of the second cluster

Definition at line 198 of file ParticleRecoveryAlgorithm.cc.

References CalculateEffectiveOverlapFractions(), f, lar_content::LArClusterHelper::GetClusterHitType(), lar_content::LArClusterHelper::GetClusterSpanX(), m_checkGaps, m_minXOverlapFraction, max, and min.

Referenced by FindOverlaps().

199 {
201  throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
202 
203  if ((0 == pCluster1->GetNCaloHits()) || (0 == pCluster2->GetNCaloHits()))
204  throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
205 
206  float xMin1(0.f), xMax1(0.f), xMin2(0.f), xMax2(0.f);
207  LArClusterHelper::GetClusterSpanX(pCluster1, xMin1, xMax1);
208  LArClusterHelper::GetClusterSpanX(pCluster2, xMin2, xMax2);
209 
210  const float xSpan1(xMax1 - xMin1), xSpan2(xMax2 - xMin2);
211 
212  if ((xSpan1 < std::numeric_limits<float>::epsilon()) || (xSpan2 < std::numeric_limits<float>::epsilon()))
213  throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
214 
215  const float xOverlap(std::min(xMax1, xMax2) - std::max(xMin1, xMin2));
216 
217  float xOverlapFraction1(xOverlap / xSpan1), xOverlapFraction2(xOverlap / xSpan2);
218 
219  if (m_checkGaps)
220  this->CalculateEffectiveOverlapFractions(pCluster1, xMin1, xMax1, pCluster2, xMin2, xMax2, xOverlapFraction1, xOverlapFraction2);
221 
222  if ((xOverlapFraction1 < m_minXOverlapFraction) || (xOverlapFraction2 < m_minXOverlapFraction))
223  return false;
224 
225  return true;
226 }
static pandora::HitType GetClusterHitType(const pandora::Cluster *const pCluster)
Get the hit type associated with a two dimensional cluster.
TFile f
Definition: plotHisto.C:6
Int_t max
Definition: plot.C:27
void CalculateEffectiveOverlapFractions(const pandora::Cluster *const pCluster1, const float xMin1, const float xMax1, const pandora::Cluster *const pCluster2, const float xMin2, const float xMax2, float &xOverlapFraction1, float &xOverlapFraction2) const
Calculate effective overlap fractions taking into account gaps.
float m_minXOverlapFraction
The min x overlap fraction required in order to id overlapping clusters.
bool m_checkGaps
Whether to check for gaps in the calculation of the overlap.
Int_t min
Definition: plot.C:26
static void GetClusterSpanX(const pandora::Cluster *const pCluster, float &xmin, float &xmax)
Get minimum and maximum X positions of the calo hits in a cluster.
StatusCode lar_content::ParticleRecoveryAlgorithm::ReadSettings ( const pandora::TiXmlHandle  xmlHandle)
private

Definition at line 472 of file ParticleRecoveryAlgorithm.cc.

References m_checkGaps, m_inputClusterListNames, m_maxVertexTransverseDistance, m_minClusterCaloHits, m_minClusterLengthSquared, m_minClusterXSpan, m_minVertexLongitudinalDistance, m_minXOverlapFraction, m_minXOverlapFractionGaps, m_outputPfoListName, m_pseudoChi2Cut, m_sampleStepSize, m_slidingFitHalfWindow, and m_vertexClusterMode.

473 {
474  PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, XmlHelper::ReadVectorOfValues(xmlHandle, "InputClusterListNames", m_inputClusterListNames));
475  PANDORA_RETURN_RESULT_IF(STATUS_CODE_SUCCESS, !=, XmlHelper::ReadValue(xmlHandle, "OutputPfoListName", m_outputPfoListName));
476 
477  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
478  "MinClusterCaloHits", m_minClusterCaloHits));
479 
480  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
481  "CheckGaps", m_checkGaps));
482 
483  float minClusterLength = std::sqrt(m_minClusterLengthSquared);
484  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
485  "MinClusterLength", minClusterLength));
486  m_minClusterLengthSquared = minClusterLength * minClusterLength;
487 
488  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
489  "MinClusterXSpan", m_minClusterXSpan));
490 
491  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
492  "VertexClusterMode", m_vertexClusterMode));
493 
494  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
495  "MinVertexLongitudinalDistance", m_minVertexLongitudinalDistance));
496 
497  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
498  "MaxVertexTransverseDistance", m_maxVertexTransverseDistance));
499 
500  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
501  "MinXOverlapFraction", m_minXOverlapFraction));
502 
503  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
504  "MinXOverlapFractionGaps", m_minXOverlapFractionGaps));
505 
506  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
507  "SampleStepSize", m_sampleStepSize));
508 
509  if (m_sampleStepSize < std::numeric_limits<float>::epsilon())
510  {
511  std::cout << "ParticleRecoveryAlgorithm: Invalid value for SampleStepSize " << m_sampleStepSize << std::endl;
512  throw StatusCodeException(STATUS_CODE_INVALID_PARAMETER);
513  }
514 
515  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
516  "SlidingFitHalfWindow", m_slidingFitHalfWindow));
517 
518  PANDORA_RETURN_RESULT_IF_AND_IF(STATUS_CODE_SUCCESS, STATUS_CODE_NOT_FOUND, !=, XmlHelper::ReadValue(xmlHandle,
519  "PseudoChi2Cut", m_pseudoChi2Cut));
520 
521  return STATUS_CODE_SUCCESS;
522 }
pandora::StringVector m_inputClusterListNames
The list of cluster list names.
unsigned int m_minClusterCaloHits
The min number of hits in base cluster selection method.
float m_minClusterLengthSquared
The min length (squared) in base cluster selection method.
float m_maxVertexTransverseDistance
Vertex association check: max transverse distance cut.
float m_minXOverlapFractionGaps
The min x overlap fraction when there are gaps involved.
float m_minVertexLongitudinalDistance
Vertex association check: min longitudinal distance cut.
float m_minClusterXSpan
The min x span required in order to consider a cluster.
float m_pseudoChi2Cut
The selection cut on the matched chi2.
unsigned int m_slidingFitHalfWindow
The half window for the fit sliding result constructor.
float m_minXOverlapFraction
The min x overlap fraction required in order to id overlapping clusters.
bool m_checkGaps
Whether to check for gaps in the calculation of the overlap.
bool m_vertexClusterMode
Whether to demand clusters are associated with vertices of existing particles.
std::string m_outputPfoListName
The output pfo list name.
float m_sampleStepSize
The sampling step size used in association checks, units cm.
StatusCode lar_content::ParticleRecoveryAlgorithm::Run ( )
private

Definition at line 44 of file ParticleRecoveryAlgorithm.cc.

References ExamineTensor(), FindOverlaps(), GetInputClusters(), and SelectInputClusters().

45 {
46  ClusterList inputClusterListU, inputClusterListV, inputClusterListW;
47  this->GetInputClusters(inputClusterListU, inputClusterListV, inputClusterListW);
48 
49  ClusterList selectedClusterListU, selectedClusterListV, selectedClusterListW;
50  this->SelectInputClusters(inputClusterListU, selectedClusterListU);
51  this->SelectInputClusters(inputClusterListV, selectedClusterListV);
52  this->SelectInputClusters(inputClusterListW, selectedClusterListW);
53 
54  SimpleOverlapTensor overlapTensor;
55  this->FindOverlaps(selectedClusterListU, selectedClusterListV, overlapTensor);
56  this->FindOverlaps(selectedClusterListV, selectedClusterListW, overlapTensor);
57  this->FindOverlaps(selectedClusterListW, selectedClusterListU, overlapTensor);
58  this->ExamineTensor(overlapTensor);
59 
60  return STATUS_CODE_SUCCESS;
61 }
void GetInputClusters(pandora::ClusterList &inputClusterListU, pandora::ClusterList &inputClusterListV, pandora::ClusterList &inputClusterListW) const
Get the input cluster lists for processing in this algorithm.
void ExamineTensor(const SimpleOverlapTensor &overlapTensor) const
Identify unambiguous cluster overlaps and resolve ambiguous overlaps, creating new track particles...
void SelectInputClusters(const pandora::ClusterList &inputClusterList, pandora::ClusterList &selectedClusterList) const
Select a subset of input clusters for processing in this algorithm.
void FindOverlaps(const pandora::ClusterList &clusterList1, const pandora::ClusterList &clusterList2, SimpleOverlapTensor &overlapTensor) const
Find cluster overlaps and record these in the overlap tensor.
void lar_content::ParticleRecoveryAlgorithm::SelectInputClusters ( const pandora::ClusterList &  inputClusterList,
pandora::ClusterList &  selectedClusterList 
) const
private

Select a subset of input clusters for processing in this algorithm.

Parameters
inputClusterListthe input cluster list
selectedClusterListto receive the selected cluster list

Definition at line 96 of file ParticleRecoveryAlgorithm.cc.

References m_vertexClusterMode, StandardClusterSelection(), and VertexClusterSelection().

Referenced by Run().

97 {
99  {
100  ClusterList vertexClusterList;
101  this->VertexClusterSelection(inputClusterList, vertexClusterList);
102  this->StandardClusterSelection(vertexClusterList, selectedClusterList);
103  }
104  else
105  {
106  this->StandardClusterSelection(inputClusterList, selectedClusterList);
107  }
108 }
bool m_vertexClusterMode
Whether to demand clusters are associated with vertices of existing particles.
void StandardClusterSelection(const pandora::ClusterList &inputClusterList, pandora::ClusterList &selectedClusterList) const
Select a subset of input clusters for processing in this algorithm.
void VertexClusterSelection(const pandora::ClusterList &inputClusterList, pandora::ClusterList &selectedClusterList) const
Select a subset of input clusters nodally associated with the vertices of existing particles...
void lar_content::ParticleRecoveryAlgorithm::StandardClusterSelection ( const pandora::ClusterList &  inputClusterList,
pandora::ClusterList &  selectedClusterList 
) const
private

Select a subset of input clusters for processing in this algorithm.

Parameters
inputClusterListthe input cluster list
selectedClusterListto receive the selected cluster list

Definition at line 112 of file ParticleRecoveryAlgorithm.cc.

References f, lar_content::LArClusterHelper::GetClusterSpanX(), lar_content::LArClusterHelper::GetLengthSquared(), m_minClusterCaloHits, m_minClusterLengthSquared, and m_minClusterXSpan.

Referenced by SelectInputClusters().

113 {
114  for (ClusterList::const_iterator iter = inputClusterList.begin(), iterEnd = inputClusterList.end(); iter != iterEnd; ++iter)
115  {
116  const Cluster *const pCluster = *iter;
117 
118  if (!pCluster->IsAvailable())
119  continue;
120 
121  if (pCluster->GetNCaloHits() < m_minClusterCaloHits)
122  continue;
123 
125  continue;
126 
127  float xMin(0.f), xMax(0.f);
128  LArClusterHelper::GetClusterSpanX(pCluster, xMin, xMax);
129 
130  if ((xMax - xMin) < m_minClusterXSpan)
131  continue;
132 
133  selectedClusterList.push_back(pCluster);
134  }
135 }
unsigned int m_minClusterCaloHits
The min number of hits in base cluster selection method.
float m_minClusterLengthSquared
The min length (squared) in base cluster selection method.
TFile f
Definition: plotHisto.C:6
intermediate_table::const_iterator const_iterator
float m_minClusterXSpan
The min x span required in order to consider a cluster.
static void GetClusterSpanX(const pandora::Cluster *const pCluster, float &xmin, float &xmax)
Get minimum and maximum X positions of the calo hits in a cluster.
static float GetLengthSquared(const pandora::Cluster *const pCluster)
Get length squared of cluster.
void lar_content::ParticleRecoveryAlgorithm::VertexClusterSelection ( const pandora::ClusterList &  inputClusterList,
pandora::ClusterList &  selectedClusterList 
) const
private

Select a subset of input clusters nodally associated with the vertices of existing particles.

Parameters
inputClusterListthe input cluster list
selectedClusterListto receive the selected cluster list

Definition at line 139 of file ParticleRecoveryAlgorithm.cc.

References lar_content::LArPointingCluster::GetInnerVertex(), lar_content::LArPointingCluster::GetOuterVertex(), lar_content::LArPointingCluster::Vertex::GetPosition(), lar_content::LArPointingClusterHelper::IsNode(), m_maxVertexTransverseDistance, and m_minVertexLongitudinalDistance.

Referenced by SelectInputClusters().

140 {
141  CartesianPointVector vertexList;
142 
143  for (ClusterList::const_iterator iter = inputClusterList.begin(), iterEnd = inputClusterList.end(); iter != iterEnd; ++iter)
144  {
145  try
146  {
147  if (!(*iter)->IsAvailable())
148  {
149  const LArPointingCluster pointingCluster(*iter);
150  vertexList.push_back(pointingCluster.GetInnerVertex().GetPosition());
151  vertexList.push_back(pointingCluster.GetOuterVertex().GetPosition());
152  }
153  }
154  catch (StatusCodeException &) {}
155  }
156 
157  for (ClusterList::const_iterator iter = inputClusterList.begin(), iterEnd = inputClusterList.end(); iter != iterEnd; ++iter)
158  {
159  try
160  {
161  const Cluster *const pCluster = *iter;
162 
163  if (!pCluster->IsAvailable())
164  continue;
165 
166  const LArPointingCluster pointingCluster(pCluster);
167 
168  for (CartesianPointVector::const_iterator vIter = vertexList.begin(), vIterEnd = vertexList.end(); vIter != vIterEnd; ++vIter)
169  {
172  {
173  selectedClusterList.push_back(pCluster);
174  break;
175  }
176  }
177  }
178  catch (StatusCodeException &) {}
179  }
180 }
float m_maxVertexTransverseDistance
Vertex association check: max transverse distance cut.
intermediate_table::const_iterator const_iterator
float m_minVertexLongitudinalDistance
Vertex association check: min longitudinal distance cut.
static bool IsNode(const pandora::CartesianVector &parentVertex, const LArPointingCluster::Vertex &daughterVertex, const float minLongitudinalDistance, const float maxTransverseDistance)
Whether pointing vertex is adjacent to a given position.

Member Data Documentation

bool lar_content::ParticleRecoveryAlgorithm::m_checkGaps
private

Whether to check for gaps in the calculation of the overlap.

Definition at line 180 of file ParticleRecoveryAlgorithm.h.

Referenced by IsOverlap(), and ReadSettings().

pandora::StringVector lar_content::ParticleRecoveryAlgorithm::m_inputClusterListNames
private

The list of cluster list names.

Definition at line 177 of file ParticleRecoveryAlgorithm.h.

Referenced by GetInputClusters(), and ReadSettings().

float lar_content::ParticleRecoveryAlgorithm::m_maxVertexTransverseDistance
private

Vertex association check: max transverse distance cut.

Definition at line 188 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings(), and VertexClusterSelection().

unsigned int lar_content::ParticleRecoveryAlgorithm::m_minClusterCaloHits
private

The min number of hits in base cluster selection method.

Definition at line 182 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings(), and StandardClusterSelection().

float lar_content::ParticleRecoveryAlgorithm::m_minClusterLengthSquared
private

The min length (squared) in base cluster selection method.

Definition at line 183 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings(), and StandardClusterSelection().

float lar_content::ParticleRecoveryAlgorithm::m_minClusterXSpan
private

The min x span required in order to consider a cluster.

Definition at line 184 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings(), and StandardClusterSelection().

float lar_content::ParticleRecoveryAlgorithm::m_minVertexLongitudinalDistance
private

Vertex association check: min longitudinal distance cut.

Definition at line 187 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings(), and VertexClusterSelection().

float lar_content::ParticleRecoveryAlgorithm::m_minXOverlapFraction
private

The min x overlap fraction required in order to id overlapping clusters.

Definition at line 190 of file ParticleRecoveryAlgorithm.h.

Referenced by IsOverlap(), and ReadSettings().

float lar_content::ParticleRecoveryAlgorithm::m_minXOverlapFractionGaps
private

The min x overlap fraction when there are gaps involved.

Definition at line 191 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings().

std::string lar_content::ParticleRecoveryAlgorithm::m_outputPfoListName
private

The output pfo list name.

Definition at line 178 of file ParticleRecoveryAlgorithm.h.

Referenced by CreateTrackParticle(), and ReadSettings().

float lar_content::ParticleRecoveryAlgorithm::m_pseudoChi2Cut
private

The selection cut on the matched chi2.

Definition at line 194 of file ParticleRecoveryAlgorithm.h.

Referenced by CheckConsistency(), and ReadSettings().

float lar_content::ParticleRecoveryAlgorithm::m_sampleStepSize
private

The sampling step size used in association checks, units cm.

Definition at line 192 of file ParticleRecoveryAlgorithm.h.

Referenced by CalculateEffectiveSpan(), and ReadSettings().

unsigned int lar_content::ParticleRecoveryAlgorithm::m_slidingFitHalfWindow
private

The half window for the fit sliding result constructor.

Definition at line 193 of file ParticleRecoveryAlgorithm.h.

Referenced by CalculateEffectiveSpan(), and ReadSettings().

bool lar_content::ParticleRecoveryAlgorithm::m_vertexClusterMode
private

Whether to demand clusters are associated with vertices of existing particles.

Definition at line 186 of file ParticleRecoveryAlgorithm.h.

Referenced by ReadSettings(), and SelectInputClusters().


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