LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ShowerDirectionCheater_tool.cc
Go to the documentation of this file.
1 //############################################################################
2 //### Name: ShowerDirectionCheater ###
3 //### Author: Ed Tyley ###
4 //### Date: 16.07.19 ###
5 //### Description: Cheating tool using truth for shower direction ###
6 //############################################################################
7 
8 //Framework Includes
10 
11 //LArSoft Includes
18 
19 //ROOT Includes
20 #include "TTree.h"
21 
22 namespace ShowerRecoTools {
23 
25 
26  public:
28 
29  //Generic Direction Finder
30  int CalculateElement(const art::Ptr<recob::PFParticle>& pfparticle,
31  art::Event& Event,
32  reco::shower::ShowerElementHolder& ShowerEleHolder) override;
33 
34  private:
35  //Algorithm functions
37 
38  //Services
40 
41  //fcl
43  const unsigned int
44  fNSegments; //Number of segement to split the shower into the perforam the RMSFlip.
45  const bool fRMSFlip; //Flip the direction by considering the rms.
46  const bool
47  fVertexFlip; //Flip the direction by considering the vertex position relative to the center position.
48 
49  //TTree Branch variables
50  TTree* Tree;
52  float rmsGradient;
53 
54  const std::string fShowerStartPositionInputLabel;
55  const std::string fTrueParticleInputLabel;
56  const std::string fShowerDirectionOutputLabel;
57  };
58 
60  : IShowerTool(pset.get<fhicl::ParameterSet>("BaseTools"))
61  , fLArPandoraShowerCheatingAlg(pset.get<fhicl::ParameterSet>("LArPandoraShowerCheatingAlg"))
62  , fPFParticleLabel(pset.get<art::InputTag>("PFParticleLabel"))
63  , fNSegments(pset.get<unsigned int>("NSegments"))
64  , fRMSFlip(pset.get<bool>("RMSFlip"))
65  , fVertexFlip(pset.get<bool>("VertexFlip"))
66  , fShowerStartPositionInputLabel(pset.get<std::string>("ShowerStartPositionInputLabel"))
67  , fTrueParticleInputLabel(pset.get<std::string>("TrueParticleInputLabel"))
68  , fShowerDirectionOutputLabel(pset.get<std::string>("ShowerDirectionOutputLabel"))
69  {
70  if (fVertexFlip || fRMSFlip) {
71  Tree = tfs->make<TTree>("DebugTreeDirCheater", "DebugTree from shower direction cheater");
72  if (fVertexFlip) Tree->Branch("vertexDotProduct", &vertexDotProduct);
73  if (fRMSFlip) Tree->Branch("rmsGradient", &rmsGradient);
74  }
75  }
76 
78  art::Event& Event,
79  reco::shower::ShowerElementHolder& ShowerEleHolder)
80  {
81 
82  const simb::MCParticle* trueParticle;
83 
84  //Get the hits from the shower:
85  auto const pfpHandle = Event.getValidHandle<std::vector<recob::PFParticle>>(fPFParticleLabel);
86 
87  auto const clockData =
89  auto const detProp =
91 
92  if (ShowerEleHolder.CheckElement(fTrueParticleInputLabel)) {
93  ShowerEleHolder.GetElement(fTrueParticleInputLabel, trueParticle);
94  }
95  else {
96 
97  //Could store these in the shower element holder and just calculate once?
98  std::map<int, const simb::MCParticle*> trueParticles =
100  std::map<int, std::vector<int>> showersMothers =
102 
103  //Get the clusters
104  auto const clusHandle = Event.getValidHandle<std::vector<recob::Cluster>>(fPFParticleLabel);
105  art::FindManyP<recob::Cluster> fmc(pfpHandle, Event, fPFParticleLabel);
106  std::vector<art::Ptr<recob::Cluster>> clusters = fmc.at(pfparticle.key());
107 
108  //Get the hit association
109  art::FindManyP<recob::Hit> fmhc(clusHandle, Event, fPFParticleLabel);
110 
111  std::vector<art::Ptr<recob::Hit>> showerHits;
112  for (auto const& cluster : clusters) {
113  //Get the hits
114  std::vector<art::Ptr<recob::Hit>> hits = fmhc.at(cluster.key());
115  showerHits.insert(showerHits.end(), hits.begin(), hits.end());
116  }
117 
118  //Get the true particle from the shower
119  std::pair<int, double> ShowerTrackInfo =
121  clockData, showersMothers, showerHits, 2);
122 
123  if (ShowerTrackInfo.first == -99999) {
124  mf::LogError("ShowerDirectionCheater") << "True shower not found, returning";
125  return 1;
126  }
127  trueParticle = trueParticles[ShowerTrackInfo.first];
128  ShowerEleHolder.SetElement(trueParticle, fTrueParticleInputLabel);
129  }
130 
131  if (!trueParticle) {
132  mf::LogError("ShowerDirectionCheater") << "True shower not found, returning";
133  return 1;
134  }
135 
136  auto trueDir = geo::Vector_t{trueParticle->Px(), trueParticle->Py(), trueParticle->Pz()}.Unit();
137 
138  geo::Vector_t trueDirErr = {-999, -999, -999};
139  ShowerEleHolder.SetElement(trueDir, trueDirErr, fShowerDirectionOutputLabel);
140 
141  if (fRMSFlip || fVertexFlip) {
142 
143  // Reset the tree values to defaults
144  rmsGradient = std::numeric_limits<float>::lowest();
145  vertexDotProduct = std::numeric_limits<float>::lowest();
146 
147  //Get the SpacePoints and hits
148  art::FindManyP<recob::SpacePoint> fmspp(pfpHandle, Event, fPFParticleLabel);
149 
150  if (!fmspp.isValid()) {
151  throw cet::exception("ShowerDirectionCheater")
152  << "Trying to get the spacepoint and failed. Something is not configured correctly. "
153  "Stopping ";
154  }
155 
156  auto const spHandle = Event.getValidHandle<std::vector<recob::SpacePoint>>(fPFParticleLabel);
157  art::FindManyP<recob::Hit> fmh(spHandle, Event, fPFParticleLabel);
158  if (!fmh.isValid()) {
159  throw cet::exception("ShowerDirectionCheater")
160  << "Spacepoint and hit association not valid. Stopping.";
161  }
162  std::vector<art::Ptr<recob::SpacePoint>> spacePoints = fmspp.at(pfparticle.key());
163 
164  if (spacePoints.size() < 3) {
165  mf::LogWarning("ShowerDirectionCheater")
166  << spacePoints.size() << " spacepoints in shower, not calculating direction" << std::endl;
167  return 1;
168  }
169 
170  //Get Shower Centre
171  float TotalCharge;
172 
173  auto const ShowerCentre = IShowerTool::GetLArPandoraShowerAlg().ShowerCentre(
174  clockData, detProp, spacePoints, fmh, TotalCharge);
175 
176  //Check if we are pointing the correct direction or not, First try the start position
177  if (ShowerEleHolder.CheckElement(fShowerStartPositionInputLabel) && fVertexFlip) {
178 
179  //Get the General direction as the vector between the start position and the centre
180  geo::Point_t StartPositionVec = {-999, -999, -999};
181  ShowerEleHolder.GetElement(fShowerStartPositionInputLabel, StartPositionVec);
182 
183  auto const GeneralDir = (ShowerCentre - StartPositionVec).Unit();
184 
185  //Dot product
186  vertexDotProduct = trueDir.Dot(GeneralDir);
187  }
188 
189  if (fRMSFlip) {
190  //Otherwise Check against the RMS of the shower. Method adapated from EMShower Thanks Mike.
192  spacePoints, ShowerCentre, trueDir, fNSegments);
193  }
194  Tree->Fill();
195  }
196 
197  return 0;
198  }
199 }
200 
art::ServiceHandle< art::TFileService > tfs
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:42
double Py(const int i=0) const
Definition: MCParticle.h:232
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
std::map< int, const simb::MCParticle * > GetTrueParticleMap() const
double Px(const int i=0) const
Definition: MCParticle.h:231
void SetElement(T &dataproduct, const std::string &Name, bool checktag=false)
STL namespace.
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
Cluster finding and building.
int CalculateElement(const art::Ptr< recob::PFParticle > &pfparticle, art::Event &Event, reco::shower::ShowerElementHolder &ShowerEleHolder) override
void hits()
Definition: readHits.C:15
parameter set interface
key_type key() const noexcept
Definition: Ptr.h:166
std::map< int, std::vector< int > > GetTrueChain(std::map< int, const simb::MCParticle * > &trueParticles) const
bool CheckElement(const std::string &Name) const
int GetElement(const std::string &Name, T &Element) const
Declaration of cluster object.
shower::LArPandoraShowerCheatingAlg fLArPandoraShowerCheatingAlg
const shower::LArPandoraShowerAlg & GetLArPandoraShowerAlg() const
Definition: IShowerTool.h:82
std::pair< int, double > TrueParticleIDFromTrueChain(detinfo::DetectorClocksData const &clockData, std::map< int, std::vector< int >> const &ShowersMothers, std::vector< art::Ptr< recob::Hit >> const &hits, int planeid) const
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
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
geo::Point_t ShowerCentre(std::vector< art::Ptr< recob::SpacePoint >> const &showersps) const
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
double Pz(const int i=0) const
Definition: MCParticle.h:233
Definition: MVAAlg.h:12
double RMSShowerGradient(std::vector< art::Ptr< recob::SpacePoint >> &sps, const geo::Point_t &ShowerCentre, const geo::Vector_t &Direction, const unsigned int nSegments) const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
ShowerDirectionCheater(const fhicl::ParameterSet &pset)