LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
trkf::CosmicTrackerAlg Class Reference

#include "CosmicTrackerAlg.h"

Public Member Functions

 CosmicTrackerAlg (fhicl::ParameterSet const &pset)
 
void reconfigure (fhicl::ParameterSet const &pset)
 
void SPTReco (std::vector< art::Ptr< recob::Hit > > &fHits)
 

Public Attributes

std::vector< TVector3 > trajPos
 
std::vector< TVector3 > trajDir
 
std::vector< std::vector< art::Ptr< recob::Hit > > > trajHit
 
std::vector< TVector3 > trkPos
 
std::vector< TVector3 > trkDir
 

Private Member Functions

void TrackTrajectory (std::vector< art::Ptr< recob::Hit > > &fHits)
 
void Track3D (std::vector< art::Ptr< recob::Hit > > &fHits)
 
void MakeSPT (std::vector< art::Ptr< recob::Hit > > &fHits)
 

Private Attributes

int fSPTAlg
 
bool fTrajOnly
 
double ftmatch
 tolerance for time matching (in ticks) More...
 
double fsmatch
 tolerance for distance matching (in cm) More...
 
TrackTrajectoryAlg fTrackTrajectoryAlg
 
std::vector< std::vector< std::vector< std::vector< double > > > > vw
 
std::vector< std::vector< std::vector< std::vector< double > > > > vt
 
std::vector< std::vector< std::vector< std::vector< unsigned int > > > > vtraj
 
art::ServiceHandle< geo::Geometrygeom
 
const detinfo::LArPropertieslarprop
 
const detinfo::DetectorPropertiesdetprop
 

Detailed Description

Definition at line 27 of file CosmicTrackerAlg.h.

Constructor & Destructor Documentation

trkf::CosmicTrackerAlg::CosmicTrackerAlg ( fhicl::ParameterSet const &  pset)

Definition at line 45 of file CosmicTrackerAlg.cxx.

45  {
46  this->reconfigure(pset);
47  }
void reconfigure(fhicl::ParameterSet const &pset)

Member Function Documentation

void trkf::CosmicTrackerAlg::MakeSPT ( std::vector< art::Ptr< recob::Hit > > &  fHits)
private

Definition at line 480 of file CosmicTrackerAlg.cxx.

References geo::CryostatID::Cryostat, e, geo::PlaneID::Plane, geo::TPCID::TPC, geo::WireID::Wire, and geo::WireID::WireID().

