LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
SpacePointAna_module.cc
Go to the documentation of this file.
1 //
2 // Name: SpacePointAna_module.cc
3 //
4 // Purpose: Module SpacePointAna.
5 //
6 //
7 // Configuration parameters.
8 //
9 // HitModuleLabel: Hit module label.
10 // UseClusterHits: If true, use clustered hits (otherwise all hits).
11 // ClusterModuleLabel: Cluster module label.
12 // UseMC: Use MC truth information.
13 // SpacePointAlgTime: SpacePointAlg configuration (loose time cuts).
14 // SpacePointAlgSep: SpacePointAlg configuration (loose separation cuts).
15 // SpacePointAlgDefault: SpacePointAlg configuration (default cuts).
16 //
17 // Created: 2-Aug-2011 H. Greenlee
18 //
19 
20 #include <iostream>
21 #include <map>
22 #include <vector>
23 #include <algorithm>
24 
32 
44 
45 #include "TH1F.h"
46 
47 
48 namespace trkf {
49 
51  {
52  public:
53 
54  // Constructors, destructor
55 
56  explicit SpacePointAna(fhicl::ParameterSet const& pset);
57  virtual ~SpacePointAna();
58 
59  // Book histograms.
60 
61  void bookHistograms(bool mc);
62 
63  // Overrides.
64 
65  void analyze(const art::Event& evt);
66 
67  private:
68 
69  // Fcl Attributes.
70 
71  const SpacePointAlg fSptalgTime; // Algorithm object (increased time cut).
72  const SpacePointAlg fSptalgSep; // Algorithm object (increased sepataion cut).
73  const SpacePointAlg fSptalgDefault; // Algorithm object (default cuts).
74  std::string fHitModuleLabel;
76  std::string fClusterModuleLabel;
77  bool fUseMC;
78 
79  // Histograms.
80 
81  bool fBooked; // Have histograms been booked yet?
82  TH1F* fHDTUE; // U-drift electrons time difference.
83  TH1F* fHDTVE; // V-drift electrons time difference.
84  TH1F* fHDTWE; // W-drift electrons time difference.
85  TH1F* fHDTUPull; // U-drift electrons time pull.
86  TH1F* fHDTVPull; // V-drift electrons time pull.
87  TH1F* fHDTWPull; // W-drift electrons time pull.
88  TH1F* fHDTUV; // U-V time difference.
89  TH1F* fHDTVW; // V-W time difference.
90  TH1F* fHDTWU; // W-U time difference.
91  TH1F* fHS; // Spatial separation.
92  TH1F* fHchisq; // Space point chisquare.
93  TH1F* fHx; // X position.
94  TH1F* fHy; // Y position.
95  TH1F* fHz; // Z position.
96  TH1F* fHAmpU; // U hit amplitude.
97  TH1F* fHAmpV; // V hit amplitude.
98  TH1F* fHAmpW; // W hit amplitude.
99  TH1F* fHAreaU; // U hit area.
100  TH1F* fHAreaV; // V hit area.
101  TH1F* fHAreaW; // W hit area.
102  TH1F* fHSumU; // U hit sum ADC.
103  TH1F* fHSumV; // V hit sum ADC.
104  TH1F* fHSumW; // W hit sum ADC.
105  TH1F* fHMCdx; // X residual (reco vs. mc truth).
106  TH1F* fHMCdy; // Y residual (reco vs. mc truth).
107  TH1F* fHMCdz; // Z residual (reco vs. mc truth).
108  TH1F* fHMCxpull; // X pull (reco vs. mc truth).
109  TH1F* fHMCypull; // Y pull (reco vs. mc truth).
110  TH1F* fHMCzpull; // Z pull (reco vs. mc truth).
111 
112  // Statistics.
113 
115  };
116 
118 
119  SpacePointAna::SpacePointAna(const fhicl::ParameterSet& pset)
120  //
121  // Purpose: Constructor.
122  //
123  // Arguments: pset - Module parameters.
124  //
125  : EDAnalyzer(pset)
126  , fSptalgTime(pset.get<fhicl::ParameterSet>("SpacePointAlgTime"))
127  , fSptalgSep(pset.get<fhicl::ParameterSet>("SpacePointAlgSep"))
128  , fSptalgDefault(pset.get<fhicl::ParameterSet>("SpacePointAlgDefault"))
129  , fHitModuleLabel(pset.get<std::string>("HitModuleLabel"))
130  , fUseClusterHits(pset.get<bool>("UseClusterHits"))
131  , fClusterModuleLabel(pset.get<std::string>("ClusterModuleLabel"))
132  , fUseMC(pset.get<bool>("UseMC"))
133  , fBooked(false)
134  , fHDTUE(0)
135  , fHDTVE(0)
136  , fHDTWE(0)
137  , fHDTUPull(0)
138  , fHDTVPull(0)
139  , fHDTWPull(0)
140  , fHDTUV(0)
141  , fHDTVW(0)
142  , fHDTWU(0)
143  , fHS(0)
144  , fHchisq(0)
145  , fHx(0)
146  , fHy(0)
147  , fHz(0)
148  , fHAmpU(0)
149  , fHAmpV(0)
150  , fHAmpW(0)
151  , fHAreaU(0)
152  , fHAreaV(0)
153  , fHAreaW(0)
154  , fHSumU(0)
155  , fHSumV(0)
156  , fHSumW(0)
157  , fHMCdx(0)
158  , fHMCdy(0)
159  , fHMCdz(0)
160  , fHMCxpull(0)
161  , fHMCypull(0)
162  , fHMCzpull(0)
163  , fNumEvent(0)
164  {
165 
166  // Report.
167 
168  mf::LogInfo("SpacePointAna")
169  << "SpacePointAna configured with the following parameters:\n"
170  << " HitModuleLabel = " << fHitModuleLabel << "\n"
171  << " UseClusterHits = " << fUseClusterHits << "\n"
172  << " ClusterModuleLabel = " << fClusterModuleLabel << "\n"
173  << " UseMC = " << fUseMC;
174  }
175 
177  //
178  // Purpose: Destructor.
179  //
180  {}
181 
183  //
184  // Purpose: Book histograms.
185  //
186  {
187  if(!fBooked) {
188  fBooked = true;
189 
192  art::TFileDirectory dir = tfs->mkdir("sptana", "SpacePointAna histograms");
193 
194  if(mc && fUseMC) {
195  fHDTUE = dir.make<TH1F>("MCDTUE", "U-Drift Electrons Time Difference", 100, -5., 5.);
196  fHDTVE = dir.make<TH1F>("MCDTVE", "V-Drift Electrons Time Difference", 100, -5., 5.);
197  fHDTWE = dir.make<TH1F>("MCDTWE", "W-Drift Electrons Time Difference", 100, -5., 5.);
198  fHDTUPull = dir.make<TH1F>("MCDTUPull", "U-Drift Electrons Time Pull", 100, -50., 50.);
199  fHDTVPull = dir.make<TH1F>("MCDTVPull", "V-Drift Electrons Time Pull", 100, -50., 50.);
200  fHDTWPull = dir.make<TH1F>("MCDTWPull", "W-Drift Electrons Time Pull", 100, -50., 50.);
201  }
202  if(!fSptalgTime.merge()) {
203  fHDTUV = dir.make<TH1F>("DTUV", "U-V time difference", 100, -20., 20.);
204  fHDTVW = dir.make<TH1F>("DTVW", "V-W time difference", 100, -20., 20.);
205  fHDTWU = dir.make<TH1F>("DTWU", "W-U time difference", 100, -20., 20.);
206  fHS = dir.make<TH1F>("DS", "Spatial Separation", 100, -2., 2.);
207  }
208  fHchisq = dir.make<TH1F>("chisq", "Chisquare", 100, 0., 20.);
209 
210  fHx = dir.make<TH1F>("xpos", "X Position",
211  100, 0., 2.*geom->DetHalfWidth());
212  fHy = dir.make<TH1F>("ypos", "Y Position",
213  100, -geom->DetHalfHeight(), geom->DetHalfHeight());
214  fHz = dir.make<TH1F>("zpos", "Z Position",
215  100, 0., geom->DetLength());
216  fHAmpU = dir.make<TH1F>("ampU", "U Hit Amplitude", 50, 0., 50.);
217  fHAmpV = dir.make<TH1F>("ampV", "V Hit Amplitude", 50, 0., 50.);
218  fHAmpW = dir.make<TH1F>("ampW", "W Hit Amplitude", 50, 0., 50.);
219  fHAreaU = dir.make<TH1F>("areaU", "U Hit Area", 100, 0., 500.);
220  fHAreaV = dir.make<TH1F>("areaV", "V Hit Area", 100, 0., 500.);
221  fHAreaW = dir.make<TH1F>("areaW", "W Hit Area", 100, 0., 500.);
222  fHSumU = dir.make<TH1F>("sumU", "U Hit Sum ADC", 100, 0., 500.);
223  fHSumV = dir.make<TH1F>("sumV", "V Hit Sum ADC", 100, 0., 500.);
224  fHSumW = dir.make<TH1F>("sumW", "W Hit Sum ADC", 100, 0., 500.);
225  if(mc && fUseMC) {
226  fHMCdx = dir.make<TH1F>("MCdx", "X MC Residual", 100, -1., 1.);
227  fHMCdy = dir.make<TH1F>("MCdy", "Y MC Residual", 100, -1., 1.);
228  fHMCdz = dir.make<TH1F>("MCdz", "Z MC Residual", 100, -1., 1.);
229  fHMCxpull = dir.make<TH1F>("MCxpull", "X MC Pull", 100, -50., 50.);
230  fHMCypull = dir.make<TH1F>("MCypull", "Y MC Pull", 100, -50., 50.);
231  fHMCzpull = dir.make<TH1F>("MCzpull", "Z MC Pull", 100, -50., 50.);
232  }
233  }
234  }
235 
237  //
238  // Purpose: Analyze method.
239  //
240  // Arguments: event - Art event.
241  //
242  {
243  ++fNumEvent;
244 
245  // Make sure histograms are booked.
246 
247  bool mc = !evt.isRealData();
248  bookHistograms(mc);
249 
250  // Get Services.
251 
253  const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
254 
255  // Get Hits.
256 
258 
259  if(fUseClusterHits) {
260 
261  // Get clusters.
262 
264  evt.getByLabel(fClusterModuleLabel, clusterh);
265 
266  // Get hits from all clusters.
268 
269  if(clusterh.isValid()) {
270  int nclus = clusterh->size();
271 
272  for(int i = 0; i < nclus; ++i) {
273  art::Ptr<recob::Cluster> pclus(clusterh, i);
274  std::vector< art::Ptr<recob::Hit> > clushits = fm.at(i);
275  int nhits = clushits.size();
276  hits.reserve(hits.size() + nhits);
277 
278  for(std::vector< art::Ptr<recob::Hit> >::const_iterator ihit = clushits.begin(); ihit != clushits.end(); ++ihit) {
279  hits.push_back(*ihit);
280  }
281  }
282  }
283  }
284  else {
285 
286  // Get unclustered hits.
287 
289  evt.getByLabel(fHitModuleLabel, hith);
290  if(hith.isValid()) {
291  int nhits = hith->size();
292  hits.reserve(nhits);
293 
294  for(int i = 0; i < nhits; ++i)
295  hits.push_back(art::Ptr<recob::Hit>(hith, i));
296  }
297  }
298 
299  // Fill histograms that don't depend on space points.
300 
301  if(mc && fUseMC) {
302 
304 
305  // Loop over hits and fill hit-electron time difference histogram.
306 
307  for(art::PtrVector<recob::Hit>::const_iterator ihit = hits.begin();
308  ihit != hits.end(); ++ihit) {
309  const recob::Hit& hit = **ihit;
310 
311  //unsigned int channel = hit.Channel();
312  //geo::View_t geo_view = geom->View(channel);
313  //geo::View_t hit_view = hit.View();
314  //assert(geo_view == hit_view);
315  double tpeak = hit.PeakTime();
316  double terr = hit.SigmaPeakTime();
317 
318  //assert(channel == hit.Channel());
319 
320  // Loop over electrons associated with this hit/channel and fill
321  // hit-electron time difference histograms.
322 
323  // loop over the map of TDC to sim::IDE to get the TDC for each energy dep
324  // Find the average time in ticks for this hit.
325 
326  //double sumw = 0.;
327  //double sumt = 0.;
328 
329  //const std::map<unsigned short, std::vector<sim::IDE> > &idemap = simchan.TDCIDEMap();
330  //std::map<unsigned short, std::vector<sim::IDE> >::const_iterator mitr = idemap.begin();
331  //for(mitr = idemap.begin(); mitr != idemap.end(); mitr++) {
332  // double tdc = double((*mitr).first);
333  // if(tdc >= tstart && tdc <= tend) {
334  // const std::vector<sim::IDE>& idevec = (*mitr).second;
335  // for(std::vector<sim::IDE>::const_iterator iide=idevec.begin();
336  // iide != idevec.end(); ++iide) {
337  // const sim::IDE& ide = *iide;
338  // double w = ide.numElectrons;
339  // sumw += w;
340  // sumt += w*tdc;
341  // }
342  // }
343  //}
344  //double tav = 0.;
345  //if(sumw != 0.)
346  // tav = sumt / sumw;
347 
348  bool tav_ok = true;
349  double tav = 0.;
350  try {
351  std::vector<double> hitxyz = bt_serv->HitToXYZ(*ihit);
352  tav = detprop->ConvertXToTicks(hitxyz[0], (*ihit)->WireID().Plane, (*ihit)->WireID().TPC, (*ihit)->WireID().Cryostat);
353  }
354  catch(cet::exception& x) {
355  tav_ok = false;
356  }
357  if(tav_ok) {
358  if((*ihit)->View() == geo::kU) {
359  fHDTUE->Fill(tpeak - tav);
360  fHDTUPull->Fill((tpeak - tav) / terr);
361  }
362  else if((*ihit)->View() == geo::kV) {
363  fHDTVE->Fill(tpeak - tav);
364  fHDTVPull->Fill((tpeak - tav) / terr);
365  }
366  else if((*ihit)->View() == geo::kZ) {
367  fHDTWE->Fill(tpeak - tav);
368  fHDTWPull->Fill((tpeak - tav) / terr);
369  }
370  else
371  throw cet::exception("SpacePointAna") << "Bad view = " << (*ihit)->View() << "\n";
372  }
373  }
374  }
375 
376  std::vector<recob::SpacePoint> spts1; // For time histograms.
377  std::vector<recob::SpacePoint> spts2; // For separation histogram.
378  std::vector<recob::SpacePoint> spts3; // Default cuts.
379 
380  // If nonzero time cut is specified, make space points using that
381  // time cut (for time histograms).
382 
383  if(!fSptalgTime.merge()) {
384  if(mc && fUseMC)
385  fSptalgTime.makeMCTruthSpacePoints(hits, spts1);
386  else
387  fSptalgTime.makeSpacePoints(hits, spts1);
388 
389  // Report number of space points.
390 
391  LOG_DEBUG("SpacePointAna") << "Found " << spts1.size()
392  << " space points using special time cut.";
393  }
394 
395  // If nonzero separation cut is specified, make space points using that
396  // separation cut (for separation histogram).
397 
398  if(!fSptalgSep.merge()) {
399  if(mc && fUseMC)
400  fSptalgSep.makeMCTruthSpacePoints(hits, spts2);
401  else
402  fSptalgSep.makeSpacePoints(hits, spts2);
403 
404  // Report number of space points.
405 
406  LOG_DEBUG("SpacePointAna") << "Found " << spts2.size()
407  << " space points using special seperation cut.";
408  }
409 
410  // Make space points using default cuts.
411 
412  if(mc && fUseMC)
414  else
415  fSptalgDefault.makeSpacePoints(hits, spts3);
416 
417  // Report number of space points.
418 
419  LOG_DEBUG("SpacePointAna") << "Found " << spts3.size()
420  << " space points using default cuts.";
421 
422  if(!fSptalgTime.merge()) {
423 
424  // Loop over space points and fill time histograms.
425 
426  for(std::vector<recob::SpacePoint>::const_iterator i = spts1.begin();
427  i != spts1.end(); ++i) {
428  const recob::SpacePoint& spt = *i;
429 
430  // Get hits associated with this SpacePoint.
431 
433 
434  // Make a double loop over hits and fill hit time difference histograms.
435 
437  ihit != spthits.end(); ++ihit) {
438  const recob::Hit& hit1 = **ihit;
439 
440  geo::WireID hit1WireID = hit1.WireID();
441  unsigned int tpc1, plane1;
442  tpc1 = hit1WireID.TPC;
443  plane1 = hit1WireID.Plane;
444 
445  geo::View_t view1 = hit1.View();
446  double t1 = fSptalgTime.correctedTime(hit1);
447 
449  jhit != spthits.end(); ++jhit) {
450  const recob::Hit& hit2 = **jhit;
451 
452  geo::WireID hit2WireID = hit2.WireID();
453  unsigned int tpc2, plane2;
454  tpc2 = hit2WireID.TPC;
455  plane2 = hit2WireID.Plane;
456 
457  // Require same tpc, different view.
458 
459  if(tpc1 == tpc2 && plane1 != plane2) {
460 
461  geo::View_t view2 = hit2.View();
462  double t2 = fSptalgTime.correctedTime(hit2);
463 
464  if(view1 == geo::kU) {
465  if(view2 == geo::kV)
466  fHDTUV->Fill(t1-t2);
467  if(view2 == geo::kZ)
468  fHDTWU->Fill(t2-t1);
469  }
470  if(view1 == geo::kV) {
471  if(view2 == geo::kZ)
472  fHDTVW->Fill(t1-t2);
473  if(view2 == geo::kU)
474  fHDTUV->Fill(t2-t1);
475  }
476  if(view1 == geo::kZ) {
477  if(view2 == geo::kU)
478  fHDTWU->Fill(t1-t2);
479  if(view2 == geo::kV)
480  fHDTVW->Fill(t2-t1);
481  }
482  }
483  }
484  }
485  }
486 
487  // Loop over space points and fill seperation histograms.
488 
489  for(std::vector<recob::SpacePoint>::const_iterator i = spts2.begin();
490  i != spts2.end(); ++i) {
491  const recob::SpacePoint& spt = *i;
492 
493  // Get hits associated with this SpacePoint.
494 
496 
497  // Fill separation histogram.
498 
499  double sep = fSptalgSep.separation(spthits);
500  fHS->Fill(sep);
501  }
502  }
503 
504  // Loop over default space points and fill histograms.
505 
507  i != spts3.end(); ++i) {
508  const recob::SpacePoint& spt = *i;
509 
510  fHchisq->Fill(spt.Chisq());
511  fHx->Fill(spt.XYZ()[0]);
512  fHy->Fill(spt.XYZ()[1]);
513  fHz->Fill(spt.XYZ()[2]);
514 
515  // Get hits associated with this SpacePoint.
516 
517  std::vector<art::Ptr<recob::Hit> > spthits;
519  for( auto const& ptr : av_spthits ){ spthits.push_back(ptr);}
520 
521  // Fill single hit histograms.
522 
523  for(art::PtrVector<recob::Hit>::const_iterator ihit = spthits.begin();
524  ihit != spthits.end(); ++ihit) {
525  const recob::Hit& hit = **ihit;
526 
527  geo::View_t view = hit.View();
528 
529  if(view == geo::kU) {
530  fHAmpU->Fill(hit.PeakAmplitude());
531  fHAreaU->Fill(hit.Integral());
532  fHSumU->Fill(hit.SummedADC());
533  }
534  if(view == geo::kV) {
535  fHAmpV->Fill(hit.PeakAmplitude());
536  fHAreaV->Fill(hit.Integral());
537  fHSumV->Fill(hit.SummedADC());
538  }
539  if(view == geo::kZ) {
540  fHAmpW->Fill(hit.PeakAmplitude());
541  fHAreaW->Fill(hit.Integral());
542  fHSumW->Fill(hit.SummedADC());
543  }
544  }
545 
546  if(mc && fUseMC) {
547 
549 
550  try {
551  std::vector<double> mcxyz = bt_serv->SpacePointHitsToWeightedXYZ(spthits);
552  fHMCdx->Fill(spt.XYZ()[0] - mcxyz[0]);
553  fHMCdy->Fill(spt.XYZ()[1] - mcxyz[1]);
554  fHMCdz->Fill(spt.XYZ()[2] - mcxyz[2]);
555  if(spt.ErrXYZ()[0] > 0.)
556  fHMCxpull->Fill((spt.XYZ()[0] - mcxyz[0]) / std::sqrt(spt.ErrXYZ()[0]));
557  if(spt.ErrXYZ()[2] > 0.)
558  fHMCypull->Fill((spt.XYZ()[1] - mcxyz[1]) / std::sqrt(spt.ErrXYZ()[2]));
559  if(spt.ErrXYZ()[5] > 0.)
560  fHMCzpull->Fill((spt.XYZ()[2] - mcxyz[2]) / std::sqrt(spt.ErrXYZ()[5]));
561  }
562  catch(cet::exception& x) {}
563  }
564  }
565  }
566 }
Float_t x
Definition: compare.C:6
double correctedTime(const recob::Hit &hit) const
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
TTree * t1
Definition: plottest35.C:26
const std::vector< double > HitToXYZ(const recob::Hit &hit)
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Planes which measure V.
Definition: geo_types.h:77
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
iterator begin()
Definition: PtrVector.h:223
Declaration of signal hit object.
STL namespace.
Planes which measure Z direction.
Definition: geo_types.h:79
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:225
geo::View_t View() const
View for the plane of the hit.
Definition: Hit.h:233
bool isRealData() const
Definition: Event.h:83
const art::PtrVector< recob::Hit > & getAssociatedHits(const recob::SpacePoint &spt) const
TFileDirectory mkdir(std::string const &dir, std::string const &descr="")
const SpacePointAlg fSptalgSep
void makeMCTruthSpacePoints(const art::PtrVector< recob::Hit > &hits, std::vector< recob::SpacePoint > &spts) const
virtual double ConvertXToTicks(double X, int p, int t, int c) const =0
float PeakAmplitude() const
The estimated amplitude of the hit at its peak, in ADC units.
Definition: Hit.h:222
double separation(const art::PtrVector< recob::Hit > &hits) const
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
Planes which measure U.
Definition: geo_types.h:76
void hits()
Definition: readHits.C:15
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
parameter set interface
iterator end()
Definition: PtrVector.h:237
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
EDAnalyzer(Table< Config > const &config)
Definition: EDAnalyzer.h:100
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
Declaration of cluster object.
const std::vector< double > SpacePointHitsToWeightedXYZ(std::vector< art::Ptr< recob::Hit >> const &hits)
Detector simulation of raw signals on wires.
TTree * t2
Definition: plottest35.C:36
const double * ErrXYZ() const
Definition: SpacePoint.h:65
data_t::const_iterator const_iterator
Definition: PtrVector.h:61
void makeSpacePoints(const art::PtrVector< recob::Hit > &hits, std::vector< recob::SpacePoint > &spts) const
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:219
double Chisq() const
Definition: SpacePoint.h:66
T * make(ARGS...args) const
SpacePointAna(fhicl::ParameterSet const &pset)
const double * XYZ() const
Definition: SpacePoint.h:64
Utility object to perform functions of association.
Encapsulate the construction of a single detector plane.
const SpacePointAlg fSptalgTime
TDirectory * dir
Definition: macro.C:5
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
void analyze(const art::Event &evt)
#define LOG_DEBUG(id)
float SummedADC() const
The sum of calibrated ADC counts of the hit (0. by default)
Definition: Hit.h:224
const SpacePointAlg fSptalgDefault
float SigmaPeakTime() const
Uncertainty for the signal peak, in tick units.
Definition: Hit.h:220
bool merge() const
Definition: SpacePointAlg.h:89
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:49
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
Algorithm for generating space points from hits.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Encapsulate the construction of a single detector plane.