LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ShowerRecoTools::ShowerPCADirection Class Reference
Inheritance diagram for ShowerRecoTools::ShowerPCADirection:
ShowerRecoTools::IShowerTool

Public Member Functions

 ShowerPCADirection (const fhicl::ParameterSet &pset)
 
int CalculateElement (const art::Ptr< recob::PFParticle > &pfparticle, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder) override
 
int RunShowerTool (const art::Ptr< recob::PFParticle > &pfparticle, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder, std::string evd_display_name_append="")
 
void SetPtr (art::ProducesCollector *collector)
 
void InitaliseProducerPtr (reco::shower::ShowerProducedPtrsHolder &uniqueproducerPtrs)
 

Protected Member Functions

const shower::LArPandoraShowerAlgGetLArPandoraShowerAlg () const
 
template<class T >
art::Ptr< T > GetProducedElementPtr (std::string Name, reco::shower::ShowerElementHolder &ShowerEleHolder, int iter=-1)
 
template<class T >
void InitialiseProduct (std::string Name, std::string InstanceName="")
 
template<class T , class A , class B >
void AddSingle (A &a, B &b, std::string Name)
 
int GetVectorPtrSize (std::string Name)
 
void PrintPtrs ()
 
void PrintPtr (std::string Name)
 

Private Member Functions

void InitialiseProducers () override
 
int AddAssociations (const art::Ptr< recob::PFParticle > &pfpPtr, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder) override
 
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, geo::Point_t &ShowerCentre)
 
geo::Vector_t GetPCAxisVector (recob::PCAxis &PCAxis)
 

Private Attributes

art::InputTag fPFParticleLabel
 
int fVerbose
 
unsigned int fNSegments
 
bool fUseStartPosition
 
bool fChargeWeighted
 
std::string fShowerStartPositionInputLabel
 
std::string fShowerDirectionOutputLabel
 
std::string fShowerCentreOutputLabel
 
std::string fShowerPCAOutputLabel
 

Detailed Description

Definition at line 28 of file ShowerPCADirection_tool.cc.

Constructor & Destructor Documentation

ShowerRecoTools::ShowerPCADirection::ShowerPCADirection ( const fhicl::ParameterSet pset)

Definition at line 71 of file ShowerPCADirection_tool.cc.

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  {}
T get(std::string const &key) const
Definition: ParameterSet.h:314
IShowerTool(const fhicl::ParameterSet &pset)
Definition: IShowerTool.h:33

Member Function Documentation

int ShowerRecoTools::ShowerPCADirection::AddAssociations ( const art::Ptr< recob::PFParticle > &  pfpPtr,
art::Event Event,
reco::shower::ShowerElementHolder ShowerEleHolder 
)
overrideprivatevirtual

Reimplemented from ShowerRecoTools::IShowerTool.

Definition at line 272 of file ShowerPCADirection_tool.cc.

References reco::shower::ShowerElementHolder::CheckElement(), DEFINE_ART_CLASS_TOOL, fShowerPCAOutputLabel, fVerbose, and ShowerRecoTools::IShowerTool::GetVectorPtrSize().

275  {
276 
277  //First check the element has been set
278  if (!ShowerEleHolder.CheckElement(fShowerPCAOutputLabel)) {
279  if (fVerbose) mf::LogError("ShowerPCADirection: Add Assns") << "PCA not set." << std::endl;
280  return 1;
281  }
282 
284 
285  const art::Ptr<recob::PCAxis> pcaPtr =
286  GetProducedElementPtr<recob::PCAxis>(fShowerPCAOutputLabel, ShowerEleHolder, ptrSize - 1);
287  const art::Ptr<recob::Shower> showerPtr =
288  GetProducedElementPtr<recob::Shower>("shower", ShowerEleHolder);
289 
290  AddSingle<art::Assns<recob::Shower, recob::PCAxis>>(showerPtr, pcaPtr, "ShowerPCAxisAssn");
291  AddSingle<art::Assns<recob::PFParticle, recob::PCAxis>>(pfpPtr, pcaPtr, "PFParticlePCAxisAssn");
292 
293  return 0;
294  }
int GetVectorPtrSize(std::string Name)
Definition: IShowerTool.h:158
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
bool CheckElement(const std::string &Name) const
Definition: fwd.h:26
template<class T , class A , class B >
void ShowerRecoTools::IShowerTool::AddSingle ( A &  a,
B &  b,
std::string  Name 
)
inlineprotectedinherited

