LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
GnocchiCalorimetry_module.cc
Go to the documentation of this file.
1 // A gnew Calorimetry module
2 //
3 // Re-write of the Calorimetry module ported into art by E. Chruch
4 // from the original (ArgoNeuT) calorimetry module by Ornella Palamara and Maddalena Antonello.
5 //
6 // Gray Putnam
7 // grayputnam@uchicago.edu
8 
9 #include <cmath>
10 #include <limits> // std::numeric_limits<>
11 #include <numeric> // std::accumulate
12 #include <optional>
13 #include <string>
14 
19 
20 #include "larcorealg/CoreUtils/NumericUtils.h" // util::absDiff()
23 #include "lardata/ArtDataHelper/TrackUtils.h" // lar::util::TrackPitchInView()
30 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
31 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
32 
35 
36 // ROOT includes
37 #include <TF1.h>
38 #include <TGraph.h>
39 #include <TMath.h>
40 #include <vector>
41 
50 #include "cetlib/pow.h" // cet::sum_of_squares()
51 #include "fhiclcpp/ParameterSet.h"
54 
55 namespace {
56  constexpr unsigned int int_max_as_unsigned_int{std::numeric_limits<int>::max()};
57  typedef std::pair<unsigned, std::vector<unsigned>> OrganizedHits;
58 }
59 
60 namespace calo {
61 
63  public:
64  struct Config {
66  using Name = fhicl::Name;
67 
68  enum ChargeMethod : unsigned {
73  };
74 
76  Comment("Module label for track producer.")};
77 
79  Comment("Module label for T0 time producer."),
80  ""};
81 
83  Name("AssocHitModuleLabel"),
84  Comment("Module label for association between tracks and hits. If not set, defaults to "
85  "TrackModuleLabel."),
86  ""};
87 
89  Name("ChargeMethod"),
90  Comment("Method used to extract charge from a hit. Options: 0==Amplitude(), 1==Integral(), "
91  "2==SummedADC(), 3==SummedIntegral(). See the ChargeMethod enum.")};
92 
94  Name("FieldDistortion"),
95  Comment("True if field distortion (i.e. from space charge) is included in the input.")};
96 
98  Name("FieldDistortionEfield"),
99  Comment("True if field distortion (i.e. from space charge) is included in the input.")};
100 
102  Name("TrackIsFieldDistortionCorrected"),
103  Comment("Whether the space-points on the input tracks have their points corrected for the "
104  "field distortions. "
105  "I.e. whether the track trajectory points represent charge as seen by wires or the "
106  "3D particle trajectory.")};
107 
109  Comment("Which cryostat number the input tracks occupy.")};
110 
112  Name("FieldDistortionCorrectionXSign"),
113  Comment("Sign of the field distortion correction to be applied in the X direction. "
114  "Positive by default."),
115  1.};
116 
118  Name("CaloAlg"),
119  Comment("Configuration for the calo::CalorimetryAlg")};
120 
122  Name("NormTools"),
123  Comment("List of INormalizeCharge tool configurations to use.")};
124  };
125 
127 
128  explicit GnocchiCalorimetry(Parameters const& config);
129 
130  void produce(art::Event& evt) override;
131 
132  private:
135  std::vector<std::unique_ptr<INormalizeCharge>> fNormTools;
136 
137  // helper functions
138  std::vector<std::vector<OrganizedHits>> OrganizeHits(
140  const std::vector<const recob::TrackHitMeta*>& thms,
141  const recob::Track& track,
142  unsigned nplanes);
143  std::vector<std::vector<OrganizedHits>> OrganizeHitsIndividual(
145  const std::vector<const recob::TrackHitMeta*>& thms,
146  const recob::Track& track,
147  unsigned nplanes);
148  std::vector<std::vector<OrganizedHits>> OrganizeHitsSnippets(
150  const std::vector<const recob::TrackHitMeta*>& thms,
151  const recob::Track& track,
152  unsigned nplanes);
153  bool HitIsValid(const recob::TrackHitMeta* thm, const recob::Track& track);
156  const recob::TrackHitMeta* meta);
159  const recob::TrackHitMeta* meta);
162  double GetPitch(const recob::Track& track,
164  const recob::TrackHitMeta* meta);
165  double GetCharge(const art::Ptr<recob::Hit> hit,
166  const std::vector<recob::Hit const*>& sharedHits);
167  double GetEfield(const detinfo::DetectorPropertiesData& dprop,
168  const recob::Track& track,
170  const recob::TrackHitMeta* meta);
171  double Normalize(double dQdx,
172  const art::Event& e,
173  const recob::Hit& h,
174  const geo::Point_t& location,
175  const geo::Vector_t& direction,
176  double t0);
177  };
178 
179 } // end namespace calo
180 
183 {
184  produces<std::vector<anab::Calorimetry>>();
185  produces<art::Assns<recob::Track, anab::Calorimetry>>();
186  std::vector<fhicl::ParameterSet> norm_tool_configs =
187  fConfig.NormTools.get<std::vector<fhicl::ParameterSet>>();
188  for (const fhicl::ParameterSet& p : norm_tool_configs) {
189  fNormTools.push_back(art::make_tool<INormalizeCharge>(p));
190  }
191 }
192 
194 {
195  // Get services
197  auto const clock_data = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
198  auto const det_prop =
200 
201  size_t nplanes = geom->Nplanes();
202 
203  // Define output collections
204  std::unique_ptr<std::vector<anab::Calorimetry>> outputCalo(new std::vector<anab::Calorimetry>);
205  std::unique_ptr<art::Assns<recob::Track, anab::Calorimetry>> outputCaloAssn(
207 
208  // collect input
209  art::Handle<std::vector<recob::Track>> trackListHandle;
210  std::vector<art::Ptr<recob::Track>> tracklist;
211  if (evt.getByLabel(fConfig.TrackModuleLabel(), trackListHandle)) {
212  art::fill_ptr_vector(tracklist, trackListHandle);
213  }
214 
215  // get the label to collect this hits
216  const std::string& hitLabel = (fConfig.AssocHitModuleLabel() == "") ?
219  art::FindManyP<recob::Hit, recob::TrackHitMeta> fmHits(trackListHandle, evt, hitLabel);
220 
221  // must be valid if the T0 module label is non-empty
222  art::FindManyP<anab::T0> fmT0s(trackListHandle, evt, fConfig.T0ModuleLabel());
223 
224  for (auto const& nt : fNormTools)
225  nt->setup(evt);
226 
227  // iterate over all the tracks
228  for (unsigned trk_i = 0; trk_i < tracklist.size(); trk_i++) {
229  const recob::Track& track = *tracklist[trk_i];
230 
231  // collect input for this track
232  const std::vector<art::Ptr<recob::Hit>>& hits = fmHits.at(trk_i);
233  const std::vector<const recob::TrackHitMeta*>& thms = fmHits.data(trk_i);
234 
235  double T0 = 0;
236  if (fConfig.T0ModuleLabel().size()) {
237  const std::vector<art::Ptr<anab::T0>>& this_t0s = fmT0s.at(trk_i);
238  if (this_t0s.size()) T0 = this_t0s.at(0)->Time();
239  }
240 
241  // organize the hits by plane
242  std::vector<std::vector<OrganizedHits>> hit_indices = OrganizeHits(hits, thms, track, nplanes);
243 
244  for (unsigned plane_i = 0; plane_i < nplanes; plane_i++) {
245 
246  float kinetic_energy = 0.;
247  std::vector<float> dEdxs;
248  std::vector<float> dQdxs;
249  std::vector<float> resranges;
250  std::vector<float> deadwireresranges;
251  float range = 0.;
252  std::vector<float> pitches;
253  std::vector<geo::Point_t> xyzs;
254  std::vector<size_t> tp_indices;
255  geo::PlaneID plane;
256 
257  // setup the plane ID
258  plane.Plane = plane_i;
259  plane.TPC = 0; // arbitrary -- tracks can cross TPC boundaries
260  plane.Cryostat = fConfig.Cryostat();
261  plane.isValid = true;
262 
263  std::vector<float> lengths;
264  for (unsigned hit_i = 0; hit_i < hit_indices[plane_i].size(); hit_i++) {
265  unsigned hit_index = hit_indices[plane_i][hit_i].first;
266 
267  std::vector<recob::Hit const*> sharedHits = {};
268  sharedHits.reserve(hit_indices[plane_i][hit_i].second.size());
269  for (const unsigned shared_hit_index : hit_indices[plane_i][hit_i].second)
270  sharedHits.push_back(hits[shared_hit_index].get());
271 
272  // Get the location of this point
273  geo::Point_t location = GetLocation(track, hits[hit_index], thms[hit_index]);
274 
275  // Get the pitch
276  double pitch = GetPitch(track, hits[hit_index], thms[hit_index]);
277 
278  // And the charge
279  double charge = GetCharge(hits[hit_index], sharedHits);
280 
281  // Get the EField
282  double EField = GetEfield(det_prop, track, hits[hit_index], thms[hit_index]);
283 
284  // Angle to the drift electric field (in x direction), in units of degrees
285  geo::Vector_t direction = track.DirectionAtPoint(thms[hit_index]->Index());
286  double phi = acos(abs(direction.x())) * 180 / M_PI;
287 
288  double dQdx = charge / pitch;
289 
290  // Normalize out the detector response
291  dQdx = Normalize(dQdx,
292  evt,
293  *hits[hit_index],
294  track.LocationAtPoint(thms[hit_index]->Index()),
295  track.DirectionAtPoint(thms[hit_index]->Index()),
296  T0);
297 
298  // turn into dEdx
300  fCaloAlg.dEdx_AMP(clock_data,
301  det_prop,
302  dQdx,
303  hits[hit_index]->PeakTime(),
304  hits[hit_index]->WireID().Plane,
305  T0,
306  EField,
307  phi) :
308  fCaloAlg.dEdx_AREA(clock_data,
309  det_prop,
310  dQdx,
311  hits[hit_index]->PeakTime(),
312  hits[hit_index]->WireID().Plane,
313  T0,
314  EField,
315  phi);
316 
317  // save the length between each pair of hits
318  if (xyzs.size() == 0) { lengths.push_back(0.); }
319  else {
320  lengths.push_back((location - xyzs.back()).r());
321  }
322 
323  // save stuff
324  dEdxs.push_back(dEdx);
325  dQdxs.push_back(dQdx);
326  pitches.push_back(pitch);
327  xyzs.push_back(location);
328  kinetic_energy += dEdx * pitch;
329 
330  // TODO: FIXME
331  // It seems weird that the "trajectory-point-index" actually is the
332  // index of the hit... is this a bug in the documentation
333  // of anab::Calorimetry?
334  //
335  // i.e. -- I think this piece of code should actually be:
336  // tp_indices.push_back(thms[hit_index]->Index());
337  tp_indices.push_back(hits[hit_index].key());
338 
339  } // end iterate over hits
340 
341  // turn the lengths vector into a residual-range vector and total length
342  if (lengths.size() > 1) {
343  range = std::accumulate(lengths.begin(), lengths.end(), 0.);
344 
345  // check the direction that the hits are going in the track:
346  // upstream (end-start) or downstream (start-end)
347  bool is_downstream =
348  (track.Trajectory().Start() - xyzs[0]).r() +
349  (track.Trajectory().End() - xyzs.back()).r() <
350  (track.Trajectory().End() - xyzs[0]).r() + (track.Trajectory().Start() - xyzs.back()).r();
351 
352  resranges.resize(lengths.size());
353  if (is_downstream) {
354  resranges[lengths.size() - 1] = lengths.back() / 2.;
355  for (int i_len = lengths.size() - 2; i_len >= 0; i_len--) {
356  resranges[i_len] = resranges[i_len + 1] + lengths[i_len + 1];
357  }
358  }
359  else {
360  resranges[0] = lengths[1] / 2.;
361  for (unsigned i_len = 1; i_len < lengths.size(); i_len++) {
362  resranges[i_len] = resranges[i_len - 1] + lengths[i_len];
363  }
364  }
365  }
366 
367  // save the Calorimetry output
368  //
369  // Bogus if less than two hits on this plane
370  if (lengths.size() > 1) {
371  outputCalo->push_back(anab::Calorimetry(kinetic_energy,
372  dEdxs,
373  dQdxs,
374  resranges,
375  deadwireresranges,
376  range,
377  pitches,
378  xyzs,
379  tp_indices,
380  plane));
381  }
382  else {
383  outputCalo->push_back(
384  anab::Calorimetry(util::kBogusD, {}, {}, {}, {}, util::kBogusD, {}, {}, {}, plane));
385  }
386 
387  util::CreateAssn(evt, *outputCalo, tracklist[trk_i], *outputCaloAssn);
388 
389  } // end iterate over planes
390 
391  } // end iterate over tracks
392 
393  evt.put(std::move(outputCalo));
394  evt.put(std::move(outputCaloAssn));
395 
396  return;
397 }
398 
399 std::vector<std::vector<OrganizedHits>> calo::GnocchiCalorimetry::OrganizeHits(
401  const std::vector<const recob::TrackHitMeta*>& thms,
402  const recob::Track& track,
403  unsigned nplanes)
404 {
405  // charge is computed per hit -- we organize hits indivudally
408  return OrganizeHitsIndividual(hits, thms, track, nplanes);
409  }
410  // charge is computed per snippet -- we organize hits by snippet
411  else {
412  return OrganizeHitsSnippets(hits, thms, track, nplanes);
413  }
414 }
415 
416 std::vector<std::vector<OrganizedHits>> calo::GnocchiCalorimetry::OrganizeHitsIndividual(
418  const std::vector<const recob::TrackHitMeta*>& thms,
419  const recob::Track& track,
420  unsigned nplanes)
421 {
422  std::vector<std::vector<OrganizedHits>> ret(nplanes);
423  for (unsigned i = 0; i < hits.size(); i++) {
424  if (HitIsValid(thms[i], track)) { ret[hits[i]->WireID().Plane].push_back({i, {}}); }
425  }
426 
427  return ret;
428 }
429 
430 std::vector<std::vector<OrganizedHits>> calo::GnocchiCalorimetry::OrganizeHitsSnippets(
432  const std::vector<const recob::TrackHitMeta*>& thms,
433  const recob::Track& track,
434  unsigned nplanes)
435 {
436  // In this case, we need to only accept one hit in each snippet
437  // Snippets are counted by the Start, End, and Wire. If all these are the same for a hit, then they are on the same snippet.
438  //
439  // If there are multiple valid hits on the same snippet, we need a way to pick the best one.
440  // (TODO: find a good way). The current method is to take the one with the highest charge integral.
441  struct HitIdentifier {
442  int startTick;
443  int endTick;
444  int wire;
445  float integral;
446 
447  // construct
448  explicit HitIdentifier(const recob::Hit& hit)
449  : startTick(hit.StartTick())
450  , endTick(hit.EndTick())
451  , wire(hit.WireID().Wire)
452  , integral(hit.Integral())
453  {}
454 
455  // Defines whether two hits are on the same snippet
456  inline bool operator==(const HitIdentifier& rhs) const
457  {
458  return startTick == rhs.startTick && endTick == rhs.endTick && wire == rhs.wire;
459  }
460 
461  // Defines which hit to pick between two both on the same snippet
462  inline bool operator>(const HitIdentifier& rhs) const { return integral > rhs.integral; }
463  };
464 
465  std::vector<std::vector<OrganizedHits>> ret(nplanes);
466  std::vector<std::vector<HitIdentifier>> hit_idents(nplanes);
467  for (unsigned i = 0; i < hits.size(); i++) {
468  if (HitIsValid(thms[i], track)) {
469  HitIdentifier this_ident(*hits[i]);
470 
471  // check if we have found a hit on this snippet before
472  bool found_snippet = false;
473  auto const plane = hits[i]->WireID().Plane;
474  for (unsigned j = 0; j < ret[plane].size(); j++) {
475  if (this_ident == hit_idents[plane][j]) {
476  found_snippet = true;
477  if (this_ident > hit_idents[plane][j]) {
478  ret[plane][j].second.push_back(ret[plane][j].first);
479  ret[plane][j].first = i;
480  hit_idents[plane][j] = this_ident;
481  }
482  else {
483  ret[plane][j].second.push_back(i);
484  }
485  break;
486  }
487  }
488  if (!found_snippet) {
489  ret[plane].push_back({i, {}});
490  hit_idents[plane].push_back(this_ident);
491  }
492  }
493  }
494  return ret;
495 }
496 
498 {
499  if (thm->Index() == int_max_as_unsigned_int) return false;
500  return track.HasValidPoint(thm->Index());
501 }
502 
505  const recob::TrackHitMeta* meta)
506 {
507  geo::Point_t loc = track.LocationAtPoint(meta->Index());
509  loc;
510 }
511 
513  const geo::TPCID& tpc)
514 {
515  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
516 
517  geo::Point_t ret = loc;
518 
519  if (sce->EnableCalSpatialSCE() && fConfig.FieldDistortion()) {
520  geo::Vector_t offset = sce->GetCalPosOffsets(ret, tpc.TPC);
521 
522  ret.SetX(ret.X() + fConfig.FieldDistortionCorrectionXSign() * offset.X());
523  ret.SetY(ret.Y() + offset.Y());
524  ret.SetZ(ret.Z() + offset.Z());
525  }
526 
527  return ret;
528 }
529 
532  const recob::TrackHitMeta* meta)
533 {
534  geo::Point_t loc = track.LocationAtPoint(meta->Index());
536  loc;
537 }
538 
540  const geo::TPCID& tpc)
541 {
542  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
544 
545  geo::Point_t ret = loc;
546 
547  if (sce->EnableCalSpatialSCE() && fConfig.FieldDistortion()) {
548  // Returned X is the drift -- multiply by the drift direction to undo this
549  int corr = geom->TPC(tpc).DriftDir().X();
550 
551  geo::Vector_t offset = sce->GetPosOffsets(ret);
552 
553  ret.SetX(ret.X() + corr * fConfig.FieldDistortionCorrectionXSign() * offset.X());
554  ret.SetY(ret.Y() + offset.Y());
555  ret.SetZ(ret.Z() + offset.Z());
556  }
557 
558  return ret;
559 }
560 
563  const recob::TrackHitMeta* meta)
564 {
566  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
567 
568  double angleToVert =
569  geom->WireAngleToVertical(hit->View(), hit->WireID().asPlaneID()) - 0.5 * ::util::pi<>();
570 
572 
573  // "dir" should be the direction that the wires see. If the track already has the field
574  // distortion corrections applied, then we need to de-apply them to get the direction as
575  // seen by the wire planes
576  if (sce->EnableCalSpatialSCE() && fConfig.FieldDistortion() &&
578  geo::Point_t loc = track.LocationAtPoint(meta->Index());
579 
580  // compute the dir of the track trajectory
581  geo::Vector_t track_dir = track.DirectionAtPoint(meta->Index());
582  geo::Point_t loc_mdx = loc - track_dir * (geom->WirePitch(hit->View()) / 2.);
583  geo::Point_t loc_pdx = loc + track_dir * (geom->WirePitch(hit->View()) / 2.);
584 
585  loc_mdx = TrajectoryToWirePosition(loc_mdx, hit->WireID());
586  loc_pdx = TrajectoryToWirePosition(loc_pdx, hit->WireID());
587 
588  // Direction at wires
589  dir = (loc_pdx - loc_mdx) / (loc_mdx - loc_pdx).r();
590  }
591  // If there is no space charge or the track is not yet corrected, then the dir
592  // is the track is what we want
593  else {
594  dir = track.DirectionAtPoint(meta->Index());
595  }
596 
597  double cosgamma = std::abs(std::sin(angleToVert) * dir.Y() + std::cos(angleToVert) * dir.Z());
598  double pitch;
599  if (cosgamma) { pitch = geom->WirePitch(hit->View()) / cosgamma; }
600  else {
601  pitch = 0.;
602  }
603 
604  // now take the pitch computed on the wires and correct it back to the particle trajectory
605  geo::Point_t loc_w = GetLocationAtWires(track, hit, meta);
606 
607  geo::Point_t locw_pdx_traj = WireToTrajectoryPosition(loc_w + pitch * dir, hit->WireID());
608  geo::Point_t loc = WireToTrajectoryPosition(loc_w, hit->WireID());
609 
610  pitch = (locw_pdx_traj - loc).R();
611 
612  return pitch;
613 }
614 
616  const std::vector<recob::Hit const*>& sharedHits)
617 {
618  switch (fConfig.ChargeMethod()) {
623  return std::accumulate(
624  sharedHits.cbegin(),
625  sharedHits.cend(),
626  hit->Integral(),
627  [](double sum, recob::Hit const* sharedHit) { return sum + sharedHit->Integral(); });
628  default: return 0.;
629  }
630  return 0.;
631 }
632 
634  const art::Event& e,
635  const recob::Hit& h,
636  const geo::Point_t& location,
637  const geo::Vector_t& direction,
638  double t0)
639 {
640  double ret = dQdx;
641  for (auto const& nt : fNormTools) {
642  ret = nt->Normalize(ret, e, h, location, direction, t0);
643  }
644 
645  return ret;
646 }
647 
649  const recob::Track& track,
651  const recob::TrackHitMeta* meta)
652 {
653  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
654 
655  double EField = dprop.Efield();
656  if (sce->EnableSimEfieldSCE() && fConfig.FieldDistortionEfield()) {
657  // Gets relative E field Distortions
658  geo::Vector_t EFieldOffsets = sce->GetEfieldOffsets(GetLocation(track, hit, meta));
659  // Add 1 in X direction as this is the direction of the drift field
660  EFieldOffsets = EFieldOffsets + geo::Vector_t{1, 0, 0};
661  // Convert to Absolute E Field from relative
662  EFieldOffsets = EField * EFieldOffsets;
663  // We only care about the magnitude for recombination
664  EField = EFieldOffsets.r();
665  }
666  return EField;
667 }
668 
TRandom r
Definition: spectrum.C:23
code to link reconstructed objects back to the MC truth information
std::vector< std::vector< OrganizedHits > > OrganizeHitsIndividual(const std::vector< art::Ptr< recob::Hit >> &hits, const std::vector< const recob::TrackHitMeta * > &thms, const recob::Track &track, unsigned nplanes)
GnocchiCalorimetry(Parameters const &config)
Functions to help with numbers.
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Vector_t
Type for representation of momenta in 3D space.
Definition: geo_vectors.h:160
void produce(art::Event &evt) override
constexpr bool operator>(Interval< Q, Cat > const a, Quantity< Args... > const b) noexcept
Definition: intervals.h:507
double GetCharge(const art::Ptr< recob::Hit > hit, const std::vector< recob::Hit const * > &sharedHits)
const recob::TrackTrajectory & Trajectory() const
Access to the stored recob::TrackTrajectory.
Definition: Track.h:132
Declaration of signal hit object.
Point_t const & LocationAtPoint(size_t i) const
Access to track position at different points.
Definition: Track.h:160
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
bool isValid
Whether this ID points to a valid element.
Definition: geo_types.h:210
bool HasValidPoint(size_t i) const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:145
constexpr auto abs(T v)
Returns the absolute value of the argument.
Class to keep data related to recob::Hit associated with recob::Track.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:244
geo::View_t View() const
View for the plane of the hit.
Definition: Hit.h:276
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
data_const_reference data(size_type i) const
Definition: FindManyP.h:446
Data related to recob::Hit associated with recob::Track.The purpose is to collect several variables t...
Definition: TrackHitMeta.h:43
double Normalize(double dQdx, const art::Event &e, const recob::Hit &h, const geo::Point_t &location, const geo::Vector_t &direction, double t0)
geo::Point_t GetLocationAtWires(const recob::Track &track, const art::Ptr< recob::Hit > hit, const recob::TrackHitMeta *meta)
float PeakAmplitude() const
The estimated amplitude of the hit at its peak, in ADC units.
Definition: Hit.h:232
double Efield(unsigned int planegap=0) const
kV/cm
fhicl::Atom< std::string > TrackModuleLabel
TPCGeo const & TPC(TPCID const &tpcid=tpc_zero) const
Returns the specified TPC.
Definition: GeometryCore.h:722
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
fhicl::Atom< std::string > T0ModuleLabel
double WireAngleToVertical(View_t view, TPCID const &tpcid) const
Returns the angle of the wires in the specified view from vertical.
void hits()
Definition: readHits.C:15
double dEdx_AMP(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, recob::Hit const &hit, double pitch, double T0=0) const
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
std::vector< std::vector< OrganizedHits > > OrganizeHitsSnippets(const std::vector< art::Ptr< recob::Hit >> &hits, const std::vector< const recob::TrackHitMeta * > &thms, const recob::Track &track, unsigned nplanes)
std::vector< std::unique_ptr< INormalizeCharge > > fNormTools
double GetEfield(const detinfo::DetectorPropertiesData &dprop, const recob::Track &track, const art::Ptr< recob::Hit > hit, const recob::TrackHitMeta *meta)
Provides recob::Track data product.
geo::Point_t WireToTrajectoryPosition(const geo::Point_t &loc, const geo::TPCID &tpc)
constexpr PlaneID const & asPlaneID() const
Conversion to PlaneID (for convenience of notation).
Definition: geo_types.h:520
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:481
raw::TDCtick_t StartTick() const
Initial tdc tick for hit.
Definition: Hit.h:212
float dEdx(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const TCSlice &slc, TP3D &tp3d)
Definition: PFPUtils.cxx:2675
Detector simulation of raw signals on wires.
Encapsulate the geometry of a wire.
raw::TDCtick_t EndTick() const
Final tdc tick for hit.
Definition: Hit.h:216
fhicl::Table< calo::CalorimetryAlg::Config > CalorimetryAlgConfig
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
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
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Utility object to perform functions of association.
Encapsulate the construction of a single detector plane.
TDirectory * dir
Definition: macro.C:5
Point_t const & End() const
Returns the position of the last valid point of the trajectory [cm].
double dEdx_AREA(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, recob::Hit const &hit, double pitch, double T0=0) const
double GetPitch(const recob::Track &track, const art::Ptr< recob::Hit > hit, const recob::TrackHitMeta *meta)
geo::Point_t TrajectoryToWirePosition(const geo::Point_t &loc, const geo::TPCID &tpc)
float SummedADC() const
The sum of calibrated ADC counts of the hit (0. by default)
Definition: Hit.h:240
unsigned int Index() const
Hit index along the track trajectory.
Definition: TrackHitMeta.h:52
unsigned int Nplanes(TPCID const &tpcid=tpc_zero) const
Returns the total number of planes in the specified TPC.
Definition: GeometryCore.h:977
Vector_t DirectionAtPoint(size_t i) const
Access to track direction at different points.
Definition: Track.h:168
constexpr double kBogusD
obviously bogus double value
geo::Point_t GetLocation(const recob::Track &track, const art::Ptr< recob::Hit > hit, const recob::TrackHitMeta *meta)
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:46
TCEvent evt
Definition: DataStructs.cxx:8
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
Collection of Physical constants used in LArSoft.
Double_t sum
Definition: plot.C:31
Float_t e
Definition: plot.C:35
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:82
fhicl::Atom< std::string > AssocHitModuleLabel
Vector_t DriftDir() const
Returns the direction of the drift (vector pointing toward the planes).
Definition: TPCGeo.h:124
Point_t const & Start() const
Returns the position of the first valid point of the trajectory [cm].
bool HitIsValid(const recob::TrackHitMeta *thm, const recob::Track &track)
Float_t track
Definition: plot.C:35
Utility functions to extract information from recob::Track
bool operator==(infinite_endcount_iterator< T > const &, count_iterator< T > const &)
Definition: counter.h:278
std::vector< std::vector< OrganizedHits > > OrganizeHits(const std::vector< art::Ptr< recob::Hit >> &hits, const std::vector< const recob::TrackHitMeta * > &thms, const recob::Track &track, unsigned nplanes)
Length_t WirePitch(PlaneID const &planeid=plane_zero) const
Returns the distance between two consecutive wires.
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:49
calorimetry
This is an interface for an art Tool which scales charge by some factor given information about its a...