All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MCTruthT0Matching_module.cc
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////////
2 /// Class: MCTruthT0Matching
3 /// Module Type: producer
4 /// File: MCTruthT0Matching_module.cc
5 ///
6 /// Author: Thomas Karl Warburton
7 /// E-mail address: k.warburton@sheffield.ac.uk
8 ///
9 /// Generated at Wed Mar 25 13:54:28 2015 by Thomas Warburton using artmod
10 /// from cetpkgsupport v1_08_04.
11 ///
12 /// This module accesses the Monte Carlo Truth information stored in the ART
13 /// event and matches that with a track. It does this by looping through the
14 /// tracks in the event and looping through each hit in the track.
15 /// For each hit it uses the backtracker service to work out the charge which
16 /// each MCTruth particle contributed to the total charge desposited for the
17 /// hit.
18 /// The MCTruth particle which is ultimately assigned to the track is simply
19 /// the particle which deposited the most charge.
20 /// It then stores an ART anab::T0 object which has the following variables;
21 /// 1) Generation time of the MCTruth particle assigned to track, in ns.
22 /// 2) The trigger type used to assign T0 (in this case 2 for MCTruth)
23 /// 3) The Geant4 TrackID of the particle (so can access all MCTrtuh info in
24 /// subsequent modules).
25 /// 4) The track number of this track in this event.
26 ///
27 /// The module has been extended to also associate an anab::T0 object with a
28 /// recob::Shower. It does this following the same algorithm, where
29 /// recob::Track has been replaced with recob::Shower.
30 ///
31 /// The module takes a reconstructed track as input.
32 /// The module outputs an anab::T0 object
33 ///
34 /// * Update (25 Oct 2017) --- wketchum@fnal.gov *
35 /// --Add option for storing hit, MCParticle associations.
36 ///
37 /// * Update (6 Nov 2017) --- yuntse@slac.stanford.edu
38 /// --Add a few variables in the metadata of hit, MCParticle associations.
39 ///
40 /////////////////////////////////////////////////////////////////////////////
41 
42 // Framework includes
43 #include "art/Framework/Core/EDProducer.h"
44 #include "art/Framework/Core/ModuleMacros.h"
45 #include "art/Framework/Principal/Event.h"
46 #include "art/Framework/Principal/Handle.h"
47 #include "art_root_io/TFileService.h"
48 #include "canvas/Persistency/Common/FindManyP.h"
49 #include "canvas/Persistency/Common/Ptr.h"
50 
51 #include "canvas/Utilities/InputTag.h"
52 #include "fhiclcpp/ParameterSet.h"
53 
54 #include <iostream>
55 #include <iterator>
56 #include <map>
57 #include <memory>
58 
59 // LArSoft
72 #include "nusimdata/SimulationBase/MCParticle.h"
73 
74 // ROOT
75 #include "TTree.h"
76 
77 namespace t0 {
78  class MCTruthT0Matching;
79 }
80 
81 class t0::MCTruthT0Matching : public art::EDProducer {
82 public:
83  explicit MCTruthT0Matching(fhicl::ParameterSet const& p);
84  // The destructor generated by the compiler is fine for classes
85  // without bare pointers or other resource use.
86 
87  // Plugins should not be copied or assigned.
88  MCTruthT0Matching(MCTruthT0Matching const&) = delete;
92 
93  // Required functions.
94  void produce(art::Event& e) override;
95 
96  // Selected optional functions.
97  void beginJob() override;
98 
99 private:
100  // Params got from fcl file
101  art::InputTag fTrackModuleLabel;
102  art::InputTag fShowerModuleLabel;
103  art::InputTag fPFParticleModuleLabel;
106 
107  art::InputTag fHitModuleLabel;
109 
111 
112  // Variable in TFS branches
113  TTree* fTree;
114  int TrackID = 0;
115  int TrueTrackID = 0;
117  double TrueTrackT0 = 0;
118 
119  int ShowerID = 0;
120  int ShowerMatchID = 0;
122  double ShowerT0 = 0;
123 };
124 
125 t0::MCTruthT0Matching::MCTruthT0Matching(fhicl::ParameterSet const& p) : EDProducer{p}
126 {
127  // Call appropriate produces<>() functions here.
128  fTrackModuleLabel = (p.get<art::InputTag>("TrackModuleLabel"));
129  fShowerModuleLabel = (p.get<art::InputTag>("ShowerModuleLabel"));
130  fPFParticleModuleLabel = (p.get<art::InputTag>("PFParticleModuleLabel", "pandoraNu"));
131  fMakeT0Assns = (p.get<bool>("makeT0Assns", true));
132  fMakePFParticleAssns = (p.get<bool>("makePFParticleAssns", false));
133 
134  fMakeHitAssns = p.get<bool>("makeHitAssns", true);
135  if (fMakeHitAssns) fHitModuleLabel = p.get<art::InputTag>("HitModuleLabel");
136  fOverrideRealData = p.get<bool>("OverrideRealData", false);
137 
138  if (
139  fMakeT0Assns) { // T0 assns are deprecated - this allows one to use deprecated funcionality. Added 2017-08-15. Should not be kept around forever
140  std::cout << "WARNING - You are using deprecated functionality\n";
141  std::cout << "MCTruthT0Matching T0 assns will be removed soon\n";
142  std::cout << "set your fcl parameter makeT0Assns to false and use MCParticle direct "
143  "associations instead"
144  << std::endl;
145  produces<std::vector<anab::T0>>();
146  produces<art::Assns<recob::Track, anab::T0>>();
147  produces<art::Assns<recob::Shower, anab::T0>>();
148  if (fMakePFParticleAssns) {
149  produces<art::Assns<recob::PFParticle, anab::T0>>();
150  } // only do PFParticles if desired by the user
151  }
152 
153  produces<art::Assns<recob::Track, simb::MCParticle, anab::BackTrackerMatchingData>>();
154  produces<art::Assns<recob::Shower, simb::MCParticle, anab::BackTrackerMatchingData>>();
155  if (fMakePFParticleAssns) { // only do PFParticles if desired by the user
156  produces<art::Assns<recob::PFParticle, simb::MCParticle, anab::BackTrackerMatchingData>>();
157  }
158 
159  if (fMakeHitAssns)
160  produces<art::Assns<recob::Hit, simb::MCParticle, anab::BackTrackerHitMatchingData>>();
161 }
162 
163 void
165 {
166  // Implementation of optional member function here.
167  art::ServiceHandle<art::TFileService const> tfs;
168  fTree = tfs->make<TTree>("MCTruthT0Matching", "MCTruthT0");
169  fTree->Branch("TrueTrackT0", &TrueTrackT0, "TrueTrackT0/D");
170  fTree->Branch("TrueTrackID", &TrueTrackID, "TrueTrackID/D");
171 }
172 
173 void
175 {
176  if (evt.isRealData() && !fOverrideRealData) return;
177 
178  // Access art services...
179  art::ServiceHandle<geo::Geometry const> geom;
180  art::ServiceHandle<cheat::BackTrackerService const> bt_serv;
181  art::ServiceHandle<cheat::ParticleInventoryService const> pi_serv;
182  auto const clockData =
183  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
184 
185  //TrackList handle
186  art::Handle<std::vector<recob::Track>> trackListHandle;
187  std::vector<art::Ptr<recob::Track>> tracklist;
188  if (evt.getByLabel(fTrackModuleLabel, trackListHandle))
189  art::fill_ptr_vector(tracklist, trackListHandle);
190 
191  //ShowerList handle
192  art::Handle<std::vector<recob::Shower>> showerListHandle;
193  std::vector<art::Ptr<recob::Shower>> showerlist;
194  if (evt.getByLabel(fShowerModuleLabel, showerListHandle))
195  art::fill_ptr_vector(showerlist, showerListHandle);
196 
197  //PFParticleList handle
198  art::Handle<std::vector<recob::PFParticle>> pfparticleListHandle;
199  std::vector<art::Ptr<recob::PFParticle>> pfparticlelist;
200  if (evt.getByLabel(fPFParticleModuleLabel, pfparticleListHandle))
201  art::fill_ptr_vector(pfparticlelist, pfparticleListHandle);
202 
203  auto mcpartHandle = evt.getValidHandle<std::vector<simb::MCParticle>>("largeant");
204  // simb::MCParticle const *firstParticle = &mcpartHandle->front();
205  // Create anab::T0 objects and make association with recob::Track, recob::Shower, and recob::PFParticle objects
206  std::unique_ptr<std::vector<anab::T0>> T0col(new std::vector<anab::T0>);
207  // if (fMakeT0Assns){ // T0 assns are deprecated - this allows one to use deprecated funcionality. Added 2017-08-15. Should not be kept around forever
208  std::unique_ptr<art::Assns<recob::Track, anab::T0>> Trackassn(
209  new art::Assns<recob::Track, anab::T0>);
210  std::unique_ptr<art::Assns<recob::Shower, anab::T0>> Showerassn(
211  new art::Assns<recob::Shower, anab::T0>);
212  std::unique_ptr<art::Assns<recob::PFParticle, anab::T0>> PFParticleassn(
213  new art::Assns<recob::PFParticle, anab::T0>);
214 
215  // Create associations directly between MCParticle and either recob::Track, recob::Shower, or recob::PFParticle
216  // These associations contain metadata summarising the quality of the match
217  std::unique_ptr<art::Assns<recob::Track, simb::MCParticle, anab::BackTrackerMatchingData>>
218  MCPartTrackassn(new art::Assns<recob::Track, simb::MCParticle, anab::BackTrackerMatchingData>);
219  std::unique_ptr<art::Assns<recob::Shower, simb::MCParticle, anab::BackTrackerMatchingData>>
220  MCPartShowerassn(
221  new art::Assns<recob::Shower, simb::MCParticle, anab::BackTrackerMatchingData>);
222  std::unique_ptr<art::Assns<recob::PFParticle, simb::MCParticle, anab::BackTrackerMatchingData>>
223  MCPartPFParticleassn(
224  new art::Assns<recob::PFParticle, simb::MCParticle, anab::BackTrackerMatchingData>);
225  // Association block for the hits<-->MCParticles
226  std::unique_ptr<art::Assns<recob::Hit, simb::MCParticle, anab::BackTrackerHitMatchingData>>
227  MCPartHitassn(new art::Assns<recob::Hit, simb::MCParticle, anab::BackTrackerHitMatchingData>);
228 
229  double maxe = -1;
230  double tote = 0;
231  // int trkid = -1;
232  int maxtrkid = -1;
233  double maxn = -1;
234  double totn = 0;
235  int maxntrkid = -1;
237 
238  std::unordered_map<int, int> trkid_lookup; //indexed by geant4trkid, delivers MC particle location
239 
240  //if we want to make per-hit assns
241  if (fMakeHitAssns) {
242 
243  art::Handle<std::vector<recob::Hit>> hitListHandle;
244  evt.getByLabel(fHitModuleLabel, hitListHandle);
245 
246  if (hitListHandle.isValid()) {
247 
248  auto const& hitList(*hitListHandle);
249  auto const& mcpartList(*mcpartHandle);
250  for (size_t i_h = 0; i_h < hitList.size(); ++i_h) {
251  art::Ptr<recob::Hit> hitPtr(hitListHandle, i_h);
252  auto trkide_list = bt_serv->HitToTrackIDEs(clockData, hitPtr);
253  struct TrackIDEinfo {
254  float E;
255  float NumElectrons;
256  };
257  std::map<int, TrackIDEinfo> trkide_collector;
258  maxe = -1;
259  tote = 0;
260  maxtrkid = -1;
261  maxn = -1;
262  totn = 0;
263  maxntrkid = -1;
264  for (auto const& t : trkide_list) {
265  trkide_collector[t.trackID].E += t.energy;
266  tote += t.energy;
267  if (trkide_collector[t.trackID].E > maxe) {
268  maxe = trkide_collector[t.trackID].E;
269  maxtrkid = t.trackID;
270  }
271  trkide_collector[t.trackID].NumElectrons += t.numElectrons;
272  totn += t.numElectrons;
273  if (trkide_collector[t.trackID].NumElectrons > maxn) {
274  maxn = trkide_collector[t.trackID].NumElectrons;
275  maxntrkid = t.trackID;
276  }
277 
278  //if not found, find mc particle...
279  if (trkid_lookup.find(t.trackID) == trkid_lookup.end()) {
280  size_t i_p = 0;
281  while (i_p < mcpartList.size()) {
282  if (mcpartList[i_p].TrackId() == t.trackID) {
283  trkid_lookup[t.trackID] = (int)i_p;
284  break;
285  }
286  ++i_p;
287  }
288  if (i_p == mcpartList.size()) trkid_lookup[t.trackID] = -1;
289  }
290 
291  } //end loop on TrackIDs
292 
293  //now find the mcparticle and loop back through ...
294  for (auto const& t : trkide_collector) {
295  int mcpart_i = trkid_lookup[t.first];
296  if (mcpart_i == -1) continue; //no mcparticle here
297  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
298  bthmd.ideFraction = t.second.E / tote;
299  bthmd.isMaxIDE = (t.first == maxtrkid);
300  bthmd.ideNFraction = t.second.NumElectrons / totn;
301  bthmd.isMaxIDEN = (t.first == maxntrkid);
302  MCPartHitassn->addSingle(hitPtr, mcpartPtr, bthmd);
303  }
304 
305  } //end loop on hits
306 
307  } //end if handle valid
308 
309  } //end if make hit/mcpart assns
310 
311  if (trackListHandle.isValid()) {
312  //Access tracks and hits
313  art::FindManyP<recob::Hit> fmtht(trackListHandle, evt, fTrackModuleLabel);
314 
315  size_t NTracks = tracklist.size();
316 
317  // Now to access MCTruth for each track...
318  for (size_t iTrk = 0; iTrk < NTracks; ++iTrk) {
319  TrueTrackT0 = 0;
320  TrackID = 0;
321  TrueTrackID = 0;
323  std::vector<art::Ptr<recob::Hit>> allHits = fmtht.at(iTrk);
324 
325  std::map<int, double> trkide;
326  for (size_t h = 0; h < allHits.size(); ++h) {
327  art::Ptr<recob::Hit> hit = allHits[h];
328  std::vector<sim::IDE> ides;
329  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
330 
331  for (size_t e = 0; e < TrackIDs.size(); ++e) {
332  trkide[TrackIDs[e].trackID] += TrackIDs[e].energy;
333  }
334  }
335  // Work out which IDE despoited the most charge in the hit if there was more than one.
336  maxe = -1;
337  tote = 0;
338  for (std::map<int, double>::iterator ii = trkide.begin(); ii != trkide.end(); ++ii) {
339  tote += ii->second;
340  if ((ii->second) > maxe) {
341  maxe = ii->second;
342  TrackID = ii->first;
343  }
344  }
345  btdata.cleanliness = maxe / tote;
346 
347  // Now have trackID, so get PdG code and T0 etc.
348  const simb::MCParticle* tmpParticle = pi_serv->TrackIdToParticle_P(TrackID);
349  if (!tmpParticle)
350  continue; // Retain this check that the BackTracker can find the right particle
351  // Now, loop through the MCParticle's myself to find the correct match
352  int mcpart_i(-1);
353  for (auto const particle : *mcpartHandle) {
354  mcpart_i++;
355  if (TrackID == particle.TrackId()) { break; }
356  }
357  const simb::MCParticle particle = mcpartHandle.product()->at(mcpart_i);
358  TrueTrackT0 = particle.T();
359  TrueTrackID = particle.TrackId();
360  TrueTriggerType = 2; // Using MCTruth as trigger, so tigger type is 2.
361 
362  T0col->push_back(anab::T0(TrueTrackT0, TrueTriggerType, TrueTrackID, (*T0col).size()));
363  //auto diff = particle - firstParticle;
364  auto diff = mcpart_i; // check I have a sensible value for this counter
365  if (diff >= (int)mcpartHandle->size()) {
366  std::cout << "Error, the backtracker is doing weird things to your pointers!" << std::endl;
367  throw std::exception();
368  }
369 
370  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
371  MCPartTrackassn->addSingle(tracklist[iTrk], mcpartPtr, btdata);
372  if (fMakeT0Assns) { util::CreateAssn(*this, evt, *T0col, tracklist[iTrk], *Trackassn); }
373  fTree->Fill();
374 
375  } // Loop over tracks
376  }
377 
378  if (showerListHandle.isValid()) {
379  art::FindManyP<recob::Hit> fmsht(showerListHandle, evt, fShowerModuleLabel);
380  // Now Loop over showers....
381  size_t NShowers = showerlist.size();
382  for (size_t Shower = 0; Shower < NShowers; ++Shower) {
383  ShowerMatchID = 0;
384  ShowerID = 0;
385  ShowerT0 = 0;
386  std::vector<art::Ptr<recob::Hit>> allHits = fmsht.at(Shower);
388 
389  std::map<int, double> showeride;
390  for (size_t h = 0; h < allHits.size(); ++h) {
391  art::Ptr<recob::Hit> hit = allHits[h];
392  std::vector<sim::IDE> ides;
393  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
394 
395  for (size_t e = 0; e < TrackIDs.size(); ++e) {
396  showeride[TrackIDs[e].trackID] += TrackIDs[e].energy;
397  }
398  }
399  // Work out which IDE despoited the most charge in the hit if there was more than one.
400  maxe = -1;
401  tote = 0;
402  for (std::map<int, double>::iterator ii = showeride.begin(); ii != showeride.end(); ++ii) {
403  tote += ii->second;
404  if ((ii->second) > maxe) {
405  maxe = ii->second;
406  ShowerID = ii->first;
407  }
408  }
409  btdata.cleanliness = maxe / tote;
410 
411  // Now have MCParticle trackID corresponding to shower, so get PdG code and T0 etc.
412  const simb::MCParticle* tmpParticle = pi_serv->TrackIdToParticle_P(ShowerID);
413  if (!tmpParticle)
414  continue; // Retain this check that the BackTracker can find the right particle
415  // Now, loop through the MCParticle's myself to find the correct match
416  int mcpart_i(-1);
417  for (auto const particle : *mcpartHandle) {
418  mcpart_i++;
419  if (ShowerID == particle.TrackId()) { break; }
420  }
421  const simb::MCParticle particle = mcpartHandle.product()->at(mcpart_i);
422  ShowerT0 = particle.T();
423  ShowerID = particle.TrackId();
424  ShowerTriggerType = 2; // Using MCTruth as trigger, so tigger type is 2.
425  T0col->push_back(anab::T0(ShowerT0, ShowerTriggerType, ShowerID, (*T0col).size()));
426 
427  auto diff = mcpart_i; // check I have a sensible value for this counter
428  if (diff >= (int)mcpartHandle->size()) {
429  std::cout << "Error, the backtracker is doing weird things to your pointers!" << std::endl;
430  throw std::exception();
431  }
432  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
433  if (fMakeT0Assns) { util::CreateAssn(*this, evt, *T0col, showerlist[Shower], *Showerassn); }
434  MCPartShowerassn->addSingle(showerlist[Shower], mcpartPtr, btdata);
435 
436  } // Loop over showers
437  }
438 
439  if (pfparticleListHandle.isValid()) {
440  //Access pfparticles and hits
441  art::FindManyP<recob::Cluster> fmcp(pfparticleListHandle, evt, fPFParticleModuleLabel);
442 
443  size_t NPfparticles = pfparticlelist.size();
444 
445  // Now to access MCTruth for each pfparticle...
446  for (size_t iPfp = 0; iPfp < NPfparticles; ++iPfp) {
447  TrueTrackT0 = 0;
448  TrackID = 0;
449  TrueTrackID = 0;
451 
452  std::vector<art::Ptr<recob::Hit>> allHits;
453  //Get all hits through associated clusters
454  std::vector<art::Ptr<recob::Cluster>> allClusters = fmcp.at(iPfp);
455  art::FindManyP<recob::Hit> fmhcl(allClusters, evt, fPFParticleModuleLabel);
456  for (size_t iclu = 0; iclu < allClusters.size(); ++iclu) {
457  std::vector<art::Ptr<recob::Hit>> hits = fmhcl.at(iclu);
458  allHits.insert(allHits.end(), hits.begin(), hits.end());
459  }
460 
461  std::map<int, double> trkide;
462  for (size_t h = 0; h < allHits.size(); ++h) {
463  art::Ptr<recob::Hit> hit = allHits[h];
464  std::vector<sim::IDE> ides;
465  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToTrackIDEs(clockData, hit);
466 
467  for (size_t e = 0; e < TrackIDs.size(); ++e) {
468  trkide[TrackIDs[e].trackID] += TrackIDs[e].energy;
469  }
470  }
471  // Work out which IDE despoited the most charge in the hit if there was more than one.
472  double maxe = -1;
473  double tote = 0;
474  for (std::map<int, double>::iterator ii = trkide.begin(); ii != trkide.end(); ++ii) {
475  tote += ii->second;
476  if ((ii->second) > maxe) {
477  maxe = ii->second;
478  TrackID = ii->first;
479  }
480  }
481  btdata.cleanliness = maxe / tote;
482 
483  // Now have trackID, so get PdG code and T0 etc.
484  const simb::MCParticle* tmpParticle = pi_serv->TrackIdToParticle_P(TrackID);
485  if (!tmpParticle)
486  continue; // Retain this check that the BackTracker can find the right particle
487  // Now, loop through the MCParticle's myself to find the correct match
488  int mcpart_i(-1);
489  for (auto const particle : *mcpartHandle) {
490  mcpart_i++;
491  if (TrackID == particle.TrackId()) { break; }
492  }
493  const simb::MCParticle particle = mcpartHandle.product()->at(mcpart_i);
494  TrueTrackT0 = particle.T();
495  TrueTrackID = particle.TrackId();
496  TrueTriggerType = 2; // Using MCTruth as trigger, so tigger type is 2.
497  //std::cout << "Got particle, PDG = " << particle.PdgCode() << std::endl;
498 
499  //std::cout << "Filling T0col with " << TrueTrackT0 << " " << TrueTriggerType << " " << TrueTrackID << " " << (*T0col).size() << std::endl;
500 
501  T0col->push_back(anab::T0(TrueTrackT0, TrueTriggerType, TrueTrackID, (*T0col).size()));
502  //auto diff = particle - firstParticle;
503  auto diff = mcpart_i; // check I have a sensible value for this counter
504  if (diff >= (int)mcpartHandle->size()) {
505  std::cout << "Error, the backtracker is doing weird things to your pointers!" << std::endl;
506  throw std::exception();
507  }
508  // art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, particle - firstParticle);
509  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
510  //std::cout << "made MCParticle Ptr" << std::endl;
511  // const std::vertex< std::pair<double, std::string> > cleanlinessCompleteness;
512  // const std::pair<double, double> tmp(0.5, 0.5);
513  // MCPartassn->addSingle(tracklist[iPfp], mcpartPtr, tmp);
514  if (fMakePFParticleAssns) {
515  if (fMakeT0Assns) {
516  util::CreateAssn(*this, evt, *T0col, pfparticlelist[iPfp], *PFParticleassn);
517  }
518  MCPartPFParticleassn->addSingle(pfparticlelist[iPfp], mcpartPtr, btdata);
519  }
520  fTree->Fill();
521  } // Loop over tracks
522  }
523 
524  if (fMakeT0Assns) {
525  evt.put(std::move(T0col));
526  evt.put(std::move(Trackassn));
527  evt.put(std::move(Showerassn));
528  if (fMakePFParticleAssns) { evt.put(std::move(PFParticleassn)); }
529  }
530  evt.put(std::move(MCPartTrackassn));
531  evt.put(std::move(MCPartShowerassn));
532  if (fMakePFParticleAssns) { evt.put(std::move(MCPartPFParticleassn)); }
533  if (fMakeHitAssns) evt.put(std::move(MCPartHitassn));
534 } // Produce
535 
536 DEFINE_ART_MODULE(t0::MCTruthT0Matching)
MCTruthT0Matching & operator=(MCTruthT0Matching const &)=delete
std::vector< TrackID > TrackIDs
Declaration of signal hit object.
pdgs p
Definition: selectors.fcl:22
process_name E
Definition: Data.h:14
Definition: T0.h:16
process_name hit
Definition: cheaterreco.fcl:51
MCTruthT0Matching(fhicl::ParameterSet const &p)
void produce(art::Event &e) override
while getopts h
Declaration of cluster object.
Provides recob::Track data product.
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
do i e
art::ServiceHandle< art::TFileService > tfs
TCEvent evt
Definition: DataStructs.cxx:8
art framework interface to geometry description
BEGIN_PROLOG could also be cout