Definition at line 152 of file IShowerTool.h.

References reco::shower::ShowerProducedPtrsHolder::AddSingle().

153  {
154  UniquePtrs->AddSingle<T>(a, b, Name);
155  }
reco::shower::ShowerProducedPtrsHolder * UniquePtrs
Definition: IShowerTool.h:85
void AddSingle(A &a, B &b, const std::string &Name)
int ShowerRecoTools::ShowerPCADirection::CalculateElement ( const art::Ptr< recob::PFParticle > &  pfparticle,
art::Event Event,
reco::shower::ShowerElementHolder ShowerEleHolder 
)
overridevirtual

Implements ShowerRecoTools::IShowerTool.

Definition at line 91 of file ShowerPCADirection_tool.cc.

References CalculateShowerPCA(), reco::shower::ShowerElementHolder::CheckElement(), fNSegments, fPFParticleLabel, fShowerCentreOutputLabel, fShowerDirectionOutputLabel, fShowerPCAOutputLabel, fShowerStartPositionInputLabel, fUseStartPosition, fVerbose, reco::shower::ShowerElementHolder::GetElement(), reco::shower::ShowerElementHolder::GetFindManyP(), ShowerRecoTools::IShowerTool::GetLArPandoraShowerAlg(), GetPCAxisVector(), art::ProductRetriever::getValidHandle(), art::Ptr< T >::key(), shower::LArPandoraShowerAlg::RMSShowerGradient(), and reco::shower::ShowerElementHolder::SetElement().

