LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
TCShowerElectronLikelihood_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: 8/3/18
6 // -------------------------------------------------
7 
8 // Framework includes
14 #include "art_root_io/TFileService.h"
16 #include "cetlib/pow.h"
17 #include "fhiclcpp/ParameterSet.h"
18 
28 
29 #include "TFile.h"
30 #include "TH1.h"
31 #include "TH3.h"
32 #include "TProfile.h"
33 #include "TProfile2D.h"
34 
35 namespace shower {
36 
38  public:
39  explicit TCShowerElectronLikelihood(fhicl::ParameterSet const& pset);
40 
41  private:
42  void beginJob() override;
43  void analyze(const art::Event& evt) override;
44 
45  void getShowerProfile(detinfo::DetectorClocksData const& clockdata,
46  detinfo::DetectorPropertiesData const& detProp,
48  TVector3 shwvtx,
49  TVector3 shwdir);
50  void findEnergyBin();
51  double getLongLikelihood();
52  double getTranLikelihood();
53 
54  void resetProfiles();
55 
56  std::string fTemplateFile;
57  std::string fROOTfile;
58 
59  TH3F* longTemplate;
60  TH3F* tranTemplate;
66  TProfile2D* longTemplateProf2D;
67  TProfile2D* tranTemplateProf2D;
68  TProfile2D* tranTemplateProf2D_1;
69  TProfile2D* tranTemplateProf2D_2;
70  TProfile2D* tranTemplateProf2D_3;
71  TProfile2D* tranTemplateProf2D_4;
72  TProfile2D* tranTemplateProf2D_5;
73 
74  //TTree* fTree;
75  TH1F* energyDist;
78  // just for printing purposes
79  TH1F* longProfHist;
85 
86  TH1F* longProfile;
87  TH1F* tranProfile;
93 
95  double energyChi2;
96  int maxt;
97 
98  const int LBINS = 20;
99  const int LMIN = 0;
100  const int LMAX = 5;
101 
102  const int TBINS = 20;
103  const int TMIN = -5;
104  const int TMAX = 5;
105 
106  const double X0 = 14;
107 
108  std::string fHitModuleLabel;
109  std::string fShowerModuleLabel;
110  std::string fTemplateModuleLabel;
111  std::string fGenieGenModuleLabel;
112  std::string fDigitModuleLabel;
113 
115 
116  }; // class TCShowerElectronLikelihood
117 
118 } // shower
119 
120 // -------------------------------------------------
121 
123  : EDAnalyzer(pset)
124  , fHitModuleLabel(pset.get<std::string>("HitModuleLabel", "trajcluster"))
125  , fShowerModuleLabel(pset.get<std::string>("ShowerModuleLabel", "tcshower"))
126  , fGenieGenModuleLabel(pset.get<std::string>("GenieGenModuleLabel", "generator"))
127  , fCalorimetryAlg(pset.get<fhicl::ParameterSet>("CalorimetryAlg"))
128 {
129  fTemplateFile = pset.get<std::string>("TemplateFile");
130  cet::search_path sp("FW_SEARCH_PATH");
131  if (!sp.find_file(fTemplateFile, fROOTfile))
132  throw cet::exception("TCShowerElectronLikelihood")
133  << "cannot find the root template file: \n"
134  << fTemplateFile << "\n bail ungracefully.\n";
135 
136  std::unique_ptr<TFile> file{TFile::Open(fROOTfile.c_str())};
137 
138  longTemplate = file->Get<TH3F>("tcshowertemplate/fLongitudinal");
139  tranTemplate = file->Get<TH3F>("tcshowertemplate/fTransverse");
140  tranTemplate_1 = file->Get<TH3F>("tcshowertemplate/fTransverse_1");
141  tranTemplate_2 = file->Get<TH3F>("tcshowertemplate/fTransverse_2");
142  tranTemplate_3 = file->Get<TH3F>("tcshowertemplate/fTransverse_3");
143  tranTemplate_4 = file->Get<TH3F>("tcshowertemplate/fTransverse_4");
144  tranTemplate_5 = file->Get<TH3F>("tcshowertemplate/fTransverse_5");
145 
146  longTemplateProf2D = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoLong2D");
147  tranTemplateProf2D = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoTrans2D");
148  tranTemplateProf2D_1 = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoTrans2D_1");
149  tranTemplateProf2D_2 = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoTrans2D_2");
150  tranTemplateProf2D_3 = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoTrans2D_3");
151  tranTemplateProf2D_4 = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoTrans2D_4");
152  tranTemplateProf2D_5 = file->Get<TProfile2D>("tcshowertemplate/fShowerProfileRecoTrans2D_5");
153 
154  longProfile = new TH1F("longProfile", "longitudinal shower profile;t;Q", LBINS, LMIN, LMAX);
155  tranProfile = new TH1F("tranProfile", "transverse shower profile;dist (cm);Q", TBINS, TMIN, TMAX);
156 
157  tranProfile_1 = new TH1F(
158  "tranProfile_1", "transverse shower profile [0 <= t < 1];dist (cm);Q", TBINS, TMIN, TMAX);
159  tranProfile_2 = new TH1F(
160  "tranProfile_2", "transverse shower profile [1 <= t < 2];dist (cm);Q", TBINS, TMIN, TMAX);
161  tranProfile_3 = new TH1F(
162  "tranProfile_3", "transverse shower profile [2 <= t < 3];dist (cm);Q", TBINS, TMIN, TMAX);
163  tranProfile_4 = new TH1F(
164  "tranProfile_4", "transverse shower profile [3 <= t < 4];dist (cm);Q", TBINS, TMIN, TMAX);
165  tranProfile_5 = new TH1F(
166  "tranProfile_5", "transverse shower profile [4 <= t < 5];dist (cm);Q", TBINS, TMIN, TMAX);
167 } // TCShowerElectronLikelihood
168 
169 // -------------------------------------------------
170 
172 {
173 
175 
176  energyDist = tfs->make<TH1F>("energyDist", "true energy - guess energy", 41, -20.5, 20.5);
177  longLikelihoodHist = tfs->make<TH1F>("longLikelihoodHist", "longitudinal likelihood", 20, 0, 2);
178  tranLikelihoodHist = tfs->make<TH1F>("tranLikelihoodHist", "transverse likelihood", 20, 0, 3);
179 
180  // just for printing purposes
181  longProfHist =
182  tfs->make<TH1F>("longProfHist", "longitudinal e- profile (reco);t;Q", LBINS, LMIN, LMAX);
183  tranProfHist_1 = tfs->make<TH1F>(
184  "tranProfHist", "transverse e- profile (reco) [0 <= t < 1];dist (cm);Q", TBINS, TMIN, TMAX);
185  tranProfHist_2 = tfs->make<TH1F>(
186  "tranProfHist", "transverse e- profile (reco) [1 <= t < 2];dist (cm);Q", TBINS, TMIN, TMAX);
187  tranProfHist_3 = tfs->make<TH1F>(
188  "tranProfHist", "transverse e- profile (reco) [2 <= t < 3];dist (cm);Q", TBINS, TMIN, TMAX);
189  tranProfHist_4 = tfs->make<TH1F>(
190  "tranProfHist", "transverse e- profile (reco) [3 <= t < 4];dist (cm);Q", TBINS, TMIN, TMAX);
191  tranProfHist_5 = tfs->make<TH1F>(
192  "tranProfHist", "transverse e- profile (reco) [4 <= t < 5];dist (cm);Q", TBINS, TMIN, TMAX);
193 
194 } // beginJob
195 
196 // -------------------------------------------------
197 
199 {
200  resetProfiles();
201 
203  std::vector<art::Ptr<recob::Hit>> hitlist;
204  if (evt.getByLabel(fHitModuleLabel, hitListHandle)) art::fill_ptr_vector(hitlist, hitListHandle);
205 
206  art::Handle<std::vector<recob::Shower>> showerListHandle;
207  std::vector<art::Ptr<recob::Shower>> showerlist;
208  if (evt.getByLabel(fShowerModuleLabel, showerListHandle))
209  art::fill_ptr_vector(showerlist, showerListHandle);
210 
211  art::Handle<std::vector<simb::MCTruth>> mctruthListHandle;
212  std::vector<art::Ptr<simb::MCTruth>> mclist;
213  if (evt.getByLabel(fGenieGenModuleLabel, mctruthListHandle))
214  art::fill_ptr_vector(mclist, mctruthListHandle);
215 
216  art::FindManyP<recob::Hit> shwfm(showerListHandle, evt, fShowerModuleLabel);
217 
218  if (showerlist.size()) {
219  auto const clock_data =
221  auto const det_prop =
223  std::vector<art::Ptr<recob::Hit>> showerhits = shwfm.at(0);
225  clock_data, det_prop, showerhits, showerlist[0]->ShowerStart(), showerlist[0]->Direction());
226 
227  maxt = std::ceil((90 - showerlist[0]->ShowerStart().Z()) / X0);
228 
229  findEnergyBin();
230 
233 
234  // check true shower energy
235  if (mclist.size()) {
236  art::Ptr<simb::MCTruth> mctruth = mclist[0];
237  if (mctruth->NeutrinoSet()) {
238  if (std::abs(mctruth->GetNeutrino().Nu().PdgCode()) == 12 &&
239  mctruth->GetNeutrino().CCNC() == 0) {
240  double elep = mctruth->GetNeutrino().Lepton().E();
241  energyDist->Fill(elep - energyGuess);
242  } // cc nue
243  } // neutrinoset
244  } // mclist
245  } // showerlist
246 
253 } // analyze
254 
255 // -------------------------------------------------
256 
258 {
259 
260  longProfile->Reset();
261  tranProfile->Reset();
262  tranProfile_1->Reset();
263  tranProfile_2->Reset();
264  tranProfile_3->Reset();
265  tranProfile_4->Reset();
266  tranProfile_5->Reset();
267 
268  energyGuess = -9999;
269  energyChi2 = -9999;
270  maxt = -9999;
271 } // resetProfiles
272 
273 // -------------------------------------------------
274 
276  detinfo::DetectorClocksData const& clockdata,
277  detinfo::DetectorPropertiesData const& detProp,
278  std::vector<art::Ptr<recob::Hit>> showerhits,
279  TVector3 shwvtx,
280  TVector3 shwdir)
281 {
283 
284  auto collectionPlane = geo::PlaneID(0, 0, 1);
285 
286  using namespace geo::vect;
287  auto const shwvtx_p = toPoint(shwvtx);
288  auto const shwvtx_p2 = shwvtx_p + toVector(shwdir);
289 
290  double shwVtxTime = detProp.ConvertXToTicks(shwvtx_p.X(), collectionPlane);
291  double shwVtxWire = geom->WireCoordinate(shwvtx_p, collectionPlane);
292 
293  double shwTwoTime = detProp.ConvertXToTicks(shwvtx_p2.X(), collectionPlane);
294  double shwTwoWire = geom->WireCoordinate(shwvtx_p2, collectionPlane);
295 
296  for (size_t i = 0; i < showerhits.size(); ++i) {
297  if (showerhits[i]->WireID().Plane != collectionPlane.Plane) continue;
298 
299  double wirePitch = geom->WirePitch(showerhits[i]->WireID());
300  double tickToDist = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
301  tickToDist *= 1.e-3 * sampling_rate(clockdata); // 1e-3 is conversion of 1/us to 1/ns
302 
303  double xvtx = shwVtxTime * tickToDist;
304  double yvtx = shwVtxWire * wirePitch;
305 
306  double xtwo = shwTwoTime * tickToDist;
307  double ytwo = shwTwoWire * wirePitch;
308 
309  double xtwoorth = (ytwo - yvtx) + xvtx;
310  double ytwoorth = -(xtwo - xvtx) + yvtx;
311 
312  double xhit = showerhits[i]->PeakTime() * tickToDist;
313  double yhit = showerhits[i]->WireID().Wire * wirePitch;
314 
315  double ldist = std::abs((ytwoorth - yvtx) * xhit - (xtwoorth - xvtx) * yhit + xtwoorth * yvtx -
316  ytwoorth * xvtx) /
317  std::hypot(ytwoorth - yvtx, xtwoorth - xvtx);
318  double tdist = ((ytwo - yvtx) * xhit - (xtwo - xvtx) * yhit + xtwo * yvtx - ytwo * xvtx) /
319  std::hypot(ytwo - yvtx, xtwo - xvtx);
320 
321  double to3D = 1. / std::hypot(xvtx - xtwo,
322  yvtx - ytwo); // distance between two points in 3D space is one
323  ldist *= to3D;
324  tdist *= to3D;
325  double t = ldist / X0;
326 
327  double Q = showerhits[i]->Integral() *
328  fCalorimetryAlg.LifetimeCorrection(clockdata, detProp, showerhits[i]->PeakTime());
329 
330  longProfile->Fill(t, Q);
331  tranProfile->Fill(tdist, Q);
332 
333  if (t < 1)
334  tranProfile_1->Fill(tdist, Q);
335  else if (t < 2)
336  tranProfile_2->Fill(tdist, Q);
337  else if (t < 3)
338  tranProfile_3->Fill(tdist, Q);
339  else if (t < 4)
340  tranProfile_4->Fill(tdist, Q);
341  else if (t < 5)
342  tranProfile_5->Fill(tdist, Q);
343 
344  } // loop through showerhits
345 
346  return;
347 
348 } // getShowerProfile
349 
350 // -------------------------------------------------
351 
353 {
354 
355  if (longProfile->GetNbinsX() != longTemplate->GetNbinsX())
356  throw cet::exception("TCShowerElectronLikelihood")
357  << "Bin mismatch in longitudinal profile template \n";
358 
359  if (tranProfile->GetNbinsX() != tranTemplate->GetNbinsX())
360  throw cet::exception("TCShowerElectronLikelihood")
361  << "Bin mismatch in transverse profile template \n";
362 
363  double chi2min = 999999;
364  double bestbin = -1;
365 
366  int ebins = longTemplate->GetNbinsY();
367  int lbins = longTemplate->GetNbinsX();
368  int tbins = tranTemplate->GetNbinsX();
369 
370  for (int i = 0; i < ebins; ++i) {
371  double thischi2 = 0;
372 
373  TProfile* ltemp = longTemplateProf2D->ProfileX("_x", i + 1, i + 1);
374  TProfile* ttemp_1 = tranTemplateProf2D_1->ProfileX("_x_1", i + 1, i + 1);
375  TProfile* ttemp_2 = tranTemplateProf2D_2->ProfileX("_x_2", i + 1, i + 1);
376  TProfile* ttemp_3 = tranTemplateProf2D_3->ProfileX("_x_3", i + 1, i + 1);
377  TProfile* ttemp_4 = tranTemplateProf2D_4->ProfileX("_x_4", i + 1, i + 1);
378  TProfile* ttemp_5 = tranTemplateProf2D_5->ProfileX("_x_5", i + 1, i + 1);
379 
380  int nlbins = 0;
381  int ntbins = 0;
382 
383  for (int j = 0; j < lbins; ++j) {
384  double obs = longProfile->GetBinContent(j + 1);
385  double exp = ltemp->GetBinContent(j + 1);
386  if (obs != 0) {
387  thischi2 += cet::square(obs - exp) / exp;
388  ++nlbins;
389  }
390  } // loop through longitudinal bins
391 
392  for (int j = 0; j < tbins; ++j) {
393  double obs = tranProfile_1->GetBinContent(j + 1);
394  double exp = ttemp_1->GetBinContent(j + 1);
395  if (obs != 0) {
396  thischi2 += cet::square(obs - exp) / exp;
397  ++ntbins;
398  }
399 
400  obs = tranProfile_2->GetBinContent(j + 1);
401  exp = ttemp_2->GetBinContent(j + 1);
402  if (obs != 0) {
403  thischi2 += cet::square(obs - exp) / exp;
404  ++ntbins;
405  }
406 
407  obs = tranProfile_3->GetBinContent(j + 1);
408  exp = ttemp_3->GetBinContent(j + 1);
409  if (obs != 0) {
410  thischi2 += cet::square(obs - exp) / exp;
411  ++ntbins;
412  }
413 
414  obs = tranProfile_4->GetBinContent(j + 1);
415  exp = ttemp_4->GetBinContent(j + 1);
416  if (obs != 0) {
417  thischi2 += cet::square(obs - exp) / exp;
418  ++ntbins;
419  }
420 
421  obs = tranProfile_5->GetBinContent(j + 1);
422  exp = ttemp_5->GetBinContent(j + 1);
423  if (obs != 0) {
424  thischi2 += cet::square(obs - exp) / exp;
425  ++ntbins;
426  }
427  } // loop through longitudinal bins
428 
429  thischi2 /= (nlbins + ntbins);
430 
431  if (thischi2 < chi2min) {
432  chi2min = thischi2;
433  bestbin = i;
434  }
435 
436  } // loop through energy bins
437 
438  energyChi2 = chi2min;
439  energyGuess = bestbin + 1;
440 } // findEnergyBin
441 
442 // -------------------------------------------------
443 
445 {
446  if (energyGuess < 0) return -9999.;
447  int energyBin = energyGuess;
448 
449  double longLikelihood = 0.;
450  int nbins = 0;
451 
452  for (int i = 0; i < LBINS; ++i) {
453  double qval = longProfile->GetBinContent(i + 1);
454  int qbin = longTemplate->GetZaxis()->FindBin(qval);
455  int binentries = longTemplate->GetBinContent(i + 1, energyBin, qbin);
456  int totentries = longTemplate->Integral(i + 1, i + 1, energyBin, energyBin, 0, 100);
457  if (qval > 0) {
458  ++nbins;
459  double prob = (double)binentries / totentries * 100;
460  if (binentries > 0) longLikelihood += log(prob);
461  }
462  } // loop through
463 
464  longLikelihood /= nbins;
465 
466  std::cout << longLikelihood << std::endl;
467  return longLikelihood;
468 } // getLongLikelihood
469 
470 // -------------------------------------------------
471 
473 {
474  if (energyGuess < 0) return -9999.;
475  int energyBin = energyGuess;
476 
477  double tranLikelihood_1 = 0;
478  double tranLikelihood_2 = 0;
479  double tranLikelihood_3 = 0;
480  double tranLikelihood_4 = 0;
481  double tranLikelihood_5 = 0;
482 
483  double qval;
484  int qbin, binentries, totentries;
485 
486  int nbins = 0;
487 
488  for (int i = 0; i < TBINS; ++i) {
489  qval = tranProfile_1->GetBinContent(i + 1);
490  qbin = tranTemplate_1->GetZaxis()->FindBin(qval);
491  binentries = tranTemplate_1->GetBinContent(i + 1, energyBin, qbin);
492  totentries = tranTemplate_1->Integral(i + 1, i + 1, energyBin, energyBin, 0, 100);
493  if (qval > 0) {
494  ++nbins;
495  double prob = (double)binentries / totentries * 100;
496  if (binentries > 0) tranLikelihood_1 += log(prob);
497  }
498 
499  qval = tranProfile_2->GetBinContent(i + 1);
500  qbin = tranTemplate_2->GetZaxis()->FindBin(qval);
501  binentries = tranTemplate_2->GetBinContent(i + 1, energyBin, qbin);
502  totentries = tranTemplate_2->Integral(i + 1, i + 1, energyBin, energyBin, 0, 100);
503  if (qval > 0) {
504  ++nbins;
505  double prob = (double)binentries / totentries * 100;
506  if (binentries > 0) tranLikelihood_2 += log(prob);
507  }
508 
509  qval = tranProfile_3->GetBinContent(i + 1);
510  qbin = tranTemplate_3->GetZaxis()->FindBin(qval);
511  binentries = tranTemplate_3->GetBinContent(i + 1, energyBin, qbin);
512  totentries = tranTemplate_3->Integral(i + 1, i + 1, energyBin, energyBin, 0, 100);
513  if (qval > 0) {
514  ++nbins;
515  double prob = (double)binentries / totentries * 100;
516  if (binentries > 0) tranLikelihood_3 += log(prob);
517  }
518 
519  qval = tranProfile_4->GetBinContent(i + 1);
520  qbin = tranTemplate_4->GetZaxis()->FindBin(qval);
521  binentries = tranTemplate_4->GetBinContent(i + 1, energyBin, qbin);
522  totentries = tranTemplate_4->Integral(i + 1, i + 1, energyBin, energyBin, 0, 100);
523  if (qval > 0) {
524  ++nbins;
525  double prob = (double)binentries / totentries * 100;
526  if (binentries > 0) tranLikelihood_4 += log(prob);
527  }
528 
529  qval = tranProfile_5->GetBinContent(i + 1);
530  qbin = tranTemplate_5->GetZaxis()->FindBin(qval);
531  binentries = tranTemplate_5->GetBinContent(i + 1, energyBin, qbin);
532  totentries = tranTemplate_5->Integral(i + 1, i + 1, energyBin, energyBin, 0, 100);
533  if (qval > 0) {
534  ++nbins;
535  double prob = (double)binentries / totentries * 100;
536  if (binentries > 0) tranLikelihood_5 += log(prob);
537  }
538 
539  } // loop through
540 
541  double const tranLikelihood =
542  (tranLikelihood_1 + tranLikelihood_2 + tranLikelihood_3 + tranLikelihood_4 + tranLikelihood_5) /
543  nbins;
544  std::cout << tranLikelihood << std::endl;
545  return tranLikelihood;
546 } // getTranLikelihood
547 
548 // -------------------------------------------------
549 
double E(const int i=0) const
Definition: MCParticle.h:234
Length_t WireCoordinate(Point_t const &pos, PlaneID const &planeid) const
Returns the index of the nearest wire to the specified position.
int PdgCode() const
Definition: MCParticle.h:213
int CCNC() const
Definition: MCNeutrino.h:148
const simb::MCNeutrino & GetNeutrino() const
Definition: MCTruth.h:77
void getShowerProfile(detinfo::DetectorClocksData const &clockdata, detinfo::DetectorPropertiesData const &detProp, std::vector< art::Ptr< recob::Hit >> showerhits, TVector3 shwvtx, TVector3 shwdir)
Declaration of signal hit object.
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
double Temperature() const
In kelvin.
::geo::Point_t toPoint(Point const &p)
Convert the specified point into a geo::Point_t.
constexpr auto abs(T v)
Returns the absolute value of the argument.
::geo::Vector_t toVector(Vector const &v)
Convert the specified vector into a geo::Vector_t.
STL namespace.
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
Particle class.
TCShowerElectronLikelihood(fhicl::ParameterSet const &pset)
double Efield(unsigned int planegap=0) const
kV/cm
IDparameter< geo::PlaneID > PlaneID
Member type of validated geo::PlaneID parameter.
Float_t Z
Definition: plot.C:37
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
Utilities to manipulate geometry vectors.The utilities include generic vector interface facilities al...
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
T get(std::string const &key) const
Definition: ParameterSet.h:314
double ConvertXToTicks(double X, int p, int t, int c) const
double DriftVelocity(double efield=0., double temperature=0.) const
cm/us
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
Contains all timing reference information for the detector.
TFile * file
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
void analyze(const art::Event &evt) override
double LifetimeCorrection(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, double time, double T0=0) const
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
Length_t WirePitch(PlaneID const &planeid=plane_zero) const
Returns the distance between two consecutive wires.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33