LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
TrackAna_module.cc
Go to the documentation of this file.
1 //
2 // Name: TrackAna_module.cc
3 //
4 // Purpose: Module TrackAna.
5 //
6 // Configuration parameters.
7 //
8 // TrackModuleLabel: Track module label.
9 // MCTrackModuleLabel: MCTrack module label.
10 // MinMCKE: Minimum MC particle kinetic energy.
11 // MatchColinearity: Minimum colinearity for mc-track matching.
12 // MatchDisp: Maximum uv displacement for mc-track matching.
13 // WMatchDisp: maximum w displacement for mc-track matching.
14 // MatchLength: Minimum length fraction for mc-track matching.
15 //
16 // Created: 2-Aug-2011 H. Greenlee
17 //
18 
19 #include <map>
20 #include <iostream>
21 #include <iomanip>
22 #include <sstream>
23 #include <cmath>
24 #include <memory>
25 
33 #include "cetlib_except/exception.h"
34 
45 
46 #include "TH2F.h"
47 #include "TFile.h"
48 
49 namespace {
50 
51  // Local functions.
52 
53  // Calculate distance to boundary.
54  //----------------------------------------------------------------------------
55  double bdist(const TVector3& pos, unsigned int /*tpc*/ = 0, unsigned int /*cstat*/ = 0)
56  {
57  // Get geometry.
58 
60 
61  double d1 = pos.X(); // Distance to right side (wires).
62  double d2 = 2.*geom->DetHalfWidth() - pos.X(); // Distance to left side (cathode).
63  double d3 = pos.Y() + geom->DetHalfHeight(); // Distance to bottom.
64  double d4 = geom->DetHalfHeight() - pos.Y(); // Distance to top.
65  double d5 = pos.Z(); // Distance to front.
66  double d6 = geom->DetLength() - pos.Z(); // Distance to back.
67 
68  double result = std::min(std::min(std::min(std::min(std::min(d1, d2), d3), d4), d5), d6);
69  return result;
70  }
71 
72  // Length of reconstructed track.
73  //----------------------------------------------------------------------------
74  double length(const recob::Track& track)
75  {
76  double result = 0.;
77  TVector3 disp = track.LocationAtPoint(0);
78  int n = track.NumberTrajectoryPoints();
79 
80  for(int i = 1; i < n; ++i) {
81  const TVector3& pos = track.LocationAtPoint(i);
82  disp -= pos;
83  result += disp.Mag();
84  disp = pos;
85  }
86 
87  return result;
88  }
89 
90  // Length of MC particle.
91  //----------------------------------------------------------------------------
92  double length(const simb::MCParticle& part, double dx,
93  TVector3& start, TVector3& end, TVector3& startmom, TVector3& endmom,
94  unsigned int /*tpc*/ = 0, unsigned int /*cstat*/ = 0)
95  {
96  // Get services.
97 
99  const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
100 
101  // Get fiducial volume boundary.
102 
103  double xmin = 0.;
104  double xmax = 2.*geom->DetHalfWidth();
105  double ymin = -geom->DetHalfHeight();
106  double ymax = geom->DetHalfHeight();
107  double zmin = 0.;
108  double zmax = geom->DetLength();
109  //double ticks_max = detprop->ReadOutWindowSize();
110  double result = 0.;
111  TVector3 disp;
112  int n = part.NumberTrajectoryPoints();
113  bool first = true;
114 
115  for(int i = 0; i < n; ++i) {
116  TVector3 pos = part.Position(i).Vect();
117 
118  // Make fiducial cuts. Require the particle to be within the physical volume of
119  // the tpc, and also require the apparent x position to be within the expanded
120  // readout frame.
121 
122  if(pos.X() >= xmin &&
123  pos.X() <= xmax &&
124  pos.Y() >= ymin &&
125  pos.Y() <= ymax &&
126  pos.Z() >= zmin &&
127  pos.Z() <= zmax) {
128  pos[0] += dx;
129  double ticks = detprop->ConvertXToTicks(pos[0], 0, 0, 0);
130  if(ticks >= 0. && ticks < detprop->ReadOutWindowSize()) {
131  if(first) {
132  start = pos;
133  startmom = part.Momentum(i).Vect();
134  }
135  else {
136  disp -= pos;
137  result += disp.Mag();
138  }
139  first = false;
140  disp = pos;
141  end = pos;
142  endmom = part.Momentum(i).Vect();
143  }
144  }
145  }
146 
147  return result;
148  }
149 
150  // Length of MCTrack.
151  // In this function, the extracted start and end momenta are converted to GeV
152  // (MCTrack stores momenta in Mev).
153  //----------------------------------------------------------------------------
154  double length(const sim::MCTrack& mctrk, double dx,
155  TVector3& start, TVector3& end, TVector3& startmom, TVector3& endmom,
156  unsigned int /*tpc*/ = 0, unsigned int /*cstat*/ = 0)
157  {
158  // Get services.
159 
161  // const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
162  auto const* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
163 
164 
165  // Get fiducial volume boundary.
166 
167  double xmin = 0.;
168  double xmax = 2.*geom->DetHalfWidth();
169  double ymin = -geom->DetHalfHeight();
170  double ymax = geom->DetHalfHeight();
171  double zmin = 0.;
172  double zmax = geom->DetLength();
173  //double ticks_max = detprop->ReadOutWindowSize();
174  double result = 0.;
175  TVector3 disp;
176  int n = mctrk.size();
177  bool first = true;
178 
179  for(int i = 0; i < n; ++i) {
180  TVector3 pos = mctrk[i].Position().Vect();
181 
182  // Make fiducial cuts. Require the particle to be within the physical volume of
183  // the tpc, and also require the apparent x position to be within the expanded
184  // readout frame.
185 
186  if(pos.X() >= xmin &&
187  pos.X() <= xmax &&
188  pos.Y() >= ymin &&
189  pos.Y() <= ymax &&
190  pos.Z() >= zmin &&
191  pos.Z() <= zmax) {
192  pos[0] += dx;
193  double ticks = detprop->ConvertXToTicks(pos[0], 0, 0, 0);
194  if(ticks >= 0. && ticks < detprop->ReadOutWindowSize()) {
195  if(first) {
196  start = pos;
197  startmom = 0.001 * mctrk[i].Momentum().Vect();
198  }
199  else {
200  disp -= pos;
201  result += disp.Mag();
202  }
203  first = false;
204  disp = pos;
205  end = pos;
206  endmom = 0.001 * mctrk[i].Momentum().Vect();
207  }
208  }
209  }
210 
211  return result;
212  }
213 
214  // Fill efficiency histogram assuming binomial errors.
215 
216  void effcalc(const TH1* hnum, const TH1* hden, TH1* heff)
217  {
218  int nbins = hnum->GetNbinsX();
219  if (nbins != hden->GetNbinsX())
220  throw cet::exception("TrackAna") << "effcalc[" __FILE__ "]: incompatible histograms (I)\n";
221  if (nbins != heff->GetNbinsX())
222  throw cet::exception("TrackAna") << "effcalc[" __FILE__ "]: incompatible histograms (II)\n";
223 
224  // Loop over bins, including underflow and overflow.
225 
226  for(int ibin = 0; ibin <= nbins+1; ++ibin) {
227  double num = hnum->GetBinContent(ibin);
228  double den = hden->GetBinContent(ibin);
229  if(den == 0.) {
230  heff->SetBinContent(ibin, 0.);
231  heff->SetBinError(ibin, 0.);
232  }
233  else {
234  double eff = num / den;
235  if(eff < 0.)
236  eff = 0.;
237  if(eff > 1.)
238  eff = 1.;
239  double err = std::sqrt(eff * (1.-eff) / den);
240  heff->SetBinContent(ibin, eff);
241  heff->SetBinError(ibin, err);
242  }
243  }
244 
245  heff->SetMinimum(0.);
246  heff->SetMaximum(1.05);
247  heff->SetMarkerStyle(20);
248  }
249 
250 
251 class flattener : public std::vector<unsigned int> {
252 
253 public:
254 
255  flattener() : std::vector<unsigned int>() {};
256 
257  flattener(const std::vector<std::vector<unsigned int> >& input)
258  { Convert(input); }
259 
260  ~flattener(){}
261 
262  void Convert(const std::vector<std::vector<unsigned int> >& input)
263  {
264  clear();
265  size_t length=0;
266  for(auto const& vec : input)
267  length += vec.size();
268  reserve(length);
269 
270  for(auto const& vec : input)
271  for(auto const& value : vec)
272  push_back(value);
273 
274  }
275 }; // end class flattener
276 
277 } // end namespace
278 
279 namespace trkf {
280 
281  class TrackAna : public art::EDAnalyzer
282  {
283  public:
284 
285  // Embedded structs.
286 
287  // Struct for histograms that depend on reco track only.
288 
289  struct RecoHists
290  {
291  // Constructors.
292 
293  RecoHists();
294  RecoHists(const std::string& subdir);
295 
296  // Pure reco track histograms.
297 
298  TH1F* fHstartx; // Starting x position.
299  TH1F* fHstarty; // Starting y position.
300  TH1F* fHstartz; // Starting z position.
301  TH1F* fHstartd; // Starting distance to boundary.
302  TH1F* fHendx; // Ending x position.
303  TH1F* fHendy; // Ending y position.
304  TH1F* fHendz; // Ending z position.
305  TH1F* fHendd; // Ending distance to boundary.
306  TH1F* fHtheta; // Theta.
307  TH1F* fHphi; // Phi.
308  TH1F* fHtheta_xz; // Theta_xz.
309  TH1F* fHtheta_yz; // Theta_yz.
310  TH1F* fHmom; // Momentum.
311  TH1F* fHmoml; // Momentum (low momentum).
312  TH1F* fHlen; // Length.
313  TH1F* fHlens; // Length (short tracks).
314 
315  // Histograms for the consituent Hits
316 
317  TH1F* fHHitChg; // hit charge
318  TH1F* fHHitWidth; // hit width
319  TH1F* fHHitPdg; // Pdg primarily responsible.
320  TH1F* fHHitTrkId; // TrkId
321  TH1F* fModeFrac; // mode fraction
322  TH1F* fNTrkIdTrks; // # of stitched tracks in which unique TrkId appears
323  TH2F* fNTrkIdTrks2;
324  TH2F* fNTrkIdTrks3;
325  };
326 
327  // Struct for mc particles and mc-matched tracks.
328 
329  struct MCHists
330  {
331  // Constructors.
332 
333  MCHists();
334  MCHists(const std::string& subdir);
335 
336  // Reco-MC matching.
337 
338  TH2F* fHduvcosth; // 2D mc vs. data matching, duv vs. cos(theta).
339  TH1F* fHcosth; // 1D direction matching, cos(theta).
340  TH1F* fHmcu; // 1D endpoint truth u.
341  TH1F* fHmcv; // 1D endpoint truth v.
342  TH1F* fHmcw; // 1D endpoint truth w.
343  TH1F* fHupull; // 1D endpoint u pull.
344  TH1F* fHvpull; // 1D endpoint v pull.
345  TH1F* fHmcdudw; // Truth du/dw.
346  TH1F* fHmcdvdw; // Truth dv/dw.
347  TH1F* fHdudwpull; // du/dw pull.
348  TH1F* fHdvdwpull; // dv/dw pull.
349  TH1F* fHHitEff; // Hit efficiency.
350  TH1F* fHHitPurity; // Hit purity.
351 
352  // Histograms for matched tracks.
353 
354  TH1F* fHstartdx; // Start dx.
355  TH1F* fHstartdy; // Start dy.
356  TH1F* fHstartdz; // Start dz.
357  TH1F* fHenddx; // End dx.
358  TH1F* fHenddy; // End dy.
359  TH1F* fHenddz; // End dz.
360  TH2F* fHlvsl; // MC vs. reco length.
361  TH1F* fHdl; // Delta(length).
362  TH2F* fHpvsp; // MC vs. reco momentum.
363  TH2F* fHpvspc; // MC vs. reco momentum (contained tracks).
364  TH1F* fHdp; // Momentum difference.
365  TH1F* fHdpc; // Momentum difference (contained tracks).
366  TH1F* fHppull; // Momentum pull.
367  TH1F* fHppullc; // Momentum pull (contained tracks).
368 
369  // Pure MC particle histograms (efficiency denominator).
370 
371  TH1F* fHmcstartx; // Starting x position.
372  TH1F* fHmcstarty; // Starting y position.
373  TH1F* fHmcstartz; // Starting z position.
374  TH1F* fHmcendx; // Ending x position.
375  TH1F* fHmcendy; // Ending y position.
376  TH1F* fHmcendz; // Ending z position.
377  TH1F* fHmctheta; // Theta.
378  TH1F* fHmcphi; // Phi.
379  TH1F* fHmctheta_xz; // Theta_xz.
380  TH1F* fHmctheta_yz; // Theta_yz.
381  TH1F* fHmcmom; // Momentum.
382  TH1F* fHmcmoml; // Momentum (low momentum).
383  TH1F* fHmcke; // Kinetic energy.
384  TH1F* fHmckel; // Kinetic energy (low energy).
385  TH1F* fHmclen; // Length.
386  TH1F* fHmclens; // Length (short tracks).
387 
388  // Histograms for well-reconstructed matched tracks (efficiency numerator).
389 
390  TH1F* fHgstartx; // Starting x position.
391  TH1F* fHgstarty; // Starting y position.
392  TH1F* fHgstartz; // Starting z position.
393  TH1F* fHgendx; // Ending x position.
394  TH1F* fHgendy; // Ending y position.
395  TH1F* fHgendz; // Ending z position.
396  TH1F* fHgtheta; // Theta.
397  TH1F* fHgphi; // Phi.
398  TH1F* fHgtheta_xz; // Theta_xz.
399  TH1F* fHgtheta_yz; // Theta_yz.
400  TH1F* fHgmom; // Momentum.
401  TH1F* fHgmoml; // Momentum (low momentum).
402  TH1F* fHgke; // Kinetic energy.
403  TH1F* fHgkel; // Kinetic energy (low momentum).
404  TH1F* fHglen; // Length.
405  TH1F* fHglens; // Length (short tracks).
406 
407  // Efficiency histograms.
408 
409  TH1F* fHestartx; // Starting x position.
410  TH1F* fHestarty; // Starting y position.
411  TH1F* fHestartz; // Starting z position.
412  TH1F* fHeendx; // Ending x position.
413  TH1F* fHeendy; // Ending y position.
414  TH1F* fHeendz; // Ending z position.
415  TH1F* fHetheta; // Theta.
416  TH1F* fHephi; // Phi.
417  TH1F* fHetheta_xz; // Theta_xz.
418  TH1F* fHetheta_yz; // Theta_yz.
419  TH1F* fHemom; // Momentum.
420  TH1F* fHemoml; // Momentum (low momentum).
421  TH1F* fHeke; // Kinetic energy.
422  TH1F* fHekel; // Kinetic energy (low momentum).
423  TH1F* fHelen; // Length.
424  TH1F* fHelens; // Length (short tracks).
425 
426 
427  };
428 
429  // Constructors, destructor
430 
431  explicit TrackAna(fhicl::ParameterSet const& pset);
432  virtual ~TrackAna();
433 
434  // Overrides.
435 
436  void analyze(const art::Event& evt);
437  void anaStitch(const art::Event& evt);
438  void endJob();
439 
440  private:
441 
442  template <typename T> std::vector<size_t> fsort_indexes(const std::vector<T> &v) ;
443 
444  // Fcl Attributes.
445 
446  std::string fTrackModuleLabel;
447  std::string fMCTrackModuleLabel;
449  std::string fStitchModuleLabel;
452  std::string fHitModuleLabel;
453 
454  int fDump; // Number of events to dump to debug message facility.
455  double fMinMCKE; // Minimum MC particle kinetic energy (GeV).
456  double fMinMCLen; // Minimum MC particle length in tpc (cm).
457  double fMatchColinearity; // Minimum matching colinearity.
458  double fMatchDisp; // Maximum matching displacement.
459  double fWMatchDisp; // Maximum matching displacement in the w direction.
460  double fMatchLength; // Minimum length fraction.
461  bool fIgnoreSign; // Ignore sign of mc particle if true.
462  bool fStitchedAnalysis; // if true, do the whole drill-down from stitched track to assd hits
463 
464  std::string fOrigin;
467  int fPrintLevel; // 0 = none, 1 = event summary, 2 = track detail
468 
469  // Histograms.
470 
471  std::map<int, MCHists> fMCHistMap; // Indexed by pdg id.
472  std::map<int, RecoHists> fRecoHistMap; // Indexed by pdg id.
473 
474  // Statistics.
475 
477  };
478 
480 
481  // RecoHists methods.
482 
484  //
485  // Purpose: Default constructor.
486  //
487  fHstartx(0),
488  fHstarty(0),
489  fHstartz(0),
490  fHstartd(0),
491  fHendx(0),
492  fHendy(0),
493  fHendz(0),
494  fHendd(0),
495  fHtheta(0),
496  fHphi(0),
497  fHtheta_xz(0),
498  fHtheta_yz(0),
499  fHmom(0),
500  fHmoml(0),
501  fHlen(0),
502  fHlens(0)
503  ,fHHitChg(0)
504  ,fHHitWidth(0)
505  ,fHHitPdg(0)
506  ,fHHitTrkId(0)
507  ,fModeFrac(0)
508  ,fNTrkIdTrks(0)
509  ,fNTrkIdTrks2(0)
510  ,fNTrkIdTrks3(0)
511  {}
512 
513  TrackAna::RecoHists::RecoHists(const std::string& subdir)
514  //
515  // Purpose: Initializing constructor.
516  //
517  {
518  // Make sure all histogram pointers are initially zero.
519 
520  *this = RecoHists();
521 
522  // Get services.
523 
526 
527  // Make histogram directory.
528 
529  art::TFileDirectory topdir = tfs->mkdir("trkana", "TrackAna histograms");
530  art::TFileDirectory dir = topdir.mkdir(subdir);
531 
532  // Book histograms.
533 
534  fHstartx = dir.make<TH1F>("xstart", "X Start Position",
535  100, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
536  fHstarty = dir.make<TH1F>("ystart", "Y Start Position",
537  100, -geom->DetHalfHeight(), geom->DetHalfHeight());
538  fHstartz = dir.make<TH1F>("zstart", "Z Start Position",
539  100, 0., geom->DetLength());
540  fHstartd = dir.make<TH1F>("dstart", "Start Position Distance to Boundary",
541  100, -10., geom->DetHalfWidth());
542  fHendx = dir.make<TH1F>("xend", "X End Position",
543  100, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
544  fHendy = dir.make<TH1F>("yend", "Y End Position",
545  100, -geom->DetHalfHeight(), geom->DetHalfHeight());
546  fHendz = dir.make<TH1F>("zend", "Z End Position",
547  100, 0., geom->DetLength());
548  fHendd = dir.make<TH1F>("dend", "End Position Distance to Boundary",
549  100, -10., geom->DetHalfWidth());
550  fHtheta = dir.make<TH1F>("theta", "Theta", 100, 0., 3.142);
551  fHphi = dir.make<TH1F>("phi", "Phi", 100, -3.142, 3.142);
552  fHtheta_xz = dir.make<TH1F>("theta_xz", "Theta_xz", 100, -3.142, 3.142);
553  fHtheta_yz = dir.make<TH1F>("theta_yz", "Theta_yz", 100, -3.142, 3.142);
554  fHmom = dir.make<TH1F>("mom", "Momentum", 100, 0., 10.);
555  fHmoml = dir.make<TH1F>("moml", "Momentum", 100, 0., 1.);
556  fHlen = dir.make<TH1F>("len", "Track Length", 100, 0., 1.1 * geom->DetLength());
557  fHlens = dir.make<TH1F>("lens", "Track Length", 100, 0., 0.1 * geom->DetLength());
558  fHHitChg = dir.make<TH1F>("hchg", "Hit Charge (ADC counts)", 100, 0., 4000.);
559  fHHitWidth = dir.make<TH1F>("hwid", "Hit Width (ticks)", 40, 0., 20.);
560  fHHitPdg = dir.make<TH1F>("hpdg", "Hit Pdg code",5001, -2500.5, +2500.5);
561  fHHitTrkId = dir.make<TH1F>("htrkid", "Hit Track ID", 401, -200.5, +200.5);
562  fModeFrac = dir.make<TH1F>("hmodefrac", "quasi-Purity: Fraction of component tracks with the Track mode value", 20, 0.01, 1.01);
563  fNTrkIdTrks = dir.make<TH1F>("hntrkids", "quasi-Efficiency: Number of stitched tracks in which TrkId appears", 20, 0., +10.0);
564  fNTrkIdTrks2 = dir.make<TH2F>("hntrkids2", "Number of stitched tracks in which TrkId appears vs KE [GeV]", 20, 0., +10.0, 20, 0.0, 1.5);
565  fNTrkIdTrks3 = dir.make<TH2F>("hntrkids3", "MC Track vs Reco Track, wtd by nhits on Collection Plane", 10, -0.5, 9.5, 10, -0.5, 9.5);
566  fNTrkIdTrks3->GetXaxis()->SetTitle("Sorted-by-Descending-CPlane-Hits outer Track Number");
567  fNTrkIdTrks3->GetYaxis()->SetTitle("Sorted-by-Descending-True-Length G4Track");
568 
569  }
570 
571  // MCHists methods.
572 
573  TrackAna::MCHists::MCHists() :
574  //
575  // Purpose: Default constructor.
576  //
577  fHduvcosth(0),
578  fHcosth(0),
579  fHmcu(0),
580  fHmcv(0),
581  fHmcw(0),
582  fHupull(0),
583  fHvpull(0),
584  fHmcdudw(0),
585  fHmcdvdw(0),
586  fHdudwpull(0),
587  fHdvdwpull(0),
588  fHHitEff(0),
589  fHHitPurity(0),
590  fHstartdx(0),
591  fHstartdy(0),
592  fHstartdz(0),
593  fHenddx(0),
594  fHenddy(0),
595  fHenddz(0),
596  fHlvsl(0),
597  fHdl(0),
598  fHpvsp(0),
599  fHpvspc(0),
600  fHdp(0),
601  fHdpc(0),
602  fHppull(0),
603  fHppullc(0),
604  fHmcstartx(0),
605  fHmcstarty(0),
606  fHmcstartz(0),
607  fHmcendx(0),
608  fHmcendy(0),
609  fHmcendz(0),
610  fHmctheta(0),
611  fHmcphi(0),
612  fHmctheta_xz(0),
613  fHmctheta_yz(0),
614  fHmcmom(0),
615  fHmcmoml(0),
616  fHmcke(0),
617  fHmckel(0),
618  fHmclen(0),
619  fHmclens(0),
620  fHgstartx(0),
621  fHgstarty(0),
622  fHgstartz(0),
623  fHgendx(0),
624  fHgendy(0),
625  fHgendz(0),
626  fHgtheta(0),
627  fHgphi(0),
628  fHgtheta_xz(0),
629  fHgtheta_yz(0),
630  fHgmom(0),
631  fHgmoml(0),
632  fHgke(0),
633  fHgkel(0),
634  fHglen(0),
635  fHglens(0),
636  fHestartx(0),
637  fHestarty(0),
638  fHestartz(0),
639  fHeendx(0),
640  fHeendy(0),
641  fHeendz(0),
642  fHetheta(0),
643  fHephi(0),
644  fHetheta_xz(0),
645  fHetheta_yz(0),
646  fHemom(0),
647  fHemoml(0),
648  fHeke(0),
649  fHekel(0),
650  fHelen(0),
651  fHelens(0)
652  {}
653 
654  TrackAna::MCHists::MCHists(const std::string& subdir)
655  //
656  // Purpose: Initializing constructor.
657  //
658  {
659  // Make sure all histogram pointers are initially zero.
660 
661  *this = MCHists();
662 
663  // Get services.
664 
667 
668  // Make histogram directory.
669 
670  art::TFileDirectory topdir = tfs->mkdir("trkana", "TrackAna histograms");
671  art::TFileDirectory dir = topdir.mkdir(subdir);
672 
673  // Book histograms.
674 
675  fHduvcosth = dir.make<TH2F>("duvcosth", "Delta(uv) vs. Colinearity",
676  100, 0.95, 1., 100, 0., 1.);
677  fHcosth = dir.make<TH1F>("colin", "Colinearity", 100, 0.95, 1.);
678  fHmcu = dir.make<TH1F>("mcu", "MC Truth U", 100, -5., 5.);
679  fHmcv = dir.make<TH1F>("mcv", "MC Truth V", 100, -5., 5.);
680  fHmcw = dir.make<TH1F>("mcw", "MC Truth W", 100, -20., 20.);
681  fHupull = dir.make<TH1F>("dupull", "U Pull", 100, -20., 20.);
682  fHvpull = dir.make<TH1F>("dvpull", "V Pull", 100, -20., 20.);
683  fHmcdudw = dir.make<TH1F>("mcdudw", "MC Truth U Slope", 100, -0.2, 0.2);
684  fHmcdvdw = dir.make<TH1F>("mcdvdw", "MV Truth V Slope", 100, -0.2, 0.2);
685  fHdudwpull = dir.make<TH1F>("dudwpull", "U Slope Pull", 100, -10., 10.);
686  fHdvdwpull = dir.make<TH1F>("dvdwpull", "V Slope Pull", 100, -10., 10.);
687  fHHitEff = dir.make<TH1F>("hiteff", "MC Hit Efficiency", 100, 0., 1.0001);
688  fHHitPurity = dir.make<TH1F>("hitpurity", "MC Hit Purity", 100, 0., 1.0001);
689  fHstartdx = dir.make<TH1F>("startdx", "Start Delta x", 100, -10., 10.);
690  fHstartdy = dir.make<TH1F>("startdy", "Start Delta y", 100, -10., 10.);
691  fHstartdz = dir.make<TH1F>("startdz", "Start Delta z", 100, -10., 10.);
692  fHenddx = dir.make<TH1F>("enddx", "End Delta x", 100, -10., 10.);
693  fHenddy = dir.make<TH1F>("enddy", "End Delta y", 100, -10., 10.);
694  fHenddz = dir.make<TH1F>("enddz", "End Delta z", 100, -10., 10.);
695  fHlvsl = dir.make<TH2F>("lvsl", "Reco Length vs. MC Truth Length",
696  100, 0., 1.1 * geom->DetLength(), 100, 0., 1.1 * geom->DetLength());
697  fHdl = dir.make<TH1F>("dl", "Track Length Minus MC Particle Length", 100, -50., 50.);
698  fHpvsp = dir.make<TH2F>("pvsp", "Reco Momentum vs. MC Truth Momentum",
699  100, 0., 5., 100, 0., 5.);
700  fHpvspc = dir.make<TH2F>("pvspc", "Reco Momentum vs. MC Truth Momentum (Contained Tracks)",
701  100, 0., 5., 100, 0., 5.);
702  fHdp = dir.make<TH1F>("dp", "Reco-MC Momentum Difference", 100, -5., 5.);
703  fHdpc = dir.make<TH1F>("dpc", "Reco-MC Momentum Difference (Contained Tracks)",
704  100, -5., 5.);
705  fHppull = dir.make<TH1F>("ppull", "Momentum Pull", 100, -10., 10.);
706  fHppullc = dir.make<TH1F>("ppullc", "Momentum Pull (Contained Tracks)", 100, -10., 10.);
707 
708  fHmcstartx = dir.make<TH1F>("mcxstart", "MC X Start Position",
709  10, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
710  fHmcstarty = dir.make<TH1F>("mcystart", "MC Y Start Position",
711  10, -geom->DetHalfHeight(), geom->DetHalfHeight());
712  fHmcstartz = dir.make<TH1F>("mczstart", "MC Z Start Position",
713  10, 0., geom->DetLength());
714  fHmcendx = dir.make<TH1F>("mcxend", "MC X End Position",
715  10, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
716  fHmcendy = dir.make<TH1F>("mcyend", "MC Y End Position",
717  10, -geom->DetHalfHeight(), geom->DetHalfHeight());
718  fHmcendz = dir.make<TH1F>("mczend", "MC Z End Position",
719  10, 0., geom->DetLength());
720  fHmctheta = dir.make<TH1F>("mctheta", "MC Theta", 20, 0., 3.142);
721  fHmcphi = dir.make<TH1F>("mcphi", "MC Phi", 10, -3.142, 3.142);
722  fHmctheta_xz = dir.make<TH1F>("mctheta_xz", "MC Theta_xz", 40, -3.142, 3.142);
723  fHmctheta_yz = dir.make<TH1F>("mctheta_yz", "MC Theta_yz", 40, -3.142, 3.142);
724  fHmcmom = dir.make<TH1F>("mcmom", "MC Momentum", 10, 0., 10.);
725  fHmcmoml = dir.make<TH1F>("mcmoml", "MC Momentum", 10, 0., 1.);
726  fHmcke = dir.make<TH1F>("mcke", "MC Kinetic Energy", 10, 0., 10.);
727  fHmckel = dir.make<TH1F>("mckel", "MC Kinetic Energy", 10, 0., 1.);
728  fHmclen = dir.make<TH1F>("mclen", "MC Particle Length", 10, 0., 1.1 * geom->DetLength());
729  fHmclens = dir.make<TH1F>("mclens", "MC Particle Length", 10, 0., 0.1 * geom->DetLength());
730 
731  fHgstartx = dir.make<TH1F>("gxstart", "Good X Start Position",
732  10, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
733  fHgstarty = dir.make<TH1F>("gystart", "Good Y Start Position",
734  10, -geom->DetHalfHeight(), geom->DetHalfHeight());
735  fHgstartz = dir.make<TH1F>("gzstart", "Good Z Start Position",
736  10, 0., geom->DetLength());
737  fHgendx = dir.make<TH1F>("gxend", "Good X End Position",
738  10, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
739  fHgendy = dir.make<TH1F>("gyend", "Good Y End Position",
740  10, -geom->DetHalfHeight(), geom->DetHalfHeight());
741  fHgendz = dir.make<TH1F>("gzend", "Good Z End Position",
742  10, 0., geom->DetLength());
743  fHgtheta = dir.make<TH1F>("gtheta", "Good Theta", 20, 0., 3.142);
744  fHgphi = dir.make<TH1F>("gphi", "Good Phi", 10, -3.142, 3.142);
745  fHgtheta_xz = dir.make<TH1F>("gtheta_xz", "Good Theta_xz", 40, -3.142, 3.142);
746  fHgtheta_yz = dir.make<TH1F>("gtheta_yz", "Good Theta_yz", 40, -3.142, 3.142);
747  fHgmom = dir.make<TH1F>("gmom", "Good Momentum", 10, 0., 10.);
748  fHgmoml = dir.make<TH1F>("gmoml", "Good Momentum", 10, 0., 1.);
749  fHgke = dir.make<TH1F>("gke", "Good Kinetic Energy", 10, 0., 10.);
750  fHgkel = dir.make<TH1F>("gkel", "Good Kinetic Energy", 10, 0., 1.);
751  fHglen = dir.make<TH1F>("glen", "Good Particle Length", 10, 0., 1.1 * geom->DetLength());
752  fHglens = dir.make<TH1F>("glens", "Good Particle Length", 10, 0., 0.1 * geom->DetLength());
753 
754  fHestartx = dir.make<TH1F>("exstart", "Efficiency vs. X Start Position",
755  10, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
756  fHestarty = dir.make<TH1F>("eystart", "Efficiency vs. Y Start Position",
757  10, -geom->DetHalfHeight(), geom->DetHalfHeight());
758  fHestartz = dir.make<TH1F>("ezstart", "Efficiency vs. Z Start Position",
759  10, 0., geom->DetLength());
760  fHeendx = dir.make<TH1F>("exend", "Efficiency vs. X End Position",
761  10, -2.*geom->DetHalfWidth(), 4.*geom->DetHalfWidth());
762  fHeendy = dir.make<TH1F>("eyend", "Efficiency vs. Y End Position",
763  10, -geom->DetHalfHeight(), geom->DetHalfHeight());
764  fHeendz = dir.make<TH1F>("ezend", "Efficiency vs. Z End Position",
765  10, 0., geom->DetLength());
766  fHetheta = dir.make<TH1F>("etheta", "Efficiency vs. Theta", 20, 0., 3.142);
767  fHephi = dir.make<TH1F>("ephi", "Efficiency vs. Phi", 10, -3.142, 3.142);
768  fHetheta_xz = dir.make<TH1F>("etheta_xz", "Efficiency vs. Theta_xz", 40, -3.142, 3.142);
769  fHetheta_yz = dir.make<TH1F>("etheta_yz", "Efficiency vs. Theta_yz", 40, -3.142, 3.142);
770  fHemom = dir.make<TH1F>("emom", "Efficiency vs. Momentum", 10, 0., 10.);
771  fHemoml = dir.make<TH1F>("emoml", "Efficiency vs. Momentum", 10, 0., 1.);
772  fHeke = dir.make<TH1F>("eke", "Efficiency vs. Kinetic Energy", 10, 0., 10.);
773  fHekel = dir.make<TH1F>("ekel", "Efficiency vs. Kinetic Energy", 10, 0., 1.);
774  fHelen = dir.make<TH1F>("elen", "Efficiency vs. Particle Length",
775  10, 0., 1.1 * geom->DetLength());
776  fHelens = dir.make<TH1F>("elens", "Efficiency vs. Particle Length",
777  10, 0., 0.1 * geom->DetLength());
778  }
779 
781  //
782  // Purpose: Constructor.
783  //
784  // Arguments: pset - Module parameters.
785  //
786  : EDAnalyzer(pset)
787  , fTrackModuleLabel(pset.get<std::string>("TrackModuleLabel"))
788  , fMCTrackModuleLabel(pset.get<std::string>("MCTrackModuleLabel"))
789  , fSpacepointModuleLabel(pset.get<std::string>("SpacepointModuleLabel"))
790  , fStitchModuleLabel(pset.get<std::string>("StitchModuleLabel"))
791  , fTrkSpptAssocModuleLabel(pset.get<std::string>("TrkSpptAssocModuleLabel"))
792  , fHitSpptAssocModuleLabel(pset.get<std::string>("HitSpptAssocModuleLabel"))
793  , fHitModuleLabel(pset.get<std::string>("HitModuleLabel"))
794  , fDump(pset.get<int>("Dump"))
795  , fMinMCKE(pset.get<double>("MinMCKE"))
796  , fMinMCLen(pset.get<double>("MinMCLen"))
797  , fMatchColinearity(pset.get<double>("MatchColinearity"))
798  , fMatchDisp(pset.get<double>("MatchDisp"))
799  , fWMatchDisp(pset.get<double>("WMatchDisp"))
800  , fMatchLength(pset.get<double>("MatchLength"))
801  , fIgnoreSign(pset.get<bool>("IgnoreSign"))
802  , fStitchedAnalysis(pset.get<bool>("StitchedAnalysis",false))
803  , fOrigin(pset.get<std::string>("MCTrackOrigin", "Any"))
804  , fPrintLevel(pset.get<int>("PrintLevel",0))
805  , fNumEvent(0)
806  {
807 
808  // Decide whether to check MCTrack origin
809  fCheckOrigin = false;
811  if(fOrigin.find("Beam") != std::string::npos) {
812  fCheckOrigin = true;
814  } else if(fOrigin.find("Cosmic") != std::string::npos) {
815  fCheckOrigin = true;
817  } else if(fOrigin.find("Super") != std::string::npos) {
818  fCheckOrigin = true;
820  } else if(fOrigin.find("Single") != std::string::npos) {
821  fCheckOrigin = true;
823  }
824 
825  // Report.
826 
827  mf::LogInfo("TrackAna")
828  << "TrackAna configured with the following parameters:\n"
829  << " TrackModuleLabel = " << fTrackModuleLabel << "\n"
830  << " MCTrackModuleLabel = " << fMCTrackModuleLabel << "\n"
831  << " StitchModuleLabel = " << fStitchModuleLabel << "\n"
832  << " TrkSpptAssocModuleLabel = " << fTrkSpptAssocModuleLabel << "\n"
833  << " HitSpptAssocModuleLabel = " << fHitSpptAssocModuleLabel << "\n"
834  << " HitModuleLabel = " << fHitModuleLabel << "\n"
835  << " Dump = " << fDump << "\n"
836  << " MinMCKE = " << fMinMCKE << "\n"
837  << " MinMCLen = " << fMinMCLen
838  << " Origin = " << fOrigin<<" Origin value "<<fOriginValue;
839  }
840 
842  //
843  // Purpose: Destructor.
844  //
845  {}
846 
847  void TrackAna::analyze(const art::Event& evt)
848  //
849  // Purpose: Analyze method.
850  //
851  // Arguments: event - Art event.
852  //
853  {
854  const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
858 
859 
860  ++fNumEvent;
861 
862  // Optional dump stream.
863 
864  std::unique_ptr<mf::LogInfo> pdump;
865  if(fDump > 0) {
866  --fDump;
867  pdump = std::unique_ptr<mf::LogInfo>(new mf::LogInfo("TrackAna"));
868  }
869 
870  // Make sure histograms are booked.
871 
872  bool mc = !evt.isRealData();
873 
874  // Get MCTracks.
875 
877  evt.getByLabel(fMCTrackModuleLabel, mctrackh);
878  // pair of MCTrack and index of matched reco track
879  std::vector<std::pair<const sim::MCTrack*, int>> selected_mctracks;
880 
881  if(mc && mctrackh.isValid()) {
882 
883  selected_mctracks.reserve(mctrackh->size());
884 
885  // Dump MCTracks.
886 
887  if(pdump) {
888  *pdump << "MC Tracks\n"
889  << " Id pdg x y z dx dy dz p\n"
890  << "-------------------------------------------------------------------------------------------\n";
891  }
892 
893  // Loop over mc tracks, and fill histograms that depend only
894  // on mc information. Also, fill a secondary list of mc tracks
895  // that pass various selection criteria.
896 
897  for(std::vector<sim::MCTrack>::const_iterator imctrk = mctrackh->begin();
898  imctrk != mctrackh->end(); ++imctrk) {
899  const sim::MCTrack& mctrk = *imctrk;
900  int pdg = mctrk.PdgCode();
901  if(fIgnoreSign)
902  pdg = std::abs(pdg);
903 
904  // Ignore everything except stable charged nonshowering particles.
905 
906  int apdg = std::abs(pdg);
907  if(apdg == 13 || // Muon
908  apdg == 211 || // Charged pion
909  apdg == 321 || // Charged kaon
910  apdg == 2212) { // (Anti)proton
911 
912  // check MC track origin?
913  if(fCheckOrigin && mctrk.Origin() != fOriginValue) continue;
914 
915  // Apply minimum energy cut.
916 
917  if(mctrk.Start().E() >= mctrk.Start().Momentum().Mag() + 1000.*fMinMCKE) {
918 
919  // Calculate the x offset due to nonzero mc particle time.
920 
921  double mctime = mctrk.Start().T(); // nsec
922  double mcdx = mctime * 1.e-3 * detprop->DriftVelocity(); // cm
923 
924  // Calculate the length of this mc particle inside the fiducial volume.
925 
926  TVector3 mcstart;
927  TVector3 mcend;
928  TVector3 mcstartmom;
929  TVector3 mcendmom;
930  double plen = length(mctrk, mcdx, mcstart, mcend, mcstartmom, mcendmom);
931 
932  // Apply minimum fiducial length cut. Always reject particles that have
933  // zero length in the tpc regardless of the configured cut.
934 
935  if(plen > 0. && plen > fMinMCLen) {
936 
937  // This is a good mc particle (capable of making a track).
938 
939  selected_mctracks.push_back(std::make_pair(&mctrk, -1));
940 
941  // Dump MC particle information here.
942 
943  if(pdump) {
944  double pstart = mcstartmom.Mag();
945  double pend = mcendmom.Mag();
946  *pdump << "\nOffset"
947  << std::setw(3) << mctrk.TrackID()
948  << std::setw(6) << mctrk.PdgCode()
949  << " "
950  << std::fixed << std::setprecision(2)
951  << std::setw(10) << mcdx
952  << "\nStart "
953  << std::setw(3) << mctrk.TrackID()
954  << std::setw(6) << mctrk.PdgCode()
955  << " "
956  << std::fixed << std::setprecision(2)
957  << std::setw(10) << mcstart[0]
958  << std::setw(10) << mcstart[1]
959  << std::setw(10) << mcstart[2];
960  if(pstart > 0.) {
961  *pdump << " "
962  << std::fixed << std::setprecision(3)
963  << std::setw(10) << mcstartmom[0]/pstart
964  << std::setw(10) << mcstartmom[1]/pstart
965  << std::setw(10) << mcstartmom[2]/pstart;
966  }
967  else
968  *pdump << std::setw(32) << " ";
969  *pdump << std::setw(12) << std::fixed << std::setprecision(3) << pstart;
970  *pdump << "\nEnd "
971  << std::setw(3) << mctrk.TrackID()
972  << std::setw(6) << mctrk.PdgCode()
973  << " "
974  << std::fixed << std::setprecision(2)
975  << std::setw(10) << mcend[0]
976  << std::setw(10) << mcend[1]
977  << std::setw(10) << mcend[2];
978  if(pend > 0.01) {
979  *pdump << " "
980  << std::fixed << std::setprecision(3)
981  << std::setw(10) << mcendmom[0]/pend
982  << std::setw(10) << mcendmom[1]/pend
983  << std::setw(10) << mcendmom[2]/pend;
984  }
985  else
986  *pdump << std::setw(32)<< " ";
987  *pdump << std::setw(12) << std::fixed << std::setprecision(3) << pend
988  << "\nLength: " << plen << "\n";
989  }
990 
991  // Fill histograms.
992 
993  if(fMCHistMap.count(pdg) == 0) {
994  std::ostringstream ostr;
995  ostr << "MC" << (fIgnoreSign ? "All" : (pdg > 0 ? "Pos" : "Neg")) << std::abs(pdg);
996  fMCHistMap[pdg] = MCHists(ostr.str());
997  }
998  const MCHists& mchists = fMCHistMap[pdg];
999 
1000  double mctheta_xz = std::atan2(mcstartmom.X(), mcstartmom.Z());
1001  double mctheta_yz = std::atan2(mcstartmom.Y(), mcstartmom.Z());
1002  double mcmom = mcstartmom.Mag();
1003  double mcmass = 0.001 * mctrk.Start().Momentum().Mag();
1004  double mcke = mcmom*mcmom / (std::sqrt(mcmom*mcmom + mcmass*mcmass) + mcmass);
1005 
1006  mchists.fHmcstartx->Fill(mcstart.X());
1007  mchists.fHmcstarty->Fill(mcstart.Y());
1008  mchists.fHmcstartz->Fill(mcstart.Z());
1009  mchists.fHmcendx->Fill(mcend.X());
1010  mchists.fHmcendy->Fill(mcend.Y());
1011  mchists.fHmcendz->Fill(mcend.Z());
1012  mchists.fHmctheta->Fill(mcstartmom.Theta());
1013  mchists.fHmcphi->Fill(mcstartmom.Phi());
1014  mchists.fHmctheta_xz->Fill(mctheta_xz);
1015  mchists.fHmctheta_yz->Fill(mctheta_yz);
1016  mchists.fHmcmom->Fill(mcmom);
1017  mchists.fHmcmoml->Fill(mcmom);
1018  mchists.fHmcke->Fill(mcke);
1019  mchists.fHmckel->Fill(mcke);
1020  mchists.fHmclen->Fill(plen);
1021  mchists.fHmclens->Fill(plen);
1022  }
1023  }
1024  }
1025  }
1026  } //mc
1027 
1028 
1029  // Get tracks and spacepoints and hits
1033 
1034  evt.getByLabel(fTrackModuleLabel, trackh);
1035  evt.getByLabel(fStitchModuleLabel,trackvh);
1036  evt.getByLabel(fHitModuleLabel, hith);
1037 
1038  // Extract all hits into a vector of art::Ptrs (the format used by back tracker).
1039 
1040  std::vector<art::Ptr<recob::Hit> > allhits;
1041  if(hith.isValid()) {
1042  allhits.reserve(hith->size());
1043  for(unsigned int i=0; i<hith->size(); ++i) {
1044  allhits.emplace_back(hith, i);
1045  }
1046  }
1047 
1048  // Construct FindManyP object to be used for finding track-hit associations.
1049 
1050  art::FindManyP<recob::Hit> tkhit_find(trackh, evt, fTrackModuleLabel);
1051 
1052  // This new top part of TrackAna between two long lines of ************s
1053  // is particular to analyzing Stitched Tracks.
1054  // *******************************************************************//
1055 
1056  if (trackvh.isValid() && fStitchedAnalysis)
1057  {
1058  mf::LogDebug("TrackAna")
1059  << "TrackAna read " << trackvh->size()
1060  << " vectors of Stitched PtrVectorsof tracks.";
1061  anaStitch(evt);
1062  }
1063 
1064  if(trackh.isValid()) {
1065 
1066  if(pdump) {
1067  *pdump << "\nReconstructed Tracks\n"
1068  << " Id MCid x y z dx dy dz p\n"
1069  << "-------------------------------------------------------------------------------------------\n";
1070  }
1071 
1072  // Loop over tracks.
1073 
1074  int ntrack = trackh->size();
1075  for(int i = 0; i < ntrack; ++i) {
1076  art::Ptr<recob::Track> ptrack(trackh, i);
1077  const recob::Track& track = *ptrack;
1078 
1079  // Extract hits associated with this track.
1080 
1081  std::vector<art::Ptr<recob::Hit> > trackhits;
1082  tkhit_find.get(i, trackhits);
1083 
1084  // Calculate the x offset due to nonzero reconstructed time.
1085 
1086  //double recotime = track.Time() * detprop->SamplingRate(); // nsec
1087 // double recotime = 0.;
1088 // double trackdx = recotime * 1.e-3 * detprop->DriftVelocity(); // cm
1089  double trackdx = 0;
1090 
1091  // Fill histograms involving reco tracks only.
1092 
1093  int ntraj = track.NumberTrajectoryPoints();
1094  if(ntraj > 0) {
1095  TVector3 pos = track.Vertex();
1096  TVector3 dir = track.VertexDirection();
1097  TVector3 end = track.End();
1098  pos[0] += trackdx;
1099  end[0] += trackdx;
1100 
1101  double dpos = bdist(pos);
1102  double dend = bdist(end);
1103  double tlen = length(track);
1104  double theta_xz = std::atan2(dir.X(), dir.Z());
1105  double theta_yz = std::atan2(dir.Y(), dir.Z());
1106 
1107  if(fRecoHistMap.count(0) == 0)
1108  fRecoHistMap[0] = RecoHists("Reco");
1109  const RecoHists& rhists = fRecoHistMap[0];
1110 
1111  rhists.fHstartx->Fill(pos.X());
1112  rhists.fHstarty->Fill(pos.Y());
1113  rhists.fHstartz->Fill(pos.Z());
1114  rhists.fHstartd->Fill(dpos);
1115  rhists.fHendx->Fill(end.X());
1116  rhists.fHendy->Fill(end.Y());
1117  rhists.fHendz->Fill(end.Z());
1118  rhists.fHendd->Fill(dend);
1119  rhists.fHtheta->Fill(dir.Theta());
1120  rhists.fHphi->Fill(dir.Phi());
1121  rhists.fHtheta_xz->Fill(theta_xz);
1122  rhists.fHtheta_yz->Fill(theta_yz);
1123 
1124  double mom = 0.;
1125  if(track.HasMomentum())
1126  mom = track.VertexMomentum();
1127  rhists.fHmom->Fill(mom);
1128  rhists.fHmoml->Fill(mom);
1129  rhists.fHlen->Fill(tlen);
1130  rhists.fHlens->Fill(tlen);
1131 
1132  // Id of matching mc particle.
1133 
1134  int mcid = -1;
1135 
1136  // Loop over direction.
1137 
1138  for(int swap=0; swap<2; ++swap) {
1139 
1140  // Analyze reversed tracks only if start momentum = end momentum.
1141 
1142  if(swap != 0 && track.HasMomentum() &&
1143  std::abs(track.VertexMomentum() - track.EndMomentum()) > 1.e-3)
1144  continue;
1145 
1146  // Calculate the global-to-local rotation matrix.
1147 
1148  TMatrixD rot(3,3);
1149  int start_point = (swap == 0 ? 0 : ntraj-1);
1150  track.GlobalToLocalRotationAtPoint(start_point, rot);
1151 
1152  // Update track data for reversed case.
1153 
1154  if(swap != 0) {
1155  rot(1, 0) = -rot(1, 0);
1156  rot(2, 0) = -rot(2, 0);
1157  rot(1, 1) = -rot(1, 1);
1158  rot(2, 1) = -rot(2, 1);
1159  rot(1, 2) = -rot(1, 2);
1160  rot(2, 2) = -rot(2, 2);
1161 
1162  pos = track.End();
1163  dir = -track.EndDirection();
1164  end = track.Vertex();
1165  pos[0] += trackdx;
1166  end[0] += trackdx;
1167 
1168  dpos = bdist(pos);
1169  dend = bdist(end);
1170  theta_xz = std::atan2(dir.X(), dir.Z());
1171  theta_yz = std::atan2(dir.Y(), dir.Z());
1172 
1173  if(track.HasMomentum()) mom = track.EndMomentum();
1174  }
1175 
1176  // Get covariance matrix.
1177 
1178  const TMatrixD& cov = (swap == 0 ? track.VertexCovariance() : track.EndCovariance());
1179 
1180  // Loop over track-like mc particles.
1181 
1182  for(unsigned int imc = 0; imc < selected_mctracks.size(); ++imc) {
1183  const sim::MCTrack& mctrk = *selected_mctracks[imc].first;
1184  int pdg = mctrk.PdgCode();
1185  if(fIgnoreSign) pdg = std::abs(pdg);
1186  auto iMCHistMap = fMCHistMap.find(pdg);
1187  if (iMCHistMap == fMCHistMap.end())
1188  throw cet::exception("TrackAna") << "no particle with ID=" << pdg << "\n";
1189  const MCHists& mchists = iMCHistMap->second;
1190 
1191  // Calculate the x offset due to nonzero mc particle time.
1192 
1193  double mctime = mctrk.Start().T(); // nsec
1194  double mcdx = mctime * 1.e-3 * detprop->DriftVelocity(); // cm
1195 
1196  // Calculate the points where this mc particle enters and leaves the
1197  // fiducial volume, and the length in the fiducial volume.
1198 
1199  TVector3 mcstart;
1200  TVector3 mcend;
1201  TVector3 mcstartmom;
1202  TVector3 mcendmom;
1203  double plen = length(mctrk, mcdx, mcstart, mcend, mcstartmom, mcendmom);
1204 
1205  // Get the displacement of this mc particle in the global coordinate system.
1206 
1207  TVector3 mcpos = mcstart - pos;
1208 
1209  // Rotate the momentum and position to the
1210  // track-local coordinate system.
1211 
1212  TVector3 mcmoml = rot * mcstartmom;
1213  TVector3 mcposl = rot * mcpos;
1214 
1215  double colinearity = mcmoml.Z() / mcmoml.Mag();
1216 
1217  double u = mcposl.X();
1218  double v = mcposl.Y();
1219  double w = mcposl.Z();
1220 
1221  double pu = mcmoml.X();
1222  double pv = mcmoml.Y();
1223  double pw = mcmoml.Z();
1224 
1225  double dudw = pu / pw;
1226  double dvdw = pv / pw;
1227 
1228  double u0 = u - w * dudw;
1229  double v0 = v - w * dvdw;
1230  double uv0 = std::sqrt(u0*u0 + v0*v0);
1231 
1232  mchists.fHduvcosth->Fill(colinearity, uv0);
1233 
1234  if(std::abs(uv0) < fMatchDisp) {
1235 
1236  // Fill slope matching histograms.
1237 
1238  mchists.fHmcdudw->Fill(dudw);
1239  mchists.fHmcdvdw->Fill(dvdw);
1240  mchists.fHdudwpull->Fill(dudw / std::sqrt(cov(2,2)));
1241  mchists.fHdvdwpull->Fill(dvdw / std::sqrt(cov(3,3)));
1242  }
1243  mchists.fHcosth->Fill(colinearity);
1244  if(colinearity > fMatchColinearity) {
1245 
1246  // Fill displacement matching histograms.
1247 
1248  mchists.fHmcu->Fill(u0);
1249  mchists.fHmcv->Fill(v0);
1250  mchists.fHmcw->Fill(w);
1251  mchists.fHupull->Fill(u0 / std::sqrt(cov(0,0)));
1252  mchists.fHvpull->Fill(v0 / std::sqrt(cov(1,1)));
1253 
1254  if(std::abs(uv0) < fMatchDisp) {
1255 
1256  // Fill matching histograms.
1257 
1258  double mctheta_xz = std::atan2(mcstartmom.X(), mcstartmom.Z());
1259  double mctheta_yz = std::atan2(mcstartmom.Y(), mcstartmom.Z());
1260  double mcmom = mcstartmom.Mag();
1261  double mcmass = 0.001 * mctrk.Start().Momentum().Mag();
1262  double mcke = mcmom*mcmom / (std::sqrt(mcmom*mcmom + mcmass*mcmass) + mcmass);
1263 
1264  mchists.fHstartdx->Fill(pos.X() - mcstart.X());
1265  mchists.fHstartdy->Fill(pos.Y() - mcstart.Y());
1266  mchists.fHstartdz->Fill(pos.Z() - mcstart.Z());
1267  mchists.fHenddx->Fill(end.X() - mcend.X());
1268  mchists.fHenddy->Fill(end.Y() - mcend.Y());
1269  mchists.fHenddz->Fill(end.Z() - mcend.Z());
1270  mchists.fHlvsl->Fill(plen, tlen);
1271  mchists.fHdl->Fill(tlen - plen);
1272  mchists.fHpvsp->Fill(mcmom, mom);
1273  double dp = mom - mcmom;
1274  mchists.fHdp->Fill(dp);
1275  mchists.fHppull->Fill(dp / std::sqrt(cov(4,4)));
1276  if(std::abs(dpos) >= 5. && std::abs(dend) >= 5.) {
1277  mchists.fHpvspc->Fill(mcmom, mom);
1278  mchists.fHdpc->Fill(dp);
1279  mchists.fHppullc->Fill(dp / std::sqrt(cov(4,4)));
1280  }
1281 
1282  // Count this track as well-reconstructed if it is matched to an
1283  // mc particle (which it is if get here), and if in addition the
1284  // starting position (w) matches and the reconstructed track length
1285  // is more than 0.5 of the mc particle trajectory length.
1286 
1287  bool good = std::abs(w) <= fWMatchDisp &&
1288  tlen > fMatchLength * plen;
1289  if(good) {
1290  mcid = mctrk.TrackID();
1291 
1292  // Calculate and fill hit efficiency and purity.
1293 
1294  std::set<int> tkidset;
1295  tkidset.insert(mcid);
1296  double hiteff =
1297  bt_serv->HitCollectionEfficiency(tkidset, trackhits, allhits, geo::k3D);
1298  double hitpurity = bt_serv->HitCollectionPurity(tkidset, trackhits);
1299  mchists.fHHitEff->Fill(hiteff);
1300  mchists.fHHitPurity->Fill(hitpurity);
1301 
1302  // Fill efficiency numerator histograms.
1303 
1304  mchists.fHgstartx->Fill(mcstart.X());
1305  mchists.fHgstarty->Fill(mcstart.Y());
1306  mchists.fHgstartz->Fill(mcstart.Z());
1307  mchists.fHgendx->Fill(mcend.X());
1308  mchists.fHgendy->Fill(mcend.Y());
1309  mchists.fHgendz->Fill(mcend.Z());
1310  mchists.fHgtheta->Fill(mcstartmom.Theta());
1311  mchists.fHgphi->Fill(mcstartmom.Phi());
1312  mchists.fHgtheta_xz->Fill(mctheta_xz);
1313  mchists.fHgtheta_yz->Fill(mctheta_yz);
1314  mchists.fHgmom->Fill(mcmom);
1315  mchists.fHgmoml->Fill(mcmom);
1316  mchists.fHgke->Fill(mcke);
1317  mchists.fHgkel->Fill(mcke);
1318  mchists.fHglen->Fill(plen);
1319  mchists.fHglens->Fill(plen);
1320 
1321  // set the match flag
1322  selected_mctracks[imc].second = i;
1323 
1324  if(fPrintLevel > 0) {
1325  const simb::MCParticle* ptkl = pi_serv->TrackIdToParticle_P(mcid);
1326  float KE = ptkl->E() - ptkl->Mass();
1327  std::string KEUnits = " GeV";
1328  if(mctrk.Origin() != simb::kCosmicRay) {
1329  // MeV for low energy particles
1330  KE *= 1000;
1331  KEUnits = " MeV";
1332  }
1333  mf::LogVerbatim("TrackAna")
1334  <<evt.run()<<"."<<evt.event()
1335  <<" Match MCTkID "<<std::setw(6)<<mctrk.TrackID()
1336  <<" Origin "<<mctrk.Origin()
1337  <<" PDG"<<std::setw(5)<<mctrk.PdgCode()
1338  <<" KE"<<std::setw(4)<<(int)KE<<KEUnits
1339  <<" RecoTrkID "<<track.ID()
1340  <<" hitEff "<<std::setprecision(2)<<hiteff<<" hitPur "<<hitpurity;
1341  int sWire, sTick, eWire, eTick;
1342  // this won't work for DUNE
1343  for(unsigned short ipl = 0; ipl < geom->Nplanes(); ++ipl) {
1344  sWire = geom->NearestWire(mcstart, ipl, 0, 0);
1345  sTick = detprop->ConvertXToTicks(mcstart[0], ipl, 0, 0);
1346  eWire = geom->NearestWire(mcend, ipl, 0, 0);
1347  eTick = detprop->ConvertXToTicks(mcend[0], ipl, 0, 0);
1348  mf::LogVerbatim("TrackAna")<<" Wire:Tick in Pln "<<ipl<<" W:T "<<sWire<<":"<<sTick<<" - "<<eWire<<":"<<eTick;
1349  } // ipl
1350  } // fPrintLevel == 2
1351  } // good
1352  }
1353  }
1354  }
1355  }
1356 
1357  // Dump track information here.
1358 
1359  if(pdump) {
1360  TVector3 pos = track.Vertex();
1361  TVector3 dir = track.VertexDirection();
1362  TVector3 end = track.End();
1363  pos[0] += trackdx;
1364  end[0] += trackdx;
1365  TVector3 enddir = track.EndDirection();
1366  double pstart = track.VertexMomentum();
1367  double pend = track.EndMomentum();
1368  *pdump << "\nOffset"
1369  << std::setw(3) << track.ID()
1370  << std::setw(6) << mcid
1371  << " "
1372  << std::fixed << std::setprecision(2)
1373  << std::setw(10) << trackdx
1374  << "\nStart "
1375  << std::setw(3) << track.ID()
1376  << std::setw(6) << mcid
1377  << " "
1378  << std::fixed << std::setprecision(2)
1379  << std::setw(10) << pos[0]
1380  << std::setw(10) << pos[1]
1381  << std::setw(10) << pos[2];
1382  if(pstart > 0.) {
1383  *pdump << " "
1384  << std::fixed << std::setprecision(3)
1385  << std::setw(10) << dir[0]
1386  << std::setw(10) << dir[1]
1387  << std::setw(10) << dir[2];
1388  }
1389  else
1390  *pdump << std::setw(32) << " ";
1391  *pdump << std::setw(12) << std::fixed << std::setprecision(3) << pstart;
1392  *pdump << "\nEnd "
1393  << std::setw(3) << track.ID()
1394  << std::setw(6) << mcid
1395  << " "
1396  << std::fixed << std::setprecision(2)
1397  << std::setw(10) << end[0]
1398  << std::setw(10) << end[1]
1399  << std::setw(10) << end[2];
1400  if(pend > 0.01) {
1401  *pdump << " "
1402  << std::fixed << std::setprecision(3)
1403  << std::setw(10) << enddir[0]
1404  << std::setw(10) << enddir[1]
1405  << std::setw(10) << enddir[2];
1406  }
1407  else
1408  *pdump << std::setw(32)<< " ";
1409  *pdump << std::setw(12) << std::fixed << std::setprecision(3) << pend
1410  << "\nLength: " << tlen << "\n";
1411  }
1412  }
1413  }
1414  } // i
1415 
1416  // print out un-matched MC tracks
1417  if(fPrintLevel > 0) {
1418  for(unsigned int imc = 0; imc < selected_mctracks.size(); ++imc) {
1419  if(selected_mctracks[imc].second >= 0) continue;
1420  const sim::MCTrack& mctrk = *selected_mctracks[imc].first;
1421  const simb::MCParticle* ptkl = pi_serv->TrackIdToParticle_P(mctrk.TrackID());
1422  float KE = ptkl->E() - ptkl->Mass();
1423  std::string KEUnits = " GeV";
1424  if(mctrk.Origin() != simb::kCosmicRay) {
1425  // MeV for low energy particles
1426  KE *= 1000;
1427  KEUnits = " MeV";
1428  }
1429  // find the start/end wire:time in each plane
1430  TVector3 mcstart, mcend, mcstartmom, mcendmom;
1431  double mcdx = mctrk.Start().T() * 1.e-3 * detprop->DriftVelocity(); // cm
1432  double plen = length(mctrk, mcdx, mcstart, mcend, mcstartmom, mcendmom);
1433  mf::LogVerbatim("TrackAna")<<evt.run()<<"."<<evt.event()
1434  <<" NoMat MCTkID "<<std::setw(6)<<mctrk.TrackID()
1435  <<" Origin "<<mctrk.Origin()
1436  <<" PDG"<<std::setw(5)<<mctrk.PdgCode()
1437  <<" KE"<<std::setw(4)<<(int)KE<<KEUnits
1438  <<" Length "<<std::fixed<<std::setprecision(1)<<plen<<" cm";
1439  if(fPrintLevel > 1) {
1440  int sWire, sTick, eWire, eTick;
1441  // this won't work for DUNE
1442  for(unsigned short ipl = 0; ipl < geom->Nplanes(); ++ipl) {
1443  sWire = geom->NearestWire(mcstart, ipl, 0, 0);
1444  sTick = detprop->ConvertXToTicks(mcstart[0], ipl, 0, 0);
1445  eWire = geom->NearestWire(mcend, ipl, 0, 0);
1446  eTick = detprop->ConvertXToTicks(mcend[0], ipl, 0, 0);
1447  mf::LogVerbatim("TrackAna")<<" Wire:Tick in Pln "<<ipl
1448  <<" W:T "<<sWire<<":"<<sTick<<" - "<<eWire<<":"<<eTick;
1449  } // ipl
1450  } // fPrintLevel > 1
1451  } // imc
1452  } // fPrintLevel > 0
1453 
1454  }
1455 
1457  {
1458 
1462  const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
1463 
1464  std::map<int, std::map<int, art::PtrVector<recob::Hit>> > hitmap; // trkID, otrk, hitvec
1465  std::map<int, int > KEmap; // length traveled in det [cm]?, trkID want to sort by KE
1466  bool mc = !evt.isRealData();
1470 
1471  evt.getByLabel(fTrackModuleLabel, trackh);
1472  evt.getByLabel(fSpacepointModuleLabel, sppth);
1473  evt.getByLabel(fStitchModuleLabel,trackvh);
1474  int ntv(trackvh->size());
1475 
1476  std::vector < art::PtrVector<recob::Track> >::const_iterator cti = trackvh->begin();
1478 
1479  if(trackh.isValid()) {
1481  int nsppts_assnwhole = fswhole.size();
1482  std::cout << "TrackAna: Number of clumps of Spacepoints from Assn for all Tracks: " << nsppts_assnwhole << std::endl;
1483  }
1484 
1485  if(fRecoHistMap.count(0) == 0)
1486  {
1487  fRecoHistMap[0] = RecoHists("Reco");
1488  std::cout << "\n" << "\t\t TrkAna: Fresh fRecoHistMap[0] ******* \n" << std::endl;
1489  }
1490  const RecoHists& rhistsStitched = fRecoHistMap[0];
1491 
1492  std::vector < std::vector <unsigned int> > NtrkIdsAll;
1493  std::vector < double > ntvsorted;
1494  hitmap.clear();
1495  KEmap.clear();
1496 
1497  // Look at the components of the stitched tracks. Grab their sppts/hits from Assns.
1498  for (int o = 0; o < ntv; ++o) // o for outer
1499  {
1500 
1501  const art::PtrVector<recob::Track> pvtrack(*(cti++));
1502  // auto it = pvtrack.begin();
1503  int ntrack = pvtrack.size();
1504  // if (ntrack>1) std::cout << "\t\t TrkAna: New Stitched Track ******* " << std::endl;
1505  std::vector< std::vector <unsigned int> > NtrkId_Hit; // hit IDs in inner tracks
1506  std::vector<unsigned int> vecMode;
1508 
1509  for(int i = 0; i < ntrack; ++i) {
1510 
1511  //const art::Ptr<recob::Track> ptrack(*(it++));
1512  // const recob::Track& track = *ptrack;
1513  // auto pcoll { ptrack };
1514  // art::FindManyP<recob::SpacePoint> fs( ptrack, evt, fTrkSpptAssocModuleLabel);
1515  // From gdb> ptype fs, the vector of Ptr<SpacePoint>s it appears is grabbed after fs.at(0)
1516  bool assns(true);
1517  try {
1518  // Got Spacepoints from this Track; now get Hits from those Spacepoints.
1519  // int nsppts = ptrack->NumberTrajectoryPoints();
1520 
1521  int nsppts_assn = fs.at(i).size();
1522  // if (ntrack>1) std::cout << "\t\tTrackAna: Number of Spacepoints from Track.NumTrajPts(): " << nsppts << std::endl;
1523  // if (ntrack>1) std::cout << "\t\tTrackAna: Number of Spacepoints from Assns for this Track: " << nsppts_assn << std::endl;
1524 
1525  const auto& sppt = fs.at(i);//.at(0);
1526  // since we're in a try and worried about failure, we won't pull the following
1527  // FindManyP out of the loop.
1529 
1530  // Importantly, loop on all sppts, though they don't all contribute to the track.
1531  // As opposed to looping on the trajectory pts, which is a lower number.
1532  // Also, important, in job in whch this runs I set TrackKal3DSPS parameter MaxPass=1,
1533  // cuz I don't want merely the sparse set of sppts as follows from the uncontained
1534  // MS-measurement in 2nd pass.
1535  std::vector <unsigned int> vecNtrkIds;
1536  for(int is = 0; is < nsppts_assn; ++is) {
1537  int nhits = fh.at(is).size(); // should be 2 or 3: number of planes.
1538  for(int ih = 0; ih < nhits; ++ih) {
1539  const auto& hit = fh.at(is).at(ih); // Our vector is after the .at(is) this time.
1540  if (hit->SignalType()!=geo::kCollection) continue;
1541  rhistsStitched.fHHitChg->Fill(hit->Integral());
1542  rhistsStitched.fHHitWidth->Fill(2. * hit->RMS());
1543  if (mc)
1544  {
1545  std::vector<sim::TrackIDE> tids = bt_serv->HitToTrackIDEs(hit);
1546  // more here.
1547  // Loop over track ids.
1548  bool justOne(true); // Only take first trk that contributed to this hit
1549  // std::cout << "\t\t TrkAna: TrkId tids.size() ******* " << tids.size() <<std::endl;
1550  for(std::vector<sim::TrackIDE>::const_iterator itid = tids.begin();itid != tids.end(); ++itid) {
1551  int trackID = std::abs(itid->trackID);
1552  hitmap[trackID][o].push_back(hit);
1553 
1554  if (justOne) { vecNtrkIds.push_back(trackID); justOne=false; }
1555  // Add hit to PtrVector corresponding to this track id.
1556  rhistsStitched.fHHitTrkId->Fill(trackID);
1557  const simb::MCParticle* part = pi_serv->TrackIdToParticle_P(trackID);
1558  if (!part) break;
1559 
1560  rhistsStitched.fHHitPdg->Fill(part->PdgCode());
1561  // This really needs to be indexed as KE deposited in volTPC, not just KE. EC, 24-July-2014.
1562 
1563  TVector3 mcstart;
1564  TVector3 mcend;
1565  TVector3 mcstartmom;
1566  TVector3 mcendmom;
1567  double mctime = part->T(); // nsec
1568  double mcdx = mctime * 1.e-3 * detprop->DriftVelocity(); // cm
1569 
1570  double plen = length(*part, mcdx, mcstart, mcend, mcstartmom, mcendmom);
1571 
1572  KEmap[(int)(1e6*plen)] = trackID; // multiple assignment but always the same, so fine.
1573  // std::cout << "\t\t TrkAna: TrkId trackID, KE [MeV] ******* " << trackID << ", " << (int)(1e3*(part->E()-part->Mass())) <<std::endl;
1574  }
1575 
1576  } // mc
1577  } // hits
1578 
1579  } // spacepoints
1580 
1581  if (mc)
1582  {
1583  NtrkId_Hit.push_back(vecNtrkIds);
1584  // Find the trkID mode for this i^th track
1585  unsigned int ii(1);
1586  int max(-12), n(1), ind(0);
1587  std::sort(vecNtrkIds.begin(),vecNtrkIds.end());
1588  std::vector<unsigned int> strkIds(vecNtrkIds);
1589  while ( ii < vecNtrkIds.size() )
1590  {
1591  if (strkIds.at(ii) != strkIds.at(ii-1))
1592  {
1593  n=1;
1594  }
1595  else
1596  {
1597  n++;
1598  }
1599  if (n>max) { max = n; ind = ii;}
1600  ii++;
1601  }
1602  // std::cout << "\t\t TrkAna: TrkId ind for this track is ******* " << ind <<std::endl;
1603  unsigned int mode(sim::NoParticleId);
1604  if (strkIds.begin()!=strkIds.end())
1605  mode = strkIds.at(ind);
1606  vecMode.push_back(mode);
1607 
1608  // if (ntrack>1) std::cout << "\t\t TrkAna: TrkId mode for this component track is ******* " << mode <<std::endl;
1609  if (strkIds.size()!=0)
1610  rhistsStitched.fModeFrac->Fill((double)max/(double)strkIds.size());
1611  else
1612  rhistsStitched.fModeFrac->Fill(-1.0);
1613  } // mc
1614 
1615  } // end try
1616  catch (cet::exception& x) {
1617  assns = false;
1618  }
1619  if (!assns) throw cet::exception("TrackAna") << "Bad Associations. \n";
1620 
1621  } // i
1622 
1623  if (mc)
1624  {
1625  // one vector per o trk, for all modes of stitched i trks
1626  NtrkIdsAll.push_back(vecMode);
1627 
1628  std::unique(NtrkIdsAll.back().begin(),NtrkIdsAll.back().end());
1629  double sum(0.0);
1630  for (auto const val : NtrkIdsAll.back())
1631  {
1632  // rhistsStitched.fNTrkIdTrks3->Fill(o,val%100,hitmap[val][o].size());
1633  sum += hitmap[val][o].size();
1634  }
1635  ntvsorted.push_back(sum);
1636 
1637  }
1638 
1639  //
1640  } // o
1641 
1642  int vtmp(0);
1643  // get KEmap indices by most energetic first, least last.
1644  for (auto it = KEmap.rbegin(); it!=KEmap.rend(); ++it)
1645  {
1646  // int tval = it->second; // grab trkIDs in order, since they're sorted by KE
1647  // int ke = it->first; // grab trkIDs in order, since they're sorted by KE
1648  // const simb::MCParticle* part = bt_serv->TrackIDToParticle(tval);
1649 
1650  // std::cout << "TrackAnaStitch: KEmap cntr vtmp, Length ke, tval, pdg : " << vtmp << ", " << ke <<", " << tval <<", " << part->PdgCode() << ", " << std::endl;
1651 
1652  vtmp++;
1653  }
1654 
1655  // get o trk indices by their hits. Most populous track first, least last.
1656  for (auto const o : fsort_indexes(ntvsorted))
1657  {
1658  int v(0);
1659  // get KEmap indices by longest trajectory first, least last.
1660  for (auto it = KEmap.rbegin(); it!=KEmap.rend(); ++it)
1661  {
1662  int val = it->second; // grab trkIDs in order, since they're sorted by KE
1663  // const simb::MCParticle* part = pi_serv->TrackIDToParticle(val);
1664  // std::cout << "TrackAnaStitch: trk o, KEmap cntr v, KE val, pdg hitmap[val][o].size(): " << o <<", " << v << ", " << val <<", " << part->PdgCode() << ", " << hitmap[val][o].size() << std::endl;
1665  rhistsStitched.fNTrkIdTrks3->Fill(o,v,hitmap[val][o].size());
1666  v++;
1667  }
1668 
1669  }
1670 
1671  // In how many o tracks did each trkId appear? Histo it. Would like it to be precisely 1.
1672  // Histo it vs. particle KE.
1673  flattener flat(NtrkIdsAll);
1674  std::vector <unsigned int> &v = flat;
1675  // auto const it ( std::unique(v.begin(),v.end()) ); // never use this it, perhaps.
1676  for (auto const val : v)
1677  {
1678  if (val != (unsigned int)sim::NoParticleId)
1679  {
1680  const simb::MCParticle* part = pi_serv->TrackIdToParticle_P( val );
1681  double T(part->E() - 0.001*part->Mass());
1682  rhistsStitched.fNTrkIdTrks->Fill(std::count(v.begin(),v.end(),val));
1683  rhistsStitched.fNTrkIdTrks2->Fill(std::count(v.begin(),v.end(),val),T);
1684  }
1685  else
1686  {
1687  rhistsStitched.fNTrkIdTrks2->Fill(-1.0,0.0);
1688  }
1689  }
1690 
1691  }
1692 
1694  //
1695  // Purpose: End of job.
1696  //
1697  {
1698  // Print summary.
1699 
1700  mf::LogInfo("TrackAna")
1701  << "TrackAna statistics:\n"
1702  << " Number of events = " << fNumEvent;
1703 
1704  // Fill efficiency histograms.
1705 
1707  i != fMCHistMap.end(); ++i) {
1708  const MCHists& mchists = i->second;
1709  effcalc(mchists.fHgstartx, mchists.fHmcstartx, mchists.fHestartx);
1710  effcalc(mchists.fHgstarty, mchists.fHmcstarty, mchists.fHestarty);
1711  effcalc(mchists.fHgstartz, mchists.fHmcstartz, mchists.fHestartz);
1712  effcalc(mchists.fHgendx, mchists.fHmcendx, mchists.fHeendx);
1713  effcalc(mchists.fHgendy, mchists.fHmcendy, mchists.fHeendy);
1714  effcalc(mchists.fHgendz, mchists.fHmcendz, mchists.fHeendz);
1715  effcalc(mchists.fHgtheta, mchists.fHmctheta, mchists.fHetheta);
1716  effcalc(mchists.fHgphi, mchists.fHmcphi, mchists.fHephi);
1717  effcalc(mchists.fHgtheta_xz, mchists.fHmctheta_xz, mchists.fHetheta_xz);
1718  effcalc(mchists.fHgtheta_yz, mchists.fHmctheta_yz, mchists.fHetheta_yz);
1719  effcalc(mchists.fHgmom, mchists.fHmcmom, mchists.fHemom);
1720  effcalc(mchists.fHgmoml, mchists.fHmcmoml, mchists.fHemoml);
1721  effcalc(mchists.fHgke, mchists.fHmcke, mchists.fHeke);
1722  effcalc(mchists.fHgkel, mchists.fHmckel, mchists.fHekel);
1723  effcalc(mchists.fHglen, mchists.fHmclen, mchists.fHelen);
1724  effcalc(mchists.fHglens, mchists.fHmclens, mchists.fHelens);
1725  }
1726  }
1727 
1728  // Stole this from online. Returns indices sorted by corresponding vector values.
1729  template <typename T>
1730  std::vector<size_t> TrackAna::fsort_indexes(const std::vector<T> &v) {
1731  // initialize original index locations
1732  std::vector<size_t> idx(v.size());
1733  for (size_t i = 0; i != idx.size(); ++i) idx[i] = i;
1734  // sort indexes based on comparing values in v
1735  std::sort(idx.begin(), idx.end(),
1736  [&v](size_t i1, size_t i2) {return v[i1] > v[i2];}); // Want most occupied trks first. was <, EC, 21-July.
1737  return idx;
1738  }
1739 
1740 
1741 }
Float_t x
Definition: compare.C:6
double E(const int i=0) const
Definition: MCParticle.h:237
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
const simb::MCParticle * TrackIdToParticle_P(int const &id)
TVector3 LocationAtPoint(unsigned int p) const
Covariance matrices are either set or not.
Definition: Track.h:241
unsigned int NumberTrajectoryPoints() const
Definition: MCParticle.h:222
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:223
simb::Origin_t Origin() const
Definition: MCTrack.h:40
double EndMomentum() const
Definition: Track.h:151
int PdgCode() const
Definition: MCParticle.h:216
double VertexMomentum() const
Definition: Track.h:149
std::vector< size_t > fsort_indexes(const std::vector< T > &v)
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
std::string fHitSpptAssocModuleLabel
TVector3 VertexDirection() const
Covariance matrices are either set or not.
Definition: Track.h:247
std::string fTrkSpptAssocModuleLabel
TMatrixD EndCovariance() const
Covariance matrices are either set or not.
Definition: Track.cxx:142
const std::vector< sim::TrackIDE > HitToTrackIDEs(recob::Hit const &hit)
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Float_t den
Definition: plot.C:37
std::string fSpacepointModuleLabel
Declaration of signal hit object.
double Mass() const
Definition: MCParticle.h:243
enum simb::_ev_origin Origin_t
event origin types
bool HasMomentum() const
Definition: Track.h:147
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:109
STL namespace.
std::string fStitchModuleLabel
TrackAna(fhicl::ParameterSet const &pset)
double T() const
Definition: MCStep.h:45
bool isRealData() const
Definition: Event.h:83
Particle class.
void anaStitch(const art::Event &evt)
TFileDirectory mkdir(std::string const &dir, std::string const &descr="")
std::map< int, MCHists > fMCHistMap
virtual double ConvertXToTicks(double X, int p, int t, int c) const =0
3-dimensional objects, potentially hits, clusters, prongs, etc.
Definition: geo_types.h:82
Rotation_t GlobalToLocalRotationAtPoint(size_t p) const
Definition: Track.h:197
bool isValid() const
Definition: Handle.h:190
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
Int_t max
Definition: plot.C:27
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
intermediate_table::const_iterator const_iterator
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
static const int NoParticleId
Definition: sim.h:28
single particles thrown at the detector
Definition: MCTruth.h:24
geo::WireID::WireID_t NearestWire(geo::Point_t const &point, geo::PlaneID const &planeid) const
Returns the index of wire closest to position in the specified TPC.
TString part[npart]
Definition: Style.C:32
double T(const int i=0) const
Definition: MCParticle.h:228
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
std::string fMCTrackModuleLabel
EventNumber_t event() const
Definition: Event.h:67
EDAnalyzer(Table< Config > const &config)
Definition: EDAnalyzer.h:100
void analyze(const art::Event &evt)
Class def header for mctrack data container.
Provides recob::Track data product.
size_type size() const
Definition: PtrVector.h:308
int PdgCode() const
Definition: MCTrack.h:41
Detector simulation of raw signals on wires.
const TLorentzVector & Momentum() const
Definition: MCStep.h:38
simb::Origin_t fOriginValue
int ID() const
Definition: Track.h:205
T * make(ARGS...args) const
std::string value(boost::any const &)
std::string fOrigin
TDirectory * dir
Definition: macro.C:5
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
Int_t min
Definition: plot.C:26
virtual double DriftVelocity(double efield=0., double temperature=0.) const =0
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
std::map< int, RecoHists > fRecoHistMap
TVector3 Vertex() const
Covariance matrices are either set or not.
Definition: Track.h:245
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:224
TMatrixD VertexCovariance() const
Covariance matrices are either set or not.
Definition: Track.cxx:131
Supernova neutrinos.
Definition: MCTruth.h:23
double E() const
Definition: MCStep.h:49
TVector3 EndDirection() const
Covariance matrices are either set or not.
Definition: Track.h:248
const MCStep & Start() const
Definition: MCTrack.h:44
Char_t n[5]
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
size_type get(size_type i, reference item, data_reference data) const
Definition: FindManyP.h:469
std::string fTrackModuleLabel
unsigned int TrackID() const
Definition: MCTrack.h:42
Float_t e
Definition: plot.C:34
RunNumber_t run() const
Definition: Event.h:77
Float_t track
Definition: plot.C:34
TVector3 End() const
Covariance matrices are either set or not.
Definition: Track.h:246
Float_t w
Definition: plot.C:23
std::string fHitModuleLabel
Tools and modules for checking out the basics of the Monte Carlo.
art framework interface to geometry description
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:51
vec_iX clear()
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Cosmic rays.
Definition: MCTruth.h:22
Signal from collection planes.
Definition: geo_types.h:93
Beam neutrinos.
Definition: MCTruth.h:21