94  {
95 
96  // Get the assocated pfParicle vertex PFParticles
97  auto const pfpHandle = Event.getValidHandle<std::vector<recob::PFParticle>>(fPFParticleLabel);
98 
100  ShowerEleHolder.GetFindManyP<recob::SpacePoint>(pfpHandle, Event, fPFParticleLabel);
101 
102  //Get the spacepoints handle and the hit assoication
103  auto const spHandle = Event.getValidHandle<std::vector<recob::SpacePoint>>(fPFParticleLabel);
104 
105  const art::FindManyP<recob::Hit>& fmh =
106  ShowerEleHolder.GetFindManyP<recob::Hit>(spHandle, Event, fPFParticleLabel);
107 
108  //SpacePoints
109  std::vector<art::Ptr<recob::SpacePoint>> spacePoints_pfp = fmspp.at(pfparticle.key());
110 
111  //We cannot progress with no spacepoints.
112  if (spacePoints_pfp.size() < 3) {
113  mf::LogWarning("ShowerPCADirection")
114  << spacePoints_pfp.size() << " spacepoints in shower, not calculating direction"
115  << std::endl;
116  return 1;
117  }
118 
119  auto const clockData =
121  auto const detProp =
123 
124  //Find the PCA vector
125  geo::Point_t ShowerCentre;
126  recob::PCAxis PCA = CalculateShowerPCA(clockData, detProp, spacePoints_pfp, fmh, ShowerCentre);
127  auto PCADirection = GetPCAxisVector(PCA);
128 
129  //Save the shower the center for downstream tools
130  geo::Point_t ShowerCentreErr = {-999, -999, -999};
131  ShowerEleHolder.SetElement(ShowerCentre, ShowerCentreErr, fShowerCentreOutputLabel);
132  ShowerEleHolder.SetElement(PCA, fShowerPCAOutputLabel);
133 
134  //Check if we are pointing the correct direction or not, First try the start position
135  if (fUseStartPosition) {
136  if (!ShowerEleHolder.CheckElement(fShowerStartPositionInputLabel)) {
137  if (fVerbose)
138  mf::LogError("ShowerPCADirection")
139  << "fUseStartPosition is set but ShowerStartPosition is not set. Bailing" << std::endl;
140  return 1;
141  }
142  //Get the General direction as the vector between the start position and the centre
143  geo::Point_t StartPositionVec = {-999, -999, -999};
144  ShowerEleHolder.GetElement(fShowerStartPositionInputLabel, StartPositionVec);
145 
146  // Calculate the general direction of the shower
147  auto const GeneralDir = (ShowerCentre - StartPositionVec).Unit();
148 
149  //Calculate the dot product between eigenvector and general direction
150  double DotProduct = PCADirection.Dot(GeneralDir);
151 
152  //If the dotproduct is negative the Direction needs Flipping
153  if (DotProduct < 0) { PCADirection *= -1; }
154 
155  //To do
156  geo::Vector_t PCADirectionErr = {-999, -999, -999};
157  ShowerEleHolder.SetElement(PCADirection, PCADirectionErr, fShowerDirectionOutputLabel);
158  return 0;
159  }
160 
161  //Otherwise Check against the RMS of the shower. Method adapated from EMShower Thanks Mike.
163  spacePoints_pfp, ShowerCentre, PCADirection, fNSegments);
164 
165  // If the alg fails to calculate the gradient it will return 0. In this case do nothing
166  // If the gradient is negative, flip the direction of the shower
167  if (RMSGradient < -std::numeric_limits<double>::epsilon()) { PCADirection *= -1; }
168 
169  //To do
170  geo::Vector_t PCADirectionErr = {-999, -999, -999};
171 
172  ShowerEleHolder.SetElement(PCADirection, PCADirectionErr, fShowerDirectionOutputLabel);
173  return 0;
174  }
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:160
geo::Vector_t GetPCAxisVector(recob::PCAxis &PCAxis)
void SetElement(T &dataproduct, const std::string &Name, bool checktag=false)
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, geo::Point_t &ShowerCentre)
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
key_type key() const noexcept
Definition: Ptr.h:166
bool CheckElement(const std::string &Name) const
int GetElement(const std::string &Name, T &Element) const
const shower::LArPandoraShowerAlg & GetLArPandoraShowerAlg() const
Definition: IShowerTool.h:82
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
double RMSShowerGradient(std::vector< art::Ptr< recob::SpacePoint >> &sps, const geo::Point_t &ShowerCentre, const geo::Vector_t &Direction, const unsigned int nSegments) const
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:46
const art::FindManyP< T1 > & GetFindManyP(const art::ValidHandle< std::vector< T2 >> &handle, const art::Event &evt, const art::InputTag &moduleTag)
recob::PCAxis ShowerRecoTools::ShowerPCADirection::CalculateShowerPCA ( const detinfo::DetectorClocksData clockData,
const detinfo::DetectorPropertiesData detProp,
const std::vector< art::Ptr< recob::SpacePoint >> &  spacePoints_pfp,
const art::FindManyP< recob::Hit > &  fmh,
geo::Point_t ShowerCentre 
)
private

Definition at line 177 of file ShowerPCADirection_tool.cc.

References detinfo::DetectorPropertiesData::ElectronLifetime(), fChargeWeighted, ShowerRecoTools::IShowerTool::GetLArPandoraShowerAlg(), detinfo::sampling_rate(), shower::LArPandoraShowerAlg::ShowerCentre(), shower::LArPandoraShowerAlg::SpacePointCharge(), shower::LArPandoraShowerAlg::SpacePointTime(), xx, and zz.

Referenced by CalculateElement().

