LArSoft  v09_93_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()
31 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
32 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
33 
36 
37 // ROOT includes
38 #include <TF1.h>
39 #include <TGraph.h>
40 #include <TMath.h>
41 #include <vector>
42 
51 #include "cetlib/pow.h" // cet::sum_of_squares()
52 #include "fhiclcpp/ParameterSet.h"
55 
56 namespace {
57  constexpr unsigned int int_max_as_unsigned_int{std::numeric_limits<int>::max()};
58  typedef std::pair<unsigned, std::vector<unsigned>> OrganizedHits;
59 }
60 
61 namespace calo {
62 
64  public:
65  struct Config {
67  using Name = fhicl::Name;
68 
69  enum ChargeMethod : unsigned {
74  };
75 
77  Comment("Module label for track producer.")};
78 
80  Comment("Module label for T0 time producer."),
81  ""};
82 
84  Name("PFPModuleLabel"),
85  Comment("Module label for PFP producer. To be used to associate T0 with tracks."),
86  ""};
87 
89  Name("AssocHitModuleLabel"),
90  Comment("Module label for association between tracks and hits. If not set, defaults to "
91  "TrackModuleLabel."),
92  ""};
93 
95  Name("ChargeMethod"),
96  Comment("Method used to extract charge from a hit. Options: 0==Amplitude(), 1==Integral(), "
97  "2==SummedADC(), 3==SummedIntegral(). See the ChargeMethod enum.")};
98 
100  Name("FieldDistortion"),
101  Comment("True if field distortion (i.e. from space charge) is included in the input.")};
102 
104  Name("FieldDistortionEfield"),
105  Comment("True if field distortion (i.e. from space charge) is included in the input.")};
106 
108  Name("TrackIsFieldDistortionCorrected"),
109  Comment("Whether the space-points on the input tracks have their points corrected for the "
110  "field distortions. "
111  "I.e. whether the track trajectory points represent charge as seen by wires or the "
112  "3D particle trajectory.")};
113 
115  Comment("Which cryostat number the input tracks occupy.")};
116 
118  Name("FieldDistortionCorrectionXSign"),
119  Comment("Sign of the field distortion correction to be applied in the X direction. "
120  "Positive by default."),
121  1.f};
122 
124  Name("CaloAlg"),
125  Comment("Configuration for the calo::CalorimetryAlg")};
126 
128  Name("NormTools"),
129  Comment("List of INormalizeCharge tool configurations to use.")};
130  };
131 
133 
134  explicit GnocchiCalorimetry(Parameters const& config);
135 
136  void produce(art::Event& evt) override;
137 
138  private:
141  std::vector<std::unique_ptr<INormalizeCharge>> fNormTools;
142 
143  // helper functions
144  std::vector<std::vector<OrganizedHits>> OrganizeHits(
146  const std::vector<const recob::TrackHitMeta*>& thms,
147  const recob::Track& track,
148  unsigned nplanes);
149  std::vector<std::vector<OrganizedHits>> OrganizeHitsIndividual(
151  const std::vector<const recob::TrackHitMeta*>& thms,
152  const recob::Track& track,
153  unsigned nplanes);
154  std::vector<std::vector<OrganizedHits>> OrganizeHitsSnippets(
156  const std::vector<const recob::TrackHitMeta*>& thms,
157  const recob::Track& track,
158  unsigned nplanes);
159  bool HitIsValid(const recob::TrackHitMeta* thm, const recob::Track& track);
162  const recob::TrackHitMeta* meta);
165  const recob::TrackHitMeta* meta);
168  double GetPitch(const recob::Track& track,
170  const recob::TrackHitMeta* meta);
171  double GetCharge(const art::Ptr<recob::Hit> hit,
172  const std::vector<recob::Hit const*>& sharedHits);
173  double GetEfield(const detinfo::DetectorPropertiesData& dprop,
174  const recob::Track& track,
176  const recob::TrackHitMeta* meta);
177  double Normalize(double dQdx,
178  const art::Event& e,
179  const recob::Hit& h,
180  const geo::Point_t& location,
181  const geo::Vector_t& direction,
182  double t0);
183  };
184 
185 } // end namespace calo
186 
189 {
190  produces<std::vector<anab::Calorimetry>>();
191  produces<art::Assns<recob::Track, anab::Calorimetry>>();
192  std::vector<fhicl::ParameterSet> norm_tool_configs =
193  fConfig.NormTools.get<std::vector<fhicl::ParameterSet>>();
194  for (const fhicl::ParameterSet& p : norm_tool_configs) {
195  fNormTools.push_back(art::make_tool<INormalizeCharge>(p));
196  }
197 }
198 
200 {
201  // Get services
203  auto const clock_data = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
204  auto const det_prop =
206 
207  size_t nplanes = geom->Nplanes();
208 
209  // Define output collections
210  std::unique_ptr<std::vector<anab::Calorimetry>> outputCalo(new std::vector<anab::Calorimetry>);
211  std::unique_ptr<art::Assns<recob::Track, anab::Calorimetry>> outputCaloAssn(
213 
214  // collect input
215  art::Handle<std::vector<recob::Track>> trackListHandle;
216  std::vector<art::Ptr<recob::Track>> tracklist;
217  if (evt.getByLabel(fConfig.TrackModuleLabel(), trackListHandle)) {
218  art::fill_ptr_vector(tracklist, trackListHandle);
219  }
220 
221  // get the label to collect this hits
222  const std::string& hitLabel = (fConfig.AssocHitModuleLabel() == "") ?
225  art::FindManyP<recob::Hit, recob::TrackHitMeta> fmHits(trackListHandle, evt, hitLabel);
226 
227  // must be valid if the T0 module label is non-empty
228  art::FindManyP<anab::T0> fmT0s(trackListHandle, evt, fConfig.T0ModuleLabel());
229  art::FindManyP<recob::PFParticle> fmPFPs(trackListHandle, evt, fConfig.TrackModuleLabel());
230  std::vector<art::Ptr<recob::PFParticle>> pfpList;
231  if (fConfig.PFPModuleLabel().size()) {
232  for (std::size_t itrk = 0, sz = fmPFPs.size(); itrk < sz; ++itrk) {
233  // assuming every recob::Track has exactly one associated recob::PFParticle
234  pfpList.push_back(fmPFPs.at(itrk).at(0));
235  }
236  }
237  art::FindManyP<anab::T0> fmPFPT0s(pfpList, evt, fConfig.T0ModuleLabel());
238 
239  for (auto const& nt : fNormTools)
240  nt->setup(evt);
241 
242  // iterate over all the tracks
243  for (unsigned trk_i = 0; trk_i < tracklist.size(); trk_i++) {
244  const recob::Track& track = *tracklist[trk_i];
245 
246  // collect input for this track
247  const std::vector<art::Ptr<recob::Hit>>& hits = fmHits.at(trk_i);
248  const std::vector<const recob::TrackHitMeta*>& thms = fmHits.data(trk_i);
249 
250  double T0 = 0;
251  if (fConfig.T0ModuleLabel().size()) {
252  if (fConfig.PFPModuleLabel().size()) {
253  const std::vector<art::Ptr<anab::T0>>& this_t0s = fmPFPT0s.at(trk_i);
254  if (this_t0s.size()) T0 = this_t0s.at(0)->Time();
255  }
256  else {
257  const std::vector<art::Ptr<anab::T0>>& this_t0s = fmT0s.at(trk_i);
258  if (this_t0s.size()) T0 = this_t0s.at(0)->Time();
259  }
260  }
261 
262  // organize the hits by plane
263  std::vector<std::vector<OrganizedHits>> hit_indices = OrganizeHits(hits, thms, track, nplanes);
264 
265  for (unsigned plane_i = 0; plane_i < nplanes; plane_i++) {
266 
267  float kinetic_energy = 0.;
268  std::vector<float> dEdxs;
269  std::vector<float> dQdxs;
270  std::vector<float> resranges;
271  std::vector<float> deadwireresranges;
272  float range = 0.;
273  std::vector<float> pitches;
274  std::vector<geo::Point_t> xyzs;
275  std::vector<size_t> tp_indices;
276  geo::PlaneID plane;
277 
278  // setup the plane ID
279  plane.Plane = plane_i;
280  plane.TPC = 0; // arbitrary -- tracks can cross TPC boundaries
281  plane.Cryostat = fConfig.Cryostat();
282  plane.isValid = true;
283 
284  std::vector<float> lengths;
285  for (unsigned hit_i = 0; hit_i < hit_indices[plane_i].size(); hit_i++) {
286  unsigned hit_index = hit_indices[plane_i][hit_i].first;
287 
288  std::vector<recob::Hit const*> sharedHits = {};
289  sharedHits.reserve(hit_indices[plane_i][hit_i].second.size());
290  for (const unsigned shared_hit_index : hit_indices[plane_i][hit_i].second)
291  sharedHits.push_back(hits[shared_hit_index].get());
292 
293  // Get the location of this point
294  geo::Point_t location = GetLocation(track, hits[hit_index], thms[hit_index]);
295 
296  // Get the pitch
297  double pitch = GetPitch(track, hits[hit_index], thms[hit_index]);
298 
299  // And the charge
300  double charge = GetCharge(hits[hit_index], sharedHits);
301 
302  // Get the EField
303  double EField = GetEfield(det_prop, track, hits[hit_index], thms[hit_index]);
304 
305  // Angle to the drift electric field (in x direction), in units of degrees
306  geo::Vector_t direction = track.DirectionAtPoint(thms[hit_index]->Index());
307  double phi = acos(abs(direction.x())) * 180 / M_PI;
308 
309  double dQdx = charge / pitch;
310 
311  // Normalize out the detector response
312  dQdx = Normalize(dQdx,
313  evt,
314  *hits[hit_index],
315  track.LocationAtPoint(thms[hit_index]->Index()),
316  track.DirectionAtPoint(thms[hit_index]->Index()),
317  T0);
318 
319  // turn into dEdx
321  fCaloAlg.dEdx_AMP(clock_data,
322  det_prop,
323  dQdx,
324  hits[hit_index]->PeakTime(),
325  hits[hit_index]->WireID().Plane,
326  T0,
327  EField,
328  phi) :
329  fCaloAlg.dEdx_AREA(clock_data,
330  det_prop,
331  dQdx,
332  hits[hit_index]->PeakTime(),
333  hits[hit_index]->WireID().Plane,
334  T0,
335  EField,
336  phi);
337 
338  // save the length between each pair of hits
339  if (xyzs.size() == 0) { lengths.push_back(0.); }
340  else {
341  lengths.push_back((location - xyzs.back()).r());
342  }
343 
344  // save stuff
345  dEdxs.push_back(dEdx);
346  dQdxs.push_back(dQdx);
347  pitches.push_back(pitch);
348  xyzs.push_back(location);
349  kinetic_energy += dEdx * pitch;
350 
351  // TODO: FIXME
352  // It seems weird that the "trajectory-point-index" actually is the
353  // index of the hit... is this a bug in the documentation
354  // of anab::Calorimetry?
355  //
356  // i.e. -- I think this piece of code should actually be:
357  // tp_indices.push_back(thms[hit_index]->Index());
358  tp_indices.push_back(hits[hit_index].key());
359 
360  } // end iterate over hits
361 
362  // turn the lengths vector into a residual-range vector and total length
363  if (lengths.size() > 1) {
364  range = std::accumulate(lengths.begin(), lengths.end(), 0.);
365 
366  // check the direction that the hits are going in the track:
367  // upstream (end-start) or downstream (start-end)
368  bool is_downstream =
369  (track.Trajectory().Start() - xyzs[0]).r() +
370  (track.Trajectory().End() - xyzs.back()).r() <
371  (track.Trajectory().End() - xyzs[0]).r() + (track.Trajectory().Start() - xyzs.back()).r();
372 
373  resranges.resize(lengths.size());
374  if (is_downstream) {
375  resranges[lengths.size() - 1] = lengths.back() / 2.;
376  for (int i_len = lengths.size() - 2; i_len >= 0; i_len--) {
377  resranges[i_len] = resranges[i_len + 1] + lengths[i_len + 1];
378  }
379  }
380  else {
381  resranges[0] = lengths[1] / 2.;
382  for (unsigned i_len = 1; i_len < lengths.size(); i_len++) {
383  resranges[i_len] = resranges[i_len - 1] + lengths[i_len];
384  }
385  }
386  }
387 
388  // save the Calorimetry output
389  //
390  // Bogus if less than two hits on this plane
391  if (lengths.size() > 1) {
392  outputCalo->push_back(anab::Calorimetry(kinetic_energy,
393  dEdxs,
394  dQdxs,
395  resranges,
396  deadwireresranges,
397  range,
398  pitches,
399  xyzs,
400  tp_indices,
401  plane));
402  }
403  else {
404  outputCalo->push_back(
405  anab::Calorimetry(util::kBogusD, {}, {}, {}, {}, util::kBogusD, {}, {}, {}, plane));
406  }
407 
408  util::CreateAssn(evt, *outputCalo, tracklist[trk_i], *outputCaloAssn);
409 
410  } // end iterate over planes
411 
412  } // end iterate over tracks
413 
414  evt.put(std::move(outputCalo));
415  evt.put(std::move(outputCaloAssn));
416 
417  return;
418 }
419 
420 std::vector<std::vector<OrganizedHits>> calo::GnocchiCalorimetry::OrganizeHits(
422  const std::vector<const recob::TrackHitMeta*>& thms,
423  const recob::Track& track,
424  unsigned nplanes)
425 {
426  // charge is computed per hit -- we organize hits indivudally
429  return OrganizeHitsIndividual(hits, thms, track, nplanes);
430  }
431  // charge is computed per snippet -- we organize hits by snippet
432  else {
433  return OrganizeHitsSnippets(hits, thms, track, nplanes);
434  }
435 }
436 
437 std::vector<std::vector<OrganizedHits>> calo::GnocchiCalorimetry::OrganizeHitsIndividual(
439  const std::vector<const recob::TrackHitMeta*>& thms,
440  const recob::Track& track,
441  unsigned nplanes)
442 {
443  std::vector<std::vector<OrganizedHits>> ret(nplanes);
444  for (unsigned i = 0; i < hits.size(); i++) {
445  if (HitIsValid(thms[i], track)) { ret[hits[i]->WireID().Plane].push_back({i, {}}); }
446  }
447 
448  return ret;
449 }
450 
451 std::vector<std::vector<OrganizedHits>> calo::GnocchiCalorimetry::OrganizeHitsSnippets(
453  const std::vector<const recob::TrackHitMeta*>& thms,
454  const recob::Track& track,
455  unsigned nplanes)
456 {
457  // In this case, we need to only accept one hit in each snippet
458  // 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.
459  //
460  // If there are multiple valid hits on the same snippet, we need a way to pick the best one.
461  // (TODO: find a good way). The current method is to take the one with the highest charge integral.
462  struct HitIdentifier {
463  int startTick;
464  int endTick;
465  int wire;
466  float integral;
467 
468  // construct
469  explicit HitIdentifier(const recob::Hit& hit)
470  : startTick(hit.StartTick())
471  , endTick(hit.EndTick())
472  , wire(hit.WireID().Wire)
473  , integral(hit.Integral())
474  {}
475 
476  // Defines whether two hits are on the same snippet
477  inline bool operator==(const HitIdentifier& rhs) const
478  {
479  return startTick == rhs.startTick && endTick == rhs.endTick && wire == rhs.wire;
480  }
481 
482  // Defines which hit to pick between two both on the same snippet
483  inline bool operator>(const HitIdentifier& rhs) const { return integral > rhs.integral; }
484  };
485 
486  std::vector<std::vector<OrganizedHits>> ret(nplanes);
487  std::vector<std::vector<HitIdentifier>> hit_idents(nplanes);
488  for (unsigned i = 0; i < hits.size(); i++) {
489  if (HitIsValid(thms[i], track)) {
490  HitIdentifier this_ident(*hits[i]);
491 
492  // check if we have found a hit on this snippet before
493  bool found_snippet = false;
494  auto const plane = hits[i]->WireID().Plane;
495  for (unsigned j = 0; j < ret[plane].size(); j++) {
496  if (this_ident == hit_idents[plane][j]) {
497  found_snippet = true;
498  if (this_ident > hit_idents[plane][j]) {
499  ret[plane][j].second.push_back(ret[plane][j].first);
500  ret[plane][j].first = i;
501  hit_idents[plane][j] = this_ident;
502  }
503  else {
504  ret[plane][j].second.push_back(i);
505  }
506  break;
507  }
508  }
509  if (!found_snippet) {
510  ret[plane].push_back({i, {}});
511  hit_idents[plane].push_back(this_ident);
512  }
513  }
514  }
515  return ret;
516 }
517 
519 {
520  if (thm->Index() == int_max_as_unsigned_int) return false;
521  return track.HasValidPoint(thm->Index());
522 }
523 
526  const recob::TrackHitMeta* meta)
527 {
528  geo::Point_t loc = track.LocationAtPoint(meta->Index());
530  loc;
531 }
532 
534  const geo::TPCID& tpc)
535 {
536  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
537 
538  geo::Point_t ret = loc;
539 
540  if (sce->EnableCalSpatialSCE() && fConfig.FieldDistortion()) {
541  geo::Vector_t offset = sce->GetCalPosOffsets(ret, tpc.TPC);
542 
543  ret.SetX(ret.X() + fConfig.FieldDistortionCorrectionXSign() * offset.X());
544  ret.SetY(ret.Y() + offset.Y());
545  ret.SetZ(ret.Z() + offset.Z());
546  }
547 
548  return ret;
549 }
550 
553  const recob::TrackHitMeta* meta)
554 {
555  geo::Point_t loc = track.LocationAtPoint(meta->Index());
557  loc;
558 }
559 
561  const geo::TPCID& tpc)
562 {
563  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
565 
566  geo::Point_t ret = loc;
567 
568  if (sce->EnableCalSpatialSCE() && fConfig.FieldDistortion()) {
569  // Returned X is the drift -- multiply by the drift direction to undo this
570  int corr = geom->TPC(tpc).DriftDir().X();
571 
572  geo::Vector_t offset = sce->GetPosOffsets(ret);
573 
574  ret.SetX(ret.X() + corr * fConfig.FieldDistortionCorrectionXSign() * offset.X());
575  ret.SetY(ret.Y() + offset.Y());
576  ret.SetZ(ret.Z() + offset.Z());
577  }
578 
579  return ret;
580 }
581 
584  const recob::TrackHitMeta* meta)
585 {
587  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
588 
589  double angleToVert =
590  geom->WireAngleToVertical(hit->View(), hit->WireID().asPlaneID()) - 0.5 * ::util::pi<>();
591 
593 
594  // "dir" should be the direction that the wires see. If the track already has the field
595  // distortion corrections applied, then we need to de-apply them to get the direction as
596  // seen by the wire planes
597  if (sce->EnableCalSpatialSCE() && fConfig.FieldDistortion() &&
599  geo::Point_t loc = track.LocationAtPoint(meta->Index());
600 
601  // compute the dir of the track trajectory
602  geo::Vector_t track_dir = track.DirectionAtPoint(meta->Index());
603  geo::Point_t loc_mdx = loc - track_dir * (geom->WirePitch(hit->View()) / 2.);
604  geo::Point_t loc_pdx = loc + track_dir * (geom->WirePitch(hit->View()) / 2.);
605 
606  loc_mdx = TrajectoryToWirePosition(loc_mdx, hit->WireID());
607  loc_pdx = TrajectoryToWirePosition(loc_pdx, hit->WireID());
608 
609  // Direction at wires
610  dir = (loc_pdx - loc_mdx) / (loc_mdx - loc_pdx).r();
611  }
612  // If there is no space charge or the track is not yet corrected, then the dir
613  // is the track is what we want
614  else {
615  dir = track.DirectionAtPoint(meta->Index());
616  }
617 
618  double cosgamma = std::abs(std::sin(angleToVert) * dir.Y() + std::cos(angleToVert) * dir.Z());
619  double pitch;
620  if (cosgamma) { pitch = geom->WirePitch(hit->View()) / cosgamma; }
621  else {
622  pitch = 0.;
623  }
624 
625  // now take the pitch computed on the wires and correct it back to the particle trajectory
626  geo::Point_t loc_w = GetLocationAtWires(track, hit, meta);
627 
628  geo::Point_t locw_pdx_traj = WireToTrajectoryPosition(loc_w + pitch * dir, hit->WireID());
629  geo::Point_t loc = WireToTrajectoryPosition(loc_w, hit->WireID());
630 
631  pitch = (locw_pdx_traj - loc).R();
632 
633  return pitch;
634 }
635 
637  const std::vector<recob::Hit const*>& sharedHits)
638 {
639  switch (fConfig.ChargeMethod()) {
644  return std::accumulate(
645  sharedHits.cbegin(),
646  sharedHits.cend(),
647  hit->Integral(),
648  [](double sum, recob::Hit const* sharedHit) { return sum + sharedHit->Integral(); });
649  default: return 0.;
650  }
651  return 0.;
652 }
653 
655  const art::Event& e,
656  const recob::Hit& h,
657  const geo::Point_t& location,
658  const geo::Vector_t& direction,
659  double t0)
660 {
661  double ret = dQdx;
662  for (auto const& nt : fNormTools) {
663  ret = nt->Normalize(ret, e, h, location, direction, t0);
664  }
665 
666  return ret;
667 }
668 
670  const recob::Track& track,
672  const recob::TrackHitMeta* meta)
673 {
674  auto const* sce = lar::providerFrom<spacecharge::SpaceChargeService>();
675 
676  double EField = dprop.Efield();
677  if (sce->EnableSimEfieldSCE() && fConfig.FieldDistortionEfield()) {
678  // Gets relative E field Distortions
679  geo::Vector_t EFieldOffsets = sce->GetEfieldOffsets(GetLocation(track, hit, meta));
680  // Add 1 in X direction as this is the direction of the drift field
681  EFieldOffsets = EFieldOffsets + geo::Vector_t{1, 0, 0};
682  // Convert to Absolute E Field from relative
683  EFieldOffsets = EField * EFieldOffsets;
684  // We only care about the magnitude for recombination
685  EField = EFieldOffsets.r();
686  }
687  return EField;
688 }
689 
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:254
geo::View_t View() const
View for the plane of the hit.
Definition: Hit.h:286
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:290
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)
Provides recob::Track data product.
float PeakAmplitude() const
The estimated amplitude of the hit at its peak, in ADC units.
Definition: Hit.h:238
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)
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:218
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:222
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.
float ROISummedADC() const
The sum of calibrated ADC counts of the ROI (0. by default)
Definition: Hit.h:246
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)
unsigned int Index() const
Hit index along the track trajectory.
Definition: TrackHitMeta.h:52
fhicl::Atom< std::string > PFPModuleLabel
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...