480  {
481 
482  larprop = lar::providerFrom<detinfo::LArPropertiesService>();
483  detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
484 
485  double timetick = detprop->SamplingRate()*1e-3; //time sample in us
486  double Efield_drift = detprop->Efield(0); // Electric Field in the drift region in kV/cm
487  double Temperature = detprop->Temperature(); // LAr Temperature in K
488  double driftvelocity = detprop->DriftVelocity(Efield_drift,Temperature); //drift velocity in the drift region (cm/us)
489  double time_pitch = driftvelocity*timetick; //time sample (cm)
490 
491  for (size_t i = 0; i<fHits.size(); ++i){
492  int ip1 = -1;
493  int ip2 = -1;
494  int ip2p = -1;
495  int ih1 = -1;
496  int ih2 = -1;
497  int ih2p = -1;
498  double mindis1 = 1e9;
499  double mindis2 = 1e9;
500  double mindis2p = 1e9;
501  geo::WireID wireid = fHits[i]->WireID();
502  unsigned int wire = wireid.Wire;
503  unsigned int plane = wireid.Plane;
504  unsigned int tpc = wireid.TPC;
505  unsigned int cstat = wireid.Cryostat;
506  double wire_pitch = geom->WirePitch(plane,tpc,cstat); //wire pitch in cm
507  //find the two trajectory points that enclose the hit
508  //find the nearest trajectory point first
509  for (size_t j = 0; j<vw[cstat][tpc][plane].size(); ++j){
510  double dis = sqrt(pow(wire_pitch*(wire-vw[cstat][tpc][plane][j]),2)+
511  pow(time_pitch*(fHits[i]->PeakTime()-vt[cstat][tpc][plane][j]),2));
512  if (dis<mindis1){
513  mindis1 = dis;
514  ip1 = vtraj[cstat][tpc][plane][j];
515  ih1 = j;
516  }
517  }
518  //find the second nearest trajectory point, prefer the point on the other side
519  for (size_t j = 0; j<vw[cstat][tpc][plane].size(); ++j){
520  if (int(j)==ih1) continue;
521  double dis = sqrt(pow(wire_pitch*(wire-vw[cstat][tpc][plane][j]),2)+
522  pow(time_pitch*(fHits[i]->PeakTime()-vt[cstat][tpc][plane][j]),2));
523  if (dis<mindis2){
524  mindis2 = dis;
525  ip2 = vtraj[cstat][tpc][plane][j];
526  ih2 = j;
527  }
528  TVector3 v1(wire_pitch*(wire-vw[cstat][tpc][plane][ih1]),
529  time_pitch*(fHits[i]->PeakTime()-vt[cstat][tpc][plane][j]),
530  0);
531  TVector3 v2(wire_pitch*(wire-vw[cstat][tpc][plane][j]),
532  time_pitch*(fHits[i]->PeakTime()-vt[cstat][tpc][plane][j]),
533  0);
534  if (v1.Dot(v2)<0){
535  if (dis<mindis2p){
536  mindis2p = dis;
537  ip2p = vtraj[cstat][tpc][plane][j];
538  ih2p = j;
539  }
540  }
541  }
542  if (ip2p != -1){
543  ip2 = ip2p;
544  ih2 = ih2p;
545  }
546  //std::cout<<i<<" "<<ip1<<" "<<ip2<<std::endl;
547  if (ip1==-1||ip2==-1){
548  return;
549  throw cet::exception("CosmicTrackerAlg")<<"Cannot find two nearest trajectory points.";
550  }
551  TVector3 v1(wire_pitch*(wire-vw[cstat][tpc][plane][ih1]),
552  time_pitch*(fHits[i]->PeakTime()-vt[cstat][tpc][plane][ih1]),
553  0);
554  TVector3 v2(wire_pitch*(vw[cstat][tpc][plane][ih2]-vw[cstat][tpc][plane][ih1]),
555  time_pitch*(vt[cstat][tpc][plane][ih2]-vt[cstat][tpc][plane][ih1]),
556  0);
557  //std::cout<<v1.X()<<" "<<v1.Y()<<" "<<v2.X()<<" "<<v2.Y()<<" "<<v1.Mag()<<" "<<v2.Mag()<<std::endl;
558  if (!v2.Mag()){
559  throw cet::exception("CosmicTrackerAlg")<<"Two identical trajectory points.";
560  }
561  double ratio = v1.Dot(v2)/v2.Mag2();
562  TVector3 pos = trajPos[ip1]+ratio*(trajPos[ip2]-trajPos[ip1]);
563  trkPos.push_back(pos);
564  if (trajDir[ip1].Mag()){
565  trkDir.push_back(trajDir[ip1]);
566  }
567  else if (trajDir[ip2].Mag()){
568  trkDir.push_back(trajDir[ip2]);
569  }
570  else{//only got two trajectory points?
571  trkDir.push_back((trajPos[ip2]-trajPos[ip1]).Unit());
572  }
573  }//i
574  }
std::vector< std::vector< std::vector< std::vector< double > > > > vw
std::vector< std::vector< std::vector< std::vector< double > > > > vt
std::vector< TVector3 > trajPos
virtual double SamplingRate() const =0
Returns the period of the TPC readout electronics clock.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:130
std::vector< TVector3 > trkDir
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:313
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
std::vector< std::vector< std::vector< std::vector< unsigned int > > > > vtraj
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
const detinfo::LArProperties * larprop
std::vector< TVector3 > trkPos
art::ServiceHandle< geo::Geometry > geom
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
Float_t e
Definition: plot.C:34
std::vector< TVector3 > trajDir
const detinfo::DetectorProperties * detprop
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void trkf::CosmicTrackerAlg::reconfigure ( fhicl::ParameterSet const &  pset)

Definition at line 50 of file CosmicTrackerAlg.cxx.

References fhicl::ParameterSet::get().

50  {
51  fSPTAlg = pset.get< int >("SPTAlg", 0 );
52  fTrajOnly = pset.get< bool >("TrajOnly", false);
53  ftmatch = pset.get< double >("TMatch");
54  fsmatch = pset.get< double >("SMatch");
55  }
double ftmatch
tolerance for time matching (in ticks)
double fsmatch
tolerance for distance matching (in cm)
void trkf::CosmicTrackerAlg::SPTReco ( std::vector< art::Ptr< recob::Hit > > &  fHits)

Definition at line 58 of file CosmicTrackerAlg.cxx.

Referenced by trkf::CosmicTracker::produce().

