All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PFPAna_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // PFPAna class
4 //
5 // Bruce Baller
6 //
7 ////////////////////////////////////////////////////////////////////////
8 
9 #include <TH1F.h>
10 #include <TProfile.h>
11 #include <array>
12 #include <iomanip>
13 #include <string>
14 
15 // Framework includes
16 #include "art/Framework/Core/EDAnalyzer.h"
17 #include "art/Framework/Core/ModuleMacros.h"
18 #include "art/Framework/Principal/Event.h"
19 #include "art/Framework/Principal/Handle.h"
20 #include "art/Framework/Services/Registry/ServiceHandle.h"
21 #include "art_root_io/TFileService.h"
22 #include "canvas/Persistency/Common/FindManyP.h"
23 #include "canvas/Persistency/Common/Ptr.h"
24 #include "canvas/Persistency/Common/PtrVector.h"
25 #include "fhiclcpp/ParameterSet.h"
26 #include "messagefacility/MessageLogger/MessageLogger.h"
27 
28 // LArSoft includes
38 #include "nug4/ParticleNavigation/ParticleList.h"
39 #include "nusimdata/SimulationBase/MCTruth.h"
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  //--------------------------------------------------------------------
107  PFPAna::PFPAna(fhicl::ParameterSet const& pset)
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  //------------------------------------------------------------------
126  void
128  {
129  art::ServiceHandle<art::TFileService const> tfs;
130 
131  fNClusters = tfs->make<TH1F>("fNoClustersInEvent", "Number of Clusters", 40, 0, 400);
132  fNHitInCluster = tfs->make<TH1F>("fNHitInCluster", "NHitInCluster", 100, 0, 100);
133 
134  if (!fSkipCosmics) {
135  // Cosmic ray histos
136  fCREP2 = tfs->make<TH1F>("CREP2", "CREP2", 50, 0, 1);
137  fCRE = tfs->make<TH1F>("CRE", "CR Efficiency", 50, 0, 1);
138  fCRP = tfs->make<TH1F>("CRP", "CR Purity", 50, 0, 1);
139  }
140  // Neutrino Int histos
141  fNuKE_elec = tfs->make<TH1F>("NuKE_elec", "NuKE electron", 100, 0, 4000);
142  fNuKE_muon = tfs->make<TH1F>("NuKE_muon", "NuKE muon", 100, 0, 4000);
143  fNuKE_pion = tfs->make<TH1F>("NuKE_pion", "NuKE pion", 100, 0, 4000);
144  fNuKE_kaon = tfs->make<TH1F>("NuKE_kaon", "NuKE kaon", 100, 0, 4000);
145  fNuKE_prot = tfs->make<TH1F>("NuKE_prot", "NuKE proton", 100, 0, 4000);
146 
147  fNuEP2_elec = tfs->make<TH1F>("NuEP2_elec", "NuEP2 electron", 50, 0, 1);
148  fNuEP2_muon = tfs->make<TH1F>("NuEP2_muon", "NuEP2 muon", 50, 0, 1);
149  fNuEP2_pion = tfs->make<TH1F>("NuEP2_pion", "NuEP2 pion", 50, 0, 1);
150  fNuEP2_kaon = tfs->make<TH1F>("NuEP2_kaon", "NuEP2 kaon", 50, 0, 1);
151  fNuEP2_prot = tfs->make<TH1F>("NuEP2_prot", "NuEP2 proton", 50, 0, 1);
152 
153  fNuE_elec = tfs->make<TH1F>("NuE_elec", "Nu Efficiency electron", 50, 0, 1);
154  fNuE_muon = tfs->make<TH1F>("NuE_muon", "Nu Efficiency muon", 50, 0, 1);
155  fNuE_pion = tfs->make<TH1F>("NuE_pion", "Nu Efficiency pion", 50, 0, 1);
156  fNuE_kaon = tfs->make<TH1F>("NuE_kaon", "Nu Efficiency kaon", 50, 0, 1);
157  fNuE_prot = tfs->make<TH1F>("NuE_prot", "Nu Efficiency proton", 50, 0, 1);
158 
159  fNuP_elec = tfs->make<TH1F>("NuP_elec", "Nu Purity electron", 50, 0, 1);
160  fNuP_muon = tfs->make<TH1F>("NuP_muon", "Nu Purity muon", 50, 0, 1);
161  fNuP_pion = tfs->make<TH1F>("NuP_pion", "Nu Purity pion", 50, 0, 1);
162  fNuP_kaon = tfs->make<TH1F>("NuP_kaon", "Nu Purity kaon", 50, 0, 1);
163  fNuP_prot = tfs->make<TH1F>("NuP_prot", "Nu Purity proton", 50, 0, 1);
164 
165  // True - Reco vertex difference
166  fNuVtx_dx = tfs->make<TH1F>("Vtx dx", "Vtx dx", 80, -10, 10);
167  fNuVtx_dy = tfs->make<TH1F>("Vtx dy", "Vtx dy", 80, -10, 10);
168  fNuVtx_dz = tfs->make<TH1F>("Vtx dz", "Vtx dz", 80, -10, 10);
169 
170  fNuEP2_KE_elec = tfs->make<TProfile>("NuEP2_KE_elec", "NuEP2 electron vs KE", 20, 0, 2000);
171  fNuEP2_KE_muon = tfs->make<TProfile>("NuEP2_KE_muon", "NuEP2 muon vs KE", 20, 0, 2000);
172  fNuEP2_KE_pion = tfs->make<TProfile>("NuEP2_KE_pion", "NuEP2 pion vs KE", 20, 0, 2000);
173  fNuEP2_KE_kaon = tfs->make<TProfile>("NuEP2_KE_kaon", "NuEP2 kaon vs KE", 20, 0, 2000);
174  fNuEP2_KE_prot = tfs->make<TProfile>("NuEP2_KE_prot", "NuEP2 proton vs KE", 20, 0, 2000);
175  }
176 
177  void
178  PFPAna::analyze(const art::Event& evt)
179  {
180  // code stolen from TrackAna_module.cc
181  art::ServiceHandle<geo::Geometry const> geom;
182  if (geom->Nplanes() > 3) return;
183 
184  // get all hits in the event
185  art::Handle<std::vector<recob::Hit>> hitListHandle;
186  evt.getByLabel(fHitsModuleLabel, hitListHandle);
187  std::vector<art::Ptr<recob::Hit>> allhits;
188  art::fill_ptr_vector(allhits, hitListHandle);
189  if (empty(allhits)) return;
190 
191  // get clusters and cluster-hit associations
192  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
193  evt.getByLabel(fClusterModuleLabel, clusterListHandle);
194  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
195  if (clusterListHandle->size() == 0) return;
196 
197  // get 3D vertices
198  art::Handle<std::vector<recob::Vertex>> vertexListHandle;
199  evt.getByLabel(fVertexModuleLabel, vertexListHandle);
200  art::PtrVector<recob::Vertex> recoVtxList;
201  double xyz[3] = {0, 0, 0};
202  for (unsigned int ii = 0; ii < vertexListHandle->size(); ++ii) {
203  art::Ptr<recob::Vertex> vertex(vertexListHandle, ii);
204  recoVtxList.push_back(vertex);
205  vertex->XYZ(xyz);
206  }
207 
208  // get PFParticles
209  art::Handle<std::vector<recob::PFParticle>> PFPListHandle;
210  evt.getByLabel(fPFParticleModuleLabel, PFPListHandle);
211  art::PtrVector<recob::PFParticle> recoPFPList;
212  for (unsigned int ii = 0; ii < PFPListHandle->size(); ++ii) {
213  art::Ptr<recob::PFParticle> pfp(PFPListHandle, ii);
214  recoPFPList.push_back(pfp);
215  mf::LogVerbatim("PFPAna") << "PFParticle PDG " << pfp->PdgCode();
216  }
217 
218  // list of all true particles
219  art::ServiceHandle<cheat::BackTrackerService const> bt_serv;
220  art::ServiceHandle<cheat::ParticleInventoryService const> pi_serv;
221  sim::ParticleList const& plist = pi_serv->ParticleList();
222  // list of all true particles that will be considered
223  std::vector<const simb::MCParticle*> plist2;
224  // true (reconstructed) hits for each particle in plist2
225  std::vector<std::vector<art::Ptr<recob::Hit>>> hlist2;
226  // index of cluster matched to each particle in plist2 in each plane
227  std::vector<std::vector<short>> truToCl;
228  // number of true hits in each plane and cluster
229  std::vector<std::vector<unsigned short>> nTruHitInCl;
230  //number of reconstructed hits in all clusters
231  std::vector<unsigned short> nRecHitInCl;
232 
233  // calculate average EP2 for every event to facilitate code development
234  // Beam Neutrinos - muons and not-muons
235 
236  float aveNuEP2mu = 0.;
237  float numNuEP2mu = 0.;
238  float aveNuEP2nm = 0.;
239  float numNuEP2nm = 0.;
240  // Cosmic Rays
241  float aveCREP2 = 0.;
242  float numCREP2 = 0.;
243 
244  // track ID of the neutrino
245  int neutTrackID = -1;
246  std::vector<int> tidlist;
247  float neutEnergy = -1.;
248  int neutIntType = -1;
249  int neutCCNC = -1;
250 
251  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
252 
253  for (sim::ParticleList::const_iterator ipart = plist.begin(); ipart != plist.end(); ++ipart) {
254  const simb::MCParticle* part = (*ipart).second;
255  assert(part != 0);
256  int pdg = abs(part->PdgCode());
257  int trackID = part->TrackId();
258  art::Ptr<simb::MCTruth> theTruth = pi_serv->TrackIdToMCTruth_P(trackID);
259  if (fSkipCosmics && theTruth->Origin() == simb::kCosmicRay) continue;
260 
261  if (fPrintLevel > 3)
262  mf::LogVerbatim("PFPAna") << "Pre-Cuts origin " << theTruth->Origin() << " trackID "
263  << trackID << " PDG " << part->PdgCode() << " E " << part->E()
264  << " mass " << part->Mass() << " Mother "
265  << part->Mother() + neutTrackID << " Proc " << part->Process();
266 
267  // Get the neutrino track ID. Assume that there is only one neutrino
268  // interaction and it is first in the list of BeamNeutrino particles
269  if (theTruth->Origin() == simb::kBeamNeutrino && neutTrackID < 0) {
270  neutTrackID = trackID;
271  simb::MCNeutrino theNeutrino = theTruth->GetNeutrino();
272  neutEnergy = 1000. * theNeutrino.Nu().E();
273  neutIntType = theNeutrino.InteractionType();
274  neutCCNC = theNeutrino.CCNC();
275  for (unsigned short iv = 0; iv < recoVtxList.size(); ++iv) {
276  recoVtxList[iv]->XYZ(xyz);
277  fNuVtx_dx->Fill(part->Vx() - xyz[0]);
278  fNuVtx_dy->Fill(part->Vy() - xyz[1]);
279  fNuVtx_dz->Fill(part->Vz() - xyz[2]);
280  } // iv
281  } // theTruth->Origin() == simb::kBeamNeutrino && neutTrackID <
282 
283  bool isCharged = (pdg == 11) || (pdg == 13) || (pdg == 211) || (pdg == 321) || (pdg == 2212);
284 
285  if (!isCharged) continue;
286 
287  float KE = 1000 * (part->E() - part->Mass());
288  // KE (MeV) cuts
289  if (pdg == 11) {
290  if (fElecKERange[0] < 0) continue;
291  // only allow primary electrons
292  if (part->Process() != "primary") continue;
293  if (KE < fElecKERange[0] || KE > fElecKERange[1]) continue;
294  }
295  if (pdg == 13) {
296  if (fMuonKERange[0] < 0) continue;
297  if (KE < fMuonKERange[0] || KE > fMuonKERange[1]) continue;
298  }
299  if (pdg == 211) {
300  if (fPionKERange[0] < 0) continue;
301  if (KE < fPionKERange[0] || KE > fPionKERange[1]) continue;
302  }
303  if (pdg == 321) {
304  if (fKaonKERange[0] < 0) continue;
305  if (KE < fKaonKERange[0] || KE > fKaonKERange[1]) continue;
306  }
307  if (pdg == 2212) {
308  if (fProtKERange[0] < 0) continue;
309  if (KE < fProtKERange[0] || KE > fProtKERange[1]) continue;
310  }
311  // ignore secondaries from neutron interactions
312  if (part->Process() == "NeutronInelastic") continue;
313  plist2.push_back(part);
314  tidlist.push_back(trackID);
315  // initialize the true->(cluster,plane) association
316  std::vector<short> temp{-1, -1, -1};
317  truToCl.push_back(temp);
318  // initialize the true hit count
319  std::vector<unsigned short> temp2(3);
320  nTruHitInCl.push_back(temp2);
321 
322  if (fPrintLevel > 2)
323  mf::LogVerbatim("PFPAna") << plist2.size() - 1 << " Origin " << theTruth->Origin()
324  << " trackID " << trackID << " PDG " << part->PdgCode() << " KE "
325  << (int)KE << " Mother " << part->Mother() + neutTrackID
326  << " Proc " << part->Process();
327  }
328 
329  if (empty(plist2)) return;
330 
331  // get the hits (in all planes) that are matched to the true tracks
332  hlist2 = bt_serv->TrackIdsToHits_Ps(clockData, tidlist, allhits);
333  if (hlist2.size() != plist2.size()) {
334  mf::LogError("PFPAna") << "MC particle list size " << plist2.size()
335  << " != size of MC particle true hits lists " << hlist2.size();
336  return;
337  }
338  tidlist.clear();
339 
340  // vector of (mother, daughter) pairs
341  std::vector<std::pair<unsigned short, unsigned short>> moda;
342  // Deal with mother-daughter tracks
343  if (fMergeDaughters && neutTrackID >= 0) {
344  // Assume that daughters appear later in the list. Step backwards
345  // to accumulate all generations of daughters
346  for (unsigned short dpl = plist2.size() - 1; dpl > 0; --dpl) {
347  // no mother
348  if (plist2[dpl]->Mother() == 0) continue;
349  // electron
350  if (abs(plist2[dpl]->PdgCode()) == 11) continue;
351  // the actual mother trackID is offset from the neutrino trackID
352  int motherID = neutTrackID + plist2[dpl]->Mother() - 1;
353  // ensure that we are only looking at BeamNeutrino daughters
354  if (motherID < 0) continue;
355  // count the number of daughters
356  int ndtr = 0;
357  for (unsigned short kpl = 0; kpl < plist2.size(); ++kpl) {
358  if (plist2[kpl]->Mother() == motherID) ++ndtr;
359  }
360  // require only one daughter
361  if (ndtr > 1) continue;
362  // find the mother in the list
363  int mpl = -1;
364  for (unsigned short jpl = dpl - 1; jpl > 0; --jpl) {
365  if (plist2[jpl]->TrackId() == motherID) {
366  mpl = jpl;
367  break;
368  }
369  } // jpl
370  // mother not found for some reason
371  if (mpl < 0) continue;
372  // ensure that PDG code for mother and daughter are the same
373  if (plist2[dpl]->PdgCode() != plist2[mpl]->PdgCode()) continue;
374  moda.push_back(std::make_pair(mpl, dpl));
375  } // dpl
376  } // MergeDaughters
377 
378  // Now match reconstructed clusters to true particles.
379  art::PtrVector<recob::Cluster> clusters;
380  for (unsigned int ii = 0; ii < clusterListHandle->size(); ++ii) {
381  art::Ptr<recob::Cluster> clusterHolder(clusterListHandle, ii);
382  clusters.push_back(clusterHolder);
383  }
384 
385  fNClusters->Fill(clusterListHandle->size());
386  nRecHitInCl.resize(clusters.size());
387 
388  // get the plane from the view. Perhaps there is a method that does
389  // this somewhere...
390  std::map<geo::View_t, unsigned int> ViewToPlane;
391  for (unsigned int plane = 0; plane < geom->Nplanes(); ++plane) {
392  geo::View_t view = geom->Plane(plane).View();
393  ViewToPlane[view] = plane;
394  }
395  for (size_t icl = 0; icl < clusters.size(); ++icl) {
396  unsigned int plane = ViewToPlane[clusters[icl]->View()];
397  std::vector<art::Ptr<recob::Hit>> cluhits = fmh.at(icl);
398  fNHitInCluster->Fill(cluhits.size());
399  nRecHitInCl[icl] = cluhits.size();
400  // count the number of hits matched to each true particle in plist2
401  std::vector<unsigned short> nHitInPl2(plist2.size());
402  for (size_t iht = 0; iht < cluhits.size(); ++iht) {
403  /*
404  mf::LogVerbatim("PFPAna")
405  <<"Clus Hit "<<cluhits[iht]->View()
406  <<":"<<cluhits[iht]->WireID().Wire
407  <<":"<<(int)cluhits[iht]->PeakTime();
408 */
409  // look for this hit in all of the truth hit lists
410  short hitInPl2 = -1;
411  for (unsigned short ipl = 0; ipl < plist2.size(); ++ipl) {
412  unsigned short imat = 0;
413  for (imat = 0; imat < hlist2[ipl].size(); ++imat) {
414  if (cluhits[iht] == hlist2[ipl][imat]) break;
415  } // imat
416  if (imat < hlist2[ipl].size()) {
417  hitInPl2 = ipl;
418  break;
419  }
420  } // ipl
421  if (hitInPl2 < 0) continue;
422  // Assign the hit count to the mother if this is a daughter.
423  // Mother-daughter pairs are entered in the moda vector in reverse
424  // order, so assign daughter hits to the highest generation mother.
425  for (unsigned short imd = 0; imd < moda.size(); ++imd) {
426  if (moda[imd].second == hitInPl2) hitInPl2 = moda[imd].first;
427  }
428  // count
429  ++nHitInPl2[hitInPl2];
430  } // iht
431  // Associate the cluster with the truth particle that has the highest
432  // number of cluster hits
433  unsigned short nhit = 0;
434  short imtru = -1;
435  for (unsigned int ipl = 0; ipl < nHitInPl2.size(); ++ipl) {
436  if (nHitInPl2[ipl] > nhit) {
437  nhit = nHitInPl2[ipl];
438  imtru = ipl;
439  }
440  } // ipl
441  // make the cluster->(true,plane) association and save the
442  // number of true hits in the cluster
443  if (imtru != -1) {
444  // clobber a previously made association?
445  if (nhit > nTruHitInCl[imtru][plane]) {
446  truToCl[imtru][plane] = icl;
447  nTruHitInCl[imtru][plane] = nhit;
448  }
449  } // imtru != 1
450  } // icl
451 
452  // ready to calculate Efficiency, Purity in each plane and EP2
453  for (unsigned short ipl = 0; ipl < plist2.size(); ++ipl) {
454  // ignore daughters
455  bool skipit = false;
456  for (unsigned short ii = 0; ii < moda.size(); ++ii) {
457  if (moda[ii].second == ipl) {
458  skipit = true;
459  break;
460  }
461  } // ii
462  if (skipit) continue;
463  // ignore true particles with few true hits. Outside the detector
464  // or not reconstructable?
465  if (hlist2[ipl].size() < 3) continue;
466 
467  int trackID = plist2[ipl]->TrackId();
468  art::Ptr<simb::MCTruth> theTruth = pi_serv->TrackIdToMCTruth_P(trackID);
469  bool isCosmic = (theTruth->Origin() == simb::kCosmicRay);
470  float KE = 1000 * (plist2[ipl]->E() - plist2[ipl]->Mass());
471  int PDG = abs(plist2[ipl]->PdgCode());
472 
473  std::vector<short> nTru(geom->Nplanes());
474  std::vector<short> nRec(geom->Nplanes());
475  std::vector<short> nTruRec(geom->Nplanes());
476  std::vector<float> eff(geom->Nplanes());
477  std::vector<float> pur(geom->Nplanes());
478  std::vector<float> ep(geom->Nplanes());
479  for (unsigned int plane = 0; plane < geom->Nplanes(); ++plane) {
480  // count the number of true hits in this plane for the true particle.
481  // First count the mother hits
482  for (unsigned short ii = 0; ii < hlist2[ipl].size(); ++ii) {
483  if (ViewToPlane[hlist2[ipl][ii]->View()] == plane) ++nTru[plane];
484  } // ii
485  // mf::LogVerbatim("PFPAna")
486  // <<"Chk mom "<<ipl<<" plane "<<plane<<" nTru "<<nTru[plane];
487  // next look for daughters and count those hits in all generations
488  unsigned short mom = ipl;
489  std::vector<std::pair<unsigned short, unsigned short>>::reverse_iterator rit =
490  moda.rbegin();
491  while (rit != moda.rend()) {
492  if ((*rit).first == mom) {
493  unsigned short dau = (*rit).second;
494  for (unsigned short jj = 0; jj < hlist2[dau].size(); ++jj) {
495  if (ViewToPlane[hlist2[dau][jj]->View()] == plane) ++nTru[plane];
496  } // jj
497  // It is likely that one hit appears in the mother list
498  // as well as the daughter list, so subtract one from the count
499  --nTru[plane];
500  mom = (*rit).second;
501  // mf::LogVerbatim("PFPAna")<<"new mom "<<mom<<" nTru
502  // "<<nTru[plane];
503  } // (*rit).first == mom
504  ++rit;
505  } // rit
506  // mf::LogVerbatim("PFPAna")<<"Chk dau "<<nTru[plane];
507  if (nTru[plane] == 0) {
508  // mf::LogVerbatim("PFPAna")<<"No true hits in plane "<<plane
509  // <<" for truth particle "<<ipl;
510  continue;
511  }
512  short icl = truToCl[ipl][plane];
513  nRec[plane] = nRecHitInCl[icl];
514  nTruRec[plane] = nTruHitInCl[ipl][plane];
515  // mf::LogVerbatim("PFPAna")<<"icl "<<icl<<" nRec "<<nRec[plane]
516  // <<" nTruRec "<<nTruRec[plane];
517  if (nTru[plane] > 0) eff[plane] = (float)nTruRec[plane] / (float)nTru[plane];
518  if (nRec[plane] > 0) pur[plane] = (float)nTruRec[plane] / (float)nRec[plane];
519  ep[plane] = eff[plane] * pur[plane];
520  } // plane
521  // sort the ep values in ascending order
522  std::vector<float> temp;
523  temp = ep;
524  std::sort(temp.begin(), temp.end());
525  // EP2 is the second highest value
526  unsigned short ii = temp.size() - 2;
527  float ep2 = temp[ii];
528  // find the plane that defined EP2
529  short ep2Plane = 0;
530  short ep2Cluster = 0;
531  for (unsigned short jj = 0; jj < temp.size(); ++jj) {
532  if (ep[jj] == ep2) {
533  ep2Plane = jj;
534  ep2Cluster = truToCl[ipl][ep2Plane];
535  break;
536  }
537  } // jj
538  // find the US and DS ends of the cluster for printing
539  std::array<double, 2> clBeg, clEnd;
540  if (ep2Cluster >= 0) {
541  clBeg[0] = clusters[ep2Cluster]->StartWire();
542  clBeg[1] = clusters[ep2Cluster]->StartTick();
543  clEnd[0] = clusters[ep2Cluster]->EndWire();
544  clEnd[1] = clusters[ep2Cluster]->EndTick();
545  }
546  else {
547  clBeg.fill(0.);
548  clEnd.fill(0.);
549  }
550  // fill histograms
551  if (isCosmic) {
552  fCREP2->Fill(ep2);
553  fCRE->Fill(eff[ep2Plane]);
554  fCRP->Fill(pur[ep2Plane]);
555  aveCREP2 += ep2;
556  numCREP2 += 1.;
557  if (fPrintLevel > 1)
558  mf::LogVerbatim("PFPAna")
559  << ">>>CREP2 " << std::fixed << std::setprecision(2) << ep2 << " E " << eff[ep2Plane]
560  << std::setprecision(2) << " P " << pur[ep2Plane] << " P:W:T " << ep2Plane << ":"
561  << (int)clBeg[0] << ":" << (int)clBeg[1] << "-" << ep2Plane << ":" << (int)clEnd[0]
562  << ":" << (int)clEnd[1] << " PDG " << PDG << " KE " << (int)KE << " MeV";
563  } // isCosmic
564  else {
565  float wght = 1.;
566  if (fTrackWeightOption == 1) wght = KE;
567  // accumulate statistics for muons and not-muons
568  if (PDG == 13) {
569  aveNuEP2mu += ep2 * wght;
570  numNuEP2mu += wght;
571  }
572  else {
573  aveNuEP2nm += ep2 * wght;
574  numNuEP2nm += wght;
575  }
576  if (PDG == 11) {
577  fNuKE_elec->Fill(KE, wght);
578  fNuE_elec->Fill(eff[ep2Plane], wght);
579  fNuP_elec->Fill(pur[ep2Plane], wght);
580  fNuEP2_elec->Fill(ep2, wght);
581  fNuEP2_KE_elec->Fill(KE, ep2, wght);
582  }
583  else if (PDG == 13) {
584  fNuKE_muon->Fill(KE, wght);
585  fNuE_muon->Fill(eff[ep2Plane], wght);
586  fNuP_muon->Fill(pur[ep2Plane], wght);
587  fNuEP2_muon->Fill(ep2, wght);
588  fNuEP2_KE_muon->Fill(KE, ep2, wght);
589  }
590  else if (PDG == 211) {
591  fNuKE_pion->Fill(KE, wght);
592  fNuE_pion->Fill(eff[ep2Plane], wght);
593  fNuP_pion->Fill(pur[ep2Plane], wght);
594  fNuEP2_pion->Fill(ep2, wght);
595  fNuEP2_KE_pion->Fill(KE, ep2, wght);
596  }
597  else if (PDG == 321) {
598  fNuKE_kaon->Fill(KE, wght);
599  fNuE_kaon->Fill(eff[ep2Plane], wght);
600  fNuP_kaon->Fill(pur[ep2Plane], wght);
601  fNuEP2_kaon->Fill(ep2, wght);
602  fNuEP2_KE_kaon->Fill(KE, ep2, wght);
603  }
604  else if (PDG == 2212) {
605  fNuKE_prot->Fill(KE, wght);
606  fNuE_prot->Fill(eff[ep2Plane], wght);
607  fNuP_prot->Fill(pur[ep2Plane], wght);
608  fNuEP2_prot->Fill(ep2, wght);
609  fNuEP2_KE_prot->Fill(KE, ep2, wght);
610  }
611  if (fPrintLevel > 1)
612  mf::LogVerbatim("PFPAna")
613  << ">>>NuEP2 " << std::fixed << std::setprecision(2) << ep2 << " E " << eff[ep2Plane]
614  << std::setprecision(2) << " P " << pur[ep2Plane] << " P:W:T " << ep2Plane << ":"
615  << (int)clBeg[0] << ":" << (int)clBeg[1] << "-" << ep2Plane << ":" << (int)clEnd[0]
616  << ":" << (int)clEnd[1] << " PDG " << PDG << " KE " << (int)KE << " MeV ";
617  if (fPrintLevel > 2) {
618  // print out the begin/end true hits
619  mf::LogVerbatim mfp("PFPAna");
620  mfp << " Truth P:W:T ";
621  for (unsigned int plane = 0; plane < geom->Nplanes(); ++plane) {
622  unsigned short loW = 9999;
623  unsigned short loT = 0;
624  unsigned short hiW = 0;
625  unsigned short hiT = 0;
626  for (unsigned short ii = 0; ii < hlist2[ipl].size(); ++ii) {
627  if (ViewToPlane[hlist2[ipl][ii]->View()] == plane) {
628  art::Ptr<recob::Hit> theHit = hlist2[ipl][ii];
629  if (theHit->WireID().Wire < loW) {
630  loW = theHit->WireID().Wire;
631  loT = theHit->PeakTime();
632  }
633  if (theHit->WireID().Wire > hiW) {
634  hiW = theHit->WireID().Wire;
635  hiT = theHit->PeakTime();
636  }
637  } // correct view
638  } // ii
639  mfp << plane << ":" << loW << ":" << loT << "-" << plane << ":" << hiW << ":" << hiT
640  << " ";
641  } // plane
642  } // fPrintLevel > 2
643  } // !isCosmic
644  } // ipl
645 
646  float ave1 = -1.;
647  if (numNuEP2mu > 0.) ave1 = aveNuEP2mu / numNuEP2mu;
648 
649  float ave2 = -1.;
650  if (numNuEP2nm > 0.) ave2 = aveNuEP2nm / numNuEP2nm;
651 
652  float ave3 = -1.;
653  if (numCREP2 > 0.) ave3 = aveCREP2 / numCREP2;
654 
655  if (fPrintLevel > 0) {
656  std::string nuType = "Other";
657  if (neutCCNC == simb::kCC) {
658  if (neutIntType == 1001) nuType = "CCQE";
659  if (neutIntType == 1091) nuType = "DIS";
660  if (neutIntType == 1097) nuType = "COH";
661  if (neutIntType > 1002 && neutIntType < 1091) nuType = "RES";
662  }
663  else if (neutCCNC == simb::kNC) {
664  nuType = "NC";
665  }
666  else {
667  nuType = "Unknown";
668  }
669  mf::LogVerbatim("PFPAna") << "EvtEP2 " << evt.id().event() << " NuType " << nuType << " Enu "
670  << std::fixed << std::setprecision(0) << neutEnergy << std::right
671  << std::fixed << std::setprecision(2) << " NuMuons " << ave1
672  << " NuPiKp " << ave2 << " CosmicRays " << ave3 << " CCNC "
673  << neutCCNC << " IntType " << neutIntType;
674  }
675  } // analyze
676 
677  DEFINE_ART_MODULE(PFPAna)
678 
679 } // end namespace
process_name vertex
Definition: cheaterreco.fcl:51
std::vector< float > fElecKERange
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
Definition: UtilFunc.cxx:42
process_name stream1 can override from command line with o or output services DetectorPropertiesService services DetectorPropertiesService services DetectorPropertiesService services DetectorPropertiesService physics analyzers pmtresponse NeutronTrackingCut services LArG4Parameters gaussian physics producers generator PDG
TH1F * fNuEP2_prot
var pdg
Definition: selectors.fcl:14
std::string fVertexModuleLabel
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Declaration of signal hit object.
walls no right
Definition: selectors.fcl:105
TH1F * fNuP_kaon
std::string fHitsModuleLabel
TH1F * fNuKE_prot
std::vector< float > fProtKERange
static constexpr bool
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
std::string fClusterModuleLabel
std::vector< float > fPionKERange
TProfile * fNuEP2_KE_kaon
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
TH1F * fNuVtx_dy
TH1F * fNuKE_elec
T abs(T value)
Charged-current interactions.
Definition: IPrediction.h:38
TProfile * fNuEP2_KE_elec
TH1F * fNuVtx_dz
TH1F * fNuEP2_muon
TH1F * fNuE_kaon
void beginJob() override
TH1F * fNuP_muon
short fTrackWeightOption
TProfile * fNuEP2_KE_pion
TH1F * fNuEP2_kaon
TH1F * fNuP_elec
const Cut kCosmicRay
std::string fPFParticleModuleLabel
TH1F * fNHitInCluster
TH1F * fNuKE_muon
Declaration of cluster object.
TH1F * fNuP_prot
TH1F * fNuVtx_dx
TH1F * fNuE_prot
TProfile * fNuEP2_KE_prot
TH1F * fNuE_muon
TH1F * fNClusters
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
Neutral-current interactions.
Definition: IPrediction.h:39
TH1F * fNuP_pion
std::vector< float > fMuonKERange
std::string fTrackModuleLabel
art::ServiceHandle< art::TFileService > tfs
TH1F * fNuKE_kaon
TCEvent evt
Definition: DataStructs.cxx:8
bool empty(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:555
void analyze(const art::Event &evt) override
art framework interface to geometry description
std::vector< float > fKaonKERange
TH1F * fNuEP2_pion