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