183  {
184 
185  float TotalCharge = 0;
186  float sumWeights = 0;
187  float xx = 0;
188  float yy = 0;
189  float zz = 0;
190  float xy = 0;
191  float xz = 0;
192  float yz = 0;
193 
194  //Get the Shower Centre
195  if (fChargeWeighted) {
197  clockData, detProp, sps, fmh, TotalCharge);
198  }
199  else {
201  }
202 
203  //Normalise the spacepoints, charge weight and add to the PCA.
204  for (auto& sp : sps) {
205 
206  float wht = 1;
207 
208  //Normalise the spacepoint position.
209  auto const sp_position = sp->position() - ShowerCentre;
210 
211  if (fChargeWeighted) {
212 
213  //Get the charge.
214  float Charge = IShowerTool::GetLArPandoraShowerAlg().SpacePointCharge(sp, fmh);
215 
216  //Get the time of the spacepoint
217  float Time = IShowerTool::GetLArPandoraShowerAlg().SpacePointTime(sp, fmh);
218 
219  //Correct for the lifetime at the moment.
220  Charge *= std::exp((sampling_rate(clockData) * Time) / (detProp.ElectronLifetime() * 1e3));
221 
222  //Charge Weight
223  wht *= std::sqrt(Charge / TotalCharge);
224  }
225 
226  xx += sp_position.X() * sp_position.X() * wht;
227  yy += sp_position.Y() * sp_position.Y() * wht;
228  zz += sp_position.Z() * sp_position.Z() * wht;
229  xy += sp_position.X() * sp_position.Y() * wht;
230  xz += sp_position.X() * sp_position.Z() * wht;
231  yz += sp_position.Y() * sp_position.Z() * wht;
232  sumWeights += wht;
233  }
234 
235  // Using Eigen package
236  Eigen::Matrix3f matrix;
237 
238  // Construct covariance matrix
239  matrix << xx, xy, xz, xy, yy, yz, xz, yz, zz;
240 
241  // Normalise from the sum of weights
242  matrix /= sumWeights;
243 
244  // Run the PCA
245  Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eigenMatrix(matrix);
246 
247  Eigen::Vector3f eigenValuesVector = eigenMatrix.eigenvalues();
248  Eigen::Matrix3f eigenVectorsMatrix = eigenMatrix.eigenvectors();
249 
250  // Put in the required form for a recob::PCAxis
251  const bool svdOk = true; //TODO: Should probably think about this a bit more
252  const int nHits = sps.size();
253  // For some reason eigen sorts the eigenvalues from smallest to largest, reverse it
254  const double eigenValues[3] = {
255  eigenValuesVector(2), eigenValuesVector(1), eigenValuesVector(0)};
256  std::vector<std::vector<double>> eigenVectors = {
257  {eigenVectorsMatrix(0, 2), eigenVectorsMatrix(1, 2), eigenVectorsMatrix(2, 2)},
258  {eigenVectorsMatrix(0, 1), eigenVectorsMatrix(1, 1), eigenVectorsMatrix(2, 1)},
259  {eigenVectorsMatrix(0, 0), eigenVectorsMatrix(1, 0), eigenVectorsMatrix(2, 0)}};
260  const double avePos[3] = {ShowerCentre.X(), ShowerCentre.Y(), ShowerCentre.Z()};
261 
262  return recob::PCAxis(svdOk, nHits, eigenValues, eigenVectors, avePos);
263  }
Double_t xx
Definition: macro.C:12
Double_t zz
Definition: plot.C:276
double SpacePointTime(art::Ptr< recob::SpacePoint > const &sp, art::FindManyP< recob::Hit > const &fmh) const
const shower::LArPandoraShowerAlg & GetLArPandoraShowerAlg() const
Definition: IShowerTool.h:82
double SpacePointCharge(art::Ptr< recob::SpacePoint > const &sp, art::FindManyP< recob::Hit > const &fmh) const
geo::Point_t ShowerCentre(std::vector< art::Ptr< recob::SpacePoint >> const &showersps) const
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
const shower::LArPandoraShowerAlg& ShowerRecoTools::IShowerTool::GetLArPandoraShowerAlg ( ) const
inlineprotectedinherited
geo::Vector_t ShowerRecoTools::ShowerPCADirection::GetPCAxisVector ( recob::PCAxis PCAxis)
private

Definition at line 265 of file ShowerPCADirection_tool.cc.

References recob::PCAxis::getEigenVectors().

Referenced by CalculateElement().

266  {
267  //Get the Eigenvectors.
268  std::vector<double> EigenVector = PCAxis.getEigenVectors()[0];
269  return {EigenVector[0], EigenVector[1], EigenVector[2]};
270  }
const EigenVectors & getEigenVectors() const
Definition: PCAxis.h:74
template<class T >
art::Ptr<T> ShowerRecoTools::IShowerTool::GetProducedElementPtr ( std::string  Name,
reco::shower::ShowerElementHolder ShowerEleHolder,
int  iter = -1 
)
inlineprotectedinherited

Definition at line 102 of file IShowerTool.h.

References reco::shower::ShowerElementHolder::CheckElement(), reco::shower::ShowerProducedPtrsHolder::CheckUniqueProduerPtr(), reco::shower::ShowerProducedPtrsHolder::GetArtPtr(), and reco::shower::ShowerElementHolder::GetShowerNumber().

105  {
106 
107  //Check the element has been set
108  bool check_element = ShowerEleHolder.CheckElement(Name);
109  if (!check_element) {
110  throw cet::exception("IShowerTool") << "tried to get a element that does not exist. Failed "
111  "at making the art ptr for Element: "
112  << Name << std::endl;
113  }
114 
115  //Check the unique ptr has been set.
116  bool check_ptr = UniquePtrs->CheckUniqueProduerPtr(Name);
117  if (!check_ptr) {
118  throw cet::exception("IShowerTool")
119  << "tried to get a ptr that does not exist. Failed at making the art ptr for Element"
120  << Name;
121  }
122 
123  //Check if the user has defined an index if not just use the current shower index/
124  int index;
125  if (iter != -1) { index = iter; }
126  else {
127  index = ShowerEleHolder.GetShowerNumber();
128  }
129 
130  //Make the ptr
131  return UniquePtrs->GetArtPtr<T>(Name, index);
132  }
reco::shower::ShowerProducedPtrsHolder * UniquePtrs
Definition: IShowerTool.h:85
bool CheckUniqueProduerPtr(const std::string &Name) const
bool CheckElement(const std::string &Name) const
art::Ptr< T > GetArtPtr(const std::string &Name, const int &iter) const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
int ShowerRecoTools::IShowerTool::GetVectorPtrSize ( std::string  Name)
inlineprotectedinherited
void ShowerRecoTools::IShowerTool::InitaliseProducerPtr ( reco::shower::ShowerProducedPtrsHolder uniqueproducerPtrs)
inlineinherited

Definition at line 68 of file IShowerTool.h.

69  {
70  UniquePtrs = &uniqueproducerPtrs;
71  }
reco::shower::ShowerProducedPtrsHolder * UniquePtrs
Definition: IShowerTool.h:85
void ShowerRecoTools::ShowerPCADirection::InitialiseProducers ( )
overrideprivatevirtual

Reimplemented from ShowerRecoTools::IShowerTool.

Definition at line 84 of file ShowerPCADirection_tool.cc.

References fShowerPCAOutputLabel.

85  {
86  InitialiseProduct<std::vector<recob::PCAxis>>(fShowerPCAOutputLabel);
87  InitialiseProduct<art::Assns<recob::Shower, recob::PCAxis>>("ShowerPCAxisAssn");
88  InitialiseProduct<art::Assns<recob::PFParticle, recob::PCAxis>>("PFParticlePCAxisAssn");
89  }
template<class T >
void ShowerRecoTools::IShowerTool::InitialiseProduct ( std::string  Name,
std::string  InstanceName = "" 
)
inlineprotectedinherited

Definition at line 137 of file IShowerTool.h.

References art::ProducesCollector::produces(), and reco::shower::ShowerProducedPtrsHolder::SetShowerUniqueProduerPtr().

138  {
139 
140  if (collectorPtr == nullptr) {
141  mf::LogWarning("IShowerTool") << "The art::ProducesCollector ptr has not been set";
142  return;
143  }
144 
145  collectorPtr->produces<T>(InstanceName);
146  UniquePtrs->SetShowerUniqueProduerPtr(type<T>(), Name, InstanceName);
147  }
reco::shower::ShowerProducedPtrsHolder * UniquePtrs
Definition: IShowerTool.h:85
art::ProducesCollector * collectorPtr
Definition: IShowerTool.h:97
void produces(std::string const &instanceName={}, Persistable const persistable=Persistable::Yes)
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
int SetShowerUniqueProduerPtr(type< T >, const std::string &Name, const std::string &Instance="")
void ShowerRecoTools::IShowerTool::PrintPtr ( std::string  Name)
inlineprotectedinherited

Definition at line 162 of file IShowerTool.h.

References reco::shower::ShowerProducedPtrsHolder::PrintPtr().

Referenced by ShowerRecoTools::ShowerExampleTool::CalculateElement().

162 { UniquePtrs->PrintPtr(Name); }
reco::shower::ShowerProducedPtrsHolder * UniquePtrs
Definition: IShowerTool.h:85
void PrintPtr(const std::string &Name) const
void ShowerRecoTools::IShowerTool::PrintPtrs ( )
inlineprotectedinherited
int ShowerRecoTools::IShowerTool::RunShowerTool ( const art::Ptr< recob::PFParticle > &  pfparticle,
art::Event Event,
reco::shower::ShowerElementHolder ShowerEleHolder,
std::string  evd_display_name_append = "" 
)
inlineinherited

Definition at line 46 of file IShowerTool.h.

50  {
51 
52  int calculation_status = CalculateElement(pfparticle, Event, ShowerEleHolder);
53  if (calculation_status != 0) return calculation_status;
54  if (fRunEventDisplay) {
56  pfparticle, Event, ShowerEleHolder, evd_display_name_append);
57  }
58  return calculation_status;
59  }
void DebugEVD(art::Ptr< recob::PFParticle > const &pfparticle, art::Event const &Event, const reco::shower::ShowerElementHolder &ShowerEleHolder, std::string const &evd_disp_name_append="") const
virtual int CalculateElement(const art::Ptr< recob::PFParticle > &pfparticle, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder)=0
const shower::LArPandoraShowerAlg & GetLArPandoraShowerAlg() const
Definition: IShowerTool.h:82
void ShowerRecoTools::IShowerTool::SetPtr ( art::ProducesCollector collector)
inlineinherited

Definition at line 65 of file IShowerTool.h.

65 { collectorPtr = collector; }
art::ProducesCollector * collectorPtr
Definition: IShowerTool.h:97

Member Data Documentation

bool ShowerRecoTools::ShowerPCADirection::fChargeWeighted
private

Definition at line 63 of file ShowerPCADirection_tool.cc.

Referenced by CalculateShowerPCA().

unsigned int ShowerRecoTools::ShowerPCADirection::fNSegments
private

Definition at line 60 of file ShowerPCADirection_tool.cc.

Referenced by CalculateElement().

art::InputTag ShowerRecoTools::ShowerPCADirection::fPFParticleLabel
private

Definition at line 57 of file ShowerPCADirection_tool.cc.

Referenced by CalculateElement().

std::string ShowerRecoTools::ShowerPCADirection::fShowerCentreOutputLabel
private

Definition at line 67 of file ShowerPCADirection_tool.cc.

Referenced by CalculateElement().

std::string ShowerRecoTools::ShowerPCADirection::fShowerDirectionOutputLabel
private

Definition at line 66 of file ShowerPCADirection_tool.cc.

Referenced by CalculateElement().

std::string ShowerRecoTools::ShowerPCADirection::fShowerPCAOutputLabel
private
std::string ShowerRecoTools::ShowerPCADirection::fShowerStartPositionInputLabel
private

Definition at line 65 of file ShowerPCADirection_tool.cc.

Referenced by CalculateElement().

bool ShowerRecoTools::ShowerPCADirection::fUseStartPosition
private

Definition at line 61 of file ShowerPCADirection_tool.cc.

Referenced by CalculateElement().

int ShowerRecoTools::ShowerPCADirection::fVerbose
private

Definition at line 58 of file ShowerPCADirection_tool.cc.

Referenced by AddAssociations(), and CalculateElement().


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