58  {
59 
60  trajPos.clear();
61  trajDir.clear();
62  trajHit.clear();
63 
64  trkPos.clear();
65  trkDir.clear();
66 
67  if (fSPTAlg==0){
68  TrackTrajectory(fHits);
69  }
70  else if (fSPTAlg==1){
71  Track3D(fHits);
72  }
73  else{
74  throw cet::exception("CosmicTrackerAlg")<<"Unknown SPTAlg "<<fSPTAlg<<", needs to be 0 or 1";
75  }
76 
77  if (!fTrajOnly&&trajPos.size()) MakeSPT(fHits);
78  else{
79  trkPos = trajPos;
80  trkDir = trajDir;
81  }
82  }
std::vector< TVector3 > trajPos
std::vector< TVector3 > trkDir
void TrackTrajectory(std::vector< art::Ptr< recob::Hit > > &fHits)
std::vector< TVector3 > trkPos
void Track3D(std::vector< art::Ptr< recob::Hit > > &fHits)
void MakeSPT(std::vector< art::Ptr< recob::Hit > > &fHits)
std::vector< TVector3 > trajDir
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
std::vector< std::vector< art::Ptr< recob::Hit > > > trajHit
void trkf::CosmicTrackerAlg::Track3D ( std::vector< art::Ptr< recob::Hit > > &  fHits)
private

Definition at line 228 of file CosmicTrackerAlg.cxx.

References evd::details::begin(), c1, c2, e, evd::details::end(), max, min, geo::PlaneID::Plane, geo::sqr(), w, geo::WireID::Wire, geo::WireID::WireID(), geo::WireIDIntersection::y, and geo::WireIDIntersection::z.

228  {
229 
230  larprop = lar::providerFrom<detinfo::LArPropertiesService>();
231  detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
232 
233  //save time/hit information along track trajectory
234  std::vector<std::map<int,double> > vtimemap(3);
235  std::vector<std::map<int,art::Ptr<recob::Hit> > > vhitmap(3);
236 
237  //find hit on each wire along the fitted line
238  for (size_t ihit = 0; ihit<fHits.size(); ++ihit){//loop over hits
239  geo::WireID hitWireID = fHits[ihit]->WireID();
240  unsigned int w = hitWireID.Wire;
241  unsigned int pl = hitWireID.Plane;
242  double time = fHits[ihit]->PeakTime();
243  time -= detprop->GetXTicksOffset(fHits[ihit]->WireID().Plane,
244  fHits[ihit]->WireID().TPC,
245  fHits[ihit]->WireID().Cryostat);
246  vtimemap[pl][w] = time;
247  vhitmap[pl][w] = fHits[ihit];
248  }
249 
250  // Find two clusters with the most numbers of hits, and time ranges
251  int iclu1 = -1;
252  int iclu2 = -1;
253  int iclu3 = -1;
254  unsigned maxnumhits0 = 0;
255  unsigned maxnumhits1 = 0;
256 
257  std::vector<double> tmin(vtimemap.size());
258  std::vector<double> tmax(vtimemap.size());
259  for (size_t iclu = 0; iclu<vtimemap.size(); ++iclu){
260  tmin[iclu] = 1e9;
261  tmax[iclu] = -1e9;
262  }
263 
264  for (size_t iclu = 0; iclu<vtimemap.size(); ++iclu){
265  for (auto itime = vtimemap[iclu].begin(); itime!=vtimemap[iclu].end(); ++itime){
266  if (itime->second>tmax[iclu]){
267  tmax[iclu] = itime->second;
268  }
269  if (itime->second<tmin[iclu]){
270  tmin[iclu] = itime->second;
271  }
272  }
273  if (vtimemap[iclu].size()>maxnumhits0){
274  if (iclu1!=-1){
275  iclu2 = iclu1;
276  maxnumhits1 = maxnumhits0;
277  }
278  iclu1 = iclu;
279  maxnumhits0 = vtimemap[iclu].size();
280  }
281  else if (vtimemap[iclu].size()>maxnumhits1){
282  iclu2 = iclu;
283  maxnumhits1 = vtimemap[iclu].size();
284  }
285  }
286 
287  std::swap(iclu1,iclu2); //now iclu1 has fewer hits than iclu2
288 
289  //find iclu3
290  for (int iclu = 0; iclu<(int)vtimemap.size(); ++iclu){
291  if (iclu!=iclu1&&iclu!=iclu2) iclu3 = iclu;
292  }
293 
294  if (iclu1!=-1&&iclu2!=-1){//at least two good clusters
295  //select hits in a common time range
296  auto ihit = vhitmap[iclu1].begin();
297  auto itime = vtimemap[iclu1].begin();
298  while (itime!=vtimemap[iclu1].end()){
299  if (itime->second<std::max(tmin[iclu1],tmin[iclu2])-ftmatch||
300  itime->second>std::min(tmax[iclu1],tmax[iclu2])+ftmatch){
301  vtimemap[iclu1].erase(itime++);
302  vhitmap[iclu1].erase(ihit++);
303  }
304  else{
305  ++itime;
306  ++ihit;
307  }
308  }
309 
310  ihit = vhitmap[iclu2].begin();
311  itime = vtimemap[iclu2].begin();
312  while (itime!=vtimemap[iclu2].end()){
313  if (itime->second<std::max(tmin[iclu1],tmin[iclu2])-ftmatch||
314  itime->second>std::min(tmax[iclu1],tmax[iclu2])+ftmatch){
315  vtimemap[iclu2].erase(itime++);
316  vhitmap[iclu2].erase(ihit++);
317  }
318  else{
319  ++itime;
320  ++ihit;
321  }
322  }
323 
324  //if one cluster is empty, replace it with iclu3
325  if (!vtimemap[iclu1].size()){
326  if (iclu3!=-1){
327  std::swap(iclu3,iclu1);
328  }
329  }
330  if (!vtimemap[iclu2].size()){
331  if (iclu3!=-1){
332  std::swap(iclu3,iclu2);
333  std::swap(iclu1,iclu2);
334  }
335  }
336  if ((!vtimemap[iclu1].size())||(!vtimemap[iclu2].size())) return;
337 
338  bool rev = false;
339  auto times1 = vtimemap[iclu1].begin();
340  auto timee1 = vtimemap[iclu1].end();
341  --timee1;
342  auto times2 = vtimemap[iclu2].begin();
343  auto timee2 = vtimemap[iclu2].end();
344  --timee2;
345 
346  double ts1 = times1->second;
347  double te1 = timee1->second;
348  double ts2 = times2->second;
349  double te2 = timee2->second;
350 
351  //find out if we need to flip ends
352  if (std::abs(ts1-ts2)+std::abs(te1-te2)>std::abs(ts1-te2)+std::abs(te1-ts2)){
353  rev = true;
354  }
355 
356  double timetick = detprop->SamplingRate()*1e-3; //time sample in us
357  double Efield_drift = detprop->Efield(0); // Electric Field in the drift region in kV/cm
358  double Temperature = detprop->Temperature(); // LAr Temperature in K
359 
360  double driftvelocity = detprop->DriftVelocity(Efield_drift,Temperature); //drift velocity in the drift region (cm/us)
361  double timepitch = driftvelocity*timetick;
362 
363  double wire_pitch = geom->WirePitch(
364  vhitmap[0].begin()->second->WireID().Plane,
365  vhitmap[0].begin()->second->WireID().TPC,
366  vhitmap[0].begin()->second->WireID().Cryostat); //wire pitch in cm
367 
368  //find out 2d track length for all clusters associated with track candidate
369  std::vector<double> vtracklength;
370  for (size_t iclu = 0; iclu<vtimemap.size(); ++iclu){
371 
372  double tracklength = 0;
373  if (vtimemap[iclu].size()==1){
374  tracklength = wire_pitch;
375  }
376  else if (!vtimemap[iclu].empty()) {
377  std::map<int,double>::const_iterator iw = vtimemap[iclu].cbegin(),
378  wend = vtimemap[iclu].cend();
379  double t0 = iw->first, w0 = iw->second;
380  while (++iw != wend) {
381  tracklength += std::sqrt(sqr((iw->first-w0)*wire_pitch)+sqr((iw->second-t0)*timepitch));
382  w0 = iw->first;
383  t0 = iw->second;
384  } // while
385  }
386  vtracklength.push_back(tracklength);
387  }
388 
389  std::map<int,int> maxhitsMatch;
390 
391  auto ihit1 = vhitmap[iclu1].begin();
392  for (auto itime1 = vtimemap[iclu1].begin();
393  itime1 != vtimemap[iclu1].end();
394  ++itime1, ++ihit1){//loop over min-hits
395  std::vector<art::Ptr<recob::Hit>> sp_hits;
396  sp_hits.push_back(ihit1->second);
397  double hitcoord[3] = {0.,0.,0.};
398  double length1 = 0;
399  if (vtimemap[iclu1].size()==1){
400  length1 = wire_pitch;
401  }
402  else{
403  for (auto iw1 = vtimemap[iclu1].begin(); iw1!=itime1; ++iw1){
404  auto iw2 = iw1;
405  ++iw2;
406  length1 += std::sqrt(std::pow((iw1->first-iw2->first)*wire_pitch,2)
407  +std::pow((iw1->second-iw2->second)*timepitch,2));
408  }
409  }
410  double difference = 1e10; //distance between two matched hits
411  auto matchedtime = vtimemap[iclu2].end();
412  auto matchedhit = vhitmap[iclu2].end();
413 
414  auto ihit2 = vhitmap[iclu2].begin();
415  for (auto itime2 = vtimemap[iclu2].begin();
416  itime2!=vtimemap[iclu2].end();
417  ++itime2, ++ihit2){//loop over max-hits
418  if (maxhitsMatch[itime2->first]) continue;
419  double length2 = 0;
420  if (vtimemap[iclu2].size()==1){
421  length2 = wire_pitch;
422  }
423  else{
424  for (auto iw1 = vtimemap[iclu2].begin(); iw1!=itime2; ++iw1){
425  auto iw2 = iw1;
426  ++iw2;
427  length2 += std::sqrt(std::pow((iw1->first-iw2->first)*wire_pitch,2)+std::pow((iw1->second-iw2->second)*timepitch,2));
428  }
429  }
430  if (rev) length2 = vtracklength[iclu2] - length2;
431  length2 = vtracklength[iclu1]/vtracklength[iclu2]*length2;
432  bool timematch = std::abs(itime1->second-itime2->second)<ftmatch;
433  if (timematch &&std::abs(length2-length1)<difference){
434  difference = std::abs(length2-length1);
435  matchedtime = itime2;
436  matchedhit = ihit2;
437  }
438  }//loop over hits2
439  if (difference<fsmatch){
440  hitcoord[0] = detprop->ConvertTicksToX(matchedtime->second+
441  detprop->GetXTicksOffset((matchedhit->second)->WireID().Plane,
442  (matchedhit->second)->WireID().TPC,
443  (matchedhit->second)->WireID().Cryostat),
444  (matchedhit->second)->WireID().Plane,
445  (matchedhit->second)->WireID().TPC,
446  (matchedhit->second)->WireID().Cryostat);
447 
448  hitcoord[1] = -1e10;
449  hitcoord[2] = -1e10;
450 
451  //WireID is the exact segment of the wire where the hit is on (1 out of 3 for the 35t)
452 
453  geo::WireID c1=(ihit1->second)->WireID();
454  geo::WireID c2=(matchedhit->second)->WireID();
455 
456  geo::WireIDIntersection tmpWIDI;
457  bool sameTpcOrNot=geom->WireIDsIntersect(c1,c2, tmpWIDI);
458 
459  if(sameTpcOrNot){
460  hitcoord[1]=tmpWIDI.y;
461  hitcoord[2]=tmpWIDI.z;
462  }
463 
464  if (hitcoord[1]>-1e9&&hitcoord[2]>-1e9){
465  maxhitsMatch[matchedtime->first] = 1;
466  sp_hits.push_back(matchedhit->second);
467  }
468  }//if (difference<fsmatch)
469  if (sp_hits.size()>1){
470  trajPos.push_back(TVector3(hitcoord));
471  trajHit.push_back(sp_hits);
472  }
473  }//loop over hits1
474  }//if (iclu1!=-1&&iclu2!=-1){//at least two good clusters
475 
476  }
code to link reconstructed objects back to the MC truth information
double z
z position of intersection
Definition: geo_types.h:494
double ftmatch
tolerance for time matching (in ticks)
std::vector< TVector3 > trajPos
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:313
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
T sqr(T v)
Int_t max
Definition: plot.C:27
TCanvas * c1
Definition: plotHisto.C:7
intermediate_table::const_iterator const_iterator
TCanvas * c2
Definition: plot_hist.C:75
std::vector< evd::details::RawDigitInfo_t >::const_iterator begin(RawDigitCacheDataClass const &cache)
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
double fsmatch
tolerance for distance matching (in cm)
void swap(art::HLTGlobalStatus &lhs, art::HLTGlobalStatus &rhs)
bool WireIDsIntersect(WireID const &wid1, WireID const &wid2, geo::Point_t &intersection) const
Computes the intersection between two wires.
const detinfo::LArProperties * larprop
Int_t min
Definition: plot.C:26
double y
y position of intersection
Definition: geo_types.h:493
art::ServiceHandle< geo::Geometry > geom
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
Float_t e
Definition: plot.C:34
recob::tracking::Plane Plane
Definition: TrackState.h:17
Float_t w
Definition: plot.C:23
const detinfo::DetectorProperties * detprop
std::vector< std::vector< art::Ptr< recob::Hit > > > trajHit
void trkf::CosmicTrackerAlg::TrackTrajectory ( std::vector< art::Ptr< recob::Hit > > &  fHits)
private

Definition at line 85 of file CosmicTrackerAlg.cxx.

References evd::details::begin(), evd::details::end(), greaterThan1(), PlnLen::index, PlnLen::length, geo::TPCGeo::Nplanes(), X, Y, and Z.

85  {
86 
87  detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
88 
89 /*
90  // Track hit X and WireIDs in each plane
91  std::array<std::vector<std::pair<double, geo::WireID>>,3> trajXW;
92  // Track hit charge ...
93  std::array<std::vector<double>,3> trajChg;
94  std::array< std::vector<art::Ptr<recob::Hit>>,3> trajHits;
95  for (size_t i = 0; i<fHits.size(); ++i){
96  trajHits[fHits[i]->WireID().Plane].push_back(fHits[i]);
97  }
98 */
99  // Track hit vectors for fitting the trajectory
100  std::array<std::vector<geo::WireID>,3> trkWID;
101  std::array<std::vector<double>,3> trkX;
102  std::array<std::vector<double>,3> trkXErr;
103 
104  std::array< std::vector<art::Ptr<recob::Hit>>,3> trajHits;
105 
106  for (size_t i = 0; i<fHits.size(); ++i){
107  trajHits[fHits[i]->WireID().Plane].push_back(fHits[i]);
108  }
109 
110  std::vector<PlnLen> spl;
111  PlnLen plnlen;
112  for(unsigned int ipl = 0; ipl < 3; ++ipl) {
113  plnlen.index = ipl;
114  plnlen.length = trajHits[ipl].size();
115  //if(plnlen.length > 0)
116  spl.push_back(plnlen);
117  }
118  std::sort (spl.begin(),spl.end(), greaterThan1);
119  // spl[0] has the most hits and spl.back() has the least
120 
121  for (size_t ipl = 0; ipl <3; ++ipl){
122  if (!trajHits[ipl].size()) continue;
123  //if (ipl == spl[0].index) continue;
124  double xbeg0 = detprop->ConvertTicksToX(trajHits[spl[0].index][0]->PeakTime(),
125  trajHits[spl[0].index][0]->WireID().Plane,
126  trajHits[spl[0].index][0]->WireID().TPC,
127  trajHits[spl[0].index][0]->WireID().Cryostat);
128  double xend0 = detprop->ConvertTicksToX(trajHits[spl[0].index].back()->PeakTime(),
129  trajHits[spl[0].index].back()->WireID().Plane,
130  trajHits[spl[0].index].back()->WireID().TPC,
131  trajHits[spl[0].index].back()->WireID().Cryostat);
132  double xbeg1 = detprop->ConvertTicksToX(trajHits[ipl][0]->PeakTime(),
133  trajHits[ipl][0]->WireID().Plane,
134  trajHits[ipl][0]->WireID().TPC,
135  trajHits[ipl][0]->WireID().Cryostat);
136  double xend1 = detprop->ConvertTicksToX(trajHits[ipl].back()->PeakTime(),
137  trajHits[ipl].back()->WireID().Plane,
138  trajHits[ipl].back()->WireID().TPC,
139  trajHits[ipl].back()->WireID().Cryostat);
140  double dx1 = std::abs(xbeg0-xbeg1)+std::abs(xend0-xend1);
141  double dx2 = std::abs(xbeg0-xend1)+std::abs(xend0-xbeg1);
142  if (std::abs(xbeg1-xend1)>0.2 // this is to make sure the track is not completely flat in t, different by ~2.5 ticks
143  &&dx2<dx1){
144  std::reverse(trajHits[ipl].begin(),trajHits[ipl].end());
145  }
146  }
147 /*
148  // make the track trajectory
149  for(size_t ipl = 0; ipl < 3; ++ipl) {
150  //if (ipl == spl.back().index) continue;
151  trajXW[ipl].resize(trajHits[ipl].size());
152  trajChg[ipl].resize(trajHits[ipl].size());
153  for(size_t iht = 0; iht < trajHits[ipl].size(); ++iht) {
154  double xx = detprop->ConvertTicksToX(trajHits[ipl][iht]->PeakTime(),
155  ipl, trajHits[ipl][iht]->WireID().TPC,
156  trajHits[ipl][iht]->WireID().Cryostat);
157  trajXW[ipl][iht] = std::make_pair(xx, trajHits[ipl][iht]->WireID());
158  trajChg[ipl][iht] = trajHits[ipl][iht]->Integral();
159  } // iht
160  } // ip
161  fTrackTrajectoryAlg.TrackTrajectory(trajXW, trajPos, trajDir, trajChg);
162 */
163  // make the track trajectory
164  for(size_t ipl = 0; ipl < 3; ++ipl) {
165  //if (ipl == spl.back().index) continue;
166  trkWID[ipl].resize(trajHits[ipl].size());
167  trkX[ipl].resize(trajHits[ipl].size());
168  trkXErr[ipl].resize(trajHits[ipl].size());
169  for(size_t iht = 0; iht < trajHits[ipl].size(); ++iht) {
170  trkWID[ipl][iht] = trajHits[ipl][iht]->WireID();
171  trkX[ipl][iht] = detprop->ConvertTicksToX(trajHits[ipl][iht]->PeakTime(),ipl, trajHits[ipl][iht]->WireID().TPC, trajHits[ipl][iht]->WireID().Cryostat);
172  trkXErr[ipl][iht] = 0.2 * trajHits[ipl][iht]->RMS() * trajHits[ipl][iht]->Multiplicity();
173  } // iht
174  } // ip
175  fTrackTrajectoryAlg.TrackTrajectory(trkWID, trkX, trkXErr, trajPos, trajDir);
176  if (!trajPos.size()||!trajDir.size()){
177  mf::LogWarning("CosmicTrackerAlg")<<"Failed to reconstruct trajectory points.";
178  return;
179  }
180  //remove duplicated points
181  for (auto ipos = trajPos.begin(), idir = trajDir.begin(); ipos != trajPos.end() - 1; ){
182  auto ipos1 = ipos +1;
183  if (*ipos1 == *ipos){
184  ipos = trajPos.erase(ipos);
185  idir = trajDir.erase(idir);
186  }
187  else{
188  ++ipos;
189  ++idir;
190  }
191  }
192  vw.clear();
193  vt.clear();
194  vtraj.clear();
195  vw.resize(geom->Ncryostats());
196  vt.resize(geom->Ncryostats());
197  vtraj.resize(geom->Ncryostats());
198  for (size_t cstat = 0; cstat < geom->Ncryostats(); ++cstat){
199  vw[cstat].resize(geom->Cryostat(cstat).NTPC());
200  vt[cstat].resize(geom->Cryostat(cstat).NTPC());
201  vtraj[cstat].resize(geom->Cryostat(cstat).NTPC());
202  for (size_t tpc = 0; tpc < geom->Cryostat(cstat).NTPC(); ++tpc){
203  const geo::TPCGeo& tpcgeom = geom->Cryostat(cstat).TPC(tpc);
204  vw[cstat][tpc].resize(tpcgeom.Nplanes());
205  vt[cstat][tpc].resize(tpcgeom.Nplanes());
206  vtraj[cstat][tpc].resize(tpcgeom.Nplanes());
207  for (size_t plane = 0; plane < tpcgeom.Nplanes(); ++plane){
208  for (size_t i = 0; i< trajPos.size(); ++i){
209  double wirecord = geom->WireCoordinate(trajPos[i].Y(),
210  trajPos[i].Z(),
211  plane,tpc,cstat);
212  double tick = detprop->ConvertXToTicks(trajPos[i].X(),
213  plane,tpc,cstat);
214  //if (int(wirecord)>=0&&int(wirecord)<int(geom->Nwires(plane,tpc,cstat))
215  //&&int(tick)>=0&&int(tick)<int(detprop->ReadOutWindowSize())){
216  vw[cstat][tpc][plane].push_back(wirecord);
217  vt[cstat][tpc][plane].push_back(tick);
218  vtraj[cstat][tpc][plane].push_back(i);
219  //}//if wire and tick make sense
220  }//trajPos.size()
221  }//plane
222  }//tpc
223  }//cstat
224 
225  }
geo::Length_t WireCoordinate(double YPos, double ZPos, geo::PlaneID const &planeid) const
Returns the index of the nearest wire to the specified position.
std::vector< std::vector< std::vector< std::vector< double > > > > vw
std::vector< std::vector< std::vector< std::vector< double > > > > vt
unsigned int Nplanes() const
Number of planes in this tpc.
Definition: TPCGeo.h:145
Float_t Y
Definition: plot.C:39
Geometry information for a single TPC.
Definition: TPCGeo.h:37
std::vector< TVector3 > trajPos
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
virtual double ConvertXToTicks(double X, int p, int t, int c) const =0
Float_t Z
Definition: plot.C:39
std::vector< std::vector< std::vector< std::vector< unsigned int > > > > vtraj
TrackTrajectoryAlg fTrackTrajectoryAlg
std::vector< evd::details::RawDigitInfo_t >::const_iterator begin(RawDigitCacheDataClass const &cache)
unsigned int NTPC() const
Number of TPCs in this cryostat.
Definition: CryostatGeo.h:155
CryostatGeo const & Cryostat(geo::CryostatID const &cryoid) const
Returns the specified cryostat.
bool greaterThan1(PlnLen p1, PlnLen p2)
void TrackTrajectory(std::array< std::vector< geo::WireID >, 3 > trkWID, std::array< std::vector< double >, 3 > trkX, std::array< std::vector< double >, 3 > trkXErr, std::vector< TVector3 > &TrajPos, std::vector< TVector3 > &TrajDir)
unsigned short index
virtual double ConvertTicksToX(double ticks, int p, int t, int c) const =0
const TPCGeo & TPC(unsigned int itpc) const
Return the itpc&#39;th TPC in the cryostat.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
art::ServiceHandle< geo::Geometry > geom
std::vector< evd::details::RawDigitInfo_t >::const_iterator end(RawDigitCacheDataClass const &cache)
recob::tracking::Plane Plane
Definition: TrackState.h:17
std::vector< TVector3 > trajDir
Float_t X
Definition: plot.C:39
const detinfo::DetectorProperties * detprop

Member Data Documentation

const detinfo::DetectorProperties* trkf::CosmicTrackerAlg::detprop
private

Definition at line 75 of file CosmicTrackerAlg.h.

double trkf::CosmicTrackerAlg::fsmatch
private

tolerance for distance matching (in cm)

Definition at line 59 of file CosmicTrackerAlg.h.

int trkf::CosmicTrackerAlg::fSPTAlg
private

Definition at line 48 of file CosmicTrackerAlg.h.

double trkf::CosmicTrackerAlg::ftmatch
private

tolerance for time matching (in ticks)

Definition at line 58 of file CosmicTrackerAlg.h.

TrackTrajectoryAlg trkf::CosmicTrackerAlg::fTrackTrajectoryAlg
private

Definition at line 65 of file CosmicTrackerAlg.h.

bool trkf::CosmicTrackerAlg::fTrajOnly
private

Definition at line 51 of file CosmicTrackerAlg.h.

art::ServiceHandle<geo::Geometry> trkf::CosmicTrackerAlg::geom
private

Definition at line 73 of file CosmicTrackerAlg.h.

const detinfo::LArProperties* trkf::CosmicTrackerAlg::larprop
private

Definition at line 74 of file CosmicTrackerAlg.h.

std::vector<TVector3> trkf::CosmicTrackerAlg::trajDir

Definition at line 39 of file CosmicTrackerAlg.h.

std::vector<std::vector<art::Ptr<recob::Hit> > > trkf::CosmicTrackerAlg::trajHit

Definition at line 40 of file CosmicTrackerAlg.h.

Referenced by trkf::CosmicTracker::produce().

std::vector<TVector3> trkf::CosmicTrackerAlg::trajPos

Definition at line 38 of file CosmicTrackerAlg.h.

Referenced by trkf::CosmicTracker::produce().

std::vector<TVector3> trkf::CosmicTrackerAlg::trkDir

Definition at line 44 of file CosmicTrackerAlg.h.

Referenced by trkf::CosmicTracker::produce().

std::vector<TVector3> trkf::CosmicTrackerAlg::trkPos

Definition at line 43 of file CosmicTrackerAlg.h.

Referenced by trkf::CosmicTracker::produce().

std::vector<std::vector<std::vector<std::vector<double> > > > trkf::CosmicTrackerAlg::vt
private

Definition at line 69 of file CosmicTrackerAlg.h.

std::vector<std::vector<std::vector<std::vector<unsigned int> > > > trkf::CosmicTrackerAlg::vtraj
private

Definition at line 70 of file CosmicTrackerAlg.h.

std::vector<std::vector<std::vector<std::vector<double> > > > trkf::CosmicTrackerAlg::vw
private

Definition at line 68 of file CosmicTrackerAlg.h.


The documentation for this class was generated from the following files: