Pandora
Pandora source code navigator
Loading...
Searching...
No Matches
PFParticleAnalysis_module.cc
Go to the documentation of this file.
1
7#include "art/Framework/Core/EDAnalyzer.h"
8#include "art/Framework/Core/ModuleMacros.h"
9
10#include "TTree.h"
11#include "TVector3.h"
12
14
15#include <string>
16
17//------------------------------------------------------------------------------------------------------------------------------------------
18
19namespace lar_pandora {
20
24 class PFParticleAnalysis : public art::EDAnalyzer {
25 public:
31 PFParticleAnalysis(fhicl::ParameterSet const& pset);
32
36 virtual ~PFParticleAnalysis();
37
38 void beginJob();
39 void endJob();
40 void analyze(const art::Event& evt);
41 void reconfigure(fhicl::ParameterSet const& pset);
42
43 private:
44 TTree* m_pRecoTree;
45
46 int m_run;
47 int m_event;
48 int m_index;
49
50 int m_self;
59 int m_track;
63
66 int m_hits;
70
71 double m_pfovtxx;
72 double m_pfovtxy;
73 double m_pfovtxz;
74
75 double m_trkvtxx;
76 double m_trkvtxy;
77 double m_trkvtxz;
78 double m_trkvtxdirx;
79 double m_trkvtxdiry;
80 double m_trkvtxdirz;
81 double m_trkendx;
82 double m_trkendy;
83 double m_trkendz;
84 double m_trkenddirx;
85 double m_trkenddiry;
86 double m_trkenddirz;
87 double m_trklength;
89
90 double m_shwvtxx;
91 double m_shwvtxy;
92 double m_shwvtxz;
93 double m_shwvtxdirx;
94 double m_shwvtxdiry;
95 double m_shwvtxdirz;
96 double m_shwlength;
99
100 double m_t0;
101
102 std::string m_particleLabel;
103 std::string m_trackLabel;
104 std::string m_showerLabel;
106 };
107
108 DEFINE_ART_MODULE(PFParticleAnalysis)
109
110} // namespace lar_pandora
111
112//------------------------------------------------------------------------------------------------------------------------------------------
113// implementation follows
114
115#include "art/Framework/Principal/Event.h"
116#include "art/Framework/Principal/Handle.h"
117#include "art/Framework/Services/Registry/ServiceHandle.h"
118#include "art_root_io/TFileDirectory.h"
119#include "art_root_io/TFileService.h"
120#include "fhiclcpp/ParameterSet.h"
121#include "messagefacility/MessageLogger/MessageLogger.h"
122
123#include "lardataobj/RecoBase/PFParticle.h"
124#include "lardataobj/RecoBase/Seed.h"
125#include "lardataobj/RecoBase/Shower.h"
126#include "lardataobj/RecoBase/Track.h"
127#include "lardataobj/RecoBase/Vertex.h"
128
129#include "lardataobj/AnalysisBase/T0.h"
130
131#include <iostream>
132
133namespace lar_pandora {
134
135 PFParticleAnalysis::PFParticleAnalysis(fhicl::ParameterSet const& pset) : art::EDAnalyzer(pset)
136 {
137 this->reconfigure(pset);
138 }
139
140 //------------------------------------------------------------------------------------------------------------------------------------------
141
143
144 //------------------------------------------------------------------------------------------------------------------------------------------
145
146 void PFParticleAnalysis::reconfigure(fhicl::ParameterSet const& pset)
147 {
148 m_particleLabel = pset.get<std::string>("PFParticleModule", "pandora");
149 m_trackLabel = pset.get<std::string>("TrackModule", "pandora");
150 m_showerLabel = pset.get<std::string>("ShowerModule", "pandora");
151 m_printDebug = pset.get<bool>("PrintDebug", false);
152 }
153
154 //------------------------------------------------------------------------------------------------------------------------------------------
155
157 {
158 mf::LogDebug("LArPandora") << " *** PFParticleAnalysis::beginJob() *** " << std::endl;
159
160 //
161 art::ServiceHandle<art::TFileService const> tfs;
162
163 m_pRecoTree = tfs->make<TTree>("pandora", "LAr PFParticles");
164 m_pRecoTree->Branch("run", &m_run, "run/I");
165 m_pRecoTree->Branch("event", &m_event, "event/I");
166 m_pRecoTree->Branch("index", &m_index, "index/I");
167 m_pRecoTree->Branch("self", &m_self, "self/I");
168 m_pRecoTree->Branch("pdgcode", &m_pdgcode, "pdgcode/I");
169 m_pRecoTree->Branch("primary", &m_primary, "primary/I");
170 m_pRecoTree->Branch("parent", &m_parent, "parent/I");
171 m_pRecoTree->Branch("daughters", &m_daughters, "daughters/I");
172 m_pRecoTree->Branch("generation", &m_generation, "generation/I");
173 m_pRecoTree->Branch("neutrino", &m_neutrino, "neutrino/I");
174 m_pRecoTree->Branch("finalstate", &m_finalstate, "finalstate/I");
175 m_pRecoTree->Branch("vertex", &m_vertex, "vertex/I");
176 m_pRecoTree->Branch("track", &m_track, "track/I");
177 m_pRecoTree->Branch("trackid", &m_trackid, "trackid/I");
178 m_pRecoTree->Branch("shower", &m_shower, "shower/I");
179 m_pRecoTree->Branch("showerid", &m_showerid, "showerid/I");
180 m_pRecoTree->Branch("clusters", &m_clusters, "clusters/I");
181 m_pRecoTree->Branch("spacepoints", &m_spacepoints, "spacepoints/I");
182 m_pRecoTree->Branch("hits", &m_hits, "hits/I");
183 m_pRecoTree->Branch("trackhits", &m_trackhits, "trackhits/I");
184 m_pRecoTree->Branch("trajectorypoints", &m_trajectorypoints, "trajectorypoints/I");
185 m_pRecoTree->Branch("showerhits", &m_showerhits, "showerhits/I");
186 m_pRecoTree->Branch("pfovtxx", &m_pfovtxx, "pfovtxx/D");
187 m_pRecoTree->Branch("pfovtxy", &m_pfovtxy, "pfovtxy/D");
188 m_pRecoTree->Branch("pfovtxz", &m_pfovtxz, "pfovtxz/D");
189 m_pRecoTree->Branch("trkvtxx", &m_trkvtxx, "trkvtxx/D");
190 m_pRecoTree->Branch("trkvtxy", &m_trkvtxy, "trkvtxy/D");
191 m_pRecoTree->Branch("trkvtxz", &m_trkvtxz, "trkvtxz/D");
192 m_pRecoTree->Branch("trkvtxdirx", &m_trkvtxdirx, "trkvtxdirx/D");
193 m_pRecoTree->Branch("trkvtxdiry", &m_trkvtxdiry, "trkvtxdiry/D");
194 m_pRecoTree->Branch("trkvtxdirz", &m_trkvtxdirz, "trkvtxdirz/D");
195 m_pRecoTree->Branch("trkendx", &m_trkendx, "trkendx/D");
196 m_pRecoTree->Branch("trkendy", &m_trkendy, "trkendy/D");
197 m_pRecoTree->Branch("trkendz", &m_trkendz, "trkendz/D");
198 m_pRecoTree->Branch("trkenddirx", &m_trkenddirx, "trkenddirx/D");
199 m_pRecoTree->Branch("trkenddiry", &m_trkenddiry, "trkenddiry/D");
200 m_pRecoTree->Branch("trkenddirz", &m_trkenddirz, "trkenddirz/D");
201 m_pRecoTree->Branch("trklength", &m_trklength, "trklength/D");
202 m_pRecoTree->Branch("trkstraightlength", &m_trkstraightlength, "trkstraightlength/D");
203 m_pRecoTree->Branch("shwvtxx", &m_shwvtxx, "shwvtxx/D");
204 m_pRecoTree->Branch("shwvtxy", &m_shwvtxy, "shwvtxy/D");
205 m_pRecoTree->Branch("shwvtxz", &m_shwvtxz, "shwvtxz/D");
206 m_pRecoTree->Branch("shwvtxdirx", &m_shwvtxdirx, "shwvtxdirx/D");
207 m_pRecoTree->Branch("shwvtxdiry", &m_shwvtxdiry, "shwvtxdiry/D");
208 m_pRecoTree->Branch("shwvtxdirz", &m_shwvtxdirz, "shwvtxdirz/D");
209 m_pRecoTree->Branch("shwlength", &m_shwlength, "shwlength/D");
210 m_pRecoTree->Branch("shwopenangle", &m_shwopenangle, "shwopenangle/D");
211 m_pRecoTree->Branch("shwbestplane", &m_shwbestplane, "shwbestplane/D");
212 m_pRecoTree->Branch("t0", &m_t0, "t0/D");
213 }
214
215 //------------------------------------------------------------------------------------------------------------------------------------------
216
218
219 //------------------------------------------------------------------------------------------------------------------------------------------
220
221 void PFParticleAnalysis::analyze(const art::Event& evt)
222 {
223 if (m_printDebug) std::cout << " *** PFParticleAnalysis::analyze(...) *** " << std::endl;
224
225 m_run = evt.run();
226 m_event = evt.id().event();
227 m_index = 0;
228
229 m_self = 0;
230 m_pdgcode = 0;
231 m_primary = 0;
232 m_parent = 0;
233 m_daughters = 0;
234 m_generation = 0;
235 m_neutrino = 0;
236 m_finalstate = 0;
237 m_vertex = 0;
238 m_track = 0;
239 m_trackid = -999;
240 m_shower = 0;
241 m_showerid = -999;
242
243 m_clusters = 0;
244 m_spacepoints = 0;
245 m_hits = 0;
247 m_trackhits = 0;
248 m_showerhits = 0;
249
250 m_pfovtxx = 0.0;
251 m_pfovtxy = 0.0;
252 m_pfovtxz = 0.0;
253
254 m_trkvtxx = 0.0;
255 m_trkvtxy = 0.0;
256 m_trkvtxz = 0.0;
257 m_trkvtxdirx = 0.0;
258 m_trkvtxdiry = 0.0;
259 m_trkvtxdirz = 0.0;
260 m_trkendx = 0.0;
261 m_trkendy = 0.0;
262 m_trkendz = 0.0;
263 m_trkenddirx = 0.0;
264 m_trkenddiry = 0.0;
265 m_trkenddirz = 0.0;
266 m_trklength = 0.0;
268
269 m_shwvtxx = 0.0;
270 m_shwvtxy = 0.0;
271 m_shwvtxz = 0.0;
272 m_shwvtxdirx = 0.0;
273 m_shwvtxdiry = 0.0;
274 m_shwvtxdirz = 0.0;
275 m_shwlength = 0.0;
276 m_shwopenangle = 0.0;
277 m_shwbestplane = 0.0;
278
279 m_t0 = 0;
280
281 if (m_printDebug) {
282 std::cout << " Run: " << m_run << std::endl;
283 std::cout << " Event: " << m_event << std::endl;
284 }
285
286 // Get the reconstructed PFParticles
287 // =================================
288 PFParticleVector particleVector;
289 PFParticleVector particles1, particles2;
290 PFParticlesToClusters particlesToClusters;
291 PFParticlesToSpacePoints particlesToSpacePoints;
292 PFParticlesToHits particlesToHits;
293 HitsToPFParticles hitsToParticles;
294
296 LArPandoraHelper::CollectPFParticles(evt, m_particleLabel, particles1, particlesToClusters);
297 LArPandoraHelper::CollectPFParticles(evt, m_particleLabel, particles2, particlesToSpacePoints);
299 evt, m_particleLabel, particlesToHits, hitsToParticles);
300
301 if (m_printDebug) std::cout << " PFParticles: " << particleVector.size() << std::endl;
302
303 if (particleVector.empty()) {
304 m_pRecoTree->Fill();
305 return;
306 }
307
308 // Get the reconstructed vertices
309 // ==============================
310 VertexVector vertexVector;
311 PFParticlesToVertices particlesToVertices;
312 LArPandoraHelper::CollectVertices(evt, m_particleLabel, vertexVector, particlesToVertices);
313
314 // Get the reconstructed tracks
315 // ============================
316 TrackVector trackVector, trackVector2;
317 PFParticlesToTracks particlesToTracks;
318 TracksToHits tracksToHits;
319 LArPandoraHelper::CollectTracks(evt, m_trackLabel, trackVector, particlesToTracks);
320 LArPandoraHelper::CollectTracks(evt, m_trackLabel, trackVector2, tracksToHits);
321
322 // Get the reconstructed showers
323 // ============================
324 ShowerVector showerVector, showerVector2;
325 PFParticlesToShowers particlesToShowers;
326 ShowersToHits showersToHits;
327 LArPandoraHelper::CollectShowers(evt, m_showerLabel, showerVector, particlesToShowers);
328 LArPandoraHelper::CollectShowers(evt, m_showerLabel, showerVector2, showersToHits);
329
330 // Get the reconstructed T0 objects
331 // ================================
332 T0Vector t0Vector;
333 PFParticlesToT0s particlesToT0s;
334 LArPandoraHelper::CollectT0s(evt, m_particleLabel, t0Vector, particlesToT0s);
335
336 // Build an indexed map of the PFParticles
337 // =======================================
338 PFParticleMap particleMap;
339 LArPandoraHelper::BuildPFParticleMap(particleVector, particleMap);
340
341 // Write PFParticle properties to ROOT file
342 // ========================================
343 for (unsigned int n = 0; n < particleVector.size(); ++n) {
344 const art::Ptr<recob::PFParticle> particle = particleVector.at(n);
345
346 m_index = n;
347 m_self = particle->Self();
348 m_pdgcode = particle->PdgCode();
349 m_primary = particle->IsPrimary();
350 m_parent = (particle->IsPrimary() ? -1 : particle->Parent());
351 m_daughters = particle->NumDaughters();
352 m_generation = LArPandoraHelper::GetGeneration(particleMap, particle);
353 m_neutrino = LArPandoraHelper::GetParentNeutrino(particleMap, particle);
354 m_finalstate = LArPandoraHelper::IsFinalState(particleMap, particle);
355 m_vertex = 0;
356 m_track = 0;
357 m_trackid = -999;
358 m_shower = 0;
359 m_showerid = -999;
360
361 m_clusters = 0;
362 m_spacepoints = 0;
363 m_hits = 0;
365 m_trackhits = 0;
366 m_showerhits = 0;
367
368 m_pfovtxx = 0.0;
369 m_pfovtxy = 0.0;
370 m_pfovtxz = 0.0;
371
372 m_trkvtxx = 0.0;
373 m_trkvtxy = 0.0;
374 m_trkvtxz = 0.0;
375 m_trkvtxdirx = 0.0;
376 m_trkvtxdiry = 0.0;
377 m_trkvtxdirz = 0.0;
378 m_trkendx = 0.0;
379 m_trkendy = 0.0;
380 m_trkendz = 0.0;
381 m_trkenddirx = 0.0;
382 m_trkenddiry = 0.0;
383 m_trkenddirz = 0.0;
384 m_trklength = 0.0;
386
387 m_shwvtxx = 0.0;
388 m_shwvtxy = 0.0;
389 m_shwvtxz = 0.0;
390 m_shwvtxdirx = 0.0;
391 m_shwvtxdiry = 0.0;
392 m_shwvtxdirz = 0.0;
393 m_shwlength = 0.0;
394 m_shwopenangle = 0.0;
395 m_shwbestplane = 0.0;
396
397 m_t0 = 0.0;
398
399 // Particles <-> Clusters
400 PFParticlesToClusters::const_iterator cIter = particlesToClusters.find(particle);
401 if (particlesToClusters.end() != cIter) m_clusters = cIter->second.size();
402
403 // Particles <-> SpacePoints
404 PFParticlesToSpacePoints::const_iterator pIter = particlesToSpacePoints.find(particle);
405 if (particlesToSpacePoints.end() != pIter) m_spacepoints = pIter->second.size();
406
407 // Particles <-> Hits
408 PFParticlesToHits::const_iterator hIter = particlesToHits.find(particle);
409 if (particlesToHits.end() != hIter) m_hits = hIter->second.size();
410
411 // Particles <-> Vertices
412 PFParticlesToVertices::const_iterator vIter = particlesToVertices.find(particle);
413 if (particlesToVertices.end() != vIter) {
414 const VertexVector& vertexVector = vIter->second;
415 if (!vertexVector.empty()) {
416 if (vertexVector.size() != 1 && m_printDebug)
417 std::cout << " Warning: Found particle with more than one associated vertex "
418 << std::endl;
419
420 const art::Ptr<recob::Vertex> vertex = *(vertexVector.begin());
421 double xyz[3] = {0.0, 0.0, 0.0};
422 vertex->XYZ(xyz);
423
424 m_vertex = 1;
425 m_pfovtxx = xyz[0];
426 m_pfovtxy = xyz[1];
427 m_pfovtxz = xyz[2];
428 }
429 }
430
431 // Particles <-> T0s
432 PFParticlesToT0s::const_iterator t0Iter = particlesToT0s.find(particle);
433 if (particlesToT0s.end() != t0Iter) {
434 const T0Vector& t0Vector = t0Iter->second;
435 if (!t0Vector.empty()) {
436 if (t0Vector.size() != 1 && m_printDebug)
437 std::cout << " Warning: Found particle with more than one associated T0 " << std::endl;
438
439 const art::Ptr<anab::T0> t0 = *(t0Vector.begin());
440 m_t0 = t0->Time();
441 }
442 }
443
444 // Particles <-> Tracks <-> Hits, T0s
445 PFParticlesToTracks::const_iterator trkIter = particlesToTracks.find(particle);
446 if (particlesToTracks.end() != trkIter) {
447 const TrackVector& trackVector = trkIter->second;
448 if (!trackVector.empty()) {
449 if (trackVector.size() != 1 && m_printDebug)
450 std::cout << " Warning: Found particle with more than one associated track "
451 << std::endl;
452
453 const art::Ptr<recob::Track> track = *(trackVector.begin());
454 const auto& trackVtxPosition = track->Vertex();
455 const auto& trackVtxDirection = track->VertexDirection();
456 const auto& trackEndPosition = track->End();
457 const auto& trackEndDirection = track->EndDirection();
458
459 m_track = 1;
460 m_trackid = track->ID();
461 m_trajectorypoints = track->NumberTrajectoryPoints();
462 m_trkvtxx = trackVtxPosition.x();
463 m_trkvtxy = trackVtxPosition.y();
464 m_trkvtxz = trackVtxPosition.z();
465 m_trkvtxdirx = trackVtxDirection.x();
466 m_trkvtxdiry = trackVtxDirection.y();
467 m_trkvtxdirz = trackVtxDirection.z();
468 m_trkendx = trackEndPosition.x();
469 m_trkendy = trackEndPosition.y();
470 m_trkendz = trackEndPosition.z();
471 m_trkenddirx = trackEndDirection.x();
472 m_trkenddiry = trackEndDirection.y();
473 m_trkenddirz = trackEndDirection.z();
474 m_trklength = track->Length();
475 m_trkstraightlength = (trackEndPosition - trackVtxPosition).R();
476
477 TracksToHits::const_iterator trkIter2 = tracksToHits.find(track);
478 if (tracksToHits.end() != trkIter2) m_trackhits = trkIter2->second.size();
479 }
480 }
481
482 // Particles <-> Showers <-> Hits
483 PFParticlesToShowers::const_iterator shwIter = particlesToShowers.find(particle);
484 if (particlesToShowers.end() != shwIter) {
485 const ShowerVector& showerVector = shwIter->second;
486 if (!showerVector.empty()) {
487 if (showerVector.size() != 1 && m_printDebug)
488 std::cout << " Warning: Found particle with more than one associated shower "
489 << std::endl;
490
491 const art::Ptr<recob::Shower> shower = *(showerVector.begin());
492 const TVector3& showerVtxPosition = shower->ShowerStart();
493 const TVector3& showerVtxDirection = shower->Direction();
494
495 m_shower = 1;
496 m_showerid = shower->ID();
497
498 m_shwvtxx = showerVtxPosition.x();
499 m_shwvtxy = showerVtxPosition.y();
500 m_shwvtxz = showerVtxPosition.z();
501 m_shwvtxdirx = showerVtxDirection.x();
502 m_shwvtxdiry = showerVtxDirection.y();
503 m_shwvtxdirz = showerVtxDirection.z();
504
505 m_shwlength = shower->Length();
506 m_shwopenangle = shower->OpenAngle();
507 m_shwbestplane = shower->best_plane();
508
509 ShowersToHits::const_iterator shwIter2 = showersToHits.find(shower);
510 if (showersToHits.end() != shwIter2) m_showerhits = shwIter2->second.size();
511 }
512 }
513
514 if (m_printDebug)
515 std::cout << " PFParticle [" << n << "] Primary=" << m_primary
516 << " FinalState=" << m_finalstate << " Pdg=" << m_pdgcode
517 << " NuPdg=" << m_neutrino << " (Self=" << m_self << ", Parent=" << m_parent
518 << ")"
519 << " (Vertex=" << m_vertex << ", Track=" << m_track << ", Shower=" << m_shower
520 << ", Clusters=" << m_clusters << ", SpacePoints=" << m_spacepoints
521 << ", Hits=" << m_hits << ") " << std::endl;
522
523 m_pRecoTree->Fill();
524 }
525 }
526
527} //namespace lar_pandora
helper function for LArPandoraInterface producer module
static void CollectShowers(const art::Event &evt, const std::string &label, ShowerVector &showerVector, PFParticlesToShowers &particlesToShowers)
Collect the reconstructed PFParticles and associated Showers from the ART event record.
static void CollectVertices(const art::Event &evt, const std::string &label, VertexVector &vertexVector, PFParticlesToVertices &particlesToVertices)
Collect the reconstructed PFParticles and associated Vertices from the ART event record.
static void BuildPFParticleMap(const PFParticleVector &particleVector, PFParticleMap &particleMap)
Build particle maps for reconstructed particles.
static void CollectT0s(const art::Event &evt, const std::string &label, T0Vector &t0Vector, PFParticlesToT0s &particlesToT0s)
Collect a vector of T0s from the ART event record.
static int GetParentNeutrino(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Return the parent neutrino PDG code (or zero for cosmics) for a given reconstructed particle.
static void CollectPFParticles(const art::Event &evt, const std::string &label, PFParticleVector &particleVector)
Collect the reconstructed PFParticles from the ART event record.
static bool IsFinalState(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Determine whether a particle has been reconstructed as a final-state particle.
static void CollectTracks(const art::Event &evt, const std::string &label, TrackVector &trackVector, PFParticlesToTracks &particlesToTracks)
Collect the reconstructed PFParticles and associated Tracks from the ART event record.
static int GetGeneration(const PFParticleMap &particleMap, const art::Ptr< recob::PFParticle > daughterParticle)
Return the generation of this particle (first generation if primary)
static void BuildPFParticleHitMaps(const PFParticleVector &particleVector, const PFParticlesToSpacePoints &particlesToSpacePoints, const SpacePointsToHits &spacePointsToHits, PFParticlesToHits &particlesToHits, HitsToPFParticles &hitsToParticles, const DaughterMode daughterMode=kUseDaughters)
Build mapping between PFParticles and Hits using PFParticle/SpacePoint/Hit maps.
PFParticleAnalysis(fhicl::ParameterSet const &pset)
Constructor.
void reconfigure(fhicl::ParameterSet const &pset)
bool m_printDebug
switch for print statements (TODO: use message service!)
std::map< art::Ptr< recob::PFParticle >, T0Vector > PFParticlesToT0s
std::vector< art::Ptr< anab::T0 > > T0Vector
std::map< int, art::Ptr< recob::PFParticle > > PFParticleMap
std::vector< art::Ptr< recob::Track > > TrackVector
std::vector< art::Ptr< recob::Shower > > ShowerVector
std::map< art::Ptr< recob::PFParticle >, ClusterVector > PFParticlesToClusters
std::vector< art::Ptr< recob::Vertex > > VertexVector
std::map< art::Ptr< recob::PFParticle >, ShowerVector > PFParticlesToShowers
std::map< art::Ptr< recob::PFParticle >, SpacePointVector > PFParticlesToSpacePoints
std::map< art::Ptr< recob::PFParticle >, HitVector > PFParticlesToHits
std::vector< art::Ptr< recob::PFParticle > > PFParticleVector
std::map< art::Ptr< recob::Shower >, HitVector > ShowersToHits
std::map< art::Ptr< recob::PFParticle >, TrackVector > PFParticlesToTracks
std::map< art::Ptr< recob::Track >, HitVector > TracksToHits
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::PFParticle > > HitsToPFParticles
std::map< art::Ptr< recob::PFParticle >, VertexVector > PFParticlesToVertices