LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
PFPAna_module.cc
Go to the documentation of this file.
1 //
3 // PFPAna class
4 //
5 // Bruce Baller
6 //
8 
9 #include <TH1F.h>
10 #include <TProfile.h>
11 #include <array>
12 #include <iomanip>
13 #include <string>
14 
15 // Framework includes
21 #include "art_root_io/TFileService.h"
25 #include "fhiclcpp/ParameterSet.h"
27 
28 // LArSoft includes
40 
41 namespace pfpf {
42 
43  class PFPAna : public art::EDAnalyzer {
44  public:
45  explicit PFPAna(fhicl::ParameterSet const& pset);
46 
47  private:
48  void analyze(const art::Event& evt) override;
49  void beginJob() override;
50 
51  TH1F* fNClusters;
53  // Cosmic Rays
54  TH1F* fCREP2;
55  TH1F* fCRE;
56  TH1F* fCRP;
57  // Neutrino interactions
58  TH1F* fNuKE_elec;
59  TH1F* fNuKE_muon;
60  TH1F* fNuKE_pion;
61  TH1F* fNuKE_kaon;
62  TH1F* fNuKE_prot;
63  TH1F* fNuEP2_elec;
64  TH1F* fNuEP2_muon;
65  TH1F* fNuEP2_pion;
66  TH1F* fNuEP2_kaon;
67  TH1F* fNuEP2_prot;
68  TH1F* fNuE_elec;
69  TH1F* fNuE_muon;
70  TH1F* fNuE_pion;
71  TH1F* fNuE_kaon;
72  TH1F* fNuE_prot;
73  TH1F* fNuP_elec;
74  TH1F* fNuP_muon;
75  TH1F* fNuP_pion;
76  TH1F* fNuP_kaon;
77  TH1F* fNuP_prot;
78 
79  TH1F* fNuVtx_dx;
80  TH1F* fNuVtx_dy;
81  TH1F* fNuVtx_dz;
82 
83  TProfile* fNuEP2_KE_elec;
84  TProfile* fNuEP2_KE_muon;
85  TProfile* fNuEP2_KE_pion;
86  TProfile* fNuEP2_KE_kaon;
87  TProfile* fNuEP2_KE_prot;
88 
89  std::string fHitsModuleLabel;
90  std::string fClusterModuleLabel;
91  std::string fTrackModuleLabel;
93  std::string fVertexModuleLabel;
94  std::vector<float> fElecKERange;
95  std::vector<float> fMuonKERange;
96  std::vector<float> fPionKERange;
97  std::vector<float> fKaonKERange;
98  std::vector<float> fProtKERange;
102  short fPrintLevel;
103 
104  }; // class PFPAna
105 
106  //--------------------------------------------------------------------
108  : EDAnalyzer(pset)
109  , fHitsModuleLabel(pset.get<std::string>("HitsModuleLabel"))
110  , fClusterModuleLabel(pset.get<std::string>("ClusterModuleLabel"))
111  , fTrackModuleLabel(pset.get<std::string>("TrackModuleLabel"))
112  , fPFParticleModuleLabel(pset.get<std::string>("PFParticleModuleLabel"))
113  , fVertexModuleLabel(pset.get<std::string>("VertexModuleLabel"))
114  , fElecKERange(pset.get<std::vector<float>>("ElecKERange"))
115  , fMuonKERange(pset.get<std::vector<float>>("MuonKERange"))
116  , fPionKERange(pset.get<std::vector<float>>("PionKERange"))
117  , fKaonKERange(pset.get<std::vector<float>>("KaonKERange"))
118  , fProtKERange(pset.get<std::vector<float>>("ProtKERange"))
119  , fTrackWeightOption(pset.get<short>("TrackWeightOption"))
120  , fMergeDaughters(pset.get<bool>("MergeDaughters"))
121  , fSkipCosmics(pset.get<bool>("SkipCosmics"))
122  , fPrintLevel(pset.get<short>("PrintLevel"))
123  {}
124 
125  //------------------------------------------------------------------
127  {
129 
130  fNClusters = tfs->make<TH1F>("fNoClustersInEvent", "Number of Clusters", 40, 0, 400);
131  fNHitInCluster = tfs->make<TH1F>("fNHitInCluster", "NHitInCluster", 100, 0, 100);
132 
133  if (!fSkipCosmics) {
134  // Cosmic ray histos
135  fCREP2 = tfs->make<TH1F>("CREP2", "CREP2", 50, 0, 1);
136  fCRE = tfs->make<TH1F>("CRE", "CR Efficiency", 50, 0, 1);
137  fCRP = tfs->make<TH1F>("CRP", "CR Purity", 50, 0, 1);
138  }
139  // Neutrino Int histos
140  fNuKE_elec = tfs->make<TH1F>("NuKE_elec", "NuKE electron", 100, 0, 4000);
141  fNuKE_muon = tfs->make<TH1F>("NuKE_muon", "NuKE muon", 100, 0, 4000);
142  fNuKE_pion = tfs->make<TH1F>("NuKE_pion", "NuKE pion", 100, 0, 4000);
143  fNuKE_kaon = tfs->make<TH1F>("NuKE_kaon", "NuKE kaon", 100, 0, 4000);
144  fNuKE_prot = tfs->make<TH1F>("NuKE_prot", "NuKE proton", 100, 0, 4000);
145 
146  fNuEP2_elec = tfs->make<TH1F>("NuEP2_elec", "NuEP2 electron", 50, 0, 1);
147  fNuEP2_muon = tfs->make<TH1F>("NuEP2_muon", "NuEP2 muon", 50, 0, 1);
148  fNuEP2_pion = tfs->make<TH1F>("NuEP2_pion", "NuEP2 pion", 50, 0, 1);
149  fNuEP2_kaon = tfs->make<TH1F>("NuEP2_kaon", "NuEP2 kaon", 50, 0, 1);
150  fNuEP2_prot = tfs->make<TH1F>("NuEP2_prot", "NuEP2 proton", 50, 0, 1);
151 
152  fNuE_elec = tfs->make<TH1F>("NuE_elec", "Nu Efficiency electron", 50, 0, 1);
153  fNuE_muon = tfs->make<TH1F>("NuE_muon", "Nu Efficiency muon", 50, 0, 1);
154  fNuE_pion = tfs->make<TH1F>("NuE_pion", "Nu Efficiency pion", 50, 0, 1);
155  fNuE_kaon = tfs->make<TH1F>("NuE_kaon", "Nu Efficiency kaon", 50, 0, 1);
156  fNuE_prot = tfs->make<TH1F>("NuE_prot", "Nu Efficiency proton", 50, 0, 1);
157 
158  fNuP_elec = tfs->make<TH1F>("NuP_elec", "Nu Purity electron", 50, 0, 1);
159  fNuP_muon = tfs->make<TH1F>("NuP_muon", "Nu Purity muon", 50, 0, 1);
160  fNuP_pion = tfs->make<TH1F>("NuP_pion", "Nu Purity pion", 50, 0, 1);
161  fNuP_kaon = tfs->make<TH1F>("NuP_kaon", "Nu Purity kaon", 50, 0, 1);
162  fNuP_prot = tfs->make<TH1F>("NuP_prot", "Nu Purity proton", 50, 0, 1);
163 
164  // True - Reco vertex difference
165  fNuVtx_dx = tfs->make<TH1F>("Vtx dx", "Vtx dx", 80, -10, 10);
166  fNuVtx_dy = tfs->make<TH1F>("Vtx dy", "Vtx dy", 80, -10, 10);
167  fNuVtx_dz = tfs->make<TH1F>("Vtx dz", "Vtx dz", 80, -10, 10);
168 
169  fNuEP2_KE_elec = tfs->make<TProfile>("NuEP2_KE_elec", "NuEP2 electron vs KE", 20, 0, 2000);
170  fNuEP2_KE_muon = tfs->make<TProfile>("NuEP2_KE_muon", "NuEP2 muon vs KE", 20, 0, 2000);
171  fNuEP2_KE_pion = tfs->make<TProfile>("NuEP2_KE_pion", "NuEP2 pion vs KE", 20, 0, 2000);
172  fNuEP2_KE_kaon = tfs->make<TProfile>("NuEP2_KE_kaon", "NuEP2 kaon vs KE", 20, 0, 2000);
173  fNuEP2_KE_prot = tfs->make<TProfile>("NuEP2_KE_prot", "NuEP2 proton vs KE", 20, 0, 2000);
174  }
175 
177  {
178  // code stolen from TrackAna_module.cc
180  if (geom->Nplanes() > 3) return;
181 
182  // get all hits in the event
184  evt.getByLabel(fHitsModuleLabel, hitListHandle);
185  std::vector<art::Ptr<recob::Hit>> allhits;
186  art::fill_ptr_vector(allhits, hitListHandle);
187  if (empty(allhits)) return;
188 
189  // get clusters and cluster-hit associations
190  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
191  evt.getByLabel(fClusterModuleLabel, clusterListHandle);
192  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
193  if (clusterListHandle->size() == 0) return;
194 
195  // get 3D vertices
196  art::Handle<std::vector<recob::Vertex>> vertexListHandle;
197  evt.getByLabel(fVertexModuleLabel, vertexListHandle);
198  art::PtrVector<recob::Vertex> recoVtxList;
199  double xyz[3] = {0, 0, 0};
200  for (unsigned int ii = 0; ii < vertexListHandle->size(); ++ii) {
201  art::Ptr<recob::Vertex> vertex(vertexListHandle, ii);
202  recoVtxList.push_back(vertex);
203  vertex->XYZ(xyz);
204  }
205 
206  // get PFParticles
208  evt.getByLabel(fPFParticleModuleLabel, PFPListHandle);
210  for (unsigned int ii = 0; ii < PFPListHandle->size(); ++ii) {
211  art::Ptr<recob::PFParticle> pfp(PFPListHandle, ii);
212  recoPFPList.push_back(pfp);
213  mf::LogVerbatim("PFPAna") << "PFParticle PDG " << pfp->PdgCode();
214  }
215 
216  // list of all true particles
219  sim::ParticleList const& plist = pi_serv->ParticleList();
220  // list of all true particles that will be considered
221  std::vector<const simb::MCParticle*> plist2;
222  // true (reconstructed) hits for each particle in plist2
223  std::vector<std::vector<art::Ptr<recob::Hit>>> hlist2;
224  // index of cluster matched to each particle in plist2 in each plane
225  std::vector<std::vector<short>> truToCl;
226  // number of true hits in each plane and cluster
227  std::vector<std::vector<unsigned short>> nTruHitInCl;
228  //number of reconstructed hits in all clusters
229  std::vector<unsigned short> nRecHitInCl;
230 
231  // calculate average EP2 for every event to facilitate code development
232  // Beam Neutrinos - muons and not-muons
233 
234  float aveNuEP2mu = 0.;
235  float numNuEP2mu = 0.;
236  float aveNuEP2nm = 0.;
237  float numNuEP2nm = 0.;
238  // Cosmic Rays
239  float aveCREP2 = 0.;
240  float numCREP2 = 0.;
241 
242  // track ID of the neutrino
243  int neutTrackID = -1;
244  std::vector<int> tidlist;
245  float neutEnergy = -1.;
246  int neutIntType = -1;
247  int neutCCNC = -1;
248 
249  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
250 
251  for (sim::ParticleList::const_iterator ipart = plist.begin(); ipart != plist.end(); ++ipart) {
252  const simb::MCParticle* part = (*ipart).second;
253  assert(part != 0);
254  int pdg = abs(part->PdgCode());
255  int trackID = part->TrackId();
256  art::Ptr<simb::MCTruth> theTruth = pi_serv->TrackIdToMCTruth_P(trackID);
257  if (fSkipCosmics && theTruth->Origin() == simb::kCosmicRay) continue;
258 
259  if (fPrintLevel > 3)
260  mf::LogVerbatim("PFPAna") << "Pre-Cuts origin " << theTruth->Origin() << " trackID "
261  << trackID << " PDG " << part->PdgCode() << " E " << part->E()
262  << " mass " << part->Mass() << " Mother "
263  << part->Mother() + neutTrackID << " Proc " << part->Process();
264 
265  // Get the neutrino track ID. Assume that there is only one neutrino
266  // interaction and it is first in the list of BeamNeutrino particles
267  if (theTruth->Origin() == simb::kBeamNeutrino && neutTrackID < 0) {
268  neutTrackID = trackID;
269  simb::MCNeutrino theNeutrino = theTruth->GetNeutrino();
270  neutEnergy = 1000. * theNeutrino.Nu().E();
271  neutIntType = theNeutrino.InteractionType();
272  neutCCNC = theNeutrino.CCNC();
273  for (unsigned short iv = 0; iv < recoVtxList.size(); ++iv) {
274  recoVtxList[iv]->XYZ(xyz);
275  fNuVtx_dx->Fill(part->Vx() - xyz[0]);
276  fNuVtx_dy->Fill(part->Vy() - xyz[1]);
277  fNuVtx_dz->Fill(part->Vz() - xyz[2]);
278  } // iv
279  } // theTruth->Origin() == simb::kBeamNeutrino && neutTrackID <
280 
281  bool isCharged = (pdg == 11) || (pdg == 13) || (pdg == 211) || (pdg == 321) || (pdg == 2212);
282 
283  if (!isCharged) continue;
284 
285  float KE = 1000 * (part->E() - part->Mass());
286  // KE (MeV) cuts
287  if (pdg == 11) {
288  if (fElecKERange[0] < 0) continue;
289  // only allow primary electrons
290  if (part->Process() != "primary") continue;
291  if (KE < fElecKERange[0] || KE > fElecKERange[1]) continue;
292  }
293  if (pdg == 13) {
294  if (fMuonKERange[0] < 0) continue;
295  if (KE < fMuonKERange[0] || KE > fMuonKERange[1]) continue;
296  }
297  if (pdg == 211) {
298  if (fPionKERange[0] < 0) continue;
299  if (KE < fPionKERange[0] || KE > fPionKERange[1]) continue;
300  }
301  if (pdg == 321) {
302  if (fKaonKERange[0] < 0) continue;
303  if (KE < fKaonKERange[0] || KE > fKaonKERange[1]) continue;
304  }
305  if (pdg == 2212) {
306  if (fProtKERange[0] < 0) continue;
307  if (KE < fProtKERange[0] || KE > fProtKERange[1]) continue;
308  }
309  // ignore secondaries from neutron interactions
310  if (part->Process() == "NeutronInelastic") continue;
311  plist2.push_back(part);
312  tidlist.push_back(trackID);
313  // initialize the true->(cluster,plane) association
314  std::vector<short> temp{-1, -1, -1};
315  truToCl.push_back(temp);
316  // initialize the true hit count
317  std::vector<unsigned short> temp2(3);
318  nTruHitInCl.push_back(temp2);
319 
320  if (fPrintLevel > 2)
321  mf::LogVerbatim("PFPAna") << plist2.size() - 1 << " Origin " << theTruth->Origin()
322  << " trackID " << trackID << " PDG " << part->PdgCode() << " KE "
323  << (int)KE << " Mother " << part->Mother() + neutTrackID
324  << " Proc " << part->Process();
325  }
326 
327  if (empty(plist2)) return;
328 
329  // get the hits (in all planes) that are matched to the true tracks
330  hlist2 = bt_serv->TrackIdsToHits_Ps(clockData, tidlist, allhits);
331  if (hlist2.size() != plist2.size()) {
332  mf::LogError("PFPAna") << "MC particle list size " << plist2.size()
333  << " != size of MC particle true hits lists " << hlist2.size();
334  return;
335  }
336  tidlist.clear();
337 
338  // vector of (mother, daughter) pairs
339  std::vector<std::pair<unsigned short, unsigned short>> moda;
340  // Deal with mother-daughter tracks
341  if (fMergeDaughters && neutTrackID >= 0) {
342  // Assume that daughters appear later in the list. Step backwards
343  // to accumulate all generations of daughters
344  for (unsigned short dpl = plist2.size() - 1; dpl > 0; --dpl) {
345  // no mother
346  if (plist2[dpl]->Mother() == 0) continue;
347  // electron
348  if (abs(plist2[dpl]->PdgCode()) == 11) continue;
349  // the actual mother trackID is offset from the neutrino trackID
350  int motherID = neutTrackID + plist2[dpl]->Mother() - 1;
351  // ensure that we are only looking at BeamNeutrino daughters
352  if (motherID < 0) continue;
353  // count the number of daughters
354  int ndtr = 0;
355  for (unsigned short kpl = 0; kpl < plist2.size(); ++kpl) {
356  if (plist2[kpl]->Mother() == motherID) ++ndtr;
357  }
358  // require only one daughter
359  if (ndtr > 1) continue;
360  // find the mother in the list
361  int mpl = -1;
362  for (unsigned short jpl = dpl - 1; jpl > 0; --jpl) {
363  if (plist2[jpl]->TrackId() == motherID) {
364  mpl = jpl;
365  break;
366  }
367  } // jpl
368  // mother not found for some reason
369  if (mpl < 0) continue;
370  // ensure that PDG code for mother and daughter are the same
371  if (plist2[dpl]->PdgCode() != plist2[mpl]->PdgCode()) continue;
372  moda.push_back(std::make_pair(mpl, dpl));
373  } // dpl
374  } // MergeDaughters
375 
376  // Now match reconstructed clusters to true particles.
378  for (unsigned int ii = 0; ii < clusterListHandle->size(); ++ii) {
379  art::Ptr<recob::Cluster> clusterHolder(clusterListHandle, ii);
380  clusters.push_back(clusterHolder);
381  }
382 
383  fNClusters->Fill(clusterListHandle->size());
384  nRecHitInCl.resize(clusters.size());
385 
386  // get the plane from the view. Perhaps there is a method that does
387  // this somewhere...
388  std::map<geo::View_t, unsigned int> ViewToPlane;
389  for (auto const& plane : geom->Iterate<geo::PlaneGeo>(geo::TPCID{0, 0})) {
390  geo::View_t view = plane.View();
391  ViewToPlane[view] = plane.ID().Plane;
392  }
393  for (size_t icl = 0; icl < clusters.size(); ++icl) {
394  unsigned int plane = ViewToPlane[clusters[icl]->View()];
395  std::vector<art::Ptr<recob::Hit>> cluhits = fmh.at(icl);
396  fNHitInCluster->Fill(cluhits.size());
397  nRecHitInCl[icl] = cluhits.size();
398  // count the number of hits matched to each true particle in plist2
399  std::vector<unsigned short> nHitInPl2(plist2.size());
400  for (size_t iht = 0; iht < cluhits.size(); ++iht) {
401  // look for this hit in all of the truth hit lists
402  short hitInPl2 = -1;
403  for (unsigned short ipl = 0; ipl < plist2.size(); ++ipl) {
404  unsigned short imat = 0;
405  for (imat = 0; imat < hlist2[ipl].size(); ++imat) {
406  if (cluhits[iht] == hlist2[ipl][imat]) break;
407  } // imat
408  if (imat < hlist2[ipl].size()) {
409  hitInPl2 = ipl;
410  break;
411  }
412  } // ipl
413  if (hitInPl2 < 0) continue;
414  // Assign the hit count to the mother if this is a daughter.
415  // Mother-daughter pairs are entered in the moda vector in reverse
416  // order, so assign daughter hits to the highest generation mother.
417  for (unsigned short imd = 0; imd < moda.size(); ++imd) {
418  if (moda[imd].second == hitInPl2) hitInPl2 = moda[imd].first;
419  }
420  // count
421  ++nHitInPl2[hitInPl2];
422  } // iht
423  // Associate the cluster with the truth particle that has the highest
424  // number of cluster hits
425  unsigned short nhit = 0;
426  short imtru = -1;
427  for (unsigned int ipl = 0; ipl < nHitInPl2.size(); ++ipl) {
428  if (nHitInPl2[ipl] > nhit) {
429  nhit = nHitInPl2[ipl];
430  imtru = ipl;
431  }
432  } // ipl
433  // make the cluster->(true,plane) association and save the
434  // number of true hits in the cluster
435  if (imtru != -1) {
436  // clobber a previously made association?
437  if (nhit > nTruHitInCl[imtru][plane]) {
438  truToCl[imtru][plane] = icl;
439  nTruHitInCl[imtru][plane] = nhit;
440  }
441  } // imtru != 1
442  } // icl
443 
444  // ready to calculate Efficiency, Purity in each plane and EP2
445  for (unsigned short ipl = 0; ipl < plist2.size(); ++ipl) {
446  // ignore daughters
447  bool skipit = false;
448  for (unsigned short ii = 0; ii < moda.size(); ++ii) {
449  if (moda[ii].second == ipl) {
450  skipit = true;
451  break;
452  }
453  } // ii
454  if (skipit) continue;
455  // ignore true particles with few true hits. Outside the detector
456  // or not reconstructable?
457  if (hlist2[ipl].size() < 3) continue;
458 
459  int trackID = plist2[ipl]->TrackId();
460  art::Ptr<simb::MCTruth> theTruth = pi_serv->TrackIdToMCTruth_P(trackID);
461  bool isCosmic = (theTruth->Origin() == simb::kCosmicRay);
462  float KE = 1000 * (plist2[ipl]->E() - plist2[ipl]->Mass());
463  int PDG = abs(plist2[ipl]->PdgCode());
464 
465  std::vector<short> nTru(geom->Nplanes());
466  std::vector<short> nRec(geom->Nplanes());
467  std::vector<short> nTruRec(geom->Nplanes());
468  std::vector<float> eff(geom->Nplanes());
469  std::vector<float> pur(geom->Nplanes());
470  std::vector<float> ep(geom->Nplanes());
471  for (unsigned int plane = 0; plane < geom->Nplanes(); ++plane) {
472  // count the number of true hits in this plane for the true particle.
473  // First count the mother hits
474  for (unsigned short ii = 0; ii < hlist2[ipl].size(); ++ii) {
475  if (ViewToPlane[hlist2[ipl][ii]->View()] == plane) ++nTru[plane];
476  } // ii
477  // next look for daughters and count those hits in all generations
478  unsigned short mom = ipl;
479  std::vector<std::pair<unsigned short, unsigned short>>::reverse_iterator rit =
480  moda.rbegin();
481  while (rit != moda.rend()) {
482  if ((*rit).first == mom) {
483  unsigned short dau = (*rit).second;
484  for (unsigned short jj = 0; jj < hlist2[dau].size(); ++jj) {
485  if (ViewToPlane[hlist2[dau][jj]->View()] == plane) ++nTru[plane];
486  } // jj
487  // It is likely that one hit appears in the mother list
488  // as well as the daughter list, so subtract one from the count
489  --nTru[plane];
490  mom = (*rit).second;
491  } // (*rit).first == mom
492  ++rit;
493  } // rit
494  if (nTru[plane] == 0) { continue; }
495  short icl = truToCl[ipl][plane];
496  nRec[plane] = nRecHitInCl[icl];
497  nTruRec[plane] = nTruHitInCl[ipl][plane];
498  if (nTru[plane] > 0) eff[plane] = (float)nTruRec[plane] / (float)nTru[plane];
499  if (nRec[plane] > 0) pur[plane] = (float)nTruRec[plane] / (float)nRec[plane];
500  ep[plane] = eff[plane] * pur[plane];
501  } // plane
502  // sort the ep values in ascending order
503  std::vector<float> temp;
504  temp = ep;
505  std::sort(temp.begin(), temp.end());
506  // EP2 is the second highest value
507  unsigned short ii = temp.size() - 2;
508  float ep2 = temp[ii];
509  // find the plane that defined EP2
510  short ep2Plane = 0;
511  short ep2Cluster = 0;
512  for (unsigned short jj = 0; jj < temp.size(); ++jj) {
513  if (ep[jj] == ep2) {
514  ep2Plane = jj;
515  ep2Cluster = truToCl[ipl][ep2Plane];
516  break;
517  }
518  } // jj
519  // find the US and DS ends of the cluster for printing
520  std::array<double, 2> clBeg, clEnd;
521  if (ep2Cluster >= 0) {
522  clBeg[0] = clusters[ep2Cluster]->StartWire();
523  clBeg[1] = clusters[ep2Cluster]->StartTick();
524  clEnd[0] = clusters[ep2Cluster]->EndWire();
525  clEnd[1] = clusters[ep2Cluster]->EndTick();
526  }
527  else {
528  clBeg.fill(0.);
529  clEnd.fill(0.);
530  }
531  // fill histograms
532  if (isCosmic) {
533  fCREP2->Fill(ep2);
534  fCRE->Fill(eff[ep2Plane]);
535  fCRP->Fill(pur[ep2Plane]);
536  aveCREP2 += ep2;
537  numCREP2 += 1.;
538  if (fPrintLevel > 1)
539  mf::LogVerbatim("PFPAna")
540  << ">>>CREP2 " << std::fixed << std::setprecision(2) << ep2 << " E " << eff[ep2Plane]
541  << std::setprecision(2) << " P " << pur[ep2Plane] << " P:W:T " << ep2Plane << ":"
542  << (int)clBeg[0] << ":" << (int)clBeg[1] << "-" << ep2Plane << ":" << (int)clEnd[0]
543  << ":" << (int)clEnd[1] << " PDG " << PDG << " KE " << (int)KE << " MeV";
544  } // isCosmic
545  else {
546  float wght = 1.;
547  if (fTrackWeightOption == 1) wght = KE;
548  // accumulate statistics for muons and not-muons
549  if (PDG == 13) {
550  aveNuEP2mu += ep2 * wght;
551  numNuEP2mu += wght;
552  }
553  else {
554  aveNuEP2nm += ep2 * wght;
555  numNuEP2nm += wght;
556  }
557  if (PDG == 11) {
558  fNuKE_elec->Fill(KE, wght);
559  fNuE_elec->Fill(eff[ep2Plane], wght);
560  fNuP_elec->Fill(pur[ep2Plane], wght);
561  fNuEP2_elec->Fill(ep2, wght);
562  fNuEP2_KE_elec->Fill(KE, ep2, wght);
563  }
564  else if (PDG == 13) {
565  fNuKE_muon->Fill(KE, wght);
566  fNuE_muon->Fill(eff[ep2Plane], wght);
567  fNuP_muon->Fill(pur[ep2Plane], wght);
568  fNuEP2_muon->Fill(ep2, wght);
569  fNuEP2_KE_muon->Fill(KE, ep2, wght);
570  }
571  else if (PDG == 211) {
572  fNuKE_pion->Fill(KE, wght);
573  fNuE_pion->Fill(eff[ep2Plane], wght);
574  fNuP_pion->Fill(pur[ep2Plane], wght);
575  fNuEP2_pion->Fill(ep2, wght);
576  fNuEP2_KE_pion->Fill(KE, ep2, wght);
577  }
578  else if (PDG == 321) {
579  fNuKE_kaon->Fill(KE, wght);
580  fNuE_kaon->Fill(eff[ep2Plane], wght);
581  fNuP_kaon->Fill(pur[ep2Plane], wght);
582  fNuEP2_kaon->Fill(ep2, wght);
583  fNuEP2_KE_kaon->Fill(KE, ep2, wght);
584  }
585  else if (PDG == 2212) {
586  fNuKE_prot->Fill(KE, wght);
587  fNuE_prot->Fill(eff[ep2Plane], wght);
588  fNuP_prot->Fill(pur[ep2Plane], wght);
589  fNuEP2_prot->Fill(ep2, wght);
590  fNuEP2_KE_prot->Fill(KE, ep2, wght);
591  }
592  if (fPrintLevel > 1)
593  mf::LogVerbatim("PFPAna")
594  << ">>>NuEP2 " << std::fixed << std::setprecision(2) << ep2 << " E " << eff[ep2Plane]
595  << std::setprecision(2) << " P " << pur[ep2Plane] << " P:W:T " << ep2Plane << ":"
596  << (int)clBeg[0] << ":" << (int)clBeg[1] << "-" << ep2Plane << ":" << (int)clEnd[0]
597  << ":" << (int)clEnd[1] << " PDG " << PDG << " KE " << (int)KE << " MeV ";
598  if (fPrintLevel > 2) {
599  // print out the begin/end true hits
600  mf::LogVerbatim mfp("PFPAna");
601  mfp << " Truth P:W:T ";
602  for (unsigned int plane = 0; plane < geom->Nplanes(); ++plane) {
603  unsigned short loW = 9999;
604  unsigned short loT = 0;
605  unsigned short hiW = 0;
606  unsigned short hiT = 0;
607  for (unsigned short ii = 0; ii < hlist2[ipl].size(); ++ii) {
608  if (ViewToPlane[hlist2[ipl][ii]->View()] == plane) {
609  art::Ptr<recob::Hit> theHit = hlist2[ipl][ii];
610  if (theHit->WireID().Wire < loW) {
611  loW = theHit->WireID().Wire;
612  loT = theHit->PeakTime();
613  }
614  if (theHit->WireID().Wire > hiW) {
615  hiW = theHit->WireID().Wire;
616  hiT = theHit->PeakTime();
617  }
618  } // correct view
619  } // ii
620  mfp << plane << ":" << loW << ":" << loT << "-" << plane << ":" << hiW << ":" << hiT
621  << " ";
622  } // plane
623  } // fPrintLevel > 2
624  } // !isCosmic
625  } // ipl
626 
627  float ave1 = -1.;
628  if (numNuEP2mu > 0.) ave1 = aveNuEP2mu / numNuEP2mu;
629 
630  float ave2 = -1.;
631  if (numNuEP2nm > 0.) ave2 = aveNuEP2nm / numNuEP2nm;
632 
633  float ave3 = -1.;
634  if (numCREP2 > 0.) ave3 = aveCREP2 / numCREP2;
635 
636  if (fPrintLevel > 0) {
637  std::string nuType = "Other";
638  if (neutCCNC == simb::kCC) {
639  if (neutIntType == 1001) nuType = "CCQE";
640  if (neutIntType == 1091) nuType = "DIS";
641  if (neutIntType == 1097) nuType = "COH";
642  if (neutIntType > 1002 && neutIntType < 1091) nuType = "RES";
643  }
644  else if (neutCCNC == simb::kNC) {
645  nuType = "NC";
646  }
647  else {
648  nuType = "Unknown";
649  }
650  mf::LogVerbatim("PFPAna") << "EvtEP2 " << evt.id().event() << " NuType " << nuType << " Enu "
651  << std::fixed << std::setprecision(0) << neutEnergy << std::right
652  << std::fixed << std::setprecision(2) << " NuMuons " << ave1
653  << " NuPiKp " << ave2 << " CosmicRays " << ave3 << " CCNC "
654  << neutCCNC << " IntType " << neutIntType;
655  }
656  } // analyze
657 
659 
660 } // end namespace
double E(const int i=0) const
Definition: MCParticle.h:234
void XYZ(double *xyz) const
Legacy method to access vertex position, preserved to avoid breaking code. Please try to use Vertex::...
Definition: Vertex.cxx:34
details::range_type< T > Iterate() const
Initializes the specified ID with the ID of the first cryostat.
Definition: GeometryCore.h:541
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
std::vector< float > fElecKERange
int PdgCode() const
Definition: MCParticle.h:213
int CCNC() const
Definition: MCNeutrino.h:148
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:77
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:102
Float_t mfp
Definition: plot.C:20
TH1F * fNuEP2_prot
Int_t ipart
Definition: Style.C:10
std::string fVertexModuleLabel
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
int Mother() const
Definition: MCParticle.h:214
Declaration of signal hit object.
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
list_type::const_iterator const_iterator
Definition: ParticleList.h:132
simb::Origin_t Origin() const
Definition: MCTruth.h:74
double Mass() const
Definition: MCParticle.h:240
TH1F * fNuP_kaon
std::string fHitsModuleLabel
TH1F * fNuKE_prot
std::vector< float > fProtKERange
constexpr auto abs(T v)
Returns the absolute value of the argument.
STL namespace.
int PdgCode() const
Return the type of particle as a PDG ID.
Definition: PFParticle.h:79
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
std::string Process() const
Definition: MCParticle.h:216
std::string fClusterModuleLabel
int TrackId() const
Definition: MCParticle.h:211
std::vector< float > fPionKERange
TProfile * fNuEP2_KE_kaon
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:101
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
TH1F * fNuVtx_dy
TH1F * fNuKE_elec
int InteractionType() const
Definition: MCNeutrino.h:150
TString part[npart]
Definition: Style.C:32
TProfile * fNuEP2_KE_elec
TH1F * fNuVtx_dz
TH1F * fNuEP2_muon
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
TH1F * fNuE_kaon
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
void beginJob() override
TH1F * fNuP_muon
short fTrackWeightOption
const art::Ptr< simb::MCTruth > & TrackIdToMCTruth_P(int id) const
TProfile * fNuEP2_KE_pion
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:78
TH1F * fNuEP2_kaon
iterator begin()
Definition: ParticleList.h:305
TH1F * fNuP_elec
std::string fPFParticleModuleLabel
TH1F * fNHitInCluster
std::vector< std::vector< art::Ptr< recob::Hit > > > TrackIdsToHits_Ps(detinfo::DetectorClocksData const &clockData, std::vector< int > const &tkIds, std::vector< art::Ptr< recob::Hit >> const &hitsIn) const
TH1F * fNuKE_muon
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
Declaration of cluster object.
size_type size() const
Definition: PtrVector.h:302
TH1F * fNuP_prot
const sim::ParticleList & ParticleList() const
TH1F * fNuVtx_dx
TH1F * fNuE_prot
TProfile * fNuEP2_KE_prot
double Vx(const int i=0) const
Definition: MCParticle.h:222
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
TH1F * fNuE_muon
TH1F * fNClusters
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Encapsulate the construction of a single detector plane.
TProfile * fNuEP2_KE_muon
TH1F * fNuE_pion
TH1F * fNuE_elec
PFPAna(fhicl::ParameterSet const &pset)
TH1F * fNuEP2_elec
TH1F * fNuKE_pion
double Vz(const int i=0) const
Definition: MCParticle.h:224
TH1F * fNuP_pion
EventNumber_t event() const
Definition: EventID.h:116
std::vector< float > fMuonKERange
unsigned int Nplanes(TPCID const &tpcid=tpc_zero) const
Returns the total number of planes in the specified TPC.
Definition: GeometryCore.h:977
std::string fTrackModuleLabel
TH1F * fNuKE_kaon
TCEvent evt
Definition: DataStructs.cxx:8
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
Particle list in DetSim contains Monte Carlo particle information.
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:82
Event generator information.
Definition: MCNeutrino.h:18
Definition: fwd.h:26
void analyze(const art::Event &evt) override
EventID id() const
Definition: Event.cc:23
double Vy(const int i=0) const
Definition: MCParticle.h:223
art framework interface to geometry description
std::vector< float > fKaonKERange
decltype(auto) constexpr empty(T &&obj)
ADL-aware version of std::empty.
Definition: StdUtils.h:109
Cosmic rays.
Definition: MCTruth.h:24
TH1F * fNuEP2_pion
Beam neutrinos.
Definition: MCTruth.h:23
vertex reconstruction