LArSoft  v10_04_05
Liquid Argon Software toolkit - https://larsoft.org/
TCShowerTemplateMaker_module.cc
Go to the documentation of this file.
1 // -------------------------------------------------
2 // makes shower profile templates
3 //
4 // Author: Rory Fitzpatrick (roryfitz@umich.edu)
5 // Created: 7/16/18
6 // -------------------------------------------------
7 
8 // Framework includes
14 #include "art_root_io/TFileService.h"
16 #include "fhiclcpp/ParameterSet.h"
17 
30 
31 #include "TH1.h"
32 #include "TH3.h"
33 #include "TProfile.h"
34 #include "TProfile2D.h"
35 
36 namespace {
37  constexpr geo::PlaneID collectionPlaneID{0, 0, 1};
38 }
39 
40 namespace shower {
41 
43  public:
44  explicit TCShowerTemplateMaker(fhicl::ParameterSet const& pset);
45 
46  private:
47  void beginJob() override;
48  void analyze(const art::Event& evt) override;
49 
50  void showerProfile(detinfo::DetectorClocksData const& clockData,
51  detinfo::DetectorPropertiesData const& detProp,
53  TVector3 shwvtx,
54  TVector3 shwdir,
55  double elep);
56  void showerProfileTrue(detinfo::DetectorClocksData const& clockData,
57  detinfo::DetectorPropertiesData const& detProp,
59  double elep);
61  simb::MCParticle electron);
62 
66 
70 
74 
78 
82 
86 
90 
94 
98 
100  TH3F* fTransverse;
106 
107  // electrons
114  // photons
121  // other
128 
129  const int LBINS = 20;
130  const int LMIN = 0;
131  const int LMAX = 5;
132 
133  const int TBINS = 20;
134  const int TMIN = -5;
135  const int TMAX = 5;
136 
137  const int EBINS = 20;
138  const double EMIN = 0.5;
139  const double EMAX = 20.5;
140 
141  const double X0 = 14;
142 
143  std::string fHitModuleLabel;
144  std::string fShowerModuleLabel;
145  std::string fGenieGenModuleLabel;
146  std::string fDigitModuleLabel;
147 
149 
150  }; // class TCShowerAnalysis
151 
152 } // shower
153 
154 // -------------------------------------------------
155 
157  : EDAnalyzer(pset)
158  , fHitModuleLabel(pset.get<std::string>("HitModuleLabel", "trajcluster"))
159  , fShowerModuleLabel(pset.get<std::string>("ShowerModuleLabel", "tcshower"))
160  , fGenieGenModuleLabel(pset.get<std::string>("GenieGenModuleLabel", "generator"))
161  , fDigitModuleLabel(pset.get<std::string>("DigitModuleLabel", "largeant"))
162  , fCalorimetryAlg(pset.get<fhicl::ParameterSet>("CalorimetryAlg"))
163 {} // TCShowerTemplateMaker
164 
165 // -------------------------------------------------
166 
168 {
169 
172  tfs->make<TProfile>("fShowerProfileSimLong",
173  "longitudinal e- profile (true, simchannel);t;E (MeV)",
174  LBINS,
175  LMIN,
176  LMAX);
177  fShowerProfileHitLong = tfs->make<TProfile>(
178  "fShowerProfileHitLong", "longitudinal e- profile (true, hit);t;E (MeV)", LBINS, LMIN, LMAX);
179  fShowerProfileRecoLong = tfs->make<TProfile>(
180  "fShowerProfileRecoLong", "longitudinal e- profile (reco);t;Q", LBINS, LMIN, LMAX);
181 
182  fShowerProfileSimLong2D = tfs->make<TProfile2D>(
183  "fShowerProfileSimLong2D",
184  "longitudinal e- profile (true, simchannel);t;electron energy (MeV);E (MeV)",
185  LBINS,
186  LMIN,
187  LMAX,
188  EBINS,
189  EMIN,
190  EMAX);
192  tfs->make<TProfile2D>("fShowerProfileHitLong2D",
193  "longitudinal e- profile (true, hit);t;electron energy (MeV);E (MeV)",
194  LBINS,
195  LMIN,
196  LMAX,
197  EBINS,
198  EMIN,
199  EMAX);
201  tfs->make<TProfile2D>("fShowerProfileRecoLong2D",
202  "longitudinal e- profile (reco);t;electron energy (MeV);Q",
203  LBINS,
204  LMIN,
205  LMAX,
206  EBINS,
207  EMIN,
208  EMAX);
209 
211  tfs->make<TProfile>("fShowerProfileSimTrans",
212  "transverse e- profile (true, simchannel);dist (cm);E (MeV)",
213  TBINS,
214  TMIN,
215  TMAX);
217  tfs->make<TProfile>("fShowerProfileHitTrans",
218  "transverse e- profile (true, hit);dist (cm);E (MeV)",
219  TBINS,
220  TMIN,
221  TMAX);
222  fShowerProfileRecoTrans = tfs->make<TProfile>(
223  "fShowerProfileRecoTrans", "transverse e- profile (reco);dist (cm);Q", TBINS, TMIN, TMAX);
224 
225  fShowerProfileSimTrans2D = tfs->make<TProfile2D>(
226  "fShowerProfileSimTrans2D",
227  "transverse e- profile (true, simchannel);t;electron energy (MeV);E (MeV)",
228  TBINS,
229  TMIN,
230  TMAX,
231  EBINS,
232  EMIN,
233  EMAX);
235  tfs->make<TProfile2D>("fShowerProfileHitTrans2D",
236  "transverse e- profile (true, hit);t;electron energy (MeV);E (MeV)",
237  TBINS,
238  TMIN,
239  TMAX,
240  EBINS,
241  EMIN,
242  EMAX);
244  tfs->make<TProfile2D>("fShowerProfileRecoTrans2D",
245  "transverse e- profile (reco);t;electron energy (MeV);Q",
246  TBINS,
247  TMIN,
248  TMAX,
249  EBINS,
250  EMIN,
251  EMAX);
252 
253  fShowerProfileSimTrans2D_1 = tfs->make<TProfile2D>(
254  "fShowerProfileSimTrans2D_1",
255  "transverse e- profile [0 <= t < 1] (true, simchannel);t;electron energy (MeV);E (MeV)",
256  TBINS,
257  TMIN,
258  TMAX,
259  EBINS,
260  EMIN,
261  EMAX);
262  fShowerProfileHitTrans2D_1 = tfs->make<TProfile2D>(
263  "fShowerProfileHitTrans2D_1",
264  "transverse e- profile [0 <= t < 1] (true, hit);t;electron energy (MeV);E (MeV)",
265  TBINS,
266  TMIN,
267  TMAX,
268  EBINS,
269  EMIN,
270  EMAX);
272  tfs->make<TProfile2D>("fShowerProfileRecoTrans2D_1",
273  "transverse e- profile [0 <= t < 1] (reco);t;electron energy (MeV);Q",
274  TBINS,
275  TMIN,
276  TMAX,
277  EBINS,
278  EMIN,
279  EMAX);
280 
281  fShowerProfileSimTrans2D_2 = tfs->make<TProfile2D>(
282  "fShowerProfileSimTrans2D_2",
283  "transverse e- profile [1 <= t < 2] (true, simchannel);t;electron energy (MeV);E (MeV)",
284  TBINS,
285  TMIN,
286  TMAX,
287  EBINS,
288  EMIN,
289  EMAX);
290  fShowerProfileHitTrans2D_2 = tfs->make<TProfile2D>(
291  "fShowerProfileHitTrans2D_2",
292  "transverse e- profile [1 <= t < 2] (true, hit);t;electron energy (MeV);E (MeV)",
293  TBINS,
294  TMIN,
295  TMAX,
296  EBINS,
297  EMIN,
298  EMAX);
300  tfs->make<TProfile2D>("fShowerProfileRecoTrans2D_2",
301  "transverse e- profile [1 <= t < 2] (reco);t;electron energy (MeV);Q",
302  TBINS,
303  TMIN,
304  TMAX,
305  EBINS,
306  EMIN,
307  EMAX);
308 
309  fShowerProfileSimTrans2D_3 = tfs->make<TProfile2D>(
310  "fShowerProfileSimTrans2D_3",
311  "transverse e- profile [2 <= t < 3] (true, simchannel);t;electron energy (MeV);E (MeV)",
312  TBINS,
313  TMIN,
314  TMAX,
315  EBINS,
316  EMIN,
317  EMAX);
318  fShowerProfileHitTrans2D_3 = tfs->make<TProfile2D>(
319  "fShowerProfileHitTrans2D_3",
320  "transverse e- profile [2 <= t < 3] (true, hit);t;electron energy (MeV);E (MeV)",
321  TBINS,
322  TMIN,
323  TMAX,
324  EBINS,
325  EMIN,
326  EMAX);
328  tfs->make<TProfile2D>("fShowerProfileRecoTrans2D_3",
329  "transverse e- profile [2 <= t < 3] (reco);t;electron energy (MeV);Q",
330  TBINS,
331  TMIN,
332  TMAX,
333  EBINS,
334  EMIN,
335  EMAX);
336 
337  fShowerProfileSimTrans2D_4 = tfs->make<TProfile2D>(
338  "fShowerProfileSimTrans2D_4",
339  "transverse e- profile [3 <= t < 4] (true, simchannel);t;electron energy (MeV);E (MeV)",
340  TBINS,
341  TMIN,
342  TMAX,
343  EBINS,
344  EMIN,
345  EMAX);
346  fShowerProfileHitTrans2D_4 = tfs->make<TProfile2D>(
347  "fShowerProfileHitTrans2D_4",
348  "transverse e- profile [3 <= t < 4] (true, hit);t;electron energy (MeV);E (MeV)",
349  TBINS,
350  TMIN,
351  TMAX,
352  EBINS,
353  EMIN,
354  EMAX);
356  tfs->make<TProfile2D>("fShowerProfileRecoTrans2D_4",
357  "transverse e- profile [3 <= t < 4] (reco);t;electron energy (MeV);Q",
358  TBINS,
359  TMIN,
360  TMAX,
361  EBINS,
362  EMIN,
363  EMAX);
364 
365  fShowerProfileSimTrans2D_5 = tfs->make<TProfile2D>(
366  "fShowerProfileSimTrans2D_5",
367  "transverse e- profile [4 <= t < 5] (true, simchannel);t;electron energy (MeV);E (MeV)",
368  TBINS,
369  TMIN,
370  TMAX,
371  EBINS,
372  EMIN,
373  EMAX);
374  fShowerProfileHitTrans2D_5 = tfs->make<TProfile2D>(
375  "fShowerProfileHitTrans2D_5",
376  "transverse e- profile [4 <= t < 5] (true, hit);t;electron energy (MeV);E (MeV)",
377  TBINS,
378  TMIN,
379  TMAX,
380  EBINS,
381  EMIN,
382  EMAX);
384  tfs->make<TProfile2D>("fShowerProfileRecoTrans2D_5",
385  "transverse e- profile [4 <= t < 5] (reco);t;electron energy (MeV);Q",
386  TBINS,
387  TMIN,
388  TMAX,
389  EBINS,
390  EMIN,
391  EMAX);
392 
393  fLongitudinal = tfs->make<TH3F>("fLongitudinal",
394  "longitudinal e- profile;t;electron energy (MeV);Q",
395  LBINS,
396  LMIN,
397  LMAX,
398  EBINS,
399  EMIN,
400  EMAX,
401  100,
402  0,
403  150000);
404  fTransverse = tfs->make<TH3F>("fTransverse",
405  "transverse e- profile;dist (cm);electron energy (MeV);Q",
406  TBINS,
407  TMIN,
408  TMAX,
409  EBINS,
410  EMIN,
411  EMAX,
412  100,
413  0,
414  150000);
415  fTransverse_1 =
416  tfs->make<TH3F>("fTransverse_1",
417  "transverse e- profile [0 <= t < 1];dist (cm);electron energy (MeV);Q",
418  TBINS,
419  TMIN,
420  TMAX,
421  EBINS,
422  EMIN,
423  EMAX,
424  100,
425  0,
426  100000);
427  fTransverse_2 =
428  tfs->make<TH3F>("fTransverse_2",
429  "transverse e- profile [1 <= t < 2];dist (cm);electron energy (MeV);Q",
430  TBINS,
431  TMIN,
432  TMAX,
433  EBINS,
434  EMIN,
435  EMAX,
436  100,
437  0,
438  100000);
439  fTransverse_3 =
440  tfs->make<TH3F>("fTransverse_3",
441  "transverse e- profile [2 <= t < 3];dist (cm);electron energy (MeV);Q",
442  TBINS,
443  TMIN,
444  TMAX,
445  EBINS,
446  EMIN,
447  EMAX,
448  100,
449  0,
450  100000);
451  fTransverse_4 =
452  tfs->make<TH3F>("fTransverse_4",
453  "transverse e- profile [3 <= t < 4];dist (cm);electron energy (MeV);Q",
454  TBINS,
455  TMIN,
456  TMAX,
457  EBINS,
458  EMIN,
459  EMAX,
460  100,
461  0,
462  100000);
463  fTransverse_5 =
464  tfs->make<TH3F>("fTransverse_5",
465  "transverse e- profile [4 <= t < 5];dist (cm);electron energy (MeV);Q",
466  TBINS,
467  TMIN,
468  TMAX,
469  EBINS,
470  EMIN,
471  EMAX,
472  100,
473  0,
474  100000);
475 
476  // electrons
477  fLongitudinal_electron = tfs->make<TProfile>(
478  "fLongitudinal_electron", "longitudinal e- profile (reco);t;Q", LBINS, LMIN, LMAX);
480  tfs->make<TProfile>("fTransverse1_electron",
481  "transverse e- profile [0 <= t < 1] (reco);dist (cm);Q",
482  TBINS,
483  TMIN,
484  TMAX);
486  tfs->make<TProfile>("fTransverse2_electron",
487  "transverse e- profile [1 <= t < 2] (reco);dist (cm);Q",
488  TBINS,
489  TMIN,
490  TMAX);
492  tfs->make<TProfile>("fTransverse3_electron",
493  "transverse e- profile [2 <= t < 3] (reco);dist (cm);Q",
494  TBINS,
495  TMIN,
496  TMAX);
498  tfs->make<TProfile>("fTransverse4_electron",
499  "transverse e- profile [3 <= t < 4] (reco);dist (cm);Q",
500  TBINS,
501  TMIN,
502  TMAX);
504  tfs->make<TProfile>("fTransverse5_electron",
505  "transverse e- profile [4 <= t < 5] (reco);dist (cm);Q",
506  TBINS,
507  TMIN,
508  TMAX);
509 
510  // photons
511  fLongitudinal_photon = tfs->make<TProfile>(
512  "fLongitudinal_photon", "longitudinal photon profile (reco);t;Q", LBINS, LMIN, LMAX);
514  tfs->make<TProfile>("fTransverse1_photon",
515  "transverse photon profile [0 <= t < 1] (reco);dist (cm);Q",
516  TBINS,
517  TMIN,
518  TMAX);
520  tfs->make<TProfile>("fTransverse2_photon",
521  "transverse photon profile [1 <= t < 2] (reco);dist (cm);Q",
522  TBINS,
523  TMIN,
524  TMAX);
526  tfs->make<TProfile>("fTransverse3_photon",
527  "transverse photon profile [2 <= t < 3] (reco);dist (cm);Q",
528  TBINS,
529  TMIN,
530  TMAX);
532  tfs->make<TProfile>("fTransverse4_photon",
533  "transverse photon profile [3 <= t < 4] (reco);dist (cm);Q",
534  TBINS,
535  TMIN,
536  TMAX);
538  tfs->make<TProfile>("fTransverse5_photon",
539  "transverse photon profile [4 <= t < 5] (reco);dist (cm);Q",
540  TBINS,
541  TMIN,
542  TMAX);
543 
544  // other
545  fLongitudinal_other = tfs->make<TProfile>(
546  "fLongitudinal_other", "longitudinal other profile (reco);t;Q", LBINS, LMIN, LMAX);
548  tfs->make<TProfile>("fTransverse1_other",
549  "transverse other profile [0 <= t < 1] (reco);dist (cm);Q",
550  TBINS,
551  TMIN,
552  TMAX);
554  tfs->make<TProfile>("fTransverse2_other",
555  "transverse other profile [1 <= t < 2] (reco);dist (cm);Q",
556  TBINS,
557  TMIN,
558  TMAX);
560  tfs->make<TProfile>("fTransverse3_other",
561  "transverse other profile [2 <= t < 3] (reco);dist (cm);Q",
562  TBINS,
563  TMIN,
564  TMAX);
566  tfs->make<TProfile>("fTransverse4_other",
567  "transverse other profile [3 <= t < 4] (reco);dist (cm);Q",
568  TBINS,
569  TMIN,
570  TMAX);
572  tfs->make<TProfile>("fTransverse5_other",
573  "transverse other profile [4 <= t < 5] (reco);dist (cm);Q",
574  TBINS,
575  TMIN,
576  TMAX);
577 
578 } // beginJob
579 
580 // -------------------------------------------------
581 
583 {
584 
586  std::vector<art::Ptr<recob::Hit>> hitlist;
587  if (evt.getByLabel(fHitModuleLabel, hitListHandle)) art::fill_ptr_vector(hitlist, hitListHandle);
588 
590  std::vector<art::Ptr<sim::SimChannel>> simchanlist;
591  if (evt.getByLabel(fDigitModuleLabel, scListHandle))
592  art::fill_ptr_vector(simchanlist, scListHandle);
593 
594  art::Handle<std::vector<recob::Shower>> showerListHandle;
595  std::vector<art::Ptr<recob::Shower>> showerlist;
596  if (evt.getByLabel(fShowerModuleLabel, showerListHandle))
597  art::fill_ptr_vector(showerlist, showerListHandle);
598 
599  art::Handle<std::vector<simb::MCTruth>> mctruthListHandle;
600  std::vector<art::Ptr<simb::MCTruth>> mclist;
601  if (evt.getByLabel(fGenieGenModuleLabel, mctruthListHandle))
602  art::fill_ptr_vector(mclist, mctruthListHandle);
603 
604  art::FindManyP<recob::Hit> shwfm(showerListHandle, evt, fShowerModuleLabel);
605 
606  if (empty(mclist)) return;
607 
608  auto const clock_data = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
609  auto const det_prop =
611 
612  art::Ptr<simb::MCTruth> mctruth = mclist[0];
613  if (mctruth->NeutrinoSet()) {
614  if (std::abs(mctruth->GetNeutrino().Nu().PdgCode()) == 12 &&
615  mctruth->GetNeutrino().CCNC() == 0) {
616  double elep = mctruth->GetNeutrino().Lepton().E();
617  if (showerlist.size()) {
618  std::vector<art::Ptr<recob::Hit>> showerhits = shwfm.at(0);
619  showerProfile(clock_data,
620  det_prop,
621  showerhits,
622  showerlist[0]->ShowerStart(),
623  showerlist[0]->Direction(),
624  elep);
625  }
626  showerProfileTrue(clock_data, det_prop, hitlist, elep);
627  showerProfileTrue(simchanlist, mctruth->GetNeutrino().Lepton());
628  }
629  }
630 } // analyze
631 
632 // -------------------------------------------------
633 
635  detinfo::DetectorPropertiesData const& detProp,
636  std::vector<art::Ptr<recob::Hit>> showerhits,
637  TVector3 shwvtx,
638  TVector3 shwdir,
639  double elep)
640 {
641  auto const& wireReadoutGeom = art::ServiceHandle<geo::WireReadout const>()->Get();
642  auto const& collectionPlane = wireReadoutGeom.Plane(collectionPlaneID);
643 
644  double shwVtxTime = detProp.ConvertXToTicks(shwvtx[0], collectionPlaneID);
645  using geo::vect::toPoint;
646  double shwVtxWire = collectionPlane.WireCoordinate(toPoint(shwvtx));
647 
648  double shwTwoTime = detProp.ConvertXToTicks(shwvtx[0] + shwdir[0], collectionPlaneID);
649  double shwTwoWire = collectionPlane.WireCoordinate(toPoint(shwvtx + shwdir));
650 
651  TH1F* ltemp = new TH1F("ltemp", "ltemp", LBINS, LMIN, LMAX);
652  TH1F* ttemp = new TH1F("ttemp", "ttemp", TBINS, TMIN, TMAX);
653 
654  TH1F* ttemp_1 = new TH1F("ttemp_1", "ttemp_1", TBINS, TMIN, TMAX);
655  TH1F* ttemp_2 = new TH1F("ttemp_2", "ttemp_2", TBINS, TMIN, TMAX);
656  TH1F* ttemp_3 = new TH1F("ttemp_3", "ttemp_3", TBINS, TMIN, TMAX);
657  TH1F* ttemp_4 = new TH1F("ttemp_4", "ttemp_4", TBINS, TMIN, TMAX);
658  TH1F* ttemp_5 = new TH1F("ttemp_5", "ttemp_5", TBINS, TMIN, TMAX);
659 
660  for (size_t i = 0; i < showerhits.size(); ++i) {
661  if (showerhits[i]->WireID().Plane != collectionPlaneID.Plane) continue;
662 
663  double wirePitch = wireReadoutGeom.Plane(showerhits[i]->WireID()).WirePitch();
664  double tickToDist = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
665  tickToDist *= 1.e-3 * sampling_rate(clockData); // 1e-3 is conversion of 1/us to 1/ns
666 
667  double xvtx = shwVtxTime * tickToDist;
668  double yvtx = shwVtxWire * wirePitch;
669 
670  double xtwo = shwTwoTime * tickToDist;
671  double ytwo = shwTwoWire * wirePitch;
672 
673  double xtwoorth = (ytwo - yvtx) + xvtx;
674  double ytwoorth = -(xtwo - xvtx) + yvtx;
675 
676  double xhit = showerhits[i]->PeakTime() * tickToDist;
677  double yhit = showerhits[i]->WireID().Wire * wirePitch;
678 
679  double ldist = std::abs((ytwoorth - yvtx) * xhit - (xtwoorth - xvtx) * yhit + xtwoorth * yvtx -
680  ytwoorth * xvtx) /
681  std::hypot(ytwoorth - yvtx, xtwoorth - xvtx);
682  double tdist = ((ytwo - yvtx) * xhit - (xtwo - xvtx) * yhit + xtwo * yvtx - ytwo * xvtx) /
683  std::hypot(ytwo - yvtx, xtwo - xvtx);
684 
685  double to3D = 1. / std::hypot(xvtx - xtwo,
686  yvtx - ytwo); // distance between two points in 3D space is one
687  ldist *= to3D;
688  tdist *= to3D;
689  double t = ldist / X0;
690 
691  double Q = showerhits[i]->Integral() *
692  fCalorimetryAlg.LifetimeCorrection(clockData, detProp, showerhits[i]->PeakTime());
693 
694  ltemp->Fill(t, Q);
695  ttemp->Fill(tdist, Q);
696 
697  if (t < 1)
698  ttemp_1->Fill(tdist, Q);
699  else if (t < 2)
700  ttemp_2->Fill(tdist, Q);
701  else if (t < 3)
702  ttemp_3->Fill(tdist, Q);
703  else if (t < 4)
704  ttemp_4->Fill(tdist, Q);
705  else if (t < 5)
706  ttemp_5->Fill(tdist, Q);
707 
708  } // loop through showerhits
709 
710  for (int i = 0; i < LBINS; ++i) {
711  if (ltemp->GetBinContent(i + 1) == 0) continue;
712  fShowerProfileRecoLong->Fill(ltemp->GetBinCenter(i + 1), ltemp->GetBinContent(i + 1));
713  fShowerProfileRecoLong2D->Fill(ltemp->GetBinCenter(i + 1), elep, ltemp->GetBinContent(i + 1));
714  fLongitudinal->Fill(ltemp->GetBinCenter(i + 1), elep, ltemp->GetBinContent(i + 1));
715  }
716 
717  for (int i = 0; i < TBINS; ++i) {
718  if (ttemp->GetBinContent(i + 1) == 0) continue;
719  fShowerProfileRecoTrans->Fill(ttemp->GetBinCenter(i + 1), ttemp->GetBinContent(i + 1));
720  fShowerProfileRecoTrans2D->Fill(ttemp->GetBinCenter(i + 1), elep, ttemp->GetBinContent(i + 1));
721  fTransverse->Fill(ttemp->GetBinCenter(i + 1), elep, ttemp->GetBinContent(i + 1));
722 
724  ttemp_1->GetBinCenter(i + 1), elep, ttemp_1->GetBinContent(i + 1));
726  ttemp_2->GetBinCenter(i + 1), elep, ttemp_2->GetBinContent(i + 1));
728  ttemp_3->GetBinCenter(i + 1), elep, ttemp_3->GetBinContent(i + 1));
730  ttemp_4->GetBinCenter(i + 1), elep, ttemp_4->GetBinContent(i + 1));
732  ttemp_5->GetBinCenter(i + 1), elep, ttemp_5->GetBinContent(i + 1));
733 
734  fTransverse_1->Fill(ttemp_1->GetBinCenter(i + 1), elep, ttemp_1->GetBinContent(i + 1));
735  fTransverse_2->Fill(ttemp_2->GetBinCenter(i + 1), elep, ttemp_2->GetBinContent(i + 1));
736  fTransverse_3->Fill(ttemp_3->GetBinCenter(i + 1), elep, ttemp_3->GetBinContent(i + 1));
737  fTransverse_4->Fill(ttemp_4->GetBinCenter(i + 1), elep, ttemp_4->GetBinContent(i + 1));
738  fTransverse_5->Fill(ttemp_5->GetBinCenter(i + 1), elep, ttemp_5->GetBinContent(i + 1));
739  }
740 } // showerProfile
741 
742 // -------------------------------------------------
743 
745  detinfo::DetectorClocksData const& clockData,
746  detinfo::DetectorPropertiesData const& detProp,
748  double elep)
749 {
750  auto const& wireReadoutGeom = art::ServiceHandle<geo::WireReadout const>()->Get();
751  auto const& collectionPlane = wireReadoutGeom.Plane(collectionPlaneID);
754  std::map<int, double> trkID_E;
755 
756  TH1F* ltemp = new TH1F("ltemp", "ltemp", LBINS, LMIN, LMAX);
757  TH1F* ttemp = new TH1F("ttemp", "ttemp", TBINS, TMIN, TMAX);
758 
759  TH1F* ttemp_1 = new TH1F("ttemp_1", "ttemp_1", TBINS, TMIN, TMAX);
760  TH1F* ttemp_2 = new TH1F("ttemp_2", "ttemp_2", TBINS, TMIN, TMAX);
761  TH1F* ttemp_3 = new TH1F("ttemp_3", "ttemp_3", TBINS, TMIN, TMAX);
762  TH1F* ttemp_4 = new TH1F("ttemp_4", "ttemp_4", TBINS, TMIN, TMAX);
763  TH1F* ttemp_5 = new TH1F("ttemp_5", "ttemp_5", TBINS, TMIN, TMAX);
764 
765  double xvtx = -999;
766  double yvtx = -999;
767  double zvtx = -999;
768  double xtwo = -999;
769  double ytwo = -999;
770  double ztwo = -999;
771  double shwvtxT = -999;
772  double shwvtxW = -999;
773  double shwtwoT = -999;
774  double shwtwoW = -999;
775 
776  double shwvtxx = -999;
777  double shwvtxy = -999;
778  double shwtwox = -999;
779  double shwtwoy = -999;
780  double xtwoorth = -999;
781  double ytwoorth = -999;
782 
783  double wirePitch = -999;
784  double tickToDist = -999;
785 
786  bool foundParent = false;
787 
788  for (auto const& hit : allhits) {
789  if (hit->WireID().Plane != collectionPlaneID.Plane) continue;
790 
791  // art::Ptr<recob::Hit> hit = allhits[i];
792  std::vector<sim::TrackIDE> trackIDs = btserv->HitToEveTrackIDEs(clockData, hit);
793 
794  for (size_t j = 0; j < trackIDs.size(); ++j) {
795  // only want energy associated with the electron and electron must have neutrino mother
796  if (std::abs((piserv->TrackIdToParticle_P(trackIDs[j].trackID))->PdgCode()) != 11) continue;
797  if (std::abs((piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Mother()) != 0) continue;
798 
799  if (!foundParent) {
800  xvtx = (piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Vx();
801  yvtx = (piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Vy();
802  zvtx = (piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Vz();
803 
804  xtwo = xvtx + (piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Px();
805  ytwo = yvtx + (piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Py();
806  ztwo = zvtx + (piserv->TrackIdToParticle_P(trackIDs[j].trackID))->Pz();
807 
808  shwvtxT = detProp.ConvertXToTicks(xvtx, collectionPlaneID);
809  geo::Point_t const vtx{xvtx, yvtx, zvtx};
810  shwvtxW = collectionPlane.WireCoordinate(vtx);
811 
812  shwtwoT = detProp.ConvertXToTicks(xtwo, collectionPlaneID);
813  geo::Point_t const vtwo{xtwo, ytwo, ztwo};
814  shwtwoW = collectionPlane.WireCoordinate(vtwo);
815 
816  wirePitch = wireReadoutGeom.Plane(hit->WireID()).WirePitch();
817  tickToDist = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
818  tickToDist *= 1.e-3 * sampling_rate(clockData); // 1e-3 is conversion of 1/us to 1/ns
819 
820  shwvtxx = shwvtxT * tickToDist;
821  shwvtxy = shwvtxW * wirePitch;
822 
823  shwtwox = shwtwoT * tickToDist;
824  shwtwoy = shwtwoW * wirePitch;
825 
826  xtwoorth = (shwtwoy - shwvtxy) + shwvtxx;
827  ytwoorth = -(shwtwox - shwvtxx) + shwvtxy;
828 
829  foundParent = true;
830  }
831  double xhit = hit->PeakTime() * tickToDist;
832  double yhit = hit->WireID().Wire * wirePitch;
833 
834  double ldist = abs((ytwoorth - shwvtxy) * xhit - (xtwoorth - shwvtxx) * yhit +
835  xtwoorth * shwvtxy - ytwoorth * shwvtxx) /
836  std::hypot(ytwoorth - shwvtxy, xtwoorth - shwvtxx);
837  double tdist = ((shwtwoy - shwvtxy) * xhit - (shwtwox - shwvtxx) * yhit + shwtwox * shwvtxy -
838  shwtwoy * shwvtxx) /
839  std::hypot(shwtwoy - shwvtxy, shwtwox - shwvtxx);
840 
841  double to3D = std::hypot(xvtx - xtwo, yvtx - ytwo, zvtx - ztwo) /
842  std::hypot(shwvtxx - shwtwox,
843  shwvtxy - shwtwoy); // distance between two points in 3D space is one
844  ldist *= to3D;
845  tdist *= to3D;
846  double t = ldist / X0;
847 
848  double energy = trackIDs[j].energy;
849 
850  ltemp->Fill(t, energy);
851  ttemp->Fill(tdist, energy);
852 
853  if (t < 1)
854  ttemp_1->Fill(tdist, energy);
855  else if (t < 2)
856  ttemp_2->Fill(tdist, energy);
857  else if (t < 3)
858  ttemp_3->Fill(tdist, energy);
859  else if (t < 4)
860  ttemp_4->Fill(tdist, energy);
861  else if (t < 5)
862  ttemp_5->Fill(tdist, energy);
863 
864  } // loop through track IDE
865 
866  } // loop through all hits
867 
868  for (int i = 0; i < LBINS; ++i) {
869  if (ltemp->GetBinContent(i + 1) == 0) continue;
870  fShowerProfileHitLong->Fill(ltemp->GetBinCenter(i + 1), ltemp->GetBinContent(i + 1));
871  fShowerProfileHitLong2D->Fill(ltemp->GetBinCenter(i + 1), elep, ltemp->GetBinContent(i + 1));
872  }
873 
874  for (int i = 0; i < TBINS; ++i) {
875  if (ttemp->GetBinContent(i + 1) == 0) continue;
876  fShowerProfileHitTrans->Fill(ttemp->GetBinCenter(i + 1), ttemp->GetBinContent(i + 1));
877  fShowerProfileHitTrans2D->Fill(ttemp->GetBinCenter(i + 1), elep, ttemp->GetBinContent(i + 1));
878 
880  ttemp_1->GetBinCenter(i + 1), elep, ttemp_1->GetBinContent(i + 1));
882  ttemp_2->GetBinCenter(i + 1), elep, ttemp_2->GetBinContent(i + 1));
884  ttemp_3->GetBinCenter(i + 1), elep, ttemp_3->GetBinContent(i + 1));
886  ttemp_4->GetBinCenter(i + 1), elep, ttemp_4->GetBinContent(i + 1));
888  ttemp_5->GetBinCenter(i + 1), elep, ttemp_5->GetBinContent(i + 1));
889  }
890 } // showerProfileTrue
891 
892 // -------------------------------------------------
893 
896  simb::MCParticle electron)
897 {
899  auto const& wireReadoutGeom = art::ServiceHandle<geo::WireReadout const>()->Get();
900 
901  std::vector<sim::MCEnDep> alledep;
902 
903  TH1F* ltemp = new TH1F("ltemp", "ltemp", LBINS, LMIN, LMAX);
904  TH1F* ttemp = new TH1F("ttemp", "ttemp", TBINS, TMIN, TMAX);
905 
906  TH1F* ttemp_1 = new TH1F("ttemp_1", "ttemp_1", TBINS, TMIN, TMAX);
907  TH1F* ttemp_2 = new TH1F("ttemp_2", "ttemp_2", TBINS, TMIN, TMAX);
908  TH1F* ttemp_3 = new TH1F("ttemp_3", "ttemp_3", TBINS, TMIN, TMAX);
909  TH1F* ttemp_4 = new TH1F("ttemp_4", "ttemp_4", TBINS, TMIN, TMAX);
910  TH1F* ttemp_5 = new TH1F("ttemp_5", "ttemp_5", TBINS, TMIN, TMAX);
911 
912  // get electron energy depositions
913  for (size_t i = 0; i < allchan.size(); ++i) {
914  art::Ptr<sim::SimChannel> simchan = allchan[i];
915  if (wireReadoutGeom.View(simchan->Channel()) != geo::kV) continue;
916  auto tdc_ide_map = simchan->TDCIDEMap();
917 
918  for (auto const& tdc_ide_pair : tdc_ide_map) {
919  for (auto const& ide : tdc_ide_pair.second) {
920  if (piserv->TrackIdToMotherParticle_P(ide.trackID) == nullptr) continue;
921  if (std::abs(piserv->TrackIdToMotherParticle_P(ide.trackID)->PdgCode()) != 11) continue;
922 
924  edep.SetVertex(ide.x, ide.y, ide.z);
925  edep.SetEnergy(ide.energy);
926  edep.SetTrackId(ide.trackID);
927 
928  alledep.push_back(edep);
929 
930  } // loop through ide
931  } // loop through tdc_ide
932 
933  } // loop through channels
934 
935  double x0 = electron.Vx();
936  double y0 = electron.Vy();
937  double z0 = electron.Vz();
938 
939  double x2 = electron.Px();
940  double y2 = electron.Py();
941  double z2 = electron.Pz();
942 
943  TVector3 v0(x2, y2, z2);
944  v0 = v0.Unit();
945 
946  for (size_t i = 0; i < alledep.size(); ++i) {
947  double x = (double)alledep[i].Vertex()[0];
948  double y = (double)alledep[i].Vertex()[1];
949  double z = (double)alledep[i].Vertex()[2];
950 
951  TVector3 v1(x - x0, y - y0, z - z0);
952 
953  double ldist = v0.Dot(v1);
954  double t = ldist / X0;
955  double tdist = (v0.Orthogonal()).Dot(v1);
956 
957  double energy = alledep[i].Energy();
958 
959  ltemp->Fill(t, energy);
960  ttemp->Fill(tdist, energy);
961 
962  if (t < 1)
963  ttemp_1->Fill(tdist, energy);
964  else if (t < 2)
965  ttemp_2->Fill(tdist, energy);
966  else if (t < 3)
967  ttemp_3->Fill(tdist, energy);
968  else if (t < 4)
969  ttemp_4->Fill(tdist, energy);
970  else if (t < 5)
971  ttemp_5->Fill(tdist, energy);
972  }
973 
974  for (int i = 0; i < LBINS; ++i) {
975  if (ltemp->GetBinContent(i + 1) == 0) continue;
976  fShowerProfileSimLong->Fill(ltemp->GetBinCenter(i + 1), ltemp->GetBinContent(i + 1));
978  ltemp->GetBinCenter(i + 1), electron.E(), ltemp->GetBinContent(i + 1));
979  }
980 
981  for (int i = 0; i < TBINS; ++i) {
982  if (ttemp->GetBinContent(i + 1) == 0) continue;
983  fShowerProfileSimTrans->Fill(ttemp->GetBinCenter(i + 1), ttemp->GetBinContent(i + 1));
985  ttemp->GetBinCenter(i + 1), electron.E(), ttemp->GetBinContent(i + 1));
986 
988  ttemp_1->GetBinCenter(i + 1), electron.E(), ttemp_1->GetBinContent(i + 1));
990  ttemp_2->GetBinCenter(i + 1), electron.E(), ttemp_2->GetBinContent(i + 1));
992  ttemp_3->GetBinCenter(i + 1), electron.E(), ttemp_3->GetBinContent(i + 1));
994  ttemp_4->GetBinCenter(i + 1), electron.E(), ttemp_4->GetBinContent(i + 1));
996  ttemp_5->GetBinCenter(i + 1), electron.E(), ttemp_5->GetBinContent(i + 1));
997  }
998 } // showerProfileTrue
999 
1000 // -------------------------------------------------
1001 
Float_t x
Definition: compare.C:6
double E(const int i=0) const
Definition: MCParticle.h:234
TCShowerTemplateMaker(fhicl::ParameterSet const &pset)
int PdgCode() const
Definition: MCParticle.h:213
int CCNC() const
Definition: MCNeutrino.h:148
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:77
double Py(const int i=0) const
Definition: MCParticle.h:232
void SetEnergy(float e)
Definition: MCDataHolder.h:36
void analyze(const art::Event &evt) override
Planes which measure V.
Definition: geo_types.h:132
const simb::MCParticle * TrackIdToParticle_P(int id) const
Declaration of signal hit object.
Float_t y
Definition: compare.C:6
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
Double_t z
Definition: plot.C:276
The data type to uniquely identify a Plane.
Definition: geo_types.h:364
double Px(const int i=0) const
Definition: MCParticle.h:231
double Temperature() const
In kelvin.
constexpr auto abs(T v)
Returns the absolute value of the argument.
STL namespace.
cout<< "Opened file "<< fin<< " ixs= "<< ixs<< endl;if(ixs==0) hhh=(TH1F *) fff-> Get("h1")
Definition: AddMC.C:8
void SetTrackId(unsigned int id)
Definition: MCDataHolder.h:38
void showerProfile(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, std::vector< art::Ptr< recob::Hit >> showerhits, TVector3 shwvtx, TVector3 shwdir, double elep)
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
Float_t y2[n_points_geant4]
Definition: compare.C:26
double Efield(unsigned int planegap=0) const
kV/cm
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
void showerProfileTrue(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, std::vector< art::Ptr< recob::Hit >> allhits, double elep)
const simb::MCParticle & Lepton() const
Definition: MCNeutrino.h:147
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
parameter set interface
double energy
Definition: plottest35.C:25
void SetVertex(float x, float y, float z)
Definition: MCDataHolder.h:29
double ConvertXToTicks(double X, int p, int t, int c) const
Point_t toPoint(Point const &p)
Convert the specified point into a geo::Point_t.
double DriftVelocity(double efield=0., double temperature=0.) const
cm/us
Double_t edep
Definition: macro.C:13
Detector simulation of raw signals on wires.
raw::ChannelID_t Channel() const
Returns the readout channel this object describes.
Definition: SimChannel.h:323
double Vx(const int i=0) const
Definition: MCParticle.h:222
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Encapsulate the construction of a single detector plane .
Contains all timing reference information for the detector.
std::vector< sim::TrackIDE > HitToEveTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
double Pz(const int i=0) const
Definition: MCParticle.h:233
object containing MC truth information necessary for making RawDigits and doing back tracking ...
double Vz(const int i=0) const
Definition: MCParticle.h:224
TDCIDEs_t const & TDCIDEMap() const
Returns all the deposited energy information as stored.
Definition: SimChannel.h:319
boost::graph_traits< ModuleGraph >::vertex_descriptor Vertex
Definition: ModuleGraph.h:25
Direction
Definition: types.h:12
bool NeutrinoSet() const
Definition: MCTruth.h:78
TCEvent evt
Definition: DataStructs.cxx:8
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
const simb::MCParticle * TrackIdToMotherParticle_P(int id) const
double LifetimeCorrection(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, double time, double T0=0) const
Float_t x2[n_points_geant4]
Definition: compare.C:26
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
double Vy(const int i=0) const
Definition: MCParticle.h:223
decltype(auto) constexpr empty(T &&obj)
ADL-aware version of std::empty.
Definition: StdUtils.h:109