LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
SeedAna_module.cc
Go to the documentation of this file.
1 //
2 // Name: SeedAna_module.cc
3 //
4 // Purpose: Module SeedAna.
5 //
6 // Configuration parameters.
7 //
8 // SeedModuleLabel: Seed module label.
9 // MCTrackModuleLabel: MCTrack module label.
10 // MinMCKE: Minimum MC particle kinetic energy.
11 // MatchColinearity: Minimum colinearity for mc-seed matching.
12 // MatchDisp: Maximum uv displacement for mc-seed matching.
13 //
14 // Created: 2-Aug-2011 H. Greenlee
15 //
16 
17 #include <cmath>
18 #include <iomanip>
19 #include <map>
20 #include <sstream>
21 
25 #include "art_root_io/TFileService.h"
27 #include "cetlib_except/exception.h"
29 
34 
35 #include "TH1F.h"
36 #include "TH2F.h"
37 #include "TMatrixD.h"
38 
39 namespace {
40 
41  // Calculate distance to boundary.
42  //----------------------------------------------------------------------------
43  double bdist(const TVector3& pos)
44  {
45  // Get geometry.
46 
48 
49  double d1 = pos.X(); // Distance to right side (wires).
50  double d2 = 2. * geom->DetHalfWidth() - pos.X(); // Distance to left side (cathode).
51  double d3 = pos.Y() + geom->DetHalfHeight(); // Distance to bottom.
52  double d4 = geom->DetHalfHeight() - pos.Y(); // Distance to top.
53  double d5 = pos.Z(); // Distance to front.
54  double d6 = geom->DetLength() - pos.Z(); // Distance to back.
55 
56  return std::min({d1, d2, d3, d4, d5, d6});
57  }
58 
59  // Find the closest matching mc trajectory point (sim::MCStep) for a given seed.
60  // Returned value is index of the trajectory point.
61  // Return -1 in case of no match.
62  int mcmatch(detinfo::DetectorPropertiesData const& detProp,
63  const sim::MCTrack& mctrk,
64  const recob::Seed& seed)
65  {
66  // Get seed point.
67 
68  double pos[3];
69  double err[3];
70  seed.GetPoint(pos, err);
71 
72  // Calculate the x offset due to nonzero mc particle time.
73  double mctime = mctrk.Start().T(); // nsec
74  double mcdx = mctime * 1.e-3 * detProp.DriftVelocity(); // cm
75 
76  // Loop over trajectory points.
77 
78  int best_traj = -1;
79  double max_dist = 0.;
80  int ntraj = mctrk.size();
81  for (int itraj = 0; itraj < ntraj; ++itraj) {
82  const TLorentzVector& vec = mctrk[itraj].Position();
83  double dx = pos[0] - vec.X() - mcdx;
84  double dy = pos[1] - vec.Y();
85  double dz = pos[2] - vec.Z();
86  double dist = std::sqrt(dx * dx + dy * dy + dz * dz);
87  if (best_traj < 0 || dist < max_dist) {
88  best_traj = itraj;
89  max_dist = dist;
90  }
91  }
92  return best_traj;
93  }
94 
95  // Length of MCTrack.
96  // In this function, the extracted start and end momenta are converted to GeV
97  // (MCTrack stores momenta in Mev).
98  //----------------------------------------------------------------------------
99  double length(detinfo::DetectorPropertiesData const& detProp,
100  const sim::MCTrack& mctrk,
101  double dx,
102  TVector3& start,
103  TVector3& end,
104  TVector3& startmom,
105  TVector3& endmom,
106  unsigned int /*tpc*/ = 0,
107  unsigned int /*cstat*/ = 0)
108  {
109  // Get services.
110 
112 
113  // Get fiducial volume boundary.
114 
115  double xmin = 0.;
116  double xmax = 2. * geom->DetHalfWidth();
117  double ymin = -geom->DetHalfHeight();
118  double ymax = geom->DetHalfHeight();
119  double zmin = 0.;
120  double zmax = geom->DetLength();
121  double result = 0.;
122  TVector3 disp;
123  int n = mctrk.size();
124  bool first = true;
125 
126  for (int i = 0; i < n; ++i) {
127  TVector3 pos = mctrk[i].Position().Vect();
128 
129  // Make fiducial cuts. Require the particle to be within the physical volume of
130  // the tpc, and also require the apparent x position to be within the expanded
131  // readout frame.
132 
133  if (pos.X() >= xmin && pos.X() <= xmax && pos.Y() >= ymin && pos.Y() <= ymax &&
134  pos.Z() >= zmin && pos.Z() <= zmax) {
135  pos[0] += dx;
136  double ticks = detProp.ConvertXToTicks(pos[0], 0, 0, 0);
137  if (ticks >= 0. && ticks < detProp.ReadOutWindowSize()) {
138  if (first) {
139  start = pos;
140  startmom = 0.001 * mctrk[i].Momentum().Vect();
141  }
142  else {
143  disp -= pos;
144  result += disp.Mag();
145  }
146  first = false;
147  disp = pos;
148  end = pos;
149  endmom = 0.001 * mctrk[i].Momentum().Vect();
150  }
151  }
152  }
153 
154  return result;
155  }
156 
157  // Fill efficiency histogram assuming binomial errors.
158 
159  void effcalc(const TH1* hnum, const TH1* hden, TH1* heff)
160  {
161  int nbins = hnum->GetNbinsX();
162  if (nbins != hden->GetNbinsX())
163  throw cet::exception("SeedAna") << "effcalc[" __FILE__ "]: incompatible histograms (I)\n";
164  if (nbins != heff->GetNbinsX())
165  throw cet::exception("SeedAna") << "effcalc[" __FILE__ "]: incompatible histograms (II)\n";
166 
167  // Loop over bins, including underflow and overflow.
168 
169  for (int ibin = 0; ibin <= nbins + 1; ++ibin) {
170  double num = hnum->GetBinContent(ibin);
171  double den = hden->GetBinContent(ibin);
172  if (den == 0.) {
173  heff->SetBinContent(ibin, 0.);
174  heff->SetBinError(ibin, 0.);
175  }
176  else {
177  double eff = num / den;
178  if (eff < 0.) eff = 0.;
179  if (eff > 1.) eff = 1.;
180  double err = std::sqrt(eff * (1. - eff) / den);
181  heff->SetBinContent(ibin, eff);
182  heff->SetBinError(ibin, err);
183  }
184  }
185 
186  heff->SetMinimum(0.);
187  heff->SetMaximum(1.05);
188  heff->SetMarkerStyle(20);
189  }
190 
191  // Fill multiplicity histogram.
192 
193  void mulcalc(const TH1* hnum, const TH1* hden, TH1* hmul)
194  {
195  int nbins = hnum->GetNbinsX();
196  if (nbins != hden->GetNbinsX())
197  throw cet::exception("SeedAna") << "mulcalc[" __FILE__ "]: incompatible histograms (I)\n";
198  if (nbins != hmul->GetNbinsX())
199  throw cet::exception("SeedAna") << "mulcalc[" __FILE__ "]: incompatible histograms (II)\n";
200 
201  // Loop over bins, including underflow and overflow.
202 
203  for (int ibin = 0; ibin <= nbins + 1; ++ibin) {
204  double num = hnum->GetBinContent(ibin);
205  double den = hden->GetBinContent(ibin);
206  if (den == 0.) {
207  hmul->SetBinContent(ibin, 0.);
208  hmul->SetBinError(ibin, 0.);
209  }
210  else {
211  double mul = num / den;
212  if (mul < 0.) mul = 0.;
213  double err = std::sqrt((1. + mul) * mul / den);
214  hmul->SetBinContent(ibin, mul);
215  hmul->SetBinError(ibin, err);
216  }
217  }
218 
219  hmul->SetMinimum(0.);
220  hmul->SetMarkerStyle(20);
221  }
222 }
223 
224 namespace trkf {
225 
226  class SeedAna : public art::EDAnalyzer {
227  public:
228  // Embedded structs.
229 
230  // Struct for histograms that depend on seeds only.
231 
232  struct RecoHists {
233  RecoHists(const std::string& subdir);
234 
235  // Pure reco seed histograms.
236 
237  TH1F* fHx{nullptr}; // Seed x position.
238  TH1F* fHy{nullptr}; // Seed y position.
239  TH1F* fHz{nullptr}; // Seed z position.
240  TH1F* fHdist{nullptr}; // Seed distance to boundary.
241  TH1F* fHtheta{nullptr}; // Theta.
242  TH1F* fHphi{nullptr}; // Phi.
243  TH1F* fHtheta_xz{nullptr}; // Theta_xz.
244  TH1F* fHtheta_yz{nullptr}; // Theta_yz.
245  };
246 
247  // Struct for mc particles and mc-matched tracks.
248 
249  struct MCHists {
250  MCHists(const std::string& subdir);
251 
252  // Reco-MC matching.
253 
254  TH2F* fHduvcosth{nullptr}; // 2D mc vs. data matching, duv vs. cos(theta).
255  TH1F* fHcosth{nullptr}; // 1D direction matching, cos(theta).
256  TH1F* fHmcu{nullptr}; // 1D endpoint truth u.
257  TH1F* fHmcv{nullptr}; // 1D endpoint truth v.
258  TH1F* fHmcw{nullptr}; // 1D endpoint truth w.
259  TH1F* fHmcdudw{nullptr}; // Truth du/dw.
260  TH1F* fHmcdvdw{nullptr}; // Truth dv/dw.
261 
262  // Pure MC particle histograms (efficiency denominator).
263 
264  TH1F* fHmcstartx{nullptr}; // Starting x position.
265  TH1F* fHmcstarty{nullptr}; // Starting y position.
266  TH1F* fHmcstartz{nullptr}; // Starting z position.
267  TH1F* fHmcendx{nullptr}; // Ending x position.
268  TH1F* fHmcendy{nullptr}; // Ending y position.
269  TH1F* fHmcendz{nullptr}; // Ending z position.
270  TH1F* fHmctheta{nullptr}; // Theta.
271  TH1F* fHmcphi{nullptr}; // Phi.
272  TH1F* fHmctheta_xz{nullptr}; // Theta_xz.
273  TH1F* fHmctheta_yz{nullptr}; // Theta_yz.
274  TH1F* fHmcmom{nullptr}; // Momentum.
275  TH1F* fHmclen{nullptr}; // Length.
276 
277  // Matched seed histograms (multiplicity numerator).
278 
279  TH1F* fHmstartx{nullptr}; // Starting x position.
280  TH1F* fHmstarty{nullptr}; // Starting y position.
281  TH1F* fHmstartz{nullptr}; // Starting z position.
282  TH1F* fHmendx{nullptr}; // Ending x position.
283  TH1F* fHmendy{nullptr}; // Ending y position.
284  TH1F* fHmendz{nullptr}; // Ending z position.
285  TH1F* fHmtheta{nullptr}; // Theta.
286  TH1F* fHmphi{nullptr}; // Phi.
287  TH1F* fHmtheta_xz{nullptr}; // Theta_xz.
288  TH1F* fHmtheta_yz{nullptr}; // Theta_yz.
289  TH1F* fHmmom{nullptr}; // Momentum.
290  TH1F* fHmlen{nullptr}; // Length.
291 
292  // Matched seed histograms (efficiency numerator).
293 
294  TH1F* fHgstartx{nullptr}; // Starting x position.
295  TH1F* fHgstarty{nullptr}; // Starting y position.
296  TH1F* fHgstartz{nullptr}; // Starting z position.
297  TH1F* fHgendx{nullptr}; // Ending x position.
298  TH1F* fHgendy{nullptr}; // Ending y position.
299  TH1F* fHgendz{nullptr}; // Ending z position.
300  TH1F* fHgtheta{nullptr}; // Theta.
301  TH1F* fHgphi{nullptr}; // Phi.
302  TH1F* fHgtheta_xz{nullptr}; // Theta_xz.
303  TH1F* fHgtheta_yz{nullptr}; // Theta_yz.
304  TH1F* fHgmom{nullptr}; // Momentum.
305  TH1F* fHglen{nullptr}; // Length.
306 
307  // Multiplicity histograms.
308 
309  TH1F* fHmulstartx{nullptr}; // Starting x position.
310  TH1F* fHmulstarty{nullptr}; // Starting y position.
311  TH1F* fHmulstartz{nullptr}; // Starting z position.
312  TH1F* fHmulendx{nullptr}; // Ending x position.
313  TH1F* fHmulendy{nullptr}; // Ending y position.
314  TH1F* fHmulendz{nullptr}; // Ending z position.
315  TH1F* fHmultheta{nullptr}; // Theta.
316  TH1F* fHmulphi{nullptr}; // Phi.
317  TH1F* fHmultheta_xz{nullptr}; // Theta_xz.
318  TH1F* fHmultheta_yz{nullptr}; // Theta_yz.
319  TH1F* fHmulmom{nullptr}; // Momentum.
320  TH1F* fHmullen{nullptr}; // Length.
321 
322  // Efficiency histograms.
323 
324  TH1F* fHestartx{nullptr}; // Starting x position.
325  TH1F* fHestarty{nullptr}; // Starting y position.
326  TH1F* fHestartz{nullptr}; // Starting z position.
327  TH1F* fHeendx{nullptr}; // Ending x position.
328  TH1F* fHeendy{nullptr}; // Ending y position.
329  TH1F* fHeendz{nullptr}; // Ending z position.
330  TH1F* fHetheta{nullptr}; // Theta.
331  TH1F* fHephi{nullptr}; // Phi.
332  TH1F* fHetheta_xz{nullptr}; // Theta_xz.
333  TH1F* fHetheta_yz{nullptr}; // Theta_yz.
334  TH1F* fHemom{nullptr}; // Momentum.
335  TH1F* fHelen{nullptr}; // Length.
336  };
337 
338  explicit SeedAna(fhicl::ParameterSet const& pset);
339 
340  private:
341  void analyze(const art::Event& evt) override;
342  void endJob() override;
343 
344  // Fcl Attributes.
345 
346  std::string fSeedModuleLabel;
347  std::string fMCTrackModuleLabel;
348  int fDump; // Number of events to dump to debug message facility.
349  double fMinMCKE; // Minimum MC particle kinetic energy (GeV).
350  double fMinMCLen; // Minimum MC particle length in tpc (cm).
351  double fMatchColinearity; // Minimum matching colinearity.
352  double fMatchDisp; // Maximum matching displacement.
353  bool fIgnoreSign; // Ignore sign of mc particle if true.
354 
355  // Histograms.
356 
357  std::map<int, MCHists> fMCHistMap; // Indexed by pdg id.
358  std::map<int, RecoHists> fRecoHistMap; // Indexed by pdg id.
359 
360  // Statistics.
361 
363  };
364 
366 
367  // RecoHists methods.
368 
369  SeedAna::RecoHists::RecoHists(const std::string& subdir)
370  //
371  // Purpose: Initializing constructor.
372  //
373  {
374  // Get services.
375 
378 
379  // Make histogram directory.
380 
381  art::TFileDirectory topdir = tfs->mkdir("seedana", "SeedAna histograms");
382  art::TFileDirectory dir = topdir.mkdir(subdir);
383 
384  // Book histograms.
385 
386  fHx =
387  dir.make<TH1F>("x", "X Position", 100, -2. * geom->DetHalfWidth(), 4. * geom->DetHalfWidth());
388  fHy = dir.make<TH1F>("y", "Y Position", 100, -geom->DetHalfHeight(), geom->DetHalfHeight());
389  fHz = dir.make<TH1F>("z", "Z Position", 100, 0., geom->DetLength());
390  fHdist =
391  dir.make<TH1F>("dist", "Position Distance to Boundary", 100, -10., geom->DetHalfWidth());
392  fHtheta = dir.make<TH1F>("theta", "Theta", 100, 0., 3.142);
393  fHphi = dir.make<TH1F>("phi", "Phi", 100, -3.142, 3.142);
394  fHtheta_xz = dir.make<TH1F>("theta_xz", "Theta_xz", 100, -3.142, 3.142);
395  fHtheta_yz = dir.make<TH1F>("theta_yz", "Theta_yz", 100, -3.142, 3.142);
396  }
397 
398  // MCHists methods.
399 
400  SeedAna::MCHists::MCHists(const std::string& subdir)
401  //
402  // Purpose: Initializing constructor.
403  //
404  {
405  // Get services.
406 
409 
410  // Make histogram directory.
411 
412  art::TFileDirectory topdir = tfs->mkdir("seedana", "SeedAna histograms");
413  art::TFileDirectory dir = topdir.mkdir(subdir);
414 
415  // Book histograms.
416 
417  fHduvcosth =
418  dir.make<TH2F>("duvcosth", "Delta(uv) vs. Colinearity", 100, 0.95, 1., 100, 0., 1.);
419  fHcosth = dir.make<TH1F>("colin", "Colinearity", 100, 0.95, 1.);
420  fHmcu = dir.make<TH1F>("mcu", "MC Truth U", 100, -5., 5.);
421  fHmcv = dir.make<TH1F>("mcv", "MC Truth V", 100, -5., 5.);
422  fHmcw = dir.make<TH1F>("mcw", "MC Truth W", 100, -20., 20.);
423  fHmcdudw = dir.make<TH1F>("mcdudw", "MC Truth U Slope", 100, -0.2, 0.2);
424  fHmcdvdw = dir.make<TH1F>("mcdvdw", "MV Truth V Slope", 100, -0.2, 0.2);
425 
426  fHmcstartx = dir.make<TH1F>(
427  "mcxstart", "MC X Start Position", 10, -2. * geom->DetHalfWidth(), 4. * geom->DetHalfWidth());
428  fHmcstarty = dir.make<TH1F>(
429  "mcystart", "MC Y Start Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
430  fHmcstartz = dir.make<TH1F>("mczstart", "MC Z Start Position", 10, 0., geom->DetLength());
431  fHmcendx = dir.make<TH1F>(
432  "mcxend", "MC X End Position", 10, -2. * geom->DetHalfWidth(), 4. * geom->DetHalfWidth());
433  fHmcendy = dir.make<TH1F>(
434  "mcyend", "MC Y End Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
435  fHmcendz = dir.make<TH1F>("mczend", "MC Z End Position", 10, 0., geom->DetLength());
436  fHmctheta = dir.make<TH1F>("mctheta", "MC Theta", 20, 0., 3.142);
437  fHmcphi = dir.make<TH1F>("mcphi", "MC Phi", 10, -3.142, 3.142);
438  fHmctheta_xz = dir.make<TH1F>("mctheta_xz", "MC Theta_xz", 40, -3.142, 3.142);
439  fHmctheta_yz = dir.make<TH1F>("mctheta_yz", "MC Theta_yz", 40, -3.142, 3.142);
440  fHmcmom = dir.make<TH1F>("mcmom", "MC Momentum", 10, 0., 10.);
441  fHmclen = dir.make<TH1F>("mclen", "MC Particle Length", 10, 0., 1.1 * geom->DetLength());
442 
443  fHmstartx = dir.make<TH1F>("mxstart",
444  "Matched X Start Position",
445  10,
446  -2. * geom->DetHalfWidth(),
447  4. * geom->DetHalfWidth());
448  fHmstarty = dir.make<TH1F>(
449  "mystart", "Matched Y Start Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
450  fHmstartz = dir.make<TH1F>("mzstart", "Matched Z Start Position", 10, 0., geom->DetLength());
451  fHmendx = dir.make<TH1F>(
452  "mxend", "Matched X End Position", 10, -2. * geom->DetHalfWidth(), 4. * geom->DetHalfWidth());
453  fHmendy = dir.make<TH1F>(
454  "myend", "Matched Y End Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
455  fHmendz = dir.make<TH1F>("mzend", "Matched Z End Position", 10, 0., geom->DetLength());
456  fHmtheta = dir.make<TH1F>("mtheta", "Matched Theta", 20, 0., 3.142);
457  fHmphi = dir.make<TH1F>("mphi", "Matched Phi", 10, -3.142, 3.142);
458  fHmtheta_xz = dir.make<TH1F>("mtheta_xz", "Matched Theta_xz", 40, -3.142, 3.142);
459  fHmtheta_yz = dir.make<TH1F>("mtheta_yz", "Matched Theta_yz", 40, -3.142, 3.142);
460  fHmmom = dir.make<TH1F>("mmom", "Matched Momentum", 10, 0., 10.);
461  fHmlen = dir.make<TH1F>("mlen", "Matched Particle Length", 10, 0., 1.1 * geom->DetLength());
462 
463  fHgstartx = dir.make<TH1F>("gxstart",
464  "Good X Start Position",
465  10,
466  -2. * geom->DetHalfWidth(),
467  4. * geom->DetHalfWidth());
468  fHgstarty = dir.make<TH1F>(
469  "gystart", "Good Y Start Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
470  fHgstartz = dir.make<TH1F>("gzstart", "Good Z Start Position", 10, 0., geom->DetLength());
471  fHgendx = dir.make<TH1F>(
472  "gxend", "Good X End Position", 10, -2. * geom->DetHalfWidth(), 4. * geom->DetHalfWidth());
473  fHgendy = dir.make<TH1F>(
474  "gyend", "Good Y End Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
475  fHgendz = dir.make<TH1F>("gzend", "Good Z End Position", 10, 0., geom->DetLength());
476  fHgtheta = dir.make<TH1F>("gtheta", "Good Theta", 20, 0., 3.142);
477  fHgphi = dir.make<TH1F>("gphi", "Good Phi", 10, -3.142, 3.142);
478  fHgtheta_xz = dir.make<TH1F>("gtheta_xz", "Good Theta_xz", 40, -3.142, 3.142);
479  fHgtheta_yz = dir.make<TH1F>("gtheta_yz", "Good Theta_yz", 40, -3.142, 3.142);
480  fHgmom = dir.make<TH1F>("gmom", "Good Momentum", 10, 0., 10.);
481  fHglen = dir.make<TH1F>("glen", "Good Particle Length", 10, 0., 1.1 * geom->DetLength());
482 
483  fHmulstartx = dir.make<TH1F>("mulxstart",
484  "Multiplicity vs. X Start Position",
485  10,
486  -2. * geom->DetHalfWidth(),
487  4. * geom->DetHalfWidth());
488  fHmulstarty = dir.make<TH1F>("mulystart",
489  "Multiplicity vs. Y Start Position",
490  10,
491  -geom->DetHalfHeight(),
492  geom->DetHalfHeight());
493  fHmulstartz =
494  dir.make<TH1F>("mulzstart", "Multiplicity vs. Z Start Position", 10, 0., geom->DetLength());
495  fHmulendx = dir.make<TH1F>("mulxend",
496  "Multiplicity vs. X End Position",
497  10,
498  -2. * geom->DetHalfWidth(),
499  4. * geom->DetHalfWidth());
500  fHmulendy = dir.make<TH1F>("mulyend",
501  "Multiplicity vs. Y End Position",
502  10,
503  -geom->DetHalfHeight(),
504  geom->DetHalfHeight());
505  fHmulendz =
506  dir.make<TH1F>("mulzend", "Multiplicity vs. Z End Position", 10, 0., geom->DetLength());
507  fHmultheta = dir.make<TH1F>("multheta", "Multiplicity vs. Theta", 20, 0., 3.142);
508  fHmulphi = dir.make<TH1F>("mulphi", "Multiplicity vs. Phi", 10, -3.142, 3.142);
509  fHmultheta_xz = dir.make<TH1F>("multheta_xz", "Multiplicity vs. Theta_xz", 40, -3.142, 3.142);
510  fHmultheta_yz = dir.make<TH1F>("multheta_yz", "Multiplicity vs. Theta_yz", 40, -3.142, 3.142);
511  fHmulmom = dir.make<TH1F>("mulmom", "Multiplicity vs. Momentum", 10, 0., 10.);
512  fHmullen =
513  dir.make<TH1F>("mullen", "Multiplicity vs. Particle Length", 10, 0., 1.1 * geom->DetLength());
514 
515  fHestartx = dir.make<TH1F>("exstart",
516  "Efficiency vs. X Start Position",
517  10,
518  -2. * geom->DetHalfWidth(),
519  4. * geom->DetHalfWidth());
520  fHestarty = dir.make<TH1F>("eystart",
521  "Efficiency vs. Y Start Position",
522  10,
523  -geom->DetHalfHeight(),
524  geom->DetHalfHeight());
525  fHestartz =
526  dir.make<TH1F>("ezstart", "Efficiency vs. Z Start Position", 10, 0., geom->DetLength());
527  fHeendx = dir.make<TH1F>("exend",
528  "Efficiency vs. X End Position",
529  10,
530  -2. * geom->DetHalfWidth(),
531  4. * geom->DetHalfWidth());
532  fHeendy = dir.make<TH1F>(
533  "eyend", "Efficiency vs. Y End Position", 10, -geom->DetHalfHeight(), geom->DetHalfHeight());
534  fHeendz = dir.make<TH1F>("ezend", "Efficiency vs. Z End Position", 10, 0., geom->DetLength());
535  fHetheta = dir.make<TH1F>("etheta", "Efficiency vs. Theta", 20, 0., 3.142);
536  fHephi = dir.make<TH1F>("ephi", "Efficiency vs. Phi", 10, -3.142, 3.142);
537  fHetheta_xz = dir.make<TH1F>("etheta_xz", "Efficiency vs. Theta_xz", 40, -3.142, 3.142);
538  fHetheta_yz = dir.make<TH1F>("etheta_yz", "Efficiency vs. Theta_yz", 40, -3.142, 3.142);
539  fHemom = dir.make<TH1F>("emom", "Efficiency vs. Momentum", 10, 0., 10.);
540  fHelen =
541  dir.make<TH1F>("elen", "Efficiency vs. Particle Length", 10, 0., 1.1 * geom->DetLength());
542  }
543 
545  //
546  // Purpose: Constructor.
547  //
548  // Arguments: pset - Module parameters.
549  //
550  : EDAnalyzer(pset)
551  , fSeedModuleLabel(pset.get<std::string>("SeedModuleLabel"))
552  , fMCTrackModuleLabel(pset.get<std::string>("MCTrackModuleLabel"))
553  , fDump(pset.get<int>("Dump"))
554  , fMinMCKE(pset.get<double>("MinMCKE"))
555  , fMinMCLen(pset.get<double>("MinMCLen"))
556  , fMatchColinearity(pset.get<double>("MatchColinearity"))
557  , fMatchDisp(pset.get<double>("MatchDisp"))
558  , fIgnoreSign(pset.get<bool>("IgnoreSign"))
559  , fNumEvent(0)
560  {
561  mf::LogInfo("SeedAna") << "SeedAna configured with the following parameters:\n"
562  << " SeedModuleLabel = " << fSeedModuleLabel << "\n"
563  << " MCTrackModuleLabel = " << fMCTrackModuleLabel << "\n"
564  << " Dump = " << fDump << "\n"
565  << " MinMCKE = " << fMinMCKE << "\n"
566  << " MinMCLen = " << fMinMCLen;
567  }
568 
570  //
571  // Purpose: Analyze method.
572  //
573  // Arguments: event - Art event.
574  //
575  {
576  ++fNumEvent;
577 
578  // Optional dump stream.
579 
580  std::unique_ptr<mf::LogInfo> pdump;
581  if (fDump > 0) {
582  --fDump;
583  pdump = std::unique_ptr<mf::LogInfo>(new mf::LogInfo("TrackAna"));
584  }
585 
586  // Make sure histograms are booked.
587 
588  bool mc = !evt.isRealData();
589 
590  // Get seed handle.
591 
593  evt.getByLabel(fSeedModuleLabel, seedh);
594 
595  // Seed->mc track matching map.
596 
597  std::map<const recob::Seed*, int> seedmap;
598 
599  if (mc) {
600 
601  // Get MCTracks.
602 
604  evt.getByLabel(fMCTrackModuleLabel, mctrackh);
605 
606  // Dump MCTracks.
607 
608  if (pdump) {
609  *pdump << "MC Tracks\n"
610  << " Id pdg x y z dx dy dz "
611  " p\n"
612  << "--------------------------------------------------------------------------------"
613  "-----------\n";
614  }
615 
616  // Loop over mc tracks, and fill histograms that depend only
617  // on mc particles.
618  auto const detProp =
620 
621  for (std::vector<sim::MCTrack>::const_iterator imctrk = mctrackh->begin();
622  imctrk != mctrackh->end();
623  ++imctrk) {
624  const sim::MCTrack& mctrk = *imctrk;
625  int pdg = mctrk.PdgCode();
626  if (fIgnoreSign) pdg = std::abs(pdg);
627 
628  // Ignore everything except stable charged nonshowering particles.
629 
630  int apdg = std::abs(pdg);
631  if (apdg == 13 || // Muon
632  apdg == 211 || // Charged pion
633  apdg == 321 || // Charged kaon
634  apdg == 2212) { // (Anti)proton
635 
636  // Apply minimum energy cut.
637 
638  if (mctrk.Start().E() >= mctrk.Start().Momentum().Mag() + 1000. * fMinMCKE) {
639 
640  // Calculate the x offset due to nonzero mc particle time.
641 
642  double mctime = mctrk.Start().T(); // nsec
643  double mcdx = mctime * 1.e-3 * detProp.DriftVelocity(); // cm
644 
645  // Calculate the length of this mc particle inside the fiducial volume.
646 
647  TVector3 mcstart;
648  TVector3 mcend;
649  TVector3 mcstartmom;
650  TVector3 mcendmom;
651  double plen = length(detProp, mctrk, mcdx, mcstart, mcend, mcstartmom, mcendmom);
652 
653  // Apply minimum fiducial length cut. Always reject particles that have
654  // zero length in the tpc regardless of the configured cut.
655 
656  if (plen > 0. && plen > fMinMCLen) {
657 
658  // Dump MC particle information here.
659 
660  if (pdump) {
661  double pstart = mcstartmom.Mag();
662  double pend = mcendmom.Mag();
663  *pdump << "\nOffset" << std::setw(3) << mctrk.TrackID() << std::setw(6)
664  << mctrk.PdgCode() << " " << std::fixed << std::setprecision(2)
665  << std::setw(10) << mcdx << "\nStart " << std::setw(3) << mctrk.TrackID()
666  << std::setw(6) << mctrk.PdgCode() << " " << std::fixed
667  << std::setprecision(2) << std::setw(10) << mcstart[0] << std::setw(10)
668  << mcstart[1] << std::setw(10) << mcstart[2];
669  if (pstart > 0.) {
670  *pdump << " " << std::fixed << std::setprecision(3) << std::setw(10)
671  << mcstartmom[0] / pstart << std::setw(10) << mcstartmom[1] / pstart
672  << std::setw(10) << mcstartmom[2] / pstart;
673  }
674  else
675  *pdump << std::setw(32) << " ";
676  *pdump << std::setw(12) << std::fixed << std::setprecision(3) << pstart;
677  *pdump << "\nEnd " << std::setw(3) << mctrk.TrackID() << std::setw(6)
678  << mctrk.PdgCode() << " " << std::fixed << std::setprecision(2)
679  << std::setw(10) << mcend[0] << std::setw(10) << mcend[1] << std::setw(10)
680  << mcend[2];
681  if (pend > 0.01) {
682  *pdump << " " << std::fixed << std::setprecision(3) << std::setw(10)
683  << mcendmom[0] / pend << std::setw(10) << mcendmom[1] / pend
684  << std::setw(10) << mcendmom[2] / pend;
685  }
686  else
687  *pdump << std::setw(32) << " ";
688  *pdump << std::setw(12) << std::fixed << std::setprecision(3) << pend << "\n";
689  }
690 
691  // Fill histograms.
692 
693  if (fMCHistMap.count(pdg) == 0) {
694  std::ostringstream ostr;
695  ostr << "MC" << (fIgnoreSign ? "All" : (pdg > 0 ? "Pos" : "Neg")) << std::abs(pdg);
696  fMCHistMap.emplace(pdg, MCHists{ostr.str()});
697  }
698  const MCHists& mchists = fMCHistMap.at(pdg);
699 
700  double mctheta_xz = std::atan2(mcstartmom.X(), mcstartmom.Z());
701  double mctheta_yz = std::atan2(mcstartmom.Y(), mcstartmom.Z());
702 
703  mchists.fHmcstartx->Fill(mcstart.X());
704  mchists.fHmcstarty->Fill(mcstart.Y());
705  mchists.fHmcstartz->Fill(mcstart.Z());
706  mchists.fHmcendx->Fill(mcend.X());
707  mchists.fHmcendy->Fill(mcend.Y());
708  mchists.fHmcendz->Fill(mcend.Z());
709  mchists.fHmctheta->Fill(mcstartmom.Theta());
710  mchists.fHmcphi->Fill(mcstartmom.Phi());
711  mchists.fHmctheta_xz->Fill(mctheta_xz);
712  mchists.fHmctheta_yz->Fill(mctheta_yz);
713  mchists.fHmcmom->Fill(mcstartmom.Mag());
714  mchists.fHmclen->Fill(plen);
715 
716  // Loop over seeds and do matching.
717 
718  int nmatch = 0;
719  if (seedh.isValid()) {
720 
721  // Loop over seeds.
722 
723  int nseed = seedh->size();
724  for (int i = 0; i < nseed; ++i) {
725  art::Ptr<recob::Seed> pseed(seedh, i);
726  const recob::Seed& seed = *pseed;
727  if (seedmap.count(&seed) == 0) seedmap[&seed] = -1;
728 
729  // Get parameters of this seed.
730 
731  TVector3 pos;
732  TVector3 dir;
733  double err[3];
734  seed.GetPoint(&pos[0], err);
735  seed.GetDirection(&dir[0], err);
736 
737  // Calculate the global-to-local rotation matrix.
738  // Copied from Track.cxx.
739 
740  TMatrixD rot(3, 3);
741  double dirmag = dir.Mag();
742  double diryz = std::sqrt(dir.Y() * dir.Y() + dir.Z() * dir.Z());
743 
744  double sinth = dir.X() / dirmag;
745  double costh = diryz / dirmag;
746  double sinphi = 0.;
747  double cosphi = 1.;
748  if (diryz != 0) {
749  sinphi = -dir.Y() / diryz;
750  cosphi = dir.Z() / diryz;
751  }
752  rot(0, 0) = costh;
753  rot(1, 0) = 0.;
754  rot(2, 0) = sinth;
755  rot(0, 1) = sinth * sinphi;
756  rot(1, 1) = cosphi;
757  rot(2, 1) = -costh * sinphi;
758  rot(0, 2) = -sinth * cosphi;
759  rot(1, 2) = sinphi;
760  rot(2, 2) = costh * cosphi;
761 
762  // Get best matching mc trajectory point.
763 
764  int itraj = mcmatch(detProp, mctrk, seed);
765  if (itraj >= 0) {
766 
767  // Get mc relative position and direction at matching trajectory point.
768 
769  TVector3 mcpos = mctrk[itraj].Position().Vect() - pos;
770  TVector3 mcmom = mctrk[itraj].Momentum().Vect();
771  mcpos[0] += mcdx;
772 
773  // Rotate the momentum and position to the
774  // seed-local coordinate system.
775 
776  TVector3 mcmoml = rot * mcmom;
777  TVector3 mcposl = rot * mcpos;
778 
779  if (mcmoml.Z() < 0.) mcmoml = -mcmoml;
780  double costh = mcmoml.Z() / mcmoml.Mag();
781 
782  double u = mcposl.X();
783  double v = mcposl.Y();
784  double w = mcposl.Z();
785 
786  double pu = mcmoml.X();
787  double pv = mcmoml.Y();
788  double pw = mcmoml.Z();
789 
790  double dudw = pu / pw;
791  double dvdw = pv / pw;
792 
793  double u0 = u - w * dudw;
794  double v0 = v - w * dvdw;
795  double uv0 = std::sqrt(u0 * u0 + v0 * v0);
796 
797  // Fill matching histograms.
798 
799  mchists.fHduvcosth->Fill(costh, uv0);
800  if (std::abs(uv0) < fMatchDisp) {
801 
802  // Fill slope matching histograms.
803 
804  mchists.fHmcdudw->Fill(dudw);
805  mchists.fHmcdvdw->Fill(dvdw);
806  }
807  mchists.fHcosth->Fill(costh);
808  if (costh > fMatchColinearity) {
809 
810  // Fill displacement matching histograms.
811 
812  mchists.fHmcu->Fill(u0);
813  mchists.fHmcv->Fill(v0);
814  mchists.fHmcw->Fill(w);
815 
816  if (std::abs(uv0) < fMatchDisp) {
817 
818  // Now we have passed all matching cuts and we have a matching
819  // mc particle + seed pair.
820 
821  ++nmatch;
822  seedmap[&seed] = mctrk.TrackID();
823 
824  // Fill matched seed histograms (seed multiplicity).
825 
826  mchists.fHmstartx->Fill(mcstart.X());
827  mchists.fHmstarty->Fill(mcstart.Y());
828  mchists.fHmstartz->Fill(mcstart.Z());
829  mchists.fHmendx->Fill(mcend.X());
830  mchists.fHmendy->Fill(mcend.Y());
831  mchists.fHmendz->Fill(mcend.Z());
832  mchists.fHmtheta->Fill(mcstartmom.Theta());
833  mchists.fHmphi->Fill(mcstartmom.Phi());
834  mchists.fHmtheta_xz->Fill(mctheta_xz);
835  mchists.fHmtheta_yz->Fill(mctheta_yz);
836  mchists.fHmmom->Fill(mcstartmom.Mag());
837  mchists.fHmlen->Fill(plen);
838  }
839  }
840  }
841  }
842 
843  // If we found at least one matched seed, fill good
844  // particle histograms.
845 
846  if (nmatch > 0) {
847  mchists.fHgstartx->Fill(mcstart.X());
848  mchists.fHgstarty->Fill(mcstart.Y());
849  mchists.fHgstartz->Fill(mcstart.Z());
850  mchists.fHgendx->Fill(mcend.X());
851  mchists.fHgendy->Fill(mcend.Y());
852  mchists.fHgendz->Fill(mcend.Z());
853  mchists.fHgtheta->Fill(mcstartmom.Theta());
854  mchists.fHgphi->Fill(mcstartmom.Phi());
855  mchists.fHgtheta_xz->Fill(mctheta_xz);
856  mchists.fHgtheta_yz->Fill(mctheta_yz);
857  mchists.fHgmom->Fill(mcstartmom.Mag());
858  mchists.fHglen->Fill(plen);
859  }
860  }
861  }
862  }
863  }
864  }
865  }
866 
867  // Loop over seeds and fill reco-only seed histograms.
868 
869  if (seedh.isValid()) {
870 
871  // Loop over seeds.
872 
873  int nseed = seedh->size();
874 
875  if (nseed > 0 && pdump != 0) {
876  *pdump << "\nReconstructed Seeds\n"
877  << " MCid x y z dx dy dz "
878  " p\n"
879  << "--------------------------------------------------------------------------------"
880  "-----------\n";
881  }
882 
883  for (int i = 0; i < nseed; ++i) {
884  art::Ptr<recob::Seed> pseed(seedh, i);
885  const recob::Seed& seed = *pseed;
886 
887  // Fill histograms involving reco seeds only.
888 
889  TVector3 pos;
890  TVector3 dir;
891  double err[3];
892  seed.GetPoint(&pos[0], err);
893  seed.GetDirection(&dir[0], err);
894  double mdir = dir.Mag();
895  if (mdir != 0.) { dir *= (1. / mdir); }
896 
897  double dpos = bdist(pos);
898  double theta_xz = std::atan2(dir.X(), dir.Z());
899  double theta_yz = std::atan2(dir.Y(), dir.Z());
900 
901  // Dump seed information here.
902 
903  if (pdump) {
904  int mcid = seedmap[&seed];
905  *pdump << std::setw(15) << mcid << " " << std::fixed << std::setprecision(2)
906  << std::setw(10) << pos[0] << std::setw(10) << pos[1] << std::setw(10) << pos[2]
907  << " " << std::fixed << std::setprecision(3) << std::setw(10) << dir[0]
908  << std::setw(10) << dir[1] << std::setw(10) << dir[2] << "\n";
909  }
910 
911  // Fill histograms.
912 
913  if (fRecoHistMap.count(0) == 0) fRecoHistMap.emplace(0, RecoHists{"Reco"});
914  const RecoHists& rhists = fRecoHistMap.at(0);
915 
916  rhists.fHx->Fill(pos.X());
917  rhists.fHy->Fill(pos.Y());
918  rhists.fHz->Fill(pos.Z());
919  rhists.fHdist->Fill(dpos);
920  rhists.fHtheta->Fill(dir.Theta());
921  rhists.fHphi->Fill(dir.Phi());
922  rhists.fHtheta_xz->Fill(theta_xz);
923  rhists.fHtheta_yz->Fill(theta_yz);
924  }
925  }
926  }
927 
929  //
930  // Purpose: End of job.
931  //
932  {
933  // Print summary.
934 
935  mf::LogInfo("SeedAna") << "SeedAna statistics:\n"
936  << " Number of events = " << fNumEvent;
937 
938  // Fill multiplicity histograms.
939 
940  for (std::map<int, MCHists>::const_iterator i = fMCHistMap.begin(); i != fMCHistMap.end();
941  ++i) {
942  const MCHists& mchists = i->second;
943  mulcalc(mchists.fHmstartx, mchists.fHmcstartx, mchists.fHmulstartx);
944  mulcalc(mchists.fHmstarty, mchists.fHmcstarty, mchists.fHmulstarty);
945  mulcalc(mchists.fHmstartz, mchists.fHmcstartz, mchists.fHmulstartz);
946  mulcalc(mchists.fHmendx, mchists.fHmcendx, mchists.fHmulendx);
947  mulcalc(mchists.fHmendy, mchists.fHmcendy, mchists.fHmulendy);
948  mulcalc(mchists.fHmendz, mchists.fHmcendz, mchists.fHmulendz);
949  mulcalc(mchists.fHmtheta, mchists.fHmctheta, mchists.fHmultheta);
950  mulcalc(mchists.fHmphi, mchists.fHmcphi, mchists.fHmulphi);
951  mulcalc(mchists.fHmtheta_xz, mchists.fHmctheta_xz, mchists.fHmultheta_xz);
952  mulcalc(mchists.fHmtheta_yz, mchists.fHmctheta_yz, mchists.fHmultheta_yz);
953  mulcalc(mchists.fHmmom, mchists.fHmcmom, mchists.fHmulmom);
954  mulcalc(mchists.fHmlen, mchists.fHmclen, mchists.fHmullen);
955  }
956  // Fill efficiency histograms.
957 
958  for (std::map<int, MCHists>::const_iterator i = fMCHistMap.begin(); i != fMCHistMap.end();
959  ++i) {
960  const MCHists& mchists = i->second;
961  effcalc(mchists.fHgstartx, mchists.fHmcstartx, mchists.fHestartx);
962  effcalc(mchists.fHgstarty, mchists.fHmcstarty, mchists.fHestarty);
963  effcalc(mchists.fHgstartz, mchists.fHmcstartz, mchists.fHestartz);
964  effcalc(mchists.fHgendx, mchists.fHmcendx, mchists.fHeendx);
965  effcalc(mchists.fHgendy, mchists.fHmcendy, mchists.fHeendy);
966  effcalc(mchists.fHgendz, mchists.fHmcendz, mchists.fHeendz);
967  effcalc(mchists.fHgtheta, mchists.fHmctheta, mchists.fHetheta);
968  effcalc(mchists.fHgphi, mchists.fHmcphi, mchists.fHephi);
969  effcalc(mchists.fHgtheta_xz, mchists.fHmctheta_xz, mchists.fHetheta_xz);
970  effcalc(mchists.fHgtheta_yz, mchists.fHmctheta_yz, mchists.fHetheta_yz);
971  effcalc(mchists.fHgmom, mchists.fHmcmom, mchists.fHemom);
972  effcalc(mchists.fHglen, mchists.fHmclen, mchists.fHelen);
973  }
974  }
975 }
std::string fMCTrackModuleLabel
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Float_t den
Definition: plot.C:35
void GetPoint(double *Pt, double *Err) const
Definition: Seed.cxx:89
Length_t DetHalfWidth(TPCID const &tpcid=tpc_zero) const
Returns the half width of the active volume of the specified TPC.
void endJob() override
constexpr auto abs(T v)
Returns the absolute value of the argument.
std::string fSeedModuleLabel
STL namespace.
intermediate_table::const_iterator const_iterator
double T() const
Definition: MCStep.h:41
bool isRealData() const
Definition: Event.cc:53
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
tick ticks
Alias for common language habits.
Definition: electronics.h:76
Length_t DetLength(TPCID const &tpcid=tpc_zero) const
Returns the length of the active volume of the specified TPC.
MCHists(const std::string &subdir)
void analyze(const art::Event &evt) override
bool isValid() const noexcept
Definition: Handle.h:203
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
double fMatchColinearity
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
long seed
Definition: chem4.cc:67
double ConvertXToTicks(double X, int p, int t, int c) const
double DriftVelocity(double efield=0., double temperature=0.) const
cm/us
Class def header for mctrack data container.
std::map< int, RecoHists > fRecoHistMap
int PdgCode() const
Definition: MCTrack.h:39
const TLorentzVector & Momentum() const
Definition: MCStep.h:34
SeedAna(fhicl::ParameterSet const &pset)
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
std::map< int, MCHists > fMCHistMap
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
TDirectory * dir
Definition: macro.C:5
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
double E() const
Definition: MCStep.h:45
const MCStep & Start() const
Definition: MCTrack.h:42
Char_t n[5]
Length_t DetHalfHeight(TPCID const &tpcid=tpc_zero) const
Returns the half height of the active volume of the specified TPC.
unsigned int TrackID() const
Definition: MCTrack.h:40
TCEvent evt
Definition: DataStructs.cxx:8
Float_t w
Definition: plot.C:20
void GetDirection(double *Dir, double *Err) const
Definition: Seed.cxx:80
Definition: fwd.h:26
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
RecoHists(const std::string &subdir)