LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
VertexFinder2D_module.cc
Go to the documentation of this file.
1 //
3 // VertexFinder2D class
4 //
5 // tjyang@fnal.gov
6 //
7 // This algorithm is designed to reconstruct the vertices using the
8 // 2D cluster information
9 //
10 // This is Preliminary Work and needs modifications
11 // ////////////////////////////////////////////////////////////////////////
12 #include <algorithm>
13 #include <cmath>
14 #include <iomanip>
15 #include <string>
16 
17 // Framework includes
23 #include "art_root_io/TFileService.h"
27 #include "fhiclcpp/ParameterSet.h"
29 
41 
42 #include "TF1.h"
43 #include "TGraph.h"
44 #include "TH1D.h"
45 #include "TMath.h"
46 
47 namespace {
48  struct CluLen {
49  int index;
50  float length;
51  };
52 
53  bool myfunction(CluLen c1, CluLen c2)
54  {
55  return (c1.length > c2.length);
56  }
57 
58  struct SortByWire {
59  bool operator()(art::Ptr<recob::Hit> const& h1, art::Ptr<recob::Hit> const& h2) const
60  {
61  return h1->Channel() < h2->Channel();
62  }
63  };
64 }
65 
67 namespace vertex {
68 
70  public:
71  explicit VertexFinder2D(fhicl::ParameterSet const& pset);
72 
73  private:
74  void beginJob() override;
75  void produce(art::Event& evt) override;
76 
77  TH1D* dtIC;
78 
79  std::string fClusterModuleLabel;
80  };
81 
82 }
83 
84 namespace vertex {
85 
86  //-----------------------------------------------------------------------------
87  VertexFinder2D::VertexFinder2D(fhicl::ParameterSet const& pset) : EDProducer{pset}
88  {
89  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
90  produces<std::vector<recob::Vertex>>();
91  produces<std::vector<recob::EndPoint2D>>();
92  produces<art::Assns<recob::EndPoint2D, recob::Hit>>();
93  produces<art::Assns<recob::Vertex, recob::Hit>>();
94  produces<art::Assns<recob::Vertex, recob::Shower>>();
95  produces<art::Assns<recob::Vertex, recob::Track>>();
96  }
97 
98  //-------------------------------------------------------------------------
100  {
101  // get access to the TFile service
103  dtIC = tfs->make<TH1D>("dtIC", "It0-Ct0", 100, -5, 5);
104  dtIC->Sumw2();
105  }
106 
107  // //-----------------------------------------------------------------------------
109  {
111  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
112  auto const detProp =
114 
115  // define TPC parameters
116  TString tpcName = geom->GetLArTPCVolumeName();
117 
118  double YC = (geom->DetHalfHeight()) * 2.;
119 
120  // wire angle with respect to the vertical direction
121  double Angle = geom->Plane(geo::PlaneID{0, 0, 1}).Wire(0).ThetaZ(false) - TMath::Pi() / 2.;
122 
123  // Parameters temporary defined here, but possibly to be retrieved somewhere
124  // in the code
125  double timetick = sampling_rate(clockData) * 1.e-3; // time sample in us
126  double presamplings = trigger_offset(clockData);
127 
128  double wire_pitch = geom->WirePitch(); //wire pitch in cm
129  double Efield_drift = detProp.Efield(); // Electric Field in the drift region in kV/cm
130  double Temperature = detProp.Temperature(); // LAr Temperature in K
131 
132  //drift velocity in the drift region (cm/us)
133  double driftvelocity = detProp.DriftVelocity(Efield_drift, Temperature);
134 
135  //time sample (cm)
136  double timepitch = driftvelocity * timetick;
137 
138  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
139  evt.getByLabel(fClusterModuleLabel, clusterListHandle);
140 
142  for (unsigned int ii = 0; ii < clusterListHandle->size(); ++ii) {
143  art::Ptr<recob::Cluster> clusterHolder(clusterListHandle, ii);
144  clusters.push_back(clusterHolder);
145  }
146 
147  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
148 
149  //Point to a collection of vertices to output.
150  std::unique_ptr<std::vector<recob::Vertex>> vcol(new std::vector<recob::Vertex>); //3D vertex
151  std::unique_ptr<std::vector<recob::EndPoint2D>> epcol(
152  new std::vector<recob::EndPoint2D>); //2D vertex
153  std::unique_ptr<art::Assns<recob::EndPoint2D, recob::Hit>> assnep(
155  std::unique_ptr<art::Assns<recob::Vertex, recob::Shower>> assnsh(
157  std::unique_ptr<art::Assns<recob::Vertex, recob::Track>> assntr(
159  std::unique_ptr<art::Assns<recob::Vertex, recob::Hit>> assnh(
161 
162  // nplanes here is really being used as a proxy for the
163  // number of views in the detector
164  int nplanes = geom->Views().size();
165 
166  std::vector<std::vector<int>> Cls(nplanes); //index to clusters in each view
167  std::vector<std::vector<CluLen>> clulens(nplanes);
168 
169  std::vector<double> dtdwstart;
170 
171  //loop over clusters
172  for (size_t iclu = 0; iclu < clusters.size(); ++iclu) {
173 
174  float w0 = clusters[iclu]->StartWire();
175  float w1 = clusters[iclu]->EndWire();
176  float t0 = clusters[iclu]->StartTick();
177  float t1 = clusters[iclu]->EndTick();
178 
179  CluLen clulen;
180  clulen.index = iclu;
181  clulen.length = std::hypot((w0 - w1) * wire_pitch,
182  detProp.ConvertTicksToX(t0, clusters[iclu]->View(), 0, 0) -
183  detProp.ConvertTicksToX(t1, clusters[iclu]->View(), 0, 0));
184 
185  switch (clusters[iclu]->View()) {
186 
187  case geo::kU: clulens[0].push_back(clulen); break;
188  case geo::kV: clulens[1].push_back(clulen); break;
189  case geo::kZ: clulens[2].push_back(clulen); break;
190  default: break;
191  }
192 
193  std::vector<double> wires;
194  std::vector<double> times;
195 
196  std::vector<art::Ptr<recob::Hit>> hit = fmh.at(iclu);
197  std::sort(hit.begin(), hit.end(), SortByWire());
198  int n = 0;
199  for (size_t i = 0; i < hit.size(); ++i) {
200  wires.push_back(hit[i]->WireID().Wire);
201  times.push_back(hit[i]->PeakTime());
202  ++n;
203  }
204  if (n >= 2) {
205  TGraph* the2Dtrack = new TGraph(std::min(10, n), &wires[0], &times[0]);
206  try {
207  the2Dtrack->Fit("pol1", "Q");
208  TF1* pol1 = (TF1*)the2Dtrack->GetFunction("pol1");
209  double par[2];
210  pol1->GetParameters(par);
211  dtdwstart.push_back(par[1]);
212  }
213  catch (...) {
214  mf::LogWarning("VertexFinder2D") << "Fitter failed";
215  delete the2Dtrack;
216  dtdwstart.push_back(std::tan(clusters[iclu]->StartAngle()));
217  continue;
218  }
219  delete the2Dtrack;
220  }
221  else
222  dtdwstart.push_back(std::tan(clusters[iclu]->StartAngle()));
223  }
224 
225  //sort clusters based on 2D length
226  for (size_t i = 0; i < clulens.size(); ++i) {
227  std::sort(clulens[i].begin(), clulens[i].end(), myfunction);
228  for (size_t j = 0; j < clulens[i].size(); ++j) {
229  Cls[i].push_back(clulens[i][j].index);
230  }
231  }
232 
233  std::vector<std::vector<int>> cluvtx(nplanes);
234  std::vector<double> vtx_w;
235  std::vector<double> vtx_t;
236 
237  for (int i = 0; i < nplanes; ++i) {
238  if (Cls[i].size() >= 1) {
239  //at least one cluster
240  //find the longest two clusters
241  int c1 = -1;
242  int c2 = -1;
243  double ww0 = -999;
244  double wb1 = -999;
245  double we1 = -999;
246  double wb2 = -999;
247  double we2 = -999;
248  double tt1 = -999;
249  double tt2 = -999;
250  double dtdw1 = -999;
251  double dtdw2 = -999;
252  double lclu1 = -999;
253  double lclu2 = -999;
254  for (unsigned j = 0; j < Cls[i].size(); ++j) {
255  double lclu = std::sqrt(
256  pow((clusters[Cls[i][j]]->StartWire() - clusters[Cls[i][j]]->EndWire()) * 13.5, 2) +
257  pow(clusters[Cls[i][j]]->StartTick() - clusters[Cls[i][j]]->EndTick(), 2));
258  bool rev = false;
259  bool deltaraylike = false;
260  bool enoughhits = false;
261  if (c1 != -1) {
262  double wb = clusters[Cls[i][j]]->StartWire();
263  double we = clusters[Cls[i][j]]->EndWire();
264  double tt = clusters[Cls[i][j]]->StartTick();
265  double dtdw = dtdwstart[Cls[i][j]];
266  int nhits = fmh.at(Cls[i][j]).size();
267  ww0 = (tt - tt1 + dtdw1 * wb1 - dtdw * wb) / (dtdw1 - dtdw);
268  if (std::abs(wb1 - ww0) > std::abs(we1 - ww0)) rev = true; //reverse cluster dir
269  if ((!rev && ww0 > wb1 + 15) || (rev && ww0 < we1 - 15)) deltaraylike = true;
270  if (((!rev && ww0 > wb1 + 10) || (rev && ww0 < we1 - 10)) && nhits < 5)
271  deltaraylike = true;
272  if (wb > wb1 + 20 && nhits < 20) deltaraylike = true;
273  if (wb > wb1 + 50 && nhits < 20) deltaraylike = true;
274  if (wb > wb1 + 8 && TMath::Abs(dtdw1 - dtdw) < 0.15) deltaraylike = true;
275  if (std::abs(wb - wb1) > 30 && std::abs(we - we1) > 30) deltaraylike = true;
276  if (std::abs(tt - tt1) > 100)
277  deltaraylike = true; //not really deltaray, but isolated cluster
278  //make sure there are enough hits in the cluster
279  //at leaset 2 hits if goes horizentally, at leaset 4 hits if goes vertically
280  double alpha = std::atan(dtdw);
281  if (nhits >= int(2 + 3 * (1 - std::abs(std::cos(alpha))))) enoughhits = true;
282  if (nhits < 5 && (ww0 < wb1 - 20 || ww0 > we1 + 20)) enoughhits = false;
283  }
284  //do not replace the second cluster if the 3rd cluster is not consistent with the existing 2
285  bool replace = true;
286  if (c1 != -1 && c2 != -1) {
287  double wb = clusters[Cls[i][j]]->StartWire();
288  double we = clusters[Cls[i][j]]->EndWire();
289  ww0 = (tt2 - tt1 + dtdw1 * wb1 - dtdw2 * wb2) / (dtdw1 - dtdw2);
290  if ((std::abs(ww0 - wb1) < 10 || std::abs(ww0 - we1) < 10) &&
291  (std::abs(ww0 - wb2) < 10 || std::abs(ww0 - we2) < 10)) {
292  if (std::abs(ww0 - wb) > 15 && std::abs(ww0 - we) > 15) replace = false;
293  }
294  }
295  if (lclu1 < lclu) {
296  if (c1 != -1 && !deltaraylike && enoughhits) {
297  lclu2 = lclu1;
298  c2 = c1;
299  wb2 = wb1;
300  we2 = we1;
301  tt2 = tt1;
302  dtdw2 = dtdw1;
303  }
304  lclu1 = lclu;
305  c1 = Cls[i][j];
306  wb1 = clusters[Cls[i][j]]->StartWire();
307  we1 = clusters[Cls[i][j]]->EndWire();
308  tt1 = clusters[Cls[i][j]]->StartTick();
309  if (wb1 > we1) {
310  wb1 = clusters[Cls[i][j]]->EndWire();
311  we1 = clusters[Cls[i][j]]->StartWire();
312  tt1 = clusters[Cls[i][j]]->EndTick();
313  }
314  dtdw1 = dtdwstart[Cls[i][j]];
315  }
316  else if (lclu2 < lclu) {
317  if (!deltaraylike && enoughhits && replace) {
318  lclu2 = lclu;
319  c2 = Cls[i][j];
320  wb2 = clusters[Cls[i][j]]->StartWire();
321  we2 = clusters[Cls[i][j]]->EndWire();
322  tt2 = clusters[Cls[i][j]]->StartTick();
323  dtdw2 = dtdwstart[Cls[i][j]];
324  }
325  }
326  }
327  if (c1 != -1 && c2 != -1) {
328  cluvtx[i].push_back(c1);
329  cluvtx[i].push_back(c2);
330 
331  double w1 = clusters[c1]->StartWire();
332  double t1 = clusters[c1]->StartTick();
333  if (clusters[c1]->StartWire() > clusters[c1]->EndWire()) {
334  w1 = clusters[c1]->EndWire();
335  t1 = clusters[c1]->EndTick();
336  }
337  double k1 = dtdwstart[c1];
338  double w2 = clusters[c2]->StartWire();
339  double t2 = clusters[c2]->StartTick();
340  if (clusters[c2]->StartWire() > clusters[c2]->EndWire()) {
341  w1 = clusters[c2]->EndWire();
342  t1 = clusters[c2]->EndTick();
343  }
344  double k2 = dtdwstart[c2];
345  // calculate the vertex
346  if (std::abs(k1 - k2) < 0.5) {
347  vtx_w.push_back(w1);
348  vtx_t.push_back(t1);
349  }
350  else {
351  double t0 = (k1 * k2 * (w1 - w2) + k1 * t2 - k2 * t1) / (k1 - k2);
352  double w0 = (t2 - t1 + k1 * w1 - k2 * w2) / (k1 - k2);
353  vtx_w.push_back(w0);
354  vtx_t.push_back(t0);
355  }
356  }
357  else if (Cls[i].size() >= 1) {
358  if (c1 != -1) {
359  cluvtx[i].push_back(c1);
360  vtx_w.push_back(wb1);
361  vtx_t.push_back(tt1);
362  }
363  else {
364  cluvtx[i].push_back(Cls[i][0]);
365  vtx_w.push_back(clusters[Cls[i][0]]->StartWire());
366  vtx_t.push_back(clusters[Cls[i][0]]->StartTick());
367  }
368  }
369  //save 2D vertex
370  // make an empty art::PtrVector of hits
373  std::vector<art::Ptr<recob::Hit>> hits = fmh.at(Cls[i][0]);
374  double totalQ = 0.;
375  for (size_t h = 0; h < hits.size(); ++h)
376  totalQ += hits[h]->Integral();
377 
378  geo::WireID wireID(hits[0]->WireID().Cryostat,
379  hits[0]->WireID().TPC,
380  hits[0]->WireID().Plane,
381  (unsigned int)vtx_w.back()); //for update to EndPoint2D ... WK 4/22/13
382 
383  recob::EndPoint2D vertex(vtx_t.back(),
384  wireID, //for update to EndPoint2D ... WK 4/22/13
385  1,
386  epcol->size(),
387  clusters[Cls[i][0]]->View(),
388  totalQ);
389  epcol->push_back(vertex);
390 
391  util::CreateAssn(evt, *epcol, hits, *assnep);
392  }
393  else {
394  //no cluster found
395  vtx_w.push_back(-1);
396  vtx_t.push_back(-1);
397  }
398  }
399  //std::cout<<vtx_w[0]<<" "<<vtx_t[0]<<" "<<vtx_w[1]<<" "<<vtx_t[1]<<std::endl;
400 
401  Double_t vtxcoord[3];
402  if (Cls[0].size() > 0 && Cls[1].size() > 0) { //ignore w view
403  double Iw0 = (vtx_w[0] + 3.95) * wire_pitch;
404  double Cw0 = (vtx_w[1] + 1.84) * wire_pitch;
405 
406  double It0 = vtx_t[0] - presamplings;
407  It0 *= timepitch;
408  double Ct0 = vtx_t[1] - presamplings;
409  Ct0 *= timepitch;
410  vtxcoord[0] = detProp.ConvertTicksToX(vtx_t[1], 1, 0, 0);
411  vtxcoord[1] = (Cw0 - Iw0) / (2. * std::sin(Angle));
412  vtxcoord[2] = (Cw0 + Iw0) / (2. * std::cos(Angle)) - YC / 2. * std::tan(Angle);
413 
414  double yy, zz;
415  if (vtx_w[0] >= 0 && vtx_w[0] <= 239 && vtx_w[1] >= 0 && vtx_w[1] <= 239) {
416  if (geom->ChannelsIntersect(
417  geom->PlaneWireToChannel(geo::WireID(0, 0, 0, (int)((Iw0 / wire_pitch) - 3.95))),
418  geom->PlaneWireToChannel(geo::WireID(0, 0, 1, (int)((Cw0 / wire_pitch) - 1.84))),
419  yy,
420  zz)) {
421  //channelsintersect provides a slightly more accurate set of y and z coordinates.
422  // use channelsintersect in case the wires in question do cross.
423  vtxcoord[1] = yy;
424  vtxcoord[2] = zz;
425  }
426  else {
427  vtxcoord[0] = -99999;
428  vtxcoord[1] = -99999;
429  vtxcoord[2] = -99999;
430  }
431  }
432  dtIC->Fill(It0 - Ct0);
433  }
434  else {
435  vtxcoord[0] = -99999;
436  vtxcoord[1] = -99999;
437  vtxcoord[2] = -99999;
438  }
439 
442  art::PtrVector<recob::Track> vTracks_vec;
443  art::PtrVector<recob::Shower> vShowers_vec;
444 
445  recob::Vertex the3Dvertex(vtxcoord, vcol->size());
446  vcol->push_back(the3Dvertex);
447 
448  if (vShowers_vec.size() > 0) { util::CreateAssn(evt, *vcol, vShowers_vec, *assnsh); }
449 
450  if (vTracks_vec.size() > 0) { util::CreateAssn(evt, *vcol, vTracks_vec, *assntr); }
451 
452  MF_LOG_VERBATIM("Summary") << std::setfill('-') << std::setw(175) << "-" << std::setfill(' ');
453  MF_LOG_VERBATIM("Summary") << "VertexFinder2D Summary:";
454  for (size_t i = 0; i < epcol->size(); ++i)
455  MF_LOG_VERBATIM("Summary") << epcol->at(i);
456  for (size_t i = 0; i < vcol->size(); ++i)
457  MF_LOG_VERBATIM("Summary") << vcol->at(i);
458 
459  evt.put(std::move(epcol));
460  evt.put(std::move(vcol));
461  evt.put(std::move(assnep));
462  evt.put(std::move(assntr));
463  evt.put(std::move(assnsh));
464  evt.put(std::move(assnh));
465 
466  } // end of produce
467 } // end of vertex namespace
468 
469 // //-----------------------------------------------------------------------------
470 
471 namespace vertex {
472 
474 
475 } // end of vertex namespace
code to link reconstructed objects back to the MC truth information
std::string GetLArTPCVolumeName(TPCID const &tpcid=tpc_zero) const
Return the name of specified LAr TPC volume.
TTree * t1
Definition: plottest35.C:26
Planes which measure V.
Definition: geo_types.h:136
Declaration of signal hit object.
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
constexpr auto abs(T v)
Returns the absolute value of the argument.
bool ChannelsIntersect(raw::ChannelID_t c1, raw::ChannelID_t c2, double &y, double &z) const
Returns an intersection point of two channels.
Planes which measure Z direction.
Definition: geo_types.h:138
Double_t zz
Definition: plot.C:276
Definition of vertex object for LArSoft.
Definition: Vertex.h:35
double ThetaZ() const
Angle of the wires from positive z axis; .
Definition: PlaneGeo.cxx:646
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:101
Definition: type_traits.h:61
Planes which measure U.
Definition: geo_types.h:135
PlaneGeo const & Plane(PlaneID const &planeid) const
Returns the specified wire.
void hits()
Definition: readHits.C:15
bool SortByWire(art::Ptr< recob::Hit > const &h1, art::Ptr< recob::Hit > const &h2)
TCanvas * c1
Definition: plotHisto.C:7
TCanvas * c2
Definition: plot_hist.C:75
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
void beginJob()
Definition: Breakpoints.cc:14
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
void produce(art::Event &evt) override
Provides recob::Track data product.
Declaration of cluster object.
size_type size() const
Definition: PtrVector.h:302
Detector simulation of raw signals on wires.
TTree * t2
Definition: plottest35.C:36
TH1F * h2
Definition: plot.C:44
Encapsulate the geometry of a wire.
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.
raw::ChannelID_t PlaneWireToChannel(WireID const &wireid) const
Returns the ID of the TPC channel connected to the specified wire.
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.
#define MF_LOG_VERBATIM(category)
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
int trigger_offset(DetectorClocksData const &data)
std::set< View_t > const & Views() const
Returns a list of possible views in the detector.
TH1F * h1
Definition: plot.C:41
Char_t n[5]
decltype(auto) constexpr begin(T &&obj)
ADL-aware version of std::begin.
Definition: StdUtils.h:69
Length_t DetHalfHeight(TPCID const &tpcid=tpc_zero) const
Returns the half height of the active volume of the specified TPC.
TCEvent evt
Definition: DataStructs.cxx:8
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
recob::tracking::Plane Plane
Definition: TrackState.h:17
raw::ChannelID_t Channel() const
ID of the readout channel the hit was extracted from.
Definition: Hit.h:268
Length_t WirePitch(PlaneID const &planeid=plane_zero) const
Returns the distance between two consecutive wires.
art framework interface to geometry description
vertex reconstruction