LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
LArPandoraShowerAlg.cxx
Go to the documentation of this file.
3 
12 
14 #include "fhiclcpp/ParameterSet.h"
15 
16 #include "TCanvas.h"
17 #include "TH3.h"
18 #include "TPolyLine3D.h"
19 #include "TPolyMarker3D.h"
20 #include "TString.h"
21 #include "TStyle.h"
22 
23 #include <memory>
24 
26  : fUseCollectionOnly(pset.get<bool>("UseCollectionOnly"))
27  , fPFParticleLabel(pset.get<art::InputTag>("PFParticleLabel"))
28  , fSCEXFlip(pset.get<bool>("SCEXFlip"))
29  , fSCE(lar::providerFrom<spacecharge::SpaceChargeService>())
30  , fInitialTrackInputLabel(pset.get<std::string>("InitialTrackInputLabel"))
31  , fShowerStartPositionInputLabel(pset.get<std::string>("ShowerStartPositionInputLabel"))
32  , fShowerDirectionInputLabel(pset.get<std::string>("ShowerDirectionInputLabel"))
33  , fInitialTrackSpacePointsInputLabel(pset.get<std::string>("InitialTrackSpacePointsInputLabel"))
34 {}
35 
36 //Order the shower hits with regards to their projected length onto
37 //the shower direction from the shower start position. This is done
38 //in the 2D coordinate system (wire direction, x)
41  geo::Point_t const& ShowerStartPosition,
42  geo::Vector_t const& ShowerDirection) const
43 {
44 
45  std::map<double, art::Ptr<recob::Hit>> OrderedHits;
46  art::Ptr<recob::Hit> startHit = hits.front();
47 
48  //Get the wireID
49  const geo::WireID startWireID = startHit->WireID();
50 
51  //Get the plane
52  const geo::PlaneID planeid = startWireID.asPlaneID();
53 
54  //Get the pitch
55  double pitch = fGeom->WirePitch(planeid);
56 
57  TVector2 Shower2DStartPosition = {
58  fGeom->WireCoordinate(ShowerStartPosition, startHit->WireID().planeID()) * pitch,
59  ShowerStartPosition.X()};
60 
61  //Vector of the plane
62  auto const PlaneDirection = fGeom->Plane(planeid).GetIncreasingWireDirection();
63 
64  //get the shower 2D direction
65  TVector2 Shower2DDirection = {ShowerDirection.Dot(PlaneDirection), ShowerDirection.X()};
66 
67  Shower2DDirection = Shower2DDirection.Unit();
68 
69  for (auto const& hit : hits) {
70 
71  //Get the wireID
72  const geo::WireID WireID = hit->WireID();
73 
74  if (WireID.asPlaneID() != startWireID.asPlaneID()) { break; }
75 
76  //Get the hit Vector.
77  TVector2 hitcoord = HitCoordinates(detProp, hit);
78 
79  //Order the hits based on the projection
80  TVector2 pos = hitcoord - Shower2DStartPosition;
81  OrderedHits[pos * Shower2DDirection] = hit;
82  }
83 
84  //Transform the shower.
85  std::vector<art::Ptr<recob::Hit>> showerHits;
86  std::transform(OrderedHits.begin(),
87  OrderedHits.end(),
88  std::back_inserter(showerHits),
89  [](std::pair<double, art::Ptr<recob::Hit>> const& hit) { return hit.second; });
90 
91  //Sometimes get the order wrong. Depends on direction compared to the plane Correct for it here:
92  art::Ptr<recob::Hit> frontHit = showerHits.front();
93  art::Ptr<recob::Hit> backHit = showerHits.back();
94 
95  //Get the hit Vector.
96  TVector2 fronthitcoord = HitCoordinates(detProp, frontHit);
97  TVector2 frontpos = fronthitcoord - Shower2DStartPosition;
98 
99  //Get the hit Vector.
100  TVector2 backhitcoord = HitCoordinates(detProp, backHit);
101  TVector2 backpos = backhitcoord - Shower2DStartPosition;
102 
103  double frontproj = frontpos * Shower2DDirection;
104  double backproj = backpos * Shower2DDirection;
105  if (std::abs(backproj) < std::abs(frontproj)) {
106  std::reverse(showerHits.begin(), showerHits.end());
107  }
108 
109  hits = showerHits;
110  return;
111 }
112 
113 //Orders the shower spacepoints with regards to there perpendicular distance from
114 //the shower axis.
117  geo::Point_t const& vertex,
118  geo::Vector_t const& direction) const
119 {
120 
121  std::map<double, art::Ptr<recob::SpacePoint>> OrderedSpacePoints;
122 
123  //Loop over the spacepoints and get the pojected distance from the vertex.
124  for (auto const& sp : showersps) {
125 
126  // Get the perpendicular distance
127  double perp = SpacePointPerpendicular(sp, vertex, direction);
128 
129  //Add to the list
130  OrderedSpacePoints[perp] = sp;
131  }
132 
133  //Return an ordered list.
134  showersps.clear();
135  for (auto const& sp : OrderedSpacePoints) {
136  showersps.push_back(sp.second);
137  }
138 }
139 
140 //Orders the shower spacepoints with regards to there prejected length from
141 //the shower start position in the shower direction.
144  geo::Point_t const& vertex,
145  geo::Vector_t const& direction) const
146 {
147 
148  std::map<double, art::Ptr<recob::SpacePoint>> OrderedSpacePoints;
149 
150  //Loop over the spacepoints and get the pojected distance from the vertex.
151  for (auto const& sp : showersps) {
152 
153  // Get the projection of the space point along the direction
154  double len = SpacePointProjection(sp, vertex, direction);
155 
156  //Add to the list
157  OrderedSpacePoints[len] = sp;
158  }
159 
160  //Return an ordered list.
161  showersps.clear();
162  for (auto const& sp : OrderedSpacePoints) {
163  showersps.push_back(sp.second);
164  }
165 }
166 
169  geo::Point_t const& vertex) const
170 {
171  std::map<double, art::Ptr<recob::SpacePoint>> OrderedSpacePoints;
172 
173  //Loop over the spacepoints and get the pojected distance from the vertex.
174  for (auto const& sp : showersps) {
175 
176  //Get the distance away from the start
177  double mag = (sp->position() - vertex).R();
178 
179  //Add to the list
180  OrderedSpacePoints[mag] = sp;
181  }
182 
183  //Return an ordered list.
184  showersps.clear();
185  for (auto const& sp : OrderedSpacePoints) {
186  showersps.push_back(sp.second);
187  }
188 }
189 
191  std::vector<art::Ptr<recob::SpacePoint>> const& showersps) const
192 {
193  if (showersps.empty()) return {};
194 
195  double x{};
196  double y{};
197  double z{};
198  for (auto const& sp : showersps) {
199  auto const pos = sp->position();
200  x += pos.X();
201  y += pos.Y();
202  z += pos.Z();
203  }
204  geo::Point_t centre_position{x, y, z};
205  centre_position *= (1. / showersps.size());
206  return centre_position;
207 }
208 
210  detinfo::DetectorClocksData const& clockData,
211  detinfo::DetectorPropertiesData const& detProp,
212  std::vector<art::Ptr<recob::SpacePoint>> const& showerspcs,
213  art::FindManyP<recob::Hit> const& fmh) const
214 {
215  float totalCharge = 0;
217  clockData, detProp, showerspcs, fmh, totalCharge);
218 }
219 
220 //Returns the vector to the shower centre and the total charge of the shower.
222  detinfo::DetectorClocksData const& clockData,
223  detinfo::DetectorPropertiesData const& detProp,
224  std::vector<art::Ptr<recob::SpacePoint>> const& showersps,
225  art::FindManyP<recob::Hit> const& fmh,
226  float& totalCharge) const
227 {
228  geo::Point_t chargePoint{};
229 
230  //Loop over the spacepoints and get the charge weighted center.
231  for (auto const& sp : showersps) {
232 
233  //Get the associated hits
234  std::vector<art::Ptr<recob::Hit>> const& hits = fmh.at(sp.key());
235 
236  //Average the charge unless sepcified.
237  float charge = 0;
238  float charge2 = 0;
239  for (auto const& hit : hits) {
240 
241  if (fUseCollectionOnly) {
242  if (hit->SignalType() == geo::kCollection) {
243  charge = hit->Integral();
244  //Correct for the lifetime: Need to do other detproperites
245  charge *= std::exp((sampling_rate(clockData) * hit->PeakTime()) /
246  (detProp.ElectronLifetime() * 1e3));
247  break;
248  }
249  }
250  else {
251 
252  //Correct for the lifetime FIX: Need to do other detproperties somehow
253  double Q = hit->Integral() * std::exp((sampling_rate(clockData) * hit->PeakTime()) /
254  (detProp.ElectronLifetime() * 1e3));
255 
256  charge += Q;
257  charge2 += Q * Q;
258  }
259  }
260 
261  if (!fUseCollectionOnly) {
262  //Calculate the unbiased standard deviation and mean.
263  float mean = charge / ((float)hits.size());
264 
265  float rms = 1;
266 
267  if (hits.size() > 1) {
268  rms = std::sqrt((charge2 - charge * charge) / ((float)(hits.size() - 1)));
269  }
270 
271  charge = 0;
272  int n = 0;
273  for (auto const& hit : hits) {
274  double lifetimecorrection = std::exp((sampling_rate(clockData) * hit->PeakTime()) /
275  (detProp.ElectronLifetime() * 1e3));
276  if (hit->Integral() * lifetimecorrection > (mean - 2 * rms) &&
277  hit->Integral() * lifetimecorrection < (mean + 2 * rms)) {
278  charge += hit->Integral() * lifetimecorrection;
279  ++n;
280  }
281  }
282 
283  if (n == 0) {
284  mf::LogWarning("LArPandoraShowerAlg") << "no points used to make the charge value. \n";
285  }
286 
287  charge /= n;
288  }
289 
290  //Get the position of the spacepoint
291  auto const pos = sp->position();
292  double const x = pos.X();
293  double const y = pos.Y();
294  double const z = pos.Z();
295 
296  chargePoint.SetXYZ(
297  chargePoint.X() + charge * x, chargePoint.Y() + charge * y, chargePoint.Z() + charge * z);
298  totalCharge += charge;
299 
300  if (charge == 0) {
301  mf::LogWarning("LArPandoraShowerAlg") << "Averaged charge, within 2 sigma, for a spacepoint "
302  "is zero, Maybe this not a good method. \n";
303  }
304  }
305 
306  double intotalcharge = 1 / totalCharge;
307  return chargePoint * intotalcharge;
308 }
309 
311  art::Ptr<recob::SpacePoint> const& sp_a,
312  art::Ptr<recob::SpacePoint> const& sp_b) const
313 {
314  return (sp_a->position() - sp_b->position()).R();
315 }
316 
317 //Return the charge of the spacepoint in ADC.
319  art::FindManyP<recob::Hit> const& fmh) const
320 {
321  double Charge = 0;
322 
323  //Average over the charge even though there is only one
324  std::vector<art::Ptr<recob::Hit>> const& hits = fmh.at(sp.key());
325  for (auto const& hit : hits) {
326  Charge += hit->Integral();
327  }
328 
329  Charge /= (float)hits.size();
330 
331  return Charge;
332 }
333 
334 //Return the spacepoint time.
336  art::FindManyP<recob::Hit> const& fmh) const
337 {
338  double Time = 0;
339 
340  //Average over the hits
341  std::vector<art::Ptr<recob::Hit>> const& hits = fmh.at(sp.key());
342  for (auto const& hit : hits) {
343  Time += hit->PeakTime();
344  }
345 
346  Time /= (float)hits.size();
347  return Time;
348 }
349 
350 //Return the cooordinates of the hit in cm in wire direction and x.
352  art::Ptr<recob::Hit> const& hit) const
353 {
354 
355  //Get the pitch
356  const geo::WireID WireID = hit->WireID();
357  const geo::PlaneID planeid = WireID.asPlaneID();
358  double pitch = fGeom->WirePitch(planeid);
359 
360  return TVector2(WireID.Wire * pitch, detProp.ConvertTicksToX(hit->PeakTime(), planeid));
361 }
362 
364  geo::Point_t const& vertex,
365  geo::Vector_t const& direction) const
366 {
367  // Get the position of the spacepoint
368  auto const pos = sp->position() - vertex;
369 
370  // Get the the projected length
371  return pos.Dot(direction);
372 }
373 
375  geo::Point_t const& vertex,
376  geo::Vector_t const& direction) const
377 {
378  // Get the projection of the spacepoint
379  double proj = shower::LArPandoraShowerAlg::SpacePointProjection(sp, vertex, direction);
380 
381  return shower::LArPandoraShowerAlg::SpacePointPerpendicular(sp, vertex, direction, proj);
382 }
383 
385  geo::Point_t const& vertex,
386  geo::Vector_t const& direction,
387  double proj) const
388 {
389  // Get the position of the spacepoint
390  auto pos = sp->position() - vertex;
391 
392  // Take away the projection * distance to find the perpendicular vector
393  pos = pos - proj * direction;
394 
395  // Get the the projected length
396  return pos.R();
397 }
398 
399 //Function to calculate the RMS at segments of the shower and calculate the gradient of this. If negative then the direction is pointing the opposite way to the correct one
401  const geo::Point_t& ShowerCentre,
402  const geo::Vector_t& Direction,
403  const unsigned int nSegments) const
404 {
405  if (nSegments == 0)
406  throw cet::exception("LArPandoraShowerAlg")
407  << "Unable to calculate RMS Shower Gradient with 0 segments" << std::endl;
408 
409  if (sps.size() < 3) return 0;
410 
411  //Order the spacepoints
412  this->OrderShowerSpacePoints(sps, ShowerCentre, Direction);
413 
414  //Get the length of the shower.
415  const double minProj = this->SpacePointProjection(sps[0], ShowerCentre, Direction);
416  const double maxProj = this->SpacePointProjection(sps[sps.size() - 1], ShowerCentre, Direction);
417 
418  const double length = (maxProj - minProj);
419  const double segmentsize = length / nSegments;
420 
421  if (segmentsize < std::numeric_limits<double>::epsilon()) return 0;
422 
423  std::map<int, std::vector<float>> len_segment_map;
424 
425  //Split the the spacepoints into segments.
426  for (auto const& sp : sps) {
427 
428  //Get the the projected length
429  const double len = this->SpacePointProjection(sp, ShowerCentre, Direction);
430 
431  //Get the length to the projection
432  const double len_perp = this->SpacePointPerpendicular(sp, ShowerCentre, Direction, len);
433 
434  int sg_len = std::round(len / segmentsize);
435  len_segment_map[sg_len].push_back(len_perp);
436  }
437 
438  int counter = 0;
439  float sumx = 0.f;
440  float sumy = 0.f;
441  float sumx2 = 0.f;
442  float sumxy = 0.f;
443 
444  //Get the rms of the segments and caclulate the gradient.
445  for (auto const& segment : len_segment_map) {
446  if (segment.second.size() < 2) continue;
447  float RMS = this->CalculateRMS(segment.second);
448 
449  // Check if the calculation failed
450  if (RMS < 0) continue;
451 
452  //Calculate the gradient using regression
453  sumx += segment.first;
454  sumy += RMS;
455  sumx2 += segment.first * segment.first;
456  sumxy += RMS * segment.first;
457  ++counter;
458  }
459 
460  const float denom = counter * sumx2 - sumx * sumx;
461 
462  return std::abs(denom) < std::numeric_limits<float>::epsilon() ?
463  0 :
464  (counter * sumxy - sumx * sumy) / denom;
465 }
466 
467 double shower::LArPandoraShowerAlg::CalculateRMS(const std::vector<float>& perps) const
468 {
469 
470  if (perps.size() < 2) return std::numeric_limits<double>::lowest();
471 
472  float sum = 0;
473  for (const auto& perp : perps) {
474  sum += perp * perp;
475  }
476 
477  return std::sqrt(sum / (perps.size() - 1));
478 }
479 
481  geo::Point_t const& pos,
482  geo::Vector_t const& dir,
483  unsigned int const& TPC) const
484 {
485 
486  if (!fSCE || !fSCE->EnableCalSpatialSCE()) {
487  throw cet::exception("LArPandoraShowerALG")
488  << "Trying to correct SCE pitch when service is not configured" << std::endl;
489  }
490  // As the input pos is sce corrected already, find uncorrected pos
491  const geo::Point_t uncorrectedPos = pos + fSCE->GetPosOffsets(pos);
492  //Get the size of the correction at pos
493  const geo::Vector_t posOffset = fSCE->GetCalPosOffsets(uncorrectedPos, TPC);
494 
495  //Get the position of next hit
496  const geo::Point_t nextPos = uncorrectedPos + pitch * dir;
497  //Get the offsets at the next pos
498  const geo::Vector_t nextPosOffset = fSCE->GetCalPosOffsets(nextPos, TPC);
499 
500  //Calculate the corrected pitch
501  const int xFlip(fSCEXFlip ? -1 : 1);
502  geo::Vector_t pitchVec{pitch * dir.X() + xFlip * (nextPosOffset.X() - posOffset.X()),
503  pitch * dir.Y() + (nextPosOffset.Y() - posOffset.Y()),
504  pitch * dir.Z() + (nextPosOffset.Z() - posOffset.Z())};
505 
506  return pitchVec.r();
507 }
508 
510  geo::Point_t const& pos) const
511 {
512 
513  // Check the space charge service is properly configured
514  if (!fSCE || !fSCE->EnableSimEfieldSCE()) {
515  throw cet::exception("LArPandoraShowerALG")
516  << "Trying to correct SCE EField when service is not configured" << std::endl;
517  }
518  // Gets relative E field Distortions
519  geo::Vector_t EFieldOffsets = fSCE->GetEfieldOffsets(pos);
520  // Add 1 in X direction as this is the direction of the drift field
521  EFieldOffsets += geo::Vector_t{1, 0, 0};
522  // Convert to Absolute E Field from relative
523  EFieldOffsets *= EField;
524  // We only care about the magnitude for recombination
525  return EFieldOffsets.r();
526 }
527 
528 std::map<art::Ptr<recob::Hit>, std::vector<art::Ptr<recob::Hit>>>
530 {
531  // In this case, we need to only accept one hit in each snippet
532  // 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.
533  //
534  // If there are multiple valid hits on the same snippet, we need a way to pick the best one.
535  // (TODO: find a good way). The current method is to take the one with the highest charge integral.
536  typedef std::pair<unsigned, std::vector<unsigned>> OrganizedHits;
537  struct HitIdentifier {
538  int startTick;
539  int endTick;
540  int wire;
541  float integral;
542 
543  // construct
544  explicit HitIdentifier(const recob::Hit& hit)
545  : startTick(hit.StartTick())
546  , endTick(hit.EndTick())
547  , wire(hit.WireID().Wire)
548  , integral(hit.Integral())
549  {}
550 
551  // Defines whether two hits are on the same snippet
552  inline bool operator==(const HitIdentifier& rhs) const
553  {
554  return startTick == rhs.startTick && endTick == rhs.endTick && wire == rhs.wire;
555  }
556 
557  // Defines which hit to pick between two both on the same snippet
558  inline bool operator>(const HitIdentifier& rhs) const { return integral > rhs.integral; }
559  };
560 
561  unsigned nplanes = 6; // TODO FIXME
562  std::vector<std::vector<OrganizedHits>> hits_org(nplanes);
563  std::vector<std::vector<HitIdentifier>> hit_idents(nplanes);
564  for (unsigned i = 0; i < hits.size(); i++) {
565  HitIdentifier this_ident(*hits[i]);
566 
567  // check if we have found a hit on this snippet before
568  bool found_snippet = false;
569  for (unsigned j = 0; j < hits_org[hits[i]->WireID().Plane].size(); j++) {
570  if (this_ident == hit_idents[hits[i]->WireID().Plane][j]) {
571  found_snippet = true;
572  if (this_ident > hit_idents[hits[i]->WireID().Plane][j]) {
573  hits_org[hits[i]->WireID().Plane][j].second.push_back(
574  hits_org[hits[i]->WireID().Plane][j].first);
575  hits_org[hits[i]->WireID().Plane][j].first = i;
576  hit_idents[hits[i]->WireID().Plane][j] = this_ident;
577  }
578  else {
579  hits_org[hits[i]->WireID().Plane][j].second.push_back(i);
580  }
581  break;
582  }
583  }
584  if (!found_snippet) {
585  hits_org[hits[i]->WireID().Plane].push_back({i, {}});
586  hit_idents[hits[i]->WireID().Plane].push_back(this_ident);
587  }
588  }
589  std::map<art::Ptr<recob::Hit>, std::vector<art::Ptr<recob::Hit>>> ret;
590  for (auto const& planeHits : hits_org) {
591  for (const OrganizedHits& hit_org : planeHits) {
592  std::vector<art::Ptr<recob::Hit>> secondaryHits{};
593  for (const unsigned secondaryID : hit_org.second) {
594  secondaryHits.push_back(hits[secondaryID]);
595  }
596  ret[hits[hit_org.first]] = std::move(secondaryHits);
597  }
598  }
599  return ret;
600 }
601 
603  art::Event const& Event,
604  reco::shower::ShowerElementHolder const& ShowerEleHolder,
605  std::string const& evd_disp_name_append) const
606 {
607 
608  std::cout << "Making Debug Event Display" << std::endl;
609 
610  //Function for drawing reco showers to check direction and initial track selection
611 
612  // Get run info to make unique canvas names
613  int run = Event.run();
614  int subRun = Event.subRun();
615  int event = Event.event();
616  int PFPID = pfparticle->Self();
617 
618  // Create the canvas
619  TString canvasName = Form("canvas_%i_%i_%i_%i", run, subRun, event, PFPID);
620  if (evd_disp_name_append.length() > 0) canvasName += "_" + evd_disp_name_append;
621  TCanvas* canvas = tfs->make<TCanvas>(canvasName, canvasName);
622 
623  // Initialise variables
624  double x = 0;
625  double y = 0;
626  double z = 0;
627 
628  double x_min = std::numeric_limits<double>::max(), x_max = -std::numeric_limits<double>::max();
629  double y_min = std::numeric_limits<double>::max(), y_max = -std::numeric_limits<double>::max();
630  double z_min = std::numeric_limits<double>::max(), z_max = -std::numeric_limits<double>::max();
631 
632  // Get a bunch of associations (again)
633  // N.B. this is a horribly inefficient way of doing things but as this is only
634  // going to be used to debug I don't care, I would rather have generality in this case
635 
636  auto const pfpHandle = Event.getValidHandle<std::vector<recob::PFParticle>>(fPFParticleLabel);
637 
638  // Get the spacepoint - PFParticle assn
639  art::FindManyP<recob::SpacePoint> fmspp(pfpHandle, Event, fPFParticleLabel);
640  if (!fmspp.isValid()) {
641  throw cet::exception("LArPandoraShowerAlg") << "Trying to get the spacepoint and failed. Somet\
642  hing is not configured correctly. Stopping ";
643  }
644 
645  // Get the SpacePoints
646  std::vector<art::Ptr<recob::SpacePoint>> const& spacePoints = fmspp.at(pfparticle.key());
647 
648  //We cannot progress with no spacepoints.
649  if (spacePoints.empty()) { return; }
650 
651  // Get info from shower property holder
652  geo::Point_t showerStartPosition = {-999, -999, -999};
653  geo::Vector_t showerDirection = {-999, -999, -999};
654  std::vector<art::Ptr<recob::SpacePoint>> trackSpacePoints;
655 
656  //######################
657  //### Start Position ###
658  //######################
659  double startXYZ[3] = {-999, -999, -999};
660  if (!ShowerEleHolder.CheckElement(fShowerStartPositionInputLabel)) {
661  mf::LogError("LArPandoraShowerAlg") << "Start position not set, returning " << std::endl;
662  // return;
663  }
664  else {
665  ShowerEleHolder.GetElement(fShowerStartPositionInputLabel, showerStartPosition);
666  // Create 3D point at vertex, chosed to be origin for ease of use of display
667  startXYZ[0] = showerStartPosition.X();
668  startXYZ[1] = showerStartPosition.Y();
669  startXYZ[2] = showerStartPosition.Z();
670  }
671  auto startPoly = std::make_unique<TPolyMarker3D>(1, startXYZ);
672 
673  //########################
674  //### Shower Direction ###
675  //########################
676 
677  double xDirPoints[2] = {-999, -999};
678  double yDirPoints[2] = {-999, -999};
679  double zDirPoints[2] = {-999, -999};
680 
681  //initialise counter point
682  int point = 0;
683 
684  // Make 3D points for each spacepoint in the shower
685  auto allPoly = std::make_unique<TPolyMarker3D>(spacePoints.size());
686 
687  if (!ShowerEleHolder.CheckElement(fShowerDirectionInputLabel) &&
688  !ShowerEleHolder.CheckElement("ShowerStartPosition")) {
689  mf::LogError("LArPandoraShowerAlg") << "Direction not set, returning " << std::endl;
690  }
691  else {
692 
693  // Get the min and max projections along the direction to know how long to draw
694  ShowerEleHolder.GetElement(fShowerDirectionInputLabel, showerDirection);
695 
696  // the direction line
697  double minProj = std::numeric_limits<double>::max();
698  double maxProj = -std::numeric_limits<double>::max();
699 
700  //initialise counter point
701  int point = 0;
702 
703  for (auto spacePoint : spacePoints) {
704  auto const pos = spacePoint->position();
705  x = pos.X();
706  y = pos.Y();
707  z = pos.Z();
708  allPoly->SetPoint(point, x, y, z);
709  ++point;
710 
711  x_min = std::min(x, x_min);
712  x_max = std::max(x, x_max);
713  y_min = std::min(y, y_min);
714  y_max = std::max(y, y_max);
715  z_min = std::min(z, z_min);
716  z_max = std::max(z, z_max);
717 
718  // Calculate the projection of (point-startpoint) along the direction
720  spacePoint, showerStartPosition, showerDirection);
721  maxProj = std::max(proj, maxProj);
722  minProj = std::min(proj, minProj);
723  } // loop over spacepoints
724 
725  xDirPoints[0] = (showerStartPosition.X() + minProj * showerDirection.X());
726  xDirPoints[1] = (showerStartPosition.X() + maxProj * showerDirection.X());
727 
728  yDirPoints[0] = (showerStartPosition.Y() + minProj * showerDirection.Y());
729  yDirPoints[1] = (showerStartPosition.Y() + maxProj * showerDirection.Y());
730 
731  zDirPoints[0] = (showerStartPosition.Z() + minProj * showerDirection.Z());
732  zDirPoints[1] = (showerStartPosition.Z() + maxProj * showerDirection.Z());
733  }
734 
735  auto dirPoly = std::make_unique<TPolyLine3D>(2, xDirPoints, yDirPoints, zDirPoints);
736 
737  //#########################
738  //### Initial Track SPs ###
739  //#########################
740 
741  auto trackPoly = std::make_unique<TPolyMarker3D>(trackSpacePoints.size());
742  if (!ShowerEleHolder.CheckElement(fInitialTrackSpacePointsInputLabel)) {
743  mf::LogError("LArPandoraShowerAlg") << "TrackSpacePoints not set, returning " << std::endl;
744  // return;
745  }
746  else {
747  ShowerEleHolder.GetElement(fInitialTrackSpacePointsInputLabel, trackSpacePoints);
748  point = 0; // re-initialise counter
749  for (auto spacePoint : trackSpacePoints) {
750  auto const pos = spacePoint->position();
751  x = pos.X();
752  y = pos.Y();
753  z = pos.Z();
754  trackPoly->SetPoint(point, x, y, z);
755  ++point;
756  x_min = std::min(x, x_min);
757  x_max = std::max(x, x_max);
758  y_min = std::min(y, y_min);
759  y_max = std::max(y, y_max);
760  z_min = std::min(z, z_min);
761  z_max = std::max(z, z_max);
762  } // loop over track spacepoints
763  }
764 
765  //#########################
766  //### Other PFParticles ###
767  //#########################
768 
769  // we want to draw all of the PFParticles in the event
770  //Get the PFParticles
771  std::vector<art::Ptr<recob::PFParticle>> pfps;
772  art::fill_ptr_vector(pfps, pfpHandle);
773 
774  // initialse counters
775  // Split into tracks and showers to make it clearer what pandora is doing
776  int pfpTrackCounter = 0;
777  int pfpShowerCounter = 0;
778 
779  // initial loop over pfps to find nuber of spacepoints for tracks and showers
780  for (auto const& pfp : pfps) {
781  std::vector<art::Ptr<recob::SpacePoint>> const& sps = fmspp.at(pfp.key());
782  // If running pandora cheating it will call photons pdg 22
783  int pdg = abs(pfp->PdgCode()); // Track or shower
784  if (pdg == 11 || pdg == 22) { pfpShowerCounter += sps.size(); }
785  else {
786  pfpTrackCounter += sps.size();
787  }
788  }
789 
790  auto pfpPolyTrack = std::make_unique<TPolyMarker3D>(pfpTrackCounter);
791  auto pfpPolyShower = std::make_unique<TPolyMarker3D>(pfpShowerCounter);
792 
793  // initialise counters
794  int trackPoints = 0;
795  int showerPoints = 0;
796 
797  for (auto const& pfp : pfps) {
798  std::vector<art::Ptr<recob::SpacePoint>> const& sps = fmspp.at(pfp.key());
799  int pdg = abs(pfp->PdgCode()); // Track or shower
800  for (auto const& sp : sps) {
801  auto const pos = sp->position();
802  x = pos.X();
803  y = pos.Y();
804  z = pos.Z();
805  x_min = std::min(x, x_min);
806  x_max = std::max(x, x_max);
807  y_min = std::min(y, y_min);
808  y_max = std::max(y, y_max);
809  z_min = std::min(z, z_min);
810  z_max = std::max(z, z_max);
811 
812  // If running pandora cheating it will call photons pdg 22
813  if (pdg == 11 || pdg == 22) {
814  pfpPolyShower->SetPoint(showerPoints, x, y, z);
815  ++showerPoints;
816  }
817  else {
818  pfpPolyTrack->SetPoint(trackPoints, x, y, z);
819  ++trackPoints;
820  }
821  } // loop over sps
822 
823  //pfpPolyTrack->Draw();
824 
825  } // if (fDrawAllPFPs)
826 
827  //#################################
828  //### Initial Track Traj Points ###
829  //#################################
830 
831  auto TrackTrajPoly = std::make_unique<TPolyMarker3D>(TPolyMarker3D(1));
832  auto TrackInitTrajPoly = std::make_unique<TPolyMarker3D>(TPolyMarker3D(1));
833 
834  if (ShowerEleHolder.CheckElement(fInitialTrackInputLabel)) {
835 
836  //Get the track
837  recob::Track InitialTrack;
838  ShowerEleHolder.GetElement(fInitialTrackInputLabel, InitialTrack);
839 
840  if (InitialTrack.NumberTrajectoryPoints() != 0) {
841 
842  point = 0;
843  // Make 3D points for each trajectory point in the track stub
844  for (unsigned int traj = 0; traj < InitialTrack.NumberTrajectoryPoints(); ++traj) {
845 
846  //ignore bogus info.
847  auto flags = InitialTrack.FlagsAtPoint(traj);
848  if (flags.isSet(recob::TrajectoryPointFlagTraits::NoPoint)) { continue; }
849 
850  geo::Point_t TrajPositionPoint = InitialTrack.LocationAtPoint(traj);
851  x = TrajPositionPoint.X();
852  y = TrajPositionPoint.Y();
853  z = TrajPositionPoint.Z();
854  TrackTrajPoly->SetPoint(point, x, y, z);
855  ++point;
856  } // loop over trajectory points
857 
858  geo::Point_t TrajInitPositionPoint = InitialTrack.LocationAtPoint(0);
859  x = TrajInitPositionPoint.X();
860  y = TrajInitPositionPoint.Y();
861  z = TrajInitPositionPoint.Z();
862  TrackInitTrajPoly->SetPoint(TrackInitTrajPoly->GetN(), x, y, z);
863  }
864  }
865 
866  gStyle->SetOptStat(0);
867  TH3F axes("axes", "", 1, x_min, x_max, 1, y_min, y_max, 1, z_min, z_max);
868  axes.SetDirectory(0);
869  axes.GetXaxis()->SetTitle("X");
870  axes.GetYaxis()->SetTitle("Y");
871  axes.GetZaxis()->SetTitle("Z");
872  axes.Draw();
873 
874  // Draw all of the things
875  pfpPolyShower->SetMarkerStyle(20);
876  pfpPolyShower->SetMarkerColor(4);
877  pfpPolyShower->Draw();
878  pfpPolyTrack->SetMarkerStyle(20);
879  pfpPolyTrack->SetMarkerColor(6);
880  pfpPolyTrack->Draw();
881  allPoly->SetMarkerStyle(20);
882  allPoly->Draw();
883  trackPoly->SetMarkerStyle(20);
884  trackPoly->SetMarkerColor(2);
885  trackPoly->Draw();
886  startPoly->SetMarkerStyle(21);
887  startPoly->SetMarkerSize(0.5);
888  startPoly->SetMarkerColor(3);
889  startPoly->Draw();
890  dirPoly->SetLineWidth(1);
891  dirPoly->SetLineColor(6);
892  dirPoly->Draw();
893  TrackTrajPoly->SetMarkerStyle(22);
894  TrackTrajPoly->SetMarkerColor(7);
895  TrackTrajPoly->Draw();
896  TrackInitTrajPoly->SetMarkerStyle(22);
897  TrackInitTrajPoly->SetMarkerColor(4);
898  TrackInitTrajPoly->Draw();
899 
900  // Save the canvas. Don't usually need this when using TFileService but this in the alg
901  // not a module and didn't work without this so im going with it.
902  canvas->Write();
903 }
Float_t x
Definition: compare.C:6
virtual bool EnableSimEfieldSCE() const =0
SubRunNumber_t subRun() const
Definition: Event.cc:35
double SCECorrectPitch(double const &pitch, geo::Point_t const &pos, geo::Vector_t const &dir, unsigned int const &TPC) const
Length_t WireCoordinate(Point_t const &pos, PlaneID const &planeid) const
Returns the index of the nearest wire to the specified position.
double DistanceBetweenSpacePoints(art::Ptr< recob::SpacePoint > const &sp_a, art::Ptr< recob::SpacePoint > const &sp_b) const
void OrderShowerSpacePointsPerpendicular(std::vector< art::Ptr< recob::SpacePoint >> &showersps, geo::Point_t const &vertex, geo::Vector_t const &direction) const
Utilities related to art service access.
size_t Self() const
Returns the index of this particle.
Definition: PFParticle.h:88
void DebugEVD(art::Ptr< recob::PFParticle > const &pfparticle, art::Event const &Event, const reco::shower::ShowerElementHolder &ShowerEleHolder, std::string const &evd_disp_name_append="") const
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
static constexpr Flag_t NoPoint
The trajectory point is not defined.
double CalculateRMS(const std::vector< float > &perps) const
constexpr bool operator>(Interval< Q, Cat > const a, Quantity< Args... > const b) noexcept
Definition: intervals.h:507
T::provider_type const * providerFrom()
Returns a constant pointer to the provider of specified service.
Definition: ServiceUtil.h:75
virtual geo::Vector_t GetCalPosOffsets(geo::Point_t const &point, int const &TPCid) const =0
art::ServiceHandle< geo::Geometry const > fGeom
Declaration of signal hit object.
Float_t y
Definition: compare.C:6
Point_t const & LocationAtPoint(size_t i) const
Access to track position at different points.
Definition: Track.h:160
Double_t z
Definition: plot.C:276
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
geo::Point_t position() const
Returns the position of the point in world coordinates [cm].
Definition: SpacePoint.h:91
const std::string fInitialTrackInputLabel
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:136
constexpr auto abs(T v)
Returns the absolute value of the argument.
STL namespace.
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:244
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
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
double SpacePointTime(art::Ptr< recob::SpacePoint > const &sp, art::FindManyP< recob::Hit > const &fmh) const
double x_min
Definition: berger.C:15
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
PlaneGeo const & Plane(PlaneID const &planeid) const
Returns the specified wire.
LArPandoraShowerAlg(const fhicl::ParameterSet &pset)
void hits()
Definition: readHits.C:15
double SpacePointPerpendicular(art::Ptr< recob::SpacePoint > const &sp, geo::Point_t const &vertex, geo::Vector_t const &direction) const
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:295
virtual geo::Vector_t GetPosOffsets(geo::Point_t const &point) const =0
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
const std::string fInitialTrackSpacePointsInputLabel
key_type key() const noexcept
Definition: Ptr.h:166
std::map< art::Ptr< recob::Hit >, std::vector< art::Ptr< recob::Hit > > > OrganizeHits(const std::vector< art::Ptr< recob::Hit >> &hits) const
bool CheckElement(const std::string &Name) const
Provides recob::Track data product.
const std::string fShowerDirectionInputLabel
constexpr PlaneID const & asPlaneID() const
Conversion to PlaneID (for convenience of notation).
Definition: geo_types.h:520
EventNumber_t event() const
Definition: Event.cc:41
int GetElement(const std::string &Name, T &Element) const
raw::TDCtick_t StartTick() const
Initial tdc tick for hit.
Definition: Hit.h:212
void OrderShowerHits(detinfo::DetectorPropertiesData const &detProp, std::vector< art::Ptr< recob::Hit >> &hits, geo::Point_t const &ShowerPosition, geo::Vector_t const &ShowerDirection) const
double SpacePointProjection(art::Ptr< recob::SpacePoint > const &sp, geo::Point_t const &vertex, geo::Vector_t const &direction) const
const std::string fShowerStartPositionInputLabel
Detector simulation of raw signals on wires.
double ConvertTicksToX(double ticks, int p, int t, int c) const
raw::TDCtick_t EndTick() const
Final tdc tick for hit.
Definition: Hit.h:216
double x_max
Definition: berger.C:16
double SpacePointCharge(art::Ptr< recob::SpacePoint > const &sp, art::FindManyP< recob::Hit > const &fmh) const
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
art::ServiceHandle< art::TFileService > tfs
virtual geo::Vector_t GetEfieldOffsets(geo::Point_t const &point) const =0
double mean(const std::vector< short > &wf, size_t start, size_t nsample)
Definition: UtilFunc.cxx:13
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
virtual bool EnableCalSpatialSCE() const =0
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
LArSoft-specific namespace.
Contains all timing reference information for the detector.
TDirectory * dir
Definition: macro.C:5
geo::Point_t ShowerCentre(std::vector< art::Ptr< recob::SpacePoint >> const &showersps) const
TVector2 HitCoordinates(detinfo::DetectorPropertiesData const &detProp, art::Ptr< recob::Hit > const &hit) const
void OrderShowerSpacePoints(std::vector< art::Ptr< recob::SpacePoint >> &showersps, geo::Point_t const &vertex, geo::Vector_t const &direction) const
Float_t proj
Definition: plot.C:35
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
Vector_t const & GetIncreasingWireDirection() const
Returns the direction of increasing wires.
Definition: PlaneGeo.h:418
Definition: MVAAlg.h:12
PointFlags_t const & FlagsAtPoint(size_t i) const
Access to i-th TrajectoryPoint or its Flags.
Definition: Track.h:152
double RMSShowerGradient(std::vector< art::Ptr< recob::SpacePoint >> &sps, const geo::Point_t &ShowerCentre, const geo::Vector_t &Direction, const unsigned int nSegments) const
Char_t n[5]
Direction
Definition: types.h:12
double SCECorrectEField(double const &EField, geo::Point_t const &pos) const
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:46
constexpr PlaneID const & planeID() const
Definition: geo_types.h:620
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
Double_t sum
Definition: plot.C:31
RunNumber_t run() const
Definition: Event.cc:29
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
recob::tracking::Plane Plane
Definition: TrackState.h:17
spacecharge::SpaceCharge const * fSCE
bool operator==(infinite_endcount_iterator< T > const &, count_iterator< T > const &)
Definition: counter.h:278
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
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Event finding and building.
Signal from collection planes.
Definition: geo_types.h:152
vertex reconstruction