All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ShowerPCADirection_tool.cc
Go to the documentation of this file.
1 //############################################################################
2 //### Name: ShowerPCADirection ###
3 //### Author: Dominic Barker and Ed Tyley ###
4 //### Date: 13.05.19 ###
5 //### Description: Tool for finding the shower direction using PCA ###
6 //### methods. Derived from LArPandoraModularShowers Method. ###
7 //############################################################################
8 
9 //Framework Includes
10 #include "art/Utilities/ToolMacros.h"
11 
12 //LArSoft Includes
22 
23 //C++ Includes
24 #include <Eigen/Dense>
25 
26 namespace ShowerRecoTools {
27 
29 
30  public:
31  ShowerPCADirection(const fhicl::ParameterSet& pset);
32 
33  //Calculate the direction of the shower.
34  int CalculateElement(const art::Ptr<recob::PFParticle>& pfparticle,
35  art::Event& Event,
36  reco::shower::ShowerElementHolder& ShowerEleHolder) override;
37 
38  private:
39  void InitialiseProducers() override;
40 
41  //Function to add the assoctions
42  int AddAssociations(const art::Ptr<recob::PFParticle>& pfpPtr,
43  art::Event& Event,
44  reco::shower::ShowerElementHolder& ShowerEleHolder) override;
45 
46  // Define standard art tool interface
48  const detinfo::DetectorClocksData& clockData,
50  const std::vector<art::Ptr<recob::SpacePoint>>& spacePoints_pfp,
51  const art::FindManyP<recob::Hit>& fmh,
52  TVector3& ShowerCentre);
53 
54  TVector3 GetPCAxisVector(recob::PCAxis& PCAxis);
55 
56  //fcl
57  art::InputTag fPFParticleLabel;
58  int fVerbose;
59  unsigned int
60  fNSegments; //Used in the RMS gradient. How many segments should we split the shower into.
61  bool fUseStartPosition; //If we use the start position the drection of the
62  //PCA vector is decided as (Shower Centre - Shower Start Position).
63  bool fChargeWeighted; //Should the PCA axis be charge weighted.
64 
68  std::string fShowerPCAOutputLabel;
69  };
70 
71  ShowerPCADirection::ShowerPCADirection(const fhicl::ParameterSet& pset)
72  : IShowerTool(pset.get<fhicl::ParameterSet>("BaseTools"))
73  , fPFParticleLabel(pset.get<art::InputTag>("PFParticleLabel"))
74  , fVerbose(pset.get<int>("Verbose"))
75  , fNSegments(pset.get<unsigned int>("NSegments"))
76  , fUseStartPosition(pset.get<bool>("UseStartPosition"))
77  , fChargeWeighted(pset.get<bool>("ChargeWeighted"))
78  , fShowerStartPositionInputLabel(pset.get<std::string>("ShowerStartPositionInputLabel"))
79  , fShowerDirectionOutputLabel(pset.get<std::string>("ShowerDirectionOutputLabel"))
80  , fShowerCentreOutputLabel(pset.get<std::string>("ShowerCentreOutputLabel"))
81  , fShowerPCAOutputLabel(pset.get<std::string>("ShowerPCAOutputLabel"))
82  {}
83 
84  void
86  {
87  InitialiseProduct<std::vector<recob::PCAxis>>(fShowerPCAOutputLabel);
88  InitialiseProduct<art::Assns<recob::Shower, recob::PCAxis>>("ShowerPCAxisAssn");
89  InitialiseProduct<art::Assns<recob::PFParticle, recob::PCAxis>>("PFParticlePCAxisAssn");
90  }
91 
92  int
93  ShowerPCADirection::CalculateElement(const art::Ptr<recob::PFParticle>& pfparticle,
94  art::Event& Event,
95  reco::shower::ShowerElementHolder& ShowerEleHolder)
96  {
97 
98  // Get the assocated pfParicle vertex PFParticles
99  auto const pfpHandle = Event.getValidHandle<std::vector<recob::PFParticle>>(fPFParticleLabel);
100 
101  const art::FindManyP<recob::SpacePoint>& fmspp =
102  ShowerEleHolder.GetFindManyP<recob::SpacePoint>(pfpHandle, Event, fPFParticleLabel);
103 
104  //Get the spacepoints handle and the hit assoication
105  auto const spHandle = Event.getValidHandle<std::vector<recob::SpacePoint>>(fPFParticleLabel);
106 
107  const art::FindManyP<recob::Hit>& fmh =
108  ShowerEleHolder.GetFindManyP<recob::Hit>(spHandle, Event, fPFParticleLabel);
109 
110  //SpacePoints
111  std::vector<art::Ptr<recob::SpacePoint>> spacePoints_pfp = fmspp.at(pfparticle.key());
112 
113  //We cannot progress with no spacepoints.
114  if (spacePoints_pfp.size() < 3) {
115  mf::LogWarning("ShowerPCADirection")
116  << spacePoints_pfp.size() << " spacepoints in shower, not calculating direction"
117  << std::endl;
118  return 1;
119  }
120 
121  auto const clockData =
122  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(Event);
123  auto const detProp =
124  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(Event, clockData);
125 
126  //Find the PCA vector
127  TVector3 ShowerCentre;
128  recob::PCAxis PCA = CalculateShowerPCA(clockData, detProp, spacePoints_pfp, fmh, ShowerCentre);
129  TVector3 PCADirection = GetPCAxisVector(PCA);
130 
131  //Save the shower the center for downstream tools
132  TVector3 ShowerCentreErr = {-999, -999, -999};
133  ShowerEleHolder.SetElement(ShowerCentre, ShowerCentreErr, fShowerCentreOutputLabel);
134  ShowerEleHolder.SetElement(PCA, fShowerPCAOutputLabel);
135 
136  //Check if we are pointing the correct direction or not, First try the start position
137  if (fUseStartPosition) {
138  if (!ShowerEleHolder.CheckElement(fShowerStartPositionInputLabel)) {
139  if (fVerbose)
140  mf::LogError("ShowerPCADirection")
141  << "fUseStartPosition is set but ShowerStartPosition is not set. Bailing" << std::endl;
142  return 1;
143  }
144  //Get the General direction as the vector between the start position and the centre
145  TVector3 StartPositionVec = {-999, -999, -999};
146  ShowerEleHolder.GetElement(fShowerStartPositionInputLabel, StartPositionVec);
147 
148  // Calculate the general direction of the shower
149  TVector3 GeneralDir = (ShowerCentre - StartPositionVec).Unit();
150 
151  //Calculate the dot product between eigenvector and general direction
152  double DotProduct = PCADirection.Dot(GeneralDir);
153 
154  //If the dotproduct is negative the Direction needs Flipping
155  if (DotProduct < 0) {
156  PCADirection[0] = -PCADirection[0];
157  PCADirection[1] = -PCADirection[1];
158  PCADirection[2] = -PCADirection[2];
159  }
160 
161  //To do
162  TVector3 PCADirectionErr = {-999, -999, -999};
163  ShowerEleHolder.SetElement(PCADirection, PCADirectionErr, fShowerDirectionOutputLabel);
164  return 0;
165  }
166 
167  //Otherwise Check against the RMS of the shower. Method adapated from EMShower Thanks Mike.
169  spacePoints_pfp, ShowerCentre, PCADirection, fNSegments);
170 
171  // If the alg fails to calculate the gradient it will return 0. In this case do nothing
172  // If the gradient is negative, flip the direction of the shower
173  if (RMSGradient < -std::numeric_limits<double>::epsilon()) {
174  PCADirection[0] = -PCADirection[0];
175  PCADirection[1] = -PCADirection[1];
176  PCADirection[2] = -PCADirection[2];
177  }
178 
179  //To do
180  TVector3 PCADirectionErr = {-999, -999, -999};
181 
182  ShowerEleHolder.SetElement(PCADirection, PCADirectionErr, fShowerDirectionOutputLabel);
183  return 0;
184  }
185 
186  //Function to calculate the shower direction using a charge weight 3D PCA calculation.
190  const std::vector<art::Ptr<recob::SpacePoint>>& sps,
191  const art::FindManyP<recob::Hit>& fmh,
192  TVector3& ShowerCentre)
193  {
194 
195  float TotalCharge = 0;
196  float sumWeights = 0;
197  float xx = 0;
198  float yy = 0;
199  float zz = 0;
200  float xy = 0;
201  float xz = 0;
202  float yz = 0;
203 
204  //Get the Shower Centre
205  if (fChargeWeighted) {
207  clockData, detProp, sps, fmh, TotalCharge);
208  }
209  else {
211  }
212 
213  //Normalise the spacepoints, charge weight and add to the PCA.
214  for (auto& sp : sps) {
215 
216  TVector3 sp_position = IShowerTool::GetLArPandoraShowerAlg().SpacePointPosition(sp);
217 
218  float wht = 1;
219 
220  //Normalise the spacepoint position.
221  sp_position = sp_position - ShowerCentre;
222 
223  if (fChargeWeighted) {
224 
225  //Get the charge.
226  float Charge = IShowerTool::GetLArPandoraShowerAlg().SpacePointCharge(sp, fmh);
227 
228  //Get the time of the spacepoint
229  float Time = IShowerTool::GetLArPandoraShowerAlg().SpacePointTime(sp, fmh);
230 
231  //Correct for the lifetime at the moment.
232  Charge *= std::exp((sampling_rate(clockData) * Time) / (detProp.ElectronLifetime() * 1e3));
233 
234  //Charge Weight
235  wht *= std::sqrt(Charge / TotalCharge);
236  }
237 
238  xx += sp_position.X() * sp_position.X() * wht;
239  yy += sp_position.Y() * sp_position.Y() * wht;
240  zz += sp_position.Z() * sp_position.Z() * wht;
241  xy += sp_position.X() * sp_position.Y() * wht;
242  xz += sp_position.X() * sp_position.Z() * wht;
243  yz += sp_position.Y() * sp_position.Z() * wht;
244  sumWeights += wht;
245  }
246 
247  // Using Eigen package
248  Eigen::Matrix3f matrix;
249 
250  // Construct covariance matrix
251  matrix << xx, xy, xz, xy, yy, yz, xz, yz, zz;
252 
253  // Normalise from the sum of weights
254  matrix /= sumWeights;
255 
256  // Run the PCA
257  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigenMatrix(matrix);
258 
259  Eigen::Vector3f eigenValuesVector = eigenMatrix.eigenvalues();
260  Eigen::Matrix3f eigenVectorsMatrix = eigenMatrix.eigenvectors();
261 
262  // Put in the required form for a recob::PCAxis
263  const bool svdOk = true; //TODO: Should probably think about this a bit more
264  const int nHits = sps.size();
265  // For some reason eigen sorts the eigenvalues from smallest to largest, reverse it
266  const double eigenValues[3] = {
267  eigenValuesVector(2), eigenValuesVector(1), eigenValuesVector(0)};
268  std::vector<std::vector<double>> eigenVectors = {
269  {eigenVectorsMatrix(0, 2), eigenVectorsMatrix(1, 2), eigenVectorsMatrix(2, 2)},
270  {eigenVectorsMatrix(0, 1), eigenVectorsMatrix(1, 1), eigenVectorsMatrix(2, 1)},
271  {eigenVectorsMatrix(0, 0), eigenVectorsMatrix(1, 0), eigenVectorsMatrix(2, 0)}};
272  const double avePos[3] = {ShowerCentre[0], ShowerCentre[1], ShowerCentre[2]};
273 
274  return recob::PCAxis(svdOk, nHits, eigenValues, eigenVectors, avePos);
275  }
276 
277  TVector3
279  {
280 
281  //Get the Eigenvectors.
282  std::vector<double> EigenVector = PCAxis.getEigenVectors()[0];
283 
284  return TVector3(EigenVector[0], EigenVector[1], EigenVector[2]);
285  }
286 
287  int
288  ShowerPCADirection::AddAssociations(const art::Ptr<recob::PFParticle>& pfpPtr,
289  art::Event& Event,
290  reco::shower::ShowerElementHolder& ShowerEleHolder)
291  {
292 
293  //First check the element has been set
294  if (!ShowerEleHolder.CheckElement(fShowerPCAOutputLabel)) {
295  if (fVerbose) mf::LogError("ShowerPCADirection: Add Assns") << "PCA not set." << std::endl;
296  return 1;
297  }
298 
300 
301  const art::Ptr<recob::PCAxis> pcaPtr =
302  GetProducedElementPtr<recob::PCAxis>(fShowerPCAOutputLabel, ShowerEleHolder, ptrSize - 1);
303  const art::Ptr<recob::Shower> showerPtr =
304  GetProducedElementPtr<recob::Shower>("shower", ShowerEleHolder);
305 
306  AddSingle<art::Assns<recob::Shower, recob::PCAxis>>(showerPtr, pcaPtr, "ShowerPCAxisAssn");
307  AddSingle<art::Assns<recob::PFParticle, recob::PCAxis>>(pfpPtr, pcaPtr, "PFParticlePCAxisAssn");
308 
309  return 0;
310  }
311 }
312 DEFINE_ART_CLASS_TOOL(ShowerRecoTools::ShowerPCADirection)
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
Definition: UtilFunc.cxx:42
const EigenVectors & getEigenVectors() const
Definition: PCAxis.h:66
int AddAssociations(const art::Ptr< recob::PFParticle > &pfpPtr, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder) override
Declaration of signal hit object.
int GetVectorPtrSize(std::string Name)
Definition: IShowerTool.h:168
void SetElement(T &dataproduct, const std::string &Name, bool checktag=false)
static constexpr bool
const art::FindManyP< T1 > & GetFindManyP(const art::ValidHandle< std::vector< T2 > > &handle, const art::Event &evt, const art::InputTag &moduleTag)
double RMSShowerGradient(std::vector< art::Ptr< recob::SpacePoint >> &sps, const TVector3 &ShowerCentre, const TVector3 &Direction, const unsigned int nSegments) const
double SpacePointTime(art::Ptr< recob::SpacePoint > const &sp, art::FindManyP< recob::Hit > const &fmh) const
ShowerPCADirection(const fhicl::ParameterSet &pset)
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
int CalculateElement(const art::Ptr< recob::PFParticle > &pfparticle, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder) override
TVector3 GetPCAxisVector(recob::PCAxis &PCAxis)
bool CheckElement(const std::string &Name) const
int GetElement(const std::string &Name, T &Element) const
const shower::LArPandoraShowerAlg & GetLArPandoraShowerAlg() const
Definition: IShowerTool.h:88
double SpacePointCharge(art::Ptr< recob::SpacePoint > const &sp, art::FindManyP< recob::Hit > const &fmh) const
Contains all timing reference information for the detector.
recob::PCAxis CalculateShowerPCA(const detinfo::DetectorClocksData &clockData, const detinfo::DetectorPropertiesData &detProp, const std::vector< art::Ptr< recob::SpacePoint >> &spacePoints_pfp, const art::FindManyP< recob::Hit > &fmh, TVector3 &ShowerCentre)
TVector3 ShowerCentre(std::vector< art::Ptr< recob::SpacePoint >> const &showersps) const
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
auto const detProp
TVector3 SpacePointPosition(art::Ptr< recob::SpacePoint > const &sp) const