6 #include "art/Framework/Principal/Event.h"
10 #include "TPolyLine3D.h"
11 #include "TPolyMarker3D.h"
18 : fLArPandoraShowerAlg(pset.
get<fhicl::ParameterSet>(
"LArPandoraShowerAlg"))
19 , fHitModuleLabel(pset.
get<art::InputTag>(
"HitModuleLabel"))
20 , fPFParticleLabel(pset.
get<art::InputTag>(
"PFParticleLabel"))
21 , fShowerStartPositionInputLabel(pset.
get<
std::string>(
"ShowerStartPositionInputLabel"))
22 , fShowerDirectionInputLabel(pset.
get<
std::string>(
"ShowerDirectionInputLabel"))
23 , fInitialTrackSpacePointsInputLabel(pset.
get<
std::string>(
"InitialTrackSpacePointsInputLabel"))
26 std::map<int, const simb::MCParticle*>
30 const sim::ParticleList& particles = particleInventory->ParticleList();
32 std::map<int, const simb::MCParticle*> trueParticles;
34 for (sim::ParticleList::const_iterator particleIt = particles.begin();
35 particleIt != particles.end();
37 const simb::MCParticle* particle = particleIt->second;
38 trueParticles[particle->TrackId()] = particle;
44 std::map<int, std::vector<int>>
46 std::map<int, const simb::MCParticle*>& trueParticles)
const
50 std::map<int, std::vector<int>> showerMothers;
53 for (
const auto& particleIt : trueParticles) {
54 const simb::MCParticle* particle = particleIt.second;
55 const simb::MCParticle* mother = particle;
57 if (
std::abs(particle->PdgCode()) != 11 &&
std::abs(particle->PdgCode()) != 22) {
continue; }
61 while (mother->Mother() != 0 && trueParticles.find(mother->Mother()) != trueParticles.end()) {
63 int motherId = mother->Mother();
64 if (
std::abs(trueParticles[motherId]->PdgCode()) != 11 &&
65 std::abs(trueParticles[motherId]->PdgCode()) != 22) {
68 mother = trueParticles[motherId];
70 showerMothers[mother->TrackId()].push_back(particle->TrackId());
78 const simb::MCParticle* trueParticle,
79 art::Event
const& Event,
81 const art::Ptr<recob::PFParticle>& pfparticle)
const
84 std::cout <<
"Making Debug Event Display" << std::endl;
89 int run = Event.run();
90 int subRun = Event.subRun();
91 int event = Event.event();
92 int PFPID = pfparticle->Self();
95 TString canvasName = Form(
"canvas_%i_%i_%i_%i", run, subRun, event, PFPID);
96 TCanvas* canvas =
tfs->make<TCanvas>(canvasName, canvasName);
103 std::vector<art::Ptr<recob::SpacePoint>> showerSpacePoints;
104 std::vector<art::Ptr<recob::SpacePoint>> otherSpacePoints;
106 auto const hitHandle = Event.getValidHandle<std::vector<recob::Hit>>(fHitModuleLabel);
107 std::vector<art::Ptr<recob::Hit>> hits;
108 art::fill_ptr_vector(hits, hitHandle);
111 const art::FindManyP<recob::SpacePoint> fmsph(hitHandle, Event, fPFParticleLabel);
112 if (!fmsph.isValid()) {
113 throw cet::exception(
"LArPandoraShowerCheatingAlg")
114 <<
"Spacepoint and hit association not valid. Stopping.";
117 std::map<art::Ptr<recob::SpacePoint>, art::Ptr<recob::Hit>> spacePointHitMap;
119 for (
auto hit : hits) {
121 std::vector<art::Ptr<recob::SpacePoint>> sps = fmsph.at(
hit.key());
122 if (sps.size() == 1) {
123 art::Ptr<recob::SpacePoint> sp = sps.front();
124 if (trueParticleID == trueParticle->TrackId()) { showerSpacePoints.push_back(sp); }
126 otherSpacePoints.push_back(sp);
131 if (!ShowerEleHolder.
CheckElement(fShowerStartPositionInputLabel)) {
132 mf::LogError(
"LArPandoraShowerCheatingAlg")
133 <<
"Start position not set, returning " << std::endl;
136 if (!ShowerEleHolder.
CheckElement(fShowerDirectionInputLabel)) {
137 mf::LogError(
"LArPandoraShowerCheatingAlg") <<
"Direction not set, returning " << std::endl;
140 if (!ShowerEleHolder.
CheckElement(fInitialTrackSpacePointsInputLabel)) {
141 mf::LogError(
"LArPandoraShowerCheatingAlg")
142 <<
"TrackSpacePoints not set, returning " << std::endl;
147 TVector3 showerStartPosition = {-999, -999, -999};
148 TVector3 showerDirection = {-999, -999, -999};
149 std::vector<art::Ptr<recob::SpacePoint>> trackSpacePoints;
151 ShowerEleHolder.
GetElement(fShowerStartPositionInputLabel, showerStartPosition);
152 ShowerEleHolder.
GetElement(fShowerDirectionInputLabel, showerDirection);
153 ShowerEleHolder.
GetElement(fInitialTrackSpacePointsInputLabel, trackSpacePoints);
156 double startXYZ[3] = {0, 0, 0};
157 auto startPoly = std::make_unique<TPolyMarker3D>(1, startXYZ);
161 double minProj = std::numeric_limits<double>::max();
162 double maxProj = -std::numeric_limits<double>::max();
164 double x_min = std::numeric_limits<double>::max(), x_max = -std::numeric_limits<double>::max();
165 double y_min = std::numeric_limits<double>::max(), y_max = -std::numeric_limits<double>::max();
166 double z_min = std::numeric_limits<double>::max(), z_max = -std::numeric_limits<double>::max();
172 auto showerPoly = std::make_unique<TPolyMarker3D>(showerSpacePoints.size());
173 for (
auto spacePoint : showerSpacePoints) {
174 TVector3 pos = fLArPandoraShowerAlg.SpacePointPosition(spacePoint) - showerStartPosition;
179 x_min = std::min(x, x_min);
180 x_max = std::max(x, x_max);
181 y_min = std::min(y, y_min);
182 y_max = std::max(y, y_max);
183 z_min = std::min(z, z_min);
184 z_max = std::max(z, z_max);
186 showerPoly->SetPoint(point, x, y, z);
191 fLArPandoraShowerAlg.SpacePointProjection(spacePoint, showerStartPosition, showerDirection);
193 maxProj = std::max(proj, maxProj);
194 minProj = std::min(proj, minProj);
198 double xDirPoints[2] = {minProj * showerDirection.X(), maxProj * showerDirection.X()};
199 double yDirPoints[2] = {minProj * showerDirection.Y(), maxProj * showerDirection.Y()};
200 double zDirPoints[2] = {minProj * showerDirection.Z(), maxProj * showerDirection.Z()};
202 auto dirPoly = std::make_unique<TPolyLine3D>(2, xDirPoints, yDirPoints, zDirPoints);
205 auto trackPoly = std::make_unique<TPolyMarker3D>(trackSpacePoints.size());
206 for (
auto spacePoint : trackSpacePoints) {
207 TVector3 pos = fLArPandoraShowerAlg.SpacePointPosition(spacePoint) - showerStartPosition;
211 trackPoly->SetPoint(point, x, y, z);
218 auto otherPoly = std::make_unique<TPolyMarker3D>(otherSpacePoints.size());
223 for (
auto const& sp : otherSpacePoints) {
224 TVector3 pos = fLArPandoraShowerAlg.SpacePointPosition(sp) - showerStartPosition;
228 x_min = std::min(x, x_min);
229 x_max = std::max(x, x_max);
230 y_min = std::min(y, y_min);
231 y_max = std::max(y, y_max);
232 z_min = std::min(z, z_min);
233 z_max = std::max(z, z_max);
234 otherPoly->SetPoint(point, x, y, z);
238 gStyle->SetOptStat(0);
239 TH3F axes(
"axes",
"", 1, x_min, x_max, 1, y_min, y_max, 1, z_min, z_max);
240 axes.SetDirectory(0);
241 axes.GetXaxis()->SetTitle(
"X");
242 axes.GetYaxis()->SetTitle(
"Y");
243 axes.GetZaxis()->SetTitle(
"Z");
246 otherPoly->SetMarkerStyle(20);
247 otherPoly->SetMarkerColor(4);
251 showerPoly->SetMarkerStyle(20);
253 trackPoly->SetMarkerStyle(20);
254 trackPoly->SetMarkerColor(2);
256 startPoly->SetMarkerStyle(21);
257 startPoly->SetMarkerSize(2);
258 startPoly->SetMarkerColor(3);
260 dirPoly->SetLineWidth(1);
261 dirPoly->SetLineColor(6);
271 const art::Ptr<recob::Hit>&
hit)
const
274 double particleEnergy = 0;
275 int likelyTrackID = 0;
276 art::ServiceHandle<cheat::BackTrackerService> bt_serv;
277 std::vector<sim::TrackIDE> trackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
278 for (
unsigned int idIt = 0; idIt < trackIDs.size(); ++idIt) {
279 if (trackIDs.at(idIt).energy > particleEnergy) {
280 particleEnergy = trackIDs.at(idIt).energy;
281 likelyTrackID = trackIDs.at(idIt).trackID;
284 return likelyTrackID;
287 std::pair<int, double>
290 std::map<
int, std::vector<int>>
const& ShowersMothers,
295 art::ServiceHandle<cheat::BackTrackerService> bt_serv;
296 art::ServiceHandle<cheat::ParticleInventoryService> particleInventory;
299 std::map<int, double> trackIDToEDepMap;
300 std::map<int, double> trackIDTo3EDepMap;
301 for (
std::vector<art::Ptr<recob::Hit>>::const_iterator hitIt = hits.begin(); hitIt != hits.end();
303 art::Ptr<recob::Hit>
hit = *hitIt;
308 std::vector<sim::TrackIDE> trackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
309 for (
unsigned int idIt = 0; idIt < trackIDs.size(); ++idIt) {
310 trackIDTo3EDepMap[
std::abs(trackIDs[idIt].trackID)] += trackIDs[idIt].energy;
311 if (PlaneID == planeid) {
312 trackIDToEDepMap[
std::abs(trackIDs[idIt].trackID)] += trackIDs[idIt].energy;
318 std::map<int, double> MotherIDtoEMap;
319 std::map<int, double> MotherIDto3EMap;
320 for (std::map<
int, std::vector<int>>::const_iterator showermother = ShowersMothers.begin();
321 showermother != ShowersMothers.end();
323 for (std::vector<int>::const_iterator daughter = (showermother->second).begin();
324 daughter != (showermother->second).
end();
326 MotherIDtoEMap[showermother->first] += trackIDToEDepMap[*daughter];
327 MotherIDto3EMap[showermother->first] += trackIDTo3EDepMap[*daughter];
332 double maxenergy = -1;
333 int objectTrack = -99999;
334 for (std::map<int, double>::iterator mapIt = MotherIDto3EMap.begin();
335 mapIt != MotherIDto3EMap.end();
337 double energy_three = mapIt->second;
338 double trackid = mapIt->first;
339 if (energy_three > maxenergy) {
340 maxenergy = energy_three;
341 objectTrack = trackid;
346 if (maxenergy == 0) {
return std::make_pair(-99999, -99999); }
348 return std::make_pair(objectTrack, MotherIDtoEMap[objectTrack]);
process_name opflash particleana ie ie ie z
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
process_name opflash particleana ie x
std::map< int, const simb::MCParticle * > GetTrueParticleMap() const
LArPandoraShowerCheatingAlg(const fhicl::ParameterSet &pset)
IDparameter< geo::PlaneID > PlaneID
Member type of validated geo::PlaneID parameter.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
process_name opflash particleana ie ie y
std::map< int, std::vector< int > > GetTrueChain(std::map< int, const simb::MCParticle * > &trueParticles) const
auto end(FixedBins< T, C > const &) noexcept
bool CheckElement(const std::string &Name) const
void CheatDebugEVD(detinfo::DetectorClocksData const &clockData, const simb::MCParticle *trueParticle, art::Event const &Event, reco::shower::ShowerElementHolder &ShowerEleHolder, const art::Ptr< recob::PFParticle > &pfparticle) const
PlaneID_t Plane
Index of the plane within its TPC.
int GetElement(const std::string &Name, T &Element) const
int TrueParticleID(detinfo::DetectorClocksData const &clockData, const art::Ptr< recob::Hit > &hit) const
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
Contains all timing reference information for the detector.
constexpr WireID()=default
Default constructor: an invalid TPC ID.
art::ServiceHandle< art::TFileService > tfs
BEGIN_PROLOG could also be cout
int TrueParticleID(detinfo::DetectorClocksData const &clockData, const art::Ptr< recob::Hit > hit, bool rollup_unsaved_ids=1)