LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
MCTruthT0Matching_module.cc
Go to the documentation of this file.
1 
42 // Framework includes
54 
56 #include "fhiclcpp/ParameterSet.h"
58 
59 #include <memory>
60 #include <iostream>
61 #include <map>
62 #include <iterator>
63 
64 // LArSoft
86 
87 // ROOT
88 #include "TTree.h"
89 #include "TFile.h"
90 
91 namespace t0 {
92  class MCTruthT0Matching;
93 }
94 
96 public:
97  explicit MCTruthT0Matching(fhicl::ParameterSet const & p);
98  // The destructor generated by the compiler is fine for classes
99  // without bare pointers or other resource use.
100 
101  // Plugins should not be copied or assigned.
102  MCTruthT0Matching(MCTruthT0Matching const &) = delete;
104  MCTruthT0Matching & operator = (MCTruthT0Matching const &) = delete;
106 
107  // Required functions.
108  void produce(art::Event & e) override;
109 
110  // Selected optional functions.
111  void beginJob() override;
112  void reconfigure(fhicl::ParameterSet const & p) ;
113 
114 
115 private:
116 
117  // Params got from fcl file
123 
126 
128 
129  // Variable in TFS branches
130  TTree* fTree;
131  int TrackID = 0;
132  int TrueTrackID = 0;
134  double TrueTrackT0 = 0;
135 
136  int ShowerID = 0;
137  int ShowerMatchID = 0;
139  double ShowerT0 = 0;
140 };
141 
142 
144 // :
145 // Initialize member data here, if know don't want to reconfigure on the fly
146 {
147  // Call appropriate produces<>() functions here.
148  reconfigure(p);
149  if (fMakeT0Assns){ // T0 assns are deprecated - this allows one to use deprecated funcionality. Added 2017-08-15. Should not be kept around forever
150  std::cout << "WARNING - You are using deprecated functionality\n";
151  std::cout << "MCTruthT0Matching T0 assns will be removed soon\n";
152  std::cout << "set your fcl parameter makeT0Assns to false and use MCParticle direct associations instead" << std::endl;
153  produces< std::vector<anab::T0> >();
154  produces< art::Assns<recob::Track , anab::T0> >();
155  produces< art::Assns<recob::Shower, anab::T0> > ();
156  if (fMakePFParticleAssns){produces< art::Assns<recob::PFParticle, anab::T0> > ();}// only do PFParticles if desired by the user
157  }
158 
159  produces< art::Assns<recob::Track , simb::MCParticle, anab::BackTrackerMatchingData > > ();
160  produces< art::Assns<recob::Shower, simb::MCParticle, anab::BackTrackerMatchingData > > ();
161  if (fMakePFParticleAssns){ // only do PFParticles if desired by the user
162  produces< art::Assns<recob::PFParticle, simb::MCParticle, anab::BackTrackerMatchingData > > ();
163  }
164 
165  if(fMakeHitAssns)
166  produces< art::Assns<recob::Hit , simb::MCParticle, anab::BackTrackerHitMatchingData > > ();
167 
168 }
169 
171 {
172  // Implementation of optional member function here.
173  fTrackModuleLabel = (p.get< art::InputTag > ("TrackModuleLabel" ) );
174  fShowerModuleLabel = (p.get< art::InputTag > ("ShowerModuleLabel") );
175  fPFParticleModuleLabel = (p.get< art::InputTag > ("PFParticleModuleLabel", "pandoraNu") );
176  fMakeT0Assns = (p.get< bool > ("makeT0Assns", true) );
177  fMakePFParticleAssns = (p.get< bool > ("makePFParticleAssns", false) );
178 
179  fMakeHitAssns = p.get<bool>("makeHitAssns",true);
180  if(fMakeHitAssns) fHitModuleLabel = p.get<art::InputTag>("HitModuleLabel");
181  fOverrideRealData = p.get<bool>("OverrideRealData",false);
182 }
183 
185 {
186  // Implementation of optional member function here.
188  fTree = tfs->make<TTree>("MCTruthT0Matching","MCTruthT0");
189  fTree->Branch("TrueTrackT0",&TrueTrackT0,"TrueTrackT0/D");
190  fTree->Branch("TrueTrackID",&TrueTrackID,"TrueTrackID/D");
191  //fTree->Branch("ShowerID" ,&ShowerID ,"ShowerID/I" );
192  //fTree->Branch("ShowerT0" ,&ShowerT0 ,"ShowerT0/D" );
193 }
194 
196 {
197  if (evt.isRealData() && !fOverrideRealData) return;
198 
199  // Access art services...
203 
204  //TrackList handle
205  art::Handle< std::vector<recob::Track> > trackListHandle;
206  std::vector<art::Ptr<recob::Track> > tracklist;
207  if (evt.getByLabel(fTrackModuleLabel,trackListHandle))
208  art::fill_ptr_vector(tracklist, trackListHandle);
209  //if (!trackListHandle.isValid()) trackListHandle.clear();
210 
211  //ShowerList handle
212  art::Handle< std::vector<recob::Shower> > showerListHandle;
213  std::vector<art::Ptr<recob::Shower> > showerlist;
214  if (evt.getByLabel(fShowerModuleLabel,showerListHandle))
215  art::fill_ptr_vector(showerlist, showerListHandle);
216  //if (!showerListHandle.isValid()) showerListHandle.clear();
217 
218  //PFParticleList handle
219  art::Handle< std::vector<recob::PFParticle> > pfparticleListHandle;
220  std::vector<art::Ptr<recob::PFParticle> > pfparticlelist;
221  if (evt.getByLabel(fPFParticleModuleLabel,pfparticleListHandle))
222  art::fill_ptr_vector(pfparticlelist, pfparticleListHandle);
223  //if (!pfparticleListHandle.isValid()) pfparticleListHandle.clear();
224 
225 
226  auto mcpartHandle = evt.getValidHandle< std::vector<simb::MCParticle> >("largeant");
227 // simb::MCParticle const *firstParticle = &mcpartHandle->front();
228  // Create anab::T0 objects and make association with recob::Track, recob::Shower, and recob::PFParticle objects
229  std::unique_ptr< std::vector<anab::T0> > T0col( new std::vector<anab::T0>);
230 // if (fMakeT0Assns){ // T0 assns are deprecated - this allows one to use deprecated funcionality. Added 2017-08-15. Should not be kept around forever
231  std::unique_ptr< art::Assns<recob::Track, anab::T0> > Trackassn( new art::Assns<recob::Track, anab::T0>);
232  std::unique_ptr< art::Assns<recob::Shower, anab::T0> > Showerassn( new art::Assns<recob::Shower, anab::T0>);
233 // if (fMakePFParticleAssns){
234  std::unique_ptr< art::Assns<recob::PFParticle, anab::T0> > PFParticleassn( new art::Assns<recob::PFParticle, anab::T0>);
235 // }
236 // }
237 
238  // Create associations directly between MCParticle and either recob::Track, recob::Shower, or recob::PFParticle
239  // These associations contain metadata summarising the quality of the match
240  std::unique_ptr< art::Assns<recob::Track, simb::MCParticle, anab::BackTrackerMatchingData > > MCPartTrackassn( new art::Assns<recob::Track, simb::MCParticle, anab::BackTrackerMatchingData >);
241  std::unique_ptr< art::Assns<recob::Shower, simb::MCParticle, anab::BackTrackerMatchingData> > MCPartShowerassn( new art::Assns<recob::Shower, simb::MCParticle, anab::BackTrackerMatchingData>);
242 // if (fMakePFParticleAssns){
243  std::unique_ptr< art::Assns<recob::PFParticle, simb::MCParticle,anab::BackTrackerMatchingData> > MCPartPFParticleassn( new art::Assns<recob::PFParticle, simb::MCParticle, anab::BackTrackerMatchingData>);
244 // }
245  // Association block for the hits<-->MCParticles
246  std::unique_ptr< art::Assns<recob::Hit, simb::MCParticle, anab::BackTrackerHitMatchingData > > MCPartHitassn( new art::Assns<recob::Hit, simb::MCParticle, anab::BackTrackerHitMatchingData >);
247 
248 
249  double maxe = -1;
250  double tote = 0;
251  // int trkid = -1;
252  int maxtrkid = -1;
253  double maxn = -1;
254  double totn = 0;
255  int maxntrkid = -1;
257 
258  std::unordered_map<int,int> trkid_lookup; //indexed by geant4trkid, delivers MC particle location
259 
260  //if we want to make per-hit assns
261  if(fMakeHitAssns){
262 
263 
264 
265  art::Handle< std::vector<recob::Hit> > hitListHandle;
266  evt.getByLabel(fHitModuleLabel,hitListHandle);
267 
268  if(hitListHandle.isValid()){
269 
270  auto const& hitList(*hitListHandle);
271  auto const& mcpartList(*mcpartHandle);
272  for(size_t i_h=0; i_h<hitList.size(); ++i_h){
273  art::Ptr<recob::Hit> hitPtr(hitListHandle, i_h);
274  auto trkide_list = bt_serv->HitToTrackIDEs(hitPtr);
275  struct TrackIDEinfo {
276  float E;
277  float NumElectrons;
278  };
279  std::map<int, TrackIDEinfo> trkide_collector;
280  maxe = -1; tote = 0; maxtrkid = -1;
281  maxn = -1; totn = 0; maxntrkid = -1;
282  for(auto const& t : trkide_list){
283  trkide_collector[t.trackID].E += t.energy;
284  tote += t.energy;
285  if(trkide_collector[t.trackID].E>maxe) { maxe = trkide_collector[t.trackID].E; maxtrkid = t.trackID; }
286  trkide_collector[t.trackID].NumElectrons += t.numElectrons;
287  totn += t.numElectrons;
288  if(trkide_collector[t.trackID].NumElectrons > maxn) {
289  maxn = trkide_collector[t.trackID].NumElectrons;
290  maxntrkid = t.trackID;
291  }
292 
293  //if not found, find mc particle...
294  if(trkid_lookup.find(t.trackID)==trkid_lookup.end()){
295  size_t i_p=0;
296  while(i_p<mcpartList.size()){
297  if(mcpartList[i_p].TrackId() == t.trackID) { trkid_lookup[t.trackID] = (int)i_p; break;}
298  ++i_p;
299  }
300  if(i_p==mcpartList.size()) trkid_lookup[t.trackID] = -1;
301  }
302 
303  }//end loop on TrackIDs
304 
305  //now find the mcparticle and loop back through ...
306  for(auto const& t : trkide_collector){
307  int mcpart_i = trkid_lookup[t.first];
308  if(mcpart_i==-1) continue; //no mcparticle here
309  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
310  bthmd.ideFraction = t.second.E / tote;
311  bthmd.isMaxIDE = (t.first==maxtrkid);
312  bthmd.ideNFraction = t.second.NumElectrons / totn;
313  bthmd.isMaxIDEN = ( t.first == maxntrkid );
314  MCPartHitassn->addSingle(hitPtr, mcpartPtr, bthmd);
315  }
316 
317 
318  }//end loop on hits
319 
320  }//end if handle valid
321 
322  }//end if make hit/mcpart assns
323 
324 
325  if (trackListHandle.isValid()){
326  //Access tracks and hits
327  art::FindManyP<recob::Hit> fmtht(trackListHandle, evt, fTrackModuleLabel);
328 
329  size_t NTracks = tracklist.size();
330 
331  // Now to access MCTruth for each track...
332  for(size_t iTrk=0; iTrk < NTracks; ++iTrk) {
333  TrueTrackT0 = 0;
334  TrackID = 0;
335  TrueTrackID = 0;
337  std::vector< art::Ptr<recob::Hit> > allHits = fmtht.at(iTrk);
338 
339  std::map<int,double> trkide;
340  for(size_t h = 0; h < allHits.size(); ++h){
341  art::Ptr<recob::Hit> hit = allHits[h];
342  std::vector<sim::IDE> ides;
343  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToTrackIDEs(hit);
344 
345  for(size_t e = 0; e < TrackIDs.size(); ++e){
346  trkide[TrackIDs[e].trackID] += TrackIDs[e].energy;
347  }
348 
349  }
350  // Work out which IDE despoited the most charge in the hit if there was more than one.
351  maxe = -1;
352  tote = 0;
353  for (std::map<int,double>::iterator ii = trkide.begin(); ii!=trkide.end(); ++ii){
354  tote += ii->second;
355  if ((ii->second)>maxe){
356  maxe = ii->second;
357  TrackID = ii->first;
358  }
359  }
360  btdata.cleanliness = maxe/tote;
361 
362  // Now have trackID, so get PdG code and T0 etc.
363  const simb::MCParticle *tmpParticle = pi_serv->TrackIdToParticle_P(TrackID);
364  if (!tmpParticle) continue; // Retain this check that the BackTracker can find the right particle
365  // Now, loop through the MCParticle's myself to find the correct match
366  int mcpart_i(-1);
367  for (auto const particle : *mcpartHandle){
368  mcpart_i++;
369  if (TrackID == particle.TrackId()){
370  break;
371  }
372  }
373  const simb::MCParticle particle = mcpartHandle.product()->at(mcpart_i);
374  TrueTrackT0 = particle.T();
375  TrueTrackID = particle.TrackId();
376  TrueTriggerType = 2; // Using MCTruth as trigger, so tigger type is 2.
377  //std::cout << "Got particle, PDG = " << particle.PdgCode() << std::endl;
378 
379  //std::cout << "Filling T0col with " << TrueTrackT0 << " " << TrueTriggerType << " " << TrueTrackID << " " << (*T0col).size() << std::endl;
380 
381  T0col->push_back(anab::T0(TrueTrackT0,
383  TrueTrackID,
384  (*T0col).size()
385  ));
386  //auto diff = particle - firstParticle;
387  auto diff = mcpart_i; // check I have a sensible value for this counter
388  if (diff >= (int)mcpartHandle->size()){
389  std::cout << "Error, the backtracker is doing weird things to your pointers!" << std::endl;
390  throw std::exception();
391  }
392 
393  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
394  MCPartTrackassn->addSingle(tracklist[iTrk], mcpartPtr, btdata);
395  if (fMakeT0Assns){
396  util::CreateAssn(*this, evt, *T0col, tracklist[iTrk], *Trackassn);
397  }
398  fTree -> Fill();
399 
400  } // Loop over tracks
401  }
402 
403  if (showerListHandle.isValid()){
404  art::FindManyP<recob::Hit> fmsht(showerListHandle,evt, fShowerModuleLabel);
405  // Now Loop over showers....
406  size_t NShowers = showerlist.size();
407  for (size_t Shower = 0; Shower < NShowers; ++Shower) {
408  ShowerMatchID = 0;
409  ShowerID = 0;
410  ShowerT0 = 0;
411  std::vector< art::Ptr<recob::Hit> > allHits = fmsht.at(Shower);
413 
414  std::map<int,double> showeride;
415  for(size_t h = 0; h < allHits.size(); ++h){
416  art::Ptr<recob::Hit> hit = allHits[h];
417  std::vector<sim::IDE> ides;
418  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToTrackIDEs(hit);
419 
420  for(size_t e = 0; e < TrackIDs.size(); ++e){
421  showeride[TrackIDs[e].trackID] += TrackIDs[e].energy;
422  }
423  }
424  // Work out which IDE despoited the most charge in the hit if there was more than one.
425  maxe = -1;
426  tote = 0;
427  for (std::map<int,double>::iterator ii = showeride.begin(); ii!=showeride.end(); ++ii){
428  tote += ii->second;
429  if ((ii->second)>maxe){
430  maxe = ii->second;
431  ShowerID = ii->first;
432  }
433  }
434  btdata.cleanliness = maxe/tote;
435 
436 
437 
438 
439  // Now have MCParticle trackID corresponding to shower, so get PdG code and T0 etc.
440  const simb::MCParticle *tmpParticle = pi_serv->TrackIdToParticle_P(ShowerID);
441  if (!tmpParticle) continue; // Retain this check that the BackTracker can find the right particle
442  // Now, loop through the MCParticle's myself to find the correct match
443  int mcpart_i(-1);
444  for (auto const particle : *mcpartHandle){
445  mcpart_i++;
446  if (ShowerID == particle.TrackId()){
447  break;
448  }
449  }
450  const simb::MCParticle particle = mcpartHandle.product()->at(mcpart_i);
451  ShowerT0 = particle.T();
452  ShowerID = particle.TrackId();
453  ShowerTriggerType = 2; // Using MCTruth as trigger, so tigger type is 2.
454  T0col->push_back(anab::T0(ShowerT0,
456  ShowerID,
457  (*T0col).size()
458  ));
459 
460  auto diff = mcpart_i; // check I have a sensible value for this counter
461  if (diff >= (int)mcpartHandle->size()){
462  std::cout << "Error, the backtracker is doing weird things to your pointers!" << std::endl;
463  throw std::exception();
464  }
465  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
466  if (fMakeT0Assns){
467  util::CreateAssn(*this, evt, *T0col, showerlist[Shower], *Showerassn);
468  }
469  MCPartShowerassn->addSingle(showerlist[Shower], mcpartPtr, btdata);
470 
471  }// Loop over showers
472  }
473 
474 
475  if (pfparticleListHandle.isValid()){
476  //Access pfparticles and hits
477  art::FindManyP<recob::Cluster> fmcp(pfparticleListHandle, evt, fPFParticleModuleLabel);
478  //art::FindManyP<recob::Hit> fmtht(pfparticleListHandle, evt, fPfparticleModuleLabel);
479 
480  size_t NPfparticles = pfparticlelist.size();
481 
482  // Now to access MCTruth for each pfparticle...
483  for(size_t iPfp=0; iPfp < NPfparticles; ++iPfp) {
484  TrueTrackT0 = 0;
485  TrackID = 0;
486  TrueTrackID = 0;
488 
489  std::vector< art::Ptr<recob::Hit> > allHits;
490  //Get all hits through associated clusters
491  std::vector< art::Ptr<recob::Cluster>> allClusters = fmcp.at(iPfp);
492  art::FindManyP<recob::Hit> fmhcl(allClusters, evt, fPFParticleModuleLabel);
493  for (size_t iclu = 0; iclu<allClusters.size(); ++iclu){
494  std::vector< art::Ptr<recob::Hit>> hits = fmhcl.at(iclu);
495  allHits.insert(allHits.end(), hits.begin(), hits.end());
496  }
497 
498  std::map<int,double> trkide;
499  for(size_t h = 0; h < allHits.size(); ++h){
500  art::Ptr<recob::Hit> hit = allHits[h];
501  std::vector<sim::IDE> ides;
502  std::vector<sim::TrackIDE> TrackIDs = bt_serv->HitToTrackIDEs(hit);
503 
504  for(size_t e = 0; e < TrackIDs.size(); ++e){
505  trkide[TrackIDs[e].trackID] += TrackIDs[e].energy;
506  }
507  }
508  // Work out which IDE despoited the most charge in the hit if there was more than one.
509  double maxe = -1;
510  double tote = 0;
511  for (std::map<int,double>::iterator ii = trkide.begin(); ii!=trkide.end(); ++ii){
512  tote += ii->second;
513  if ((ii->second)>maxe){
514  maxe = ii->second;
515  TrackID = ii->first;
516  }
517  }
518  btdata.cleanliness = maxe/tote;
519 
520  // Now have trackID, so get PdG code and T0 etc.
521  const simb::MCParticle *tmpParticle = pi_serv->TrackIdToParticle_P(TrackID);
522  if (!tmpParticle) continue; // Retain this check that the BackTracker can find the right particle
523  // Now, loop through the MCParticle's myself to find the correct match
524  int mcpart_i(-1);
525  for (auto const particle : *mcpartHandle){
526  mcpart_i++;
527  if (TrackID == particle.TrackId()){
528  break;
529  }
530  }
531  const simb::MCParticle particle = mcpartHandle.product()->at(mcpart_i);
532  TrueTrackT0 = particle.T();
533  TrueTrackID = particle.TrackId();
534  TrueTriggerType = 2; // Using MCTruth as trigger, so tigger type is 2.
535  //std::cout << "Got particle, PDG = " << particle.PdgCode() << std::endl;
536 
537  //std::cout << "Filling T0col with " << TrueTrackT0 << " " << TrueTriggerType << " " << TrueTrackID << " " << (*T0col).size() << std::endl;
538 
539  T0col->push_back(anab::T0(TrueTrackT0,
541  TrueTrackID,
542  (*T0col).size()
543  ));
544  //auto diff = particle - firstParticle;
545  auto diff = mcpart_i; // check I have a sensible value for this counter
546  if (diff >= (int)mcpartHandle->size()){
547  std::cout << "Error, the backtracker is doing weird things to your pointers!" << std::endl;
548  throw std::exception();
549  }
550 // art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, particle - firstParticle);
551  art::Ptr<simb::MCParticle> mcpartPtr(mcpartHandle, mcpart_i);
552  //std::cout << "made MCParticle Ptr" << std::endl;
553 // const std::vertex< std::pair<double, std::string> > cleanlinessCompleteness;
554 // const std::pair<double, double> tmp(0.5, 0.5);
555 // MCPartassn->addSingle(tracklist[iPfp], mcpartPtr, tmp);
557  if (fMakeT0Assns){
558  util::CreateAssn(*this, evt, *T0col, pfparticlelist[iPfp], *PFParticleassn);
559  }
560  MCPartPFParticleassn->addSingle(pfparticlelist[iPfp], mcpartPtr, btdata);
561  }
562  fTree -> Fill();
563  } // Loop over tracks
564  }
565 
566  if (fMakeT0Assns){
567  evt.put(std::move(T0col));
568  evt.put(std::move(Trackassn));
569  evt.put(std::move(Showerassn));
571  evt.put(std::move(PFParticleassn));
572  }
573  }
574  evt.put(std::move(MCPartTrackassn));
575  evt.put(std::move(MCPartShowerassn));
577  evt.put(std::move(MCPartPFParticleassn));
578  }
579  if(fMakeHitAssns)
580  evt.put(std::move(MCPartHitassn));
581 } // Produce
582 
584 
code to link reconstructed objects back to the MC truth information
const simb::MCParticle * TrackIdToParticle_P(int const &id)
Float_t E
Definition: plot.C:23
MCTruthT0Matching & operator=(MCTruthT0Matching const &)=delete
intermediate_table::iterator iterator
const std::vector< sim::TrackIDE > HitToTrackIDEs(recob::Hit const &hit)
std::vector< TrackID > TrackIDs
Declaration of signal hit object.
TNtupleSim Fill(f1, f2, f3, f4)
Definition: T0.h:19
bool isRealData() const
Definition: Event.h:83
Particle class.
MCTruthT0Matching(fhicl::ParameterSet const &p)
int TrackId() const
Definition: MCParticle.h:214
void produce(art::Event &e) override
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
bool isValid() const
Definition: Handle.h:190
void hits()
Definition: readHits.C:15
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
Provides recob::Track data product.
T get(std::string const &key) const
Definition: ParameterSet.h:231
double T(const int i=0) const
Definition: MCParticle.h:228
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
Declaration of cluster object.
Detector simulation of raw signals on wires.
Encapsulate the geometry of a wire.
T * make(ARGS...args) const
Utility object to perform functions of association.
Encapsulate the construction of a single detector plane.
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
TCEvent evt
Definition: DataStructs.cxx:5
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:464
Float_t e
Definition: plot.C:34
void reconfigure(fhicl::ParameterSet const &p)
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33