LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ShowerReco_module.cc
Go to the documentation of this file.
1 //
3 // \file ShowerReco_module.cc
4 //
5 // biagio.rossi@lhep.unibe.ch (FWMK : argoslope.resize(fNPlanes);neut specific)
6 // thomas.strauss@lhep.unibe.ch (ART : general detector)
7 //
8 // andrzej.szelc@yale.edu (port to detector agnostic version)
9 // jasaadi@fnal.gov (Try to rewrite the code to make it more readable and to be able to handle
10 // multiple TPC's and cryostats and more than one cluster per plane.)
11 //
12 // This algorithm is designed to reconstruct showers
13 //
15 
16 // ### Generic C++ includes ###
17 #include <cmath> // std::tan() ...
18 #include <string>
19 #include <vector>
20 
21 // ### Framework includes ###
27 #include "art_root_io/TFileService.h"
31 #include "fhiclcpp/ParameterSet.h"
33 
34 // ### ROOT includes ###
35 #include "TMath.h"
36 #include "TTree.h"
37 
38 // ### LArSoft includes ###
53 
54 #include "range/v3/view.hpp"
55 
56 #include <cmath>
57 
58 namespace shwf {
59 
60  class ShowerReco : public art::EDProducer {
61  public:
62  explicit ShowerReco(fhicl::ParameterSet const& pset);
63 
64  private:
65  void beginJob();
66  void beginRun(art::Run& run);
67  void produce(art::Event& evt);
70  unsigned int plane); // Get shower vertex and slopes.
71 
72  void LongTransEnergy(geo::GeometryCore const* geom,
73  detinfo::DetectorClocksData const& clockData,
74  detinfo::DetectorPropertiesData const& detProp,
75  unsigned int set,
76  std::vector<art::Ptr<recob::Hit>> hitlist); // Longitudinal
77 
78  void ClearandResizeVectors(unsigned int nPlanes);
79 
81  // input labels:
82  float slope[3]; // in cm, cm
83  float angle[3]; // in radians
84 
85  std::string fClusterModuleLabel;
86 
87  float ftimetick; // time sample in us
88 
89  double fMean_wire_pitch; // wire pitch in cm
91 
92  std::vector<double> fRMS_2cm;
93  std::vector<int> fNpoints_2cm;
94  std::vector<double> fCorr_MeV_2cm;
95  std::vector<double> fCorr_Charge_2cm;
96 
97  std::vector<int> fNpoints_corr_ADC_2cm;
98  std::vector<int> fNpoints_corr_MeV_2cm;
99 
100  std::vector<double> fTotChargeADC; //Total charge in ADC/cm for each plane
101  std::vector<double> fTotChargeMeV; //Total charge in MeV/cm for each plane
102  std::vector<double> fTotChargeMeV_MIPs; //Total charge in MeV/cm for each plane
103 
104  std::vector<double> fChargeADC_2cm; //Initial charge in ADC/cm for each plane first 4cm
105  std::vector<double> fChargeMeV_2cm; //initial charge in MeV/cm for each plane first 4cm
106 
107  std::vector<double> fChargeMeV_2cm_refined;
108  std::vector<double> fChargeMeV_2cm_axsum;
109 
110  std::vector<std::vector<double>> fDistribChargeADC; //vector with the first De/Dx points ADC
111  std::vector<std::vector<double>>
112  fDistribChargeMeV; //vector with the first De/Dx points converted energy
113  std::vector<std::vector<double>> fDistribHalfChargeMeV;
114  std::vector<std::vector<double>>
115  fDistribChargeposition; //vector with the first De/Dx points' positions
116 
117  std::vector<std::vector<double>> fSingleEvtAngle; //vector with the first De/Dx points
118  std::vector<std::vector<double>> fSingleEvtAngleVal; //vector with the first De/Dx points
119 
120  std::vector<unsigned int> fWire_vertex; // wire coordinate of vertex for each plane
121  std::vector<double> fTime_vertex; // time coordinate of vertex for each plane
122 
123  std::vector<double> fWire_vertexError; // wire coordinate of vertex for each plane
124  std::vector<double> fTime_vertexError; // time coordinate of vertex for each plane
125 
126  std::vector<unsigned int> fWire_last; // wire coordinate of last point for each plane
127  std::vector<double> fTime_last; // time coordinate of last point for each plane
128 
129  // reconstructed 3D position of shower start
130  std::vector<double> xyz_vertex_fit;
131 
132  std::vector<std::vector<double>>
133  fNPitch; // double array, to use each plane for each set of angles
134 
135  // calorimetry variables
136  float Kin_En;
137  std::vector<float> vdEdx;
138  std::vector<float> vresRange;
139  std::vector<float> vdQdx;
140  std::vector<float> deadwire; // residual range for dead wires
141  float Trk_Length;
142  float fTrkPitchC;
143  float fdEdxlength; // distance that gets used to determine e/gamma
144  // separation
145  float fcalodEdxlength; // cutoff distance for hits saved to the calo object.
146  bool fUseArea;
147 
148  double xphi, xtheta; // new calculated angles.
149  unsigned int fNPlanes; // number of planes
150  unsigned int fNAngles;
151  TTree* ftree_shwf;
152 
153  // conversion and useful constants
154  double fWirePitch; // wire pitch in cm
155  double fTimeTick;
158 
159  std::vector<int> fNhitsperplane;
160  std::vector<double> fTotADCperplane;
161  }; // class ShowerReco
162 
163  //------------------------------------------------------------------------------
165  {
166  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
167  fCaloPSet = pset.get<fhicl::ParameterSet>("CaloAlg");
168 
169  fdEdxlength =
170  pset.get<double>("dEdxlength"); // distance that gets used to determine e/gamma separation
172  pset.get<double>("calodEdxlength"); // cutoff distance for hits saved to the calo object.
173  fUseArea = pset.get<bool>("UseArea");
174 
175  produces<std::vector<recob::Shower>>();
176  produces<art::Assns<recob::Shower, recob::Cluster>>();
177  produces<art::Assns<recob::Shower, recob::Hit>>();
178  produces<std::vector<anab::Calorimetry>>();
179  produces<art::Assns<recob::Shower, anab::Calorimetry>>();
180  }
181 
182  // ***************** //
184  {
186 
190  fNPlanes = geo->Nplanes();
191  fMean_wire_pitch = geo->WirePitch(); // wire pitch in cm
192 
193  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
194  ftimetick = sampling_rate(clockData) / 1000.;
195 
198 
199  ftree_shwf = tfs->make<TTree>("ShowerReco",
200  "Results");
201  ftree_shwf->Branch("run", &fRun, "run/I");
202  ftree_shwf->Branch("subrun", &fSubRun, "subrun/I");
203  ftree_shwf->Branch("event", &fEvent, "event/I");
204  ftree_shwf->Branch("nplanes", &fNPlanes, "nplanes/I");
205  ftree_shwf->Branch("nangles", &fNAngles, "nangles/I");
206  ftree_shwf->Branch("xtheta", &xtheta, "xtheta/D");
207  ftree_shwf->Branch("xphi", &xphi, "xphi/D");
208  ftree_shwf->Branch("ftotChargeADC", "std::vector<double>", &fTotChargeADC);
209  ftree_shwf->Branch("ftotChargeMeV", "std::vector<double>", &fTotChargeMeV);
210  ftree_shwf->Branch("fTotChargeMeV_MIPs", "std::vector<double>", &fTotChargeMeV_MIPs);
211  ftree_shwf->Branch("NPitch", "std::vector< std::vector<double> >", &fNPitch);
212 
213  // this should be temporary - until the omega is sorted out.
214  ftree_shwf->Branch("RMS_2cm", "std::vector<double>", &fRMS_2cm);
215  ftree_shwf->Branch("Npoints_2cm", "std::vector<int>", &fNpoints_2cm);
216  ftree_shwf->Branch("ChargeADC_2cm", "std::vector<double>", &fChargeADC_2cm);
217  ftree_shwf->Branch("ChargeMeV_2cm", "std::vector<double>", &fChargeMeV_2cm);
218  ftree_shwf->Branch("ChargeMeV_2cm_refined", "std::vector<double>", &fChargeMeV_2cm_refined);
219  ftree_shwf->Branch("ChargeMeV_2cm_axsum", "std::vector<double>", &fChargeMeV_2cm_axsum);
220  ftree_shwf->Branch("fNhitsperplane", "std::vector<int>", &fNhitsperplane);
221  ftree_shwf->Branch("fTotADCperplane", "std::vector<double>", &fTotADCperplane);
222  ftree_shwf->Branch(
223  "ChargedistributionADC", "std::vector<std::vector<double>>", &fDistribChargeADC);
224  ftree_shwf->Branch(
225  "ChargedistributionMeV", "std::vector<std::vector<double>>", &fDistribChargeMeV);
226  ftree_shwf->Branch(
227  "DistribHalfChargeMeV", "std::vector<std::vector<double>>", &fDistribHalfChargeMeV);
228  ftree_shwf->Branch(
229  "ChargedistributionPosition", "std::vector<std::vector<double>>", &fDistribChargeposition);
230  ftree_shwf->Branch("xyz_vertex_fit", "std::vector<double>", &xyz_vertex_fit);
231  }
232 
234  {
235  auto const* geom = lar::providerFrom<geo::Geometry>();
236  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
237  auto const detProp =
239 
240  fWirePitch = geom->WirePitch(); // wire pitch in cm
241  fTimeTick = sampling_rate(clockData) / 1000.;
242  fDriftVelocity = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
244  }
245 
246  void ShowerReco::ShowerReco::ClearandResizeVectors(unsigned int /*nPlanes*/)
247  {
248  //calculate factorial for number of angles
249  int fact = 1;
250  for (unsigned int i = 1; i <= fNPlanes; ++i)
251  fact *= i;
252 
253  fNAngles = fact / 2;
254 
255  fDistribChargeADC.clear();
256  fDistribChargeMeV.clear();
257  fDistribChargeposition.clear();
258 
259  fDistribChargeADC.resize(fNPlanes);
260  fDistribChargeMeV.resize(fNPlanes);
262 
263  fNPitch.clear();
264  fDistribChargeADC.clear();
265  fDistribChargeMeV.clear();
266  fDistribHalfChargeMeV.clear();
267  fDistribChargeposition.clear();
268 
269  fNPitch.resize(fNAngles);
270  fDistribChargeADC.resize(fNPlanes);
271  fDistribChargeMeV.resize(fNPlanes);
274 
275  for (unsigned int ii = 0; ii < fNAngles; ii++) {
276  fNPitch[ii].resize(fNPlanes, -1);
277  }
278 
279  fNpoints_corr_ADC_2cm.clear();
280  fNpoints_corr_MeV_2cm.clear();
281 
282  fNpoints_corr_ADC_2cm.resize(fNAngles, -1);
283  fNpoints_corr_MeV_2cm.resize(fNAngles, -1);
284 
285  for (unsigned int ii = 0; ii < fNPlanes; ii++) {
286  fDistribChargeADC[ii].resize(0); //vector with the first De/Dx points
287  fDistribChargeMeV[ii].resize(0); //vector with the first De/Dx points
288  fDistribHalfChargeMeV[ii].resize(0);
289  fDistribChargeposition[ii].resize(0); //vector with the first De/Dx points' positions
290  }
291 
292  fWire_vertex.clear();
293  fTime_vertex.clear();
294  fWire_vertexError.clear();
295  fTime_vertexError.clear();
296  fWire_last.clear();
297  fTime_last.clear();
298  fTotChargeADC.clear();
299  fTotChargeMeV.clear();
300  fTotChargeMeV_MIPs.clear();
301  fRMS_2cm.clear();
302  fNpoints_2cm.clear();
303 
304  fNhitsperplane.clear();
305  fTotADCperplane.clear();
306 
307  fWire_vertex.resize(fNAngles, -1);
308  fTime_vertex.resize(fNAngles, -1);
309  fWire_vertexError.resize(fNPlanes, -1);
310  fTime_vertexError.resize(fNPlanes, -1);
311  fWire_last.resize(fNAngles, -1);
312  fTime_last.resize(fNAngles, -1);
313  fTotChargeADC.resize(fNAngles, 0);
314  fTotChargeMeV.resize(fNAngles, 0);
315  fTotChargeMeV_MIPs.resize(fNAngles, 0);
316  fRMS_2cm.resize(fNAngles, 0);
317  fNpoints_2cm.resize(fNAngles, 0);
318 
319  fCorr_MeV_2cm.clear();
320  fCorr_Charge_2cm.clear();
321  xyz_vertex_fit.clear();
322 
323  fChargeADC_2cm.clear();
324  fChargeMeV_2cm.clear();
325  fChargeMeV_2cm_refined.clear();
326  fChargeMeV_2cm_refined.clear();
327  fChargeMeV_2cm_axsum.clear();
328 
329  fCorr_MeV_2cm.resize(fNAngles, 0);
330  fCorr_Charge_2cm.resize(fNAngles, 0);
331  xyz_vertex_fit.resize(3);
332 
333  fChargeADC_2cm.resize(fNAngles,
334  0); //Initial charge in ADC/cm for each plane angle calculation 4cm
335  fChargeMeV_2cm.resize(fNAngles,
336  0); //initial charge in MeV/cm for each angle calculation first 4cm
337  fChargeMeV_2cm_refined.resize(0);
338  fChargeMeV_2cm_refined.resize(fNAngles, 0);
339  fChargeMeV_2cm_axsum.resize(fNAngles, 0);
340 
341  fNhitsperplane.resize(fNPlanes, -1);
342  fTotADCperplane.resize(fNPlanes, -1);
343 
344  vdEdx.clear();
345  vresRange.clear();
346  vdQdx.clear();
347  }
348 
349  // ***************** //
351  {
352  auto const* geom = lar::providerFrom<geo::Geometry>();
353  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
354  auto const detProp =
356 
357  util::GeometryUtilities const gser{*geom, clockData, detProp};
358  constexpr geo::TPCID tpcid{0, 0};
359  fNPlanes = geom->Nplanes(tpcid);
360  auto Shower3DVector = std::make_unique<std::vector<recob::Shower>>();
361  auto cassn = std::make_unique<art::Assns<recob::Shower, recob::Cluster>>();
362  auto hassn = std::make_unique<art::Assns<recob::Shower, recob::Hit>>();
363  auto calorimetrycol = std::make_unique<std::vector<anab::Calorimetry>>();
364  auto calassn = std::make_unique<art::Assns<anab::Calorimetry, recob::Shower>>();
365 
368  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
369  evt.getByLabel(fClusterModuleLabel, clusterListHandle);
370 
372  evt.getByLabel(fClusterModuleLabel, clusterAssociationHandle);
373 
374  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
375 
376  fRun = evt.id().run();
377  fSubRun = evt.id().subRun();
378  fEvent = evt.id().event();
379 
380  // find all the hits associated to all the clusters (once and for all);
381  // the index of the query matches the index of the cluster in the collection
382  // (conveniently carried around in its art pointer)
383  art::FindManyP<recob::Hit> ClusterHits(clusterListHandle, evt, fClusterModuleLabel);
384 
385  std::vector<art::PtrVector<recob::Cluster>>::const_iterator clusterSet =
386  clusterAssociationHandle->begin();
387  // loop over vector of vectors (each size of NPlanes) and reconstruct showers from each of those
388  for (size_t iClustSet = 0; iClustSet < clusterAssociationHandle->size(); iClustSet++) {
389 
390  const art::PtrVector<recob::Cluster> CurrentClusters = (*(clusterSet++));
391 
392  // do some error checking - i.e. are the clusters themselves present.
393  if (clusterListHandle->size() < 2 || CurrentClusters.size() < 2) {
394  ftree_shwf->Fill();
395  return;
396  }
397 
399 
400  std::vector<std::vector<art::Ptr<recob::Hit>>> hitlist_all;
401  hitlist_all.resize(fNPlanes);
402 
403  for (size_t iClust = 0; iClust < CurrentClusters.size(); iClust++) {
404  art::Ptr<recob::Cluster> const& pclust = CurrentClusters[iClust];
405 
406  // get all the hits for this cluster;
407  // pclust is a art::Ptr to the original cluster collection stored in the event;
408  // its key corresponds to its index in the collection
409  // (and therefore to the index in the query)
410  std::vector<art::Ptr<recob::Hit>> const& hitlist = ClusterHits.at(pclust.key());
411 
412  unsigned int p(0); //c=channel, p=plane, w=wire
413 
414  if (hitlist.size() == 0) continue;
415 
416  p = (*hitlist.begin())->WireID().Plane;
417  // get vertex position and slope information to start with - ii is the
418  // posistion of the correct cluster:
420 
421  double ADCcharge = 0;
422  //loop over cluster hits
423  for (art::Ptr<recob::Hit> const& hit : hitlist) {
424  p = hit->WireID().Plane;
425  hitlist_all[p].push_back(hit);
426  ADCcharge += hit->PeakAmplitude();
427  }
428  fNhitsperplane[p] = hitlist_all[p].size();
429  fTotADCperplane[p] = ADCcharge;
430  } // End loop on clusters.
431  // Now I have the Hitlists and the relevent clusters parameters saved.
432 
433  // find best set:
434  unsigned int bp1 = 0, bp2 = 0;
435  double minerror1 = 99999999, minerror2 = 9999999;
436  for (unsigned int ii = 0; ii < fNPlanes; ++ii) {
437  double locerror =
439  fTime_vertexError[ii] * fTime_vertexError[ii]; // time coordinate of vertex for each plane
440 
441  if (minerror1 >= locerror) // the >= sign is to favorize collection
442  {
443  minerror1 = locerror;
444  bp1 = ii;
445  }
446  }
447  for (unsigned int ij = 0; ij < fNPlanes; ++ij) {
448  double locerror =
450  fTime_vertexError[ij] * fTime_vertexError[ij]; // time coordinate of vertex for each plane
451 
452  if (minerror2 >= locerror && ij != bp1) {
453  minerror2 = locerror;
454  bp2 = ij;
455  }
456  }
457 
458  gser.Get3DaxisN(bp1, bp2, angle[bp1], angle[bp2], xphi, xtheta);
459 
461  std::vector<geo::Point_t> position;
462  position.reserve(fNPlanes);
463  // get starting positions for all planes -- FIXME: only position[0] is used.
464  for (auto const& plane : geom->Iterate<geo::PlaneGeo>(tpcid)) {
465  position.push_back(plane.GetBoxCenter());
466  }
467 
468  // Assuming there is no problem ( and we found the best pair that comes
469  // close in time ) we try to get the Y and Z coordinates for the start of
470  // the shower.
471  double fTimeTick = sampling_rate(clockData) / 1000.;
472  double fDriftVelocity = detProp.DriftVelocity(detProp.Efield(), detProp.Temperature());
473  try {
474  int chan1 = geom->PlaneWireToChannel({0, 0, bp1, fWire_vertex[bp1]});
475  int chan2 = geom->PlaneWireToChannel({0, 0, bp2, fWire_vertex[bp2]});
476 
477  double y, z;
478  geom->ChannelsIntersect(chan1, chan2, y, z);
479 
480  xyz_vertex_fit[1] = y;
481  xyz_vertex_fit[2] = z;
482  xyz_vertex_fit[0] =
483  (fTime_vertex[bp1] - trigger_offset(clockData)) * fDriftVelocity * fTimeTick +
484  position[0].X();
485  }
486  catch (cet::exception const& e) {
487  mf::LogWarning("ShowerReco") << "caught exception \n" << e;
488  xyz_vertex_fit[1] = 0;
489  xyz_vertex_fit[2] = 0;
490  xyz_vertex_fit[0] = 0;
491  }
492 
493  // if collection is not best plane, project starting point from that
494  if (bp1 != fNPlanes - 1 && bp2 != fNPlanes - 1) {
495  geo::PlaneID const lastPlaneID{0, 0, fNPlanes - 1};
496  auto pos = geom->Plane(lastPlaneID).GetBoxCenter();
497  pos.SetY(xyz_vertex_fit[1]);
498  pos.SetZ(xyz_vertex_fit[2]);
499  auto const wirevertex = geom->NearestWireID(pos, lastPlaneID).Wire;
500 
501  double drifttick =
502  (xyz_vertex_fit[0] / detProp.DriftVelocity(detProp.Efield(), detProp.Temperature())) *
503  (1. / fTimeTick);
504  fWire_vertex[fNPlanes - 1] = wirevertex; // wire coordinate of vertex for each plane
505  fTime_vertex[fNPlanes - 1] =
506  drifttick -
507  (pos.X() / detProp.DriftVelocity(detProp.Efield(), detProp.Temperature())) *
508  (1. / fTimeTick) +
509  trigger_offset(clockData);
510  }
511 
512  if (fabs(xphi) < 5.) {
513  xtheta = gser.Get3DSpecialCaseTheta(
514  bp1, bp2, fWire_last[bp1] - fWire_vertex[bp1], fWire_last[bp2] - fWire_vertex[bp2]);
515  }
516 
517  // zero the arrays just to make sure
518  for (unsigned int i = 0; i < fNAngles; ++i) {
519  fTotChargeADC[i] = 0;
520  fTotChargeMeV[i] = 0;
521  fTotChargeMeV_MIPs[i] = 0;
522  fNpoints_corr_ADC_2cm[i] = 0;
523  fNpoints_corr_MeV_2cm[i] = 0;
524 
525  fRMS_2cm[i] = 0;
526  fNpoints_2cm[i] = 0;
527 
528  fCorr_MeV_2cm[i] = 0;
529  fCorr_Charge_2cm[i] = 0;
530 
531  fChargeADC_2cm[i] = 0; //Initial charge in ADC/cm for each plane angle calculation 4cm
532  fChargeMeV_2cm[i] = 0; //initial charge in MeV/cm for each angle calculation first 4cm
533  }
534 
535  // do loop and choose Collection. With ful calorimetry can do all.
536  if (!(fabs(xphi) > 89 && fabs(xphi) < 91)) // do not calculate pitch for extreme angles
537  LongTransEnergy(geom,
538  clockData,
539  detProp,
540  0,
541  hitlist_all[fNPlanes - 1]); // temporary only plane 2.
542 
544 
545  // make an art::PtrVector of the clusters
547  for (unsigned int i = 0; i < clusterListHandle->size(); ++i) {
548  art::Ptr<recob::Cluster> prod(clusterListHandle, i);
549  prodvec.push_back(prod);
550  }
551 
552  //create a singleSpacePoint at vertex.
553  std::vector<recob::SpacePoint> spcpts;
554 
555  //get direction cosines and set them for the shower
556  // TBD determine which angle to use for the actual shower
557  double fPhi = xphi;
558  double fTheta = xtheta;
559 
560  TVector3 dcosVtx(std::cos(fPhi * TMath::Pi() / 180) * std::sin(fTheta * TMath::Pi() / 180),
561  std::cos(fTheta * TMath::Pi() / 180),
562  std::sin(fPhi * TMath::Pi() / 180) * std::sin(fTheta * TMath::Pi() / 180));
565  // fill with bogus values for now
566  TVector3 dcosVtxErr(util::kBogusD, util::kBogusD, util::kBogusD);
567  recob::Shower singShower;
568  singShower.set_direction(dcosVtx);
569  singShower.set_direction_err(dcosVtxErr);
570 
571  Shower3DVector->push_back(singShower);
572  // associate the shower with its clusters
573  util::CreateAssn(evt, *Shower3DVector, prodvec, *cassn);
574 
575  // get the hits associated with each cluster and associate those with the shower
576  for (size_t p = 0; p < prodvec.size(); ++p) {
577  std::vector<art::Ptr<recob::Hit>> hits = fmh.at(p);
578  util::CreateAssn(evt, *Shower3DVector, hits, *hassn);
579  }
580 
581  geo::PlaneID planeID(0, 0, fNPlanes - 1);
582  calorimetrycol->emplace_back(
584 
586 
587  art::ProductID aid = evt.getProductID<std::vector<recob::Shower>>();
588  art::Ptr<recob::Shower> aptr(aid, 0, evt.productGetter(aid));
589  ssvec.push_back(aptr);
590 
591  util::CreateAssn(evt, *calorimetrycol, ssvec, *calassn);
592  ftree_shwf->Fill();
593  } // end loop on Vectors of "Associated clusters"
594 
595  evt.put(std::move(Shower3DVector));
596  evt.put(std::move(cassn));
597  evt.put(std::move(hassn));
598  evt.put(std::move(calorimetrycol));
599  evt.put(std::move(calassn));
600  }
601 
602  //------------------------------------------------------------------------------
604  detinfo::DetectorClocksData const& clockData,
605  detinfo::DetectorPropertiesData const& detProp,
606  unsigned int set,
608  {
609  // alogorithm for energy vs dx of the shower (roto-translation) COLLECTION
610  // VIEW
611 
613 
614  double totCnrg = 0,
615  totCnrg_corr = 0; // tot enegry of the shower in collection
616 
617  double time;
618  unsigned int wire = 0, plane = fNPlanes - 1;
619 
620  double mevav2cm = 0.;
621  double npoints_calo = 0;
622 
623  int direction = -1;
624 
625  //override direction if phi (XZ angle) is less than 90 degrees
626  if (fabs(xphi) < 90) direction = 1;
627 
628  //variables to check whether a hit is close to the shower axis.
629  double ortdist, linedist;
630  double wire_on_line, time_on_line;
631 
632  //get effective pitch using 3D angles
633  util::GeometryUtilities const gser{*geom, clockData, detProp};
634  double newpitch = gser.PitchInView(plane, xphi, xtheta);
635 
636  using lar::to_element;
637  using ranges::views::transform;
638  for (auto const& hit : hitlist | transform(to_element)) {
639  time = hit.PeakTime();
640  wire = hit.WireID().Wire;
641  plane = hit.WireID().Plane;
642 
643  double dEdx_new;
644 
645  if (fUseArea) { dEdx_new = calalg.dEdx_AREA(clockData, detProp, hit, newpitch); }
646  else // this will hopefully go away, once all of the calibration factors
647  // are calculated.
648  {
649  dEdx_new = calalg.dEdx_AMP(clockData, detProp, hit, newpitch);
650  }
651 
652  //calculate total energy.
653  totCnrg_corr += dEdx_new;
654 
655  // calculate the wire,time coordinates of the hit projection on to the 2D shower axis
656  gser.GetPointOnLine(slope[plane] / fWireTimetoCmCm,
657  fWire_vertex[plane],
658  fTime_vertex[plane],
659  wire,
660  time,
661  wire_on_line,
662  time_on_line);
663  linedist =
664  gser.Get2DDistance(wire_on_line, time_on_line, fWire_vertex[plane], fTime_vertex[plane]);
665  ortdist = gser.Get2DDistance(wire_on_line, time_on_line, wire, time);
666 
667  //calculate the distance from the vertex using the effective pitch metric
668  double wdist = (((double)wire - (double)fWire_vertex[plane]) * newpitch) *
669  direction; //wdist is always positive
670 
671  if ((wdist < fcalodEdxlength) && (wdist > 0.2)) {
672 
673  vdEdx.push_back(dEdx_new);
674  vresRange.push_back(fabs(wdist));
675  vdQdx.push_back(hit.PeakAmplitude() / newpitch);
676  Trk_Length = wdist;
677  fTrkPitchC = fNPitch[set][plane];
678  Kin_En += dEdx_new * newpitch;
679  npoints_calo++;
680 
681  if (wdist < fdEdxlength &&
682  ((direction == 1 && wire > fWire_vertex[plane]) // take no hits before vertex
683  // (depending on direction)
684  || (direction == -1 && wire < fWire_vertex[plane])) &&
685  ortdist < 4.5 && linedist < fdEdxlength) {
686  fChargeMeV_2cm[set] += dEdx_new;
687  fNpoints_2cm[set]++;
688  }
689 
690  // fill out for 4cm preshower
691 
692  fDistribChargeMeV[set].push_back(dEdx_new); // vector with the first De/Dx points
693  fDistribChargeposition[set].push_back(
694  wdist); // vector with the first De/Dx points' positions
695 
696  } // end inside range if statement
697 
698  } // end first loop on hits.
699 
700  auto const signalType =
701  hitlist.empty() ? geo::kMysteryType : geom->SignalType(hitlist.front()->WireID());
702 
703  if (signalType == geo::kCollection) {
704  fTotChargeADC[set] = totCnrg * newpitch;
705  fTotChargeMeV[set] = totCnrg_corr * newpitch;
706  }
707 
708  // calculate average dE/dx
709  if (fNpoints_2cm[set] > 0) { mevav2cm = fChargeMeV_2cm[set] / fNpoints_2cm[set]; }
710 
711  // second loop to calculate RMS
712  for (auto const& hit : hitlist | transform(to_element)) {
713  time = hit.PeakTime();
714  wire = hit.WireID().Wire;
715  plane = hit.WireID().Plane;
716  double dEdx = 0;
717 
718  if (fUseArea) { dEdx = calalg.dEdx_AREA(clockData, detProp, hit, newpitch); }
719  else // this will hopefully go away, once all of the calibration factors
720  // are calculated.
721  {
722  dEdx = calalg.dEdx_AMP(clockData, detProp, hit, newpitch);
723  }
724 
725  gser.GetPointOnLine(slope[plane] / fWireTimetoCmCm,
726  fWire_vertex[plane],
727  fTime_vertex[plane],
728  wire,
729  time,
730  wire_on_line,
731  time_on_line);
732  linedist =
733  gser.Get2DDistance(wire_on_line, time_on_line, fWire_vertex[plane], fTime_vertex[plane]);
734  ortdist = gser.Get2DDistance(wire_on_line, time_on_line, wire, time);
735 
736  double wdist = (((double)wire - (double)fWire_vertex[plane]) * newpitch) * direction;
737 
738  if ((wdist < fcalodEdxlength) && (wdist > 0.2)) {
739  if (wdist < fdEdxlength &&
740  ((direction == 1 && wire > fWire_vertex[plane]) ||
741  (direction == -1 && wire < fWire_vertex[plane])) &&
742  ortdist < 4.5 && linedist < fdEdxlength) {
743  fRMS_2cm[set] += (dEdx - mevav2cm) * (dEdx - mevav2cm);
744  }
745 
746  } // end if on correct hits.
747  } // end RMS_calculating loop.
748 
749  if (fNpoints_2cm[set] > 0) { fRMS_2cm[set] = TMath::Sqrt(fRMS_2cm[set] / fNpoints_2cm[set]); }
750 
752 
753  for (auto const& hit : hitlist | transform(to_element)) {
754  time = hit.PeakTime();
755  wire = hit.WireID().Wire;
756  plane = hit.WireID().Plane;
757 
758  double dEdx = 0;
759  if (fUseArea) { dEdx = calalg.dEdx_AREA(clockData, detProp, hit, newpitch); }
760  else // this will hopefully go away, once all of the calibration factors
761  // are calculated.
762  {
763  dEdx = calalg.dEdx_AMP(clockData, detProp, hit, newpitch);
764  }
765 
766  gser.GetPointOnLine(slope[plane] / fWireTimetoCmCm,
767  fWire_vertex[plane],
768  fTime_vertex[plane],
769  wire,
770  time,
771  wire_on_line,
772  time_on_line);
773  linedist =
774  gser.Get2DDistance(wire_on_line, time_on_line, fWire_vertex[plane], fTime_vertex[plane]);
775  ortdist = gser.Get2DDistance(wire_on_line, time_on_line, wire, time);
776 
777  double wdist = (((double)wire - (double)fWire_vertex[plane]) * newpitch) * direction;
778 
779  if ((wdist < fcalodEdxlength) && (wdist > 0.2 &&
780  ((direction == 1 && wire > fWire_vertex[plane]) ||
781  (direction == -1 && wire < fWire_vertex[plane])) &&
782  ortdist < 4.5 && linedist < fdEdxlength)) {
783  if (wdist < fdEdxlength) {
784  if (((dEdx > (mevav2cm - fRMS_2cm[set])) && (dEdx < (mevav2cm + fRMS_2cm[set]))) ||
785  (newpitch > 0.3 * fdEdxlength)) {
786  fCorr_MeV_2cm[set] += dEdx;
787  fNpoints_corr_MeV_2cm[set]++;
788  }
789 
790  } // end if on good hits
791  }
792  } // end of third loop on hits
793 
794  if (fNpoints_corr_MeV_2cm[set] > 0) {
797  }
798  }
799 
800  //------------------------------------------------------------------------------------//
802  // Get shower vertex and slopes.
803  {
804  // convert to cm/cm units needed in the calculation
805  angle[plane] = clust->StartAngle();
806  slope[plane] = std::tan(clust->StartAngle());
807  fWire_vertex[plane] = clust->StartWire();
808  fTime_vertex[plane] = clust->StartTick();
809 
810  fWire_vertexError[plane] = clust->SigmaStartWire(); // wire coordinate of vertex for each plane
811  fTime_vertexError[plane] = clust->SigmaStartTick(); // time coordinate of vertex for each plane
812 
813  fWire_last[plane] = clust->EndWire(); // wire coordinate of last point for each plane
814  fTime_last[plane] = clust->EndTick();
815  }
816 
818 }
std::vector< double > fTime_last
fhicl::ParameterSet fCaloPSet
std::vector< double > fWire_vertexError
ShowerReco(fhicl::ParameterSet const &pset)
std::vector< double > fTotChargeADC
std::vector< double > fCorr_Charge_2cm
void set_direction_err(const TVector3 &dir_e)
Definition: Shower.h:135
Utilities related to art service access.
constexpr to_element_t to_element
Definition: ToElement.h:25
Who knows?
Definition: geo_types.h:153
std::vector< int > fNhitsperplane
ProductID getProductID(std::string const &instance_name="") const
std::vector< std::vector< double > > fNPitch
Declaration of signal hit object.
Float_t y
Definition: compare.C:6
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
Double_t z
Definition: plot.C:276
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
std::vector< std::vector< double > > fDistribChargeMeV
float StartWire() const
Returns the wire coordinate of the start of the cluster.
Definition: Cluster.h:276
std::vector< int > fNpoints_corr_MeV_2cm
std::vector< std::vector< double > > fSingleEvtAngle
shower finding
std::vector< double > fTime_vertex
float EndTick() const
Returns the tick coordinate of the end of the cluster.
Definition: Cluster.h:331
std::vector< double > fTotADCperplane
std::vector< unsigned int > fWire_last
std::vector< std::vector< double > > fDistribHalfChargeMeV
float StartAngle() const
Returns the starting angle of the cluster.
Definition: Cluster.h:461
std::vector< double > fChargeMeV_2cm_refined
std::vector< float > vdEdx
RunNumber_t run() const
Definition: EventID.h:98
void GetVertexAndAnglesFromCluster(art::Ptr< recob::Cluster > clust, unsigned int plane)
Definition: Run.h:37
std::vector< float > vdQdx
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
std::vector< int > fNpoints_2cm
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
std::vector< double > fTime_vertexError
std::vector< std::vector< double > > fSingleEvtAngleVal
void hits()
Definition: readHits.C:15
void set_direction(const TVector3 &dir)
Definition: Shower.h:134
EDProductGetter const * productGetter(ProductID const pid) const
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
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
std::vector< float > vresRange
key_type key() const noexcept
Definition: Ptr.h:166
T get(std::string const &key) const
Definition: ParameterSet.h:314
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:78
std::vector< double > fRMS_2cm
void produce(art::Event &evt)
float SigmaStartWire() const
Returns the uncertainty on wire coordinate of the start of the cluster.
Definition: Cluster.h:296
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
Description of geometry of one entire detector.
Definition: GeometryCore.h:119
Declaration of cluster object.
float dEdx(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const TCSlice &slc, TP3D &tp3d)
Definition: PFPUtils.cxx:2675
size_type size() const
Definition: PtrVector.h:302
Detector simulation of raw signals on wires.
std::vector< double > fCorr_MeV_2cm
std::string fClusterModuleLabel
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.
std::vector< double > xyz_vertex_fit
SigType_t SignalType(PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
void beginRun(art::Run &run)
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.
Contains all timing reference information for the detector.
std::vector< double > fChargeMeV_2cm
std::vector< double > fTotChargeMeV_MIPs
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
float SigmaStartTick() const
Returns the uncertainty on tick coordinate of the start of the cluster.
Definition: Cluster.h:305
std::vector< double > fChargeADC_2cm
double dEdx_AREA(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, recob::Hit const &hit, double pitch, double T0=0) const
std::vector< float > deadwire
int trigger_offset(DetectorClocksData const &data)
std::vector< std::vector< double > > fDistribChargeADC
EventNumber_t event() const
Definition: EventID.h:116
void LongTransEnergy(geo::GeometryCore const *geom, detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, unsigned int set, std::vector< art::Ptr< recob::Hit >> hitlist)
unsigned int Nplanes(TPCID const &tpcid=tpc_zero) const
Returns the total number of planes in the specified TPC.
Definition: GeometryCore.h:977
constexpr double kBogusD
obviously bogus double value
std::vector< double > fChargeMeV_2cm_axsum
void ClearandResizeVectors(unsigned int nPlanes)
TCEvent evt
Definition: DataStructs.cxx:8
Float_t e
Definition: plot.C:35
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
std::vector< int > fNpoints_corr_ADC_2cm
Float_t X
Definition: plot.C:37
Namespace collecting geometry-related classes utilities.
std::vector< std::vector< double > > fDistribChargeposition
float StartTick() const
Returns the tick coordinate of the start of the cluster.
Definition: Cluster.h:287
Length_t WirePitch(PlaneID const &planeid=plane_zero) const
Returns the distance between two consecutive wires.
std::vector< unsigned int > fWire_vertex
SubRunNumber_t subRun() const
Definition: EventID.h:110
EventID id() const
Definition: Event.cc:23
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
std::vector< double > fTotChargeMeV
Signal from collection planes.
Definition: geo_types.h:152
float EndWire() const
Returns the wire coordinate of the end of the cluster.
Definition: Cluster.h:318