All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Private Attributes | List of all members
sbnd::CRTTrackRecoAlg Class Reference

#include <CRTTrackRecoAlg.h>

Classes

struct  Config
 

Public Member Functions

 CRTTrackRecoAlg (const Config &config)
 
 CRTTrackRecoAlg (const fhicl::ParameterSet &pset)
 
 ~CRTTrackRecoAlg ()
 
void reconfigure (const Config &config)
 
std::vector< std::vector
< art::Ptr< sbn::crt::CRTHit > > > 
CreateCRTTzeros (std::vector< art::Ptr< sbn::crt::CRTHit >>)
 
sbn::crt::CRTTrack FillCrtTrack (sbn::crt::CRTHit hit1, sbn::crt::CRTHit hit2, bool complete)
 
sbn::crt::CRTTrack FillCrtTrack (sbn::crt::CRTHit hit1, sbn::crt::CRTHit hit2, size_t nhits)
 
std::vector< std::pair
< sbn::crt::CRTHit,
std::vector< int > > > 
AverageHits (std::vector< art::Ptr< sbn::crt::CRTHit >> hits, std::map< art::Ptr< sbn::crt::CRTHit >, int > hitIds)
 
std::vector< sbn::crt::CRTHitAverageHits (std::vector< art::Ptr< sbn::crt::CRTHit >> hits)
 
sbn::crt::CRTHit DoAverage (std::vector< art::Ptr< sbn::crt::CRTHit >> hits)
 
std::vector< std::pair
< sbn::crt::CRTTrack,
std::vector< int > > > 
CreateTracks (std::vector< std::pair< sbn::crt::CRTHit, std::vector< int >>> hits)
 
std::vector< sbn::crt::CRTTrackCreateTracks (std::vector< sbn::crt::CRTHit > hits)
 

Private Attributes

double fTimeLimit
 
double fAverageHitDistance
 
double fDistanceLimit
 
CRTHitRecoAlg hitAlg
 
CRTGeoAlg fCrtGeo
 

Detailed Description

Definition at line 59 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.h.

Constructor & Destructor Documentation

sbnd::CRTTrackRecoAlg::CRTTrackRecoAlg ( const Config config)

Definition at line 5 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

6  : hitAlg() {
7 
8  this->reconfigure(config);
9 
10 }
sbnd::CRTTrackRecoAlg::CRTTrackRecoAlg ( const fhicl::ParameterSet &  pset)
inline

Definition at line 85 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.h.

85  :
86  CRTTrackRecoAlg(fhicl::Table<Config>(pset, {})()) {}
sbnd::CRTTrackRecoAlg::~CRTTrackRecoAlg ( )

Definition at line 21 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

21  {
22 
23 }

Member Function Documentation

std::vector< std::pair< sbn::crt::CRTHit, std::vector< int > > > sbnd::CRTTrackRecoAlg::AverageHits ( std::vector< art::Ptr< sbn::crt::CRTHit >>  hits,
std::map< art::Ptr< sbn::crt::CRTHit >, int >  hitIds 
)

Definition at line 160 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

161 {
162 
163  std::vector<std::pair<sbn::crt::CRTHit, std::vector<int>>> returnHits;
164  std::vector<art::Ptr<sbn::crt::CRTHit>> aveHits;
165  std::vector<art::Ptr<sbn::crt::CRTHit>> spareHits;
166 
167  if (hits.size()>0){
168  // loop over size of tx
169  bool first = true;
170  TVector3 middle(0., 0., 0.);
171  for (size_t i = 0; i < hits.size(); i++){
172  // Get the position of the hit
173  TVector3 pos(hits[i]->x_pos, hits[i]->y_pos, hits[i]->z_pos);
174  // If first then set average = hit pos
175  if(first){
176  middle = pos;
177  first = false;
178  }
179  // If distance from average < limit then add to average
180  if((pos-middle).Mag() < fAverageHitDistance){
181  aveHits.push_back(hits[i]);
182  }
183  // Else add to another vector
184  else{
185  spareHits.push_back(hits[i]);
186  }
187  }
188 
189  sbn::crt::CRTHit aveHit = DoAverage(aveHits);
190  std::vector<int> ids;
191  for(size_t i = 0; i < aveHits.size(); i++){
192  ids.push_back(hitIds[aveHits[i]]);
193  }
194  returnHits.push_back(std::make_pair(aveHit, ids));
195 
196  //Do this recursively
197  std::vector<std::pair<sbn::crt::CRTHit, std::vector<int>>> moreHits = AverageHits(spareHits, hitIds);
198  returnHits.insert(returnHits.end(), moreHits.begin(), moreHits.end());
199  return returnHits;
200  }
201  else {
202  return returnHits;
203  }
204 
205 } // CRTTrackRecoAlg::AverageHits()
std::vector< std::pair< sbn::crt::CRTHit, std::vector< int > > > AverageHits(std::vector< art::Ptr< sbn::crt::CRTHit >> hits, std::map< art::Ptr< sbn::crt::CRTHit >, int > hitIds)
sbn::crt::CRTHit DoAverage(std::vector< art::Ptr< sbn::crt::CRTHit >> hits)
std::vector< sbn::crt::CRTHit > sbnd::CRTTrackRecoAlg::AverageHits ( std::vector< art::Ptr< sbn::crt::CRTHit >>  hits)

Definition at line 208 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

209 {
210 
211  std::map<art::Ptr<sbn::crt::CRTHit>, int> dummy;
212  for(size_t i = 0; i < hits.size(); i++){
213  dummy[hits[i]] = 0;
214  }
215 
216  std::vector<std::pair<sbn::crt::CRTHit, std::vector<int>>> output = AverageHits(hits, dummy);
217 
218  std::vector<sbn::crt::CRTHit> returnHits;
219  for(auto const& out : output){
220  returnHits.push_back(out.first);
221  }
222 
223  return returnHits;
224 
225 } // CRTTrackRecoAlg::AverageHits()
std::vector< std::pair< sbn::crt::CRTHit, std::vector< int > > > AverageHits(std::vector< art::Ptr< sbn::crt::CRTHit >> hits, std::map< art::Ptr< sbn::crt::CRTHit >, int > hitIds)
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG simSlidingORM6O6 effSlidingORW output
std::vector< std::vector< art::Ptr< sbn::crt::CRTHit > > > sbnd::CRTTrackRecoAlg::CreateCRTTzeros ( std::vector< art::Ptr< sbn::crt::CRTHit >>  hits)

Definition at line 37 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

38 {
39 
40  std::vector<std::vector<art::Ptr<sbn::crt::CRTHit>>> crtTzeroVect;
41  std::vector<int> iflag(hits.size(), 0);
42 
43  // Sort CRTHits by time
44  std::sort(hits.begin(), hits.end(), [](auto& left, auto& right)->bool{
45  return left->ts1_ns < right->ts1_ns;});
46 
47  // Loop over crt hits
48  for(size_t i = 0; i<hits.size(); i++){
49  if(iflag[i] == 0){
50  std::vector<art::Ptr<sbn::crt::CRTHit>> crtTzero;
51  double time_ns_A = hits[i]->ts1_ns;
52  iflag[i]=1;
53  crtTzero.push_back(hits[i]);
54 
55  // Sort into a Tzero collection
56  // Loop over all the other CRT hits
57  for(size_t j = i+1; j<hits.size(); j++){
58  if(iflag[j] == 0){
59  // If ts1_ns - ts1_ns < diff then put them in a vector
60  double time_ns_B = hits[j]->ts1_ns;
61  double diff = std::abs(time_ns_B - time_ns_A) * 1e-3; // [us]
62  if(diff < fTimeLimit){
63  iflag[j] = 1;
64  crtTzero.push_back(hits[j]);
65  }
66  }
67  }
68 
69  crtTzeroVect.push_back(crtTzero);
70  }
71  }
72  return crtTzeroVect;
73 }
walls no right
Definition: selectors.fcl:105
T abs(T value)
walls no left
Definition: selectors.fcl:105
do i e
std::vector< std::pair< sbn::crt::CRTTrack, std::vector< int > > > sbnd::CRTTrackRecoAlg::CreateTracks ( std::vector< std::pair< sbn::crt::CRTHit, std::vector< int >>>  hits)

Definition at line 273 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

274 {
275 
276  std::vector<std::pair<sbn::crt::CRTTrack, std::vector<int>>> returnTracks;
277 
278  std::vector<std::vector<size_t>> trackCandidates;
279  // Loop over all hits
280  for(size_t i = 0; i < hits.size(); i++){
281 
282  // Loop over all unique pairs
283  for(size_t j = i+1; j < hits.size(); j++){
284  if(hits[i].first.tagger == hits[j].first.tagger) continue;
285 
286  // Draw a track between the two hits
287  TVector3 start (hits[i].first.x_pos, hits[i].first.y_pos, hits[i].first.z_pos);
288  TVector3 end (hits[j].first.x_pos, hits[j].first.y_pos, hits[j].first.z_pos);
289  TVector3 diff = end - start;
290 
291  std::vector<size_t> candidate {i, j};
292 
293  // Loop over all other hits on different taggers and calculate DCA with variations
294  for(size_t k = 0; k < hits.size(); k++){
295  if(k == i || k == j || hits[k].first.tagger == hits[i].first.tagger
296  || hits[k].first.tagger == hits[j].first.tagger) continue;
297 
298  // If hit within certain distance then add it to the track candidate
300  candidate.push_back(k);
301  }
302  }
303  trackCandidates.push_back(candidate);
304  }
305  }
306 
307  // Sort track candidates by number of hits
308  std::sort(trackCandidates.begin(), trackCandidates.end(), [](auto& left, auto& right){
309  return left.size() > right.size();});
310 
311  // Loop over track candidates
312  std::vector<size_t> usedHits;
313  for(auto const& candidate : trackCandidates){
314  // Check if hits have been used
315  bool used = false;
316  for(size_t i = 0; i < candidate.size(); i++){
317  //Check if any of the hits have been used
318  if(std::find(usedHits.begin(), usedHits.end(), candidate[i]) != usedHits.end()) used = true;
319  }
320  if(used) continue;
321 
322  // Create track
323  if(candidate.size() < 2) continue;
324  sbn::crt::CRTHit ihit = hits[candidate[0]].first;
325  sbn::crt::CRTHit jhit = hits[candidate[1]].first;
326  sbn::crt::CRTTrack crtTrack = FillCrtTrack(ihit, jhit, candidate.size());
327 
328  std::vector<int> ids;
329  //TODO: Add charge matching for ambiguous cases
330  // If nhits > 2 then record used hits
331  for(size_t i = 0; i < candidate.size(); i++){
332  ids.insert(ids.end(), hits[candidate[i]].second.begin(), hits[candidate[i]].second.end());
333  if(candidate.size()>2) usedHits.push_back(candidate[i]);
334  }
335 
336  returnTracks.push_back(std::make_pair(crtTrack, ids));
337 
338  }
339 
340  return returnTracks;
341 
342 } // CRTTrackRecoAlg::CreateTracks()
walls no right
Definition: selectors.fcl:105
double DistToCrtHit(sbn::crt::CRTHit hit, TVector3 start, TVector3 end)
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
walls no left
Definition: selectors.fcl:105
sbn::crt::CRTTrack FillCrtTrack(sbn::crt::CRTHit hit1, sbn::crt::CRTHit hit2, bool complete)
pdgs k
Definition: selectors.fcl:22
std::vector< sbn::crt::CRTTrack > sbnd::CRTTrackRecoAlg::CreateTracks ( std::vector< sbn::crt::CRTHit hits)

Definition at line 345 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

346 {
347 
348  std::vector<std::pair<sbn::crt::CRTHit, std::vector<int>>> input;
349  for(auto const& hit : hits){
350  std::vector<int> dummy;
351  input.push_back(std::make_pair(hit, dummy));
352  }
353 
354  std::vector<std::pair<sbn::crt::CRTTrack, std::vector<int>>> output = CreateTracks(input);
355 
356  std::vector<sbn::crt::CRTTrack> tracks;
357  for(auto const& out : output){
358  tracks.push_back(out.first);
359  }
360 
361  return tracks;
362 
363 
364 } // CRTTrackRecoAlg::CreateTracks()
ClusterModuleLabel join with tracks
std::vector< std::pair< sbn::crt::CRTTrack, std::vector< int > > > CreateTracks(std::vector< std::pair< sbn::crt::CRTHit, std::vector< int >>> hits)
process_name hit
Definition: cheaterreco.fcl:51
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG simSlidingORM6O6 effSlidingORW output
sbn::crt::CRTHit sbnd::CRTTrackRecoAlg::DoAverage ( std::vector< art::Ptr< sbn::crt::CRTHit >>  hits)

Definition at line 230 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

231 {
232 
233  // Initialize values
234  std::string tagger = hits[0]->tagger;
235  double xpos = 0.;
236  double ypos = 0.;
237  double zpos = 0.;
238  double xmax = -99999; double xmin = 99999;
239  double ymax = -99999; double ymin = 99999;
240  double zmax = -99999; double zmin = 99999;
241  double ts1_ns = 0.;
242  int nhits = 0;
243 
244  // Loop over hits
245  for( auto& hit : hits ){
246  // Get the mean x,y,z and times
247  xpos += hit->x_pos;
248  ypos += hit->y_pos;
249  zpos += hit->z_pos;
250  ts1_ns += (double)(int)hit->ts1_ns;
251  // For the errors get the maximum limits
252  if(hit->x_pos + hit->x_err > xmax) xmax = hit->x_pos + hit->x_err;
253  if(hit->x_pos - hit->x_err < xmin) xmin = hit->x_pos - hit->x_err;
254  if(hit->y_pos + hit->y_err > ymax) ymax = hit->y_pos + hit->y_err;
255  if(hit->y_pos - hit->y_err < ymin) ymin = hit->y_pos - hit->y_err;
256  if(hit->z_pos + hit->z_err > zmax) zmax = hit->z_pos + hit->z_err;
257  if(hit->z_pos - hit->z_err < zmin) zmin = hit->z_pos - hit->z_err;
258  // Add all the unique IDs in the vector
259  nhits++;
260  }
261 
262  // Create a hit
263  sbn::crt::CRTHit crtHit = hitAlg.FillCrtHit(hits[0]->feb_id, hits[0]->pesmap, hits[0]->peshit,
264  (ts1_ns/nhits)*1e-3, 0, xpos/nhits, (xmax-xmin)/2,
265  ypos/nhits, (ymax-ymin)/2., zpos/nhits, (zmax-zmin)/2., tagger);
266 
267  return crtHit;
268 
269 } // CRTTrackRecoAlg::DoAverage()
process_name hit
Definition: cheaterreco.fcl:51
process_name pandoraGausCryo1 vertexChargeCryo1 vertexStubCryo1 xmin
sbn::crt::CRTHit FillCrtHit(std::vector< uint8_t > tfeb_id, std::map< uint8_t, std::vector< std::pair< int, float >>> tpesmap, float peshit, double time, int plane, double x, double ex, double y, double ey, double z, double ez, std::string tagger)
do i e
sbn::crt::CRTTrack sbnd::CRTTrackRecoAlg::FillCrtTrack ( sbn::crt::CRTHit  hit1,
sbn::crt::CRTHit  hit2,
bool  complete 
)

Definition at line 77 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

78 {
79 
80  sbn::crt::CRTTrack newtr;
81 
82  newtr.ts0_s = (hit1.ts0_s + hit2.ts0_s)/2.;
83  newtr.ts0_s_err = (hit1.ts0_s - hit2.ts0_s)/2.;
84  newtr.ts0_ns_h1 = hit1.ts0_ns;
85  newtr.ts0_ns_err_h1 = hit1.ts0_ns_corr;
86  newtr.ts0_ns_h2 = hit2.ts0_ns;
87  newtr.ts0_ns_err_h2 = hit2.ts0_ns_corr;
88  newtr.ts0_ns = (hit1.ts0_ns + hit2.ts0_ns)/2.;
89  newtr.ts0_ns_err = sqrt(hit1.ts0_ns_corr*hit1.ts0_ns_corr + hit2.ts0_ns_corr*hit2.ts0_ns_corr)/2.;
90  newtr.ts1_ns = (hit1.ts1_ns + hit2.ts1_ns)/2.;
91  newtr.ts1_ns_err = sqrt(hit1.ts0_ns_corr*hit1.ts0_ns_corr + hit2.ts0_ns_corr*hit2.ts0_ns_corr)/2.;
92  newtr.peshit = hit1.peshit+hit2.peshit;
93  newtr.x1_pos = hit1.x_pos;
94  newtr.x1_err = hit1.x_err;
95  newtr.y1_pos = hit1.y_pos;
96  newtr.y1_err = hit1.y_err;
97  newtr.z1_pos = hit1.z_pos;
98  newtr.z1_err = hit1.z_err;
99  newtr.x2_pos = hit2.x_pos;
100  newtr.x2_err = hit2.x_err;
101  newtr.y2_pos = hit2.y_pos;
102  newtr.y2_err = hit2.y_err;
103  newtr.z2_pos = hit2.z_pos;
104  newtr.z2_err = hit2.z_err;
105  float deltax = hit1.x_pos - hit2.x_pos;
106  float deltay = hit1.y_pos - hit2.y_pos;
107  float deltaz = hit1.z_pos - hit2.z_pos;
108  newtr.length = sqrt(deltax*deltax + deltay*deltay+deltaz*deltaz);
109  newtr.thetaxy = atan2(deltax,deltay);
110  newtr.phizy = atan2(deltaz,deltay);
111  newtr.plane1 = hit1.plane;
112  newtr.plane2 = hit2.plane;
113  newtr.complete = complete;
114 
115  return(newtr);
116 
117 } // CRTTrackRecoAlg::FillCrtTrack()
float z_err
position uncertainty in z-direction (cm).
Definition: CRTHit.hh:43
float x_err
position uncertainty in x-direction (cm).
Definition: CRTHit.hh:39
double ts0_ns_err_h2
T0 time error of second CRTHit.
Definition: CRTTrack.hh:51
double ts1_ns_err
Error on average T1 (nanosecond) of the two hits making the track.
Definition: CRTTrack.hh:29
double ts0_s
Average time (second) of the two hits making the track.
Definition: CRTTrack.hh:24
double ts1_ns
Timestamp T1 ([signal time w.r.t. Trigger time]), in UTC absolute time scale in nanoseconds from the ...
Definition: CRTHit.hh:34
double ts0_ns_err
Error on average T0 (nanosecond) of the two hits making the track.
Definition: CRTTrack.hh:27
int plane
Name of the CRT wall (in the form of numbers).
Definition: CRTHit.hh:36
float peshit
Total photo-electron (PE) in a crt hit.
Definition: CRTHit.hh:27
float y_err
position uncertainty in y-direction (cm).
Definition: CRTHit.hh:41
double ts0_s_err
Average time (second) spread of the two hits making the track.
Definition: CRTTrack.hh:25
float x1_pos
X position of first CRTHit.
Definition: CRTTrack.hh:33
float y1_err
Y position error of first CRTHit.
Definition: CRTTrack.hh:36
double ts0_ns_corr
[Honestly, not sure at this point, it was there since long time (BB)]
Definition: CRTHit.hh:33
double ts0_ns_h1
T0 time of first CRTHit.
Definition: CRTTrack.hh:48
uint64_t ts0_s
Second-only part of timestamp T0.
Definition: CRTHit.hh:29
float z_pos
position in z-direction (cm).
Definition: CRTHit.hh:42
double ts0_ns
Timestamp T0 (from White Rabbit), in UTC absolute time scale in nanoseconds from the Epoch...
Definition: CRTHit.hh:32
float z1_err
Z position error of first CRTHit.
Definition: CRTTrack.hh:38
float y2_err
Y position error of second CRTHit.
Definition: CRTTrack.hh:42
double ts1_ns
Average T1 (nanosecond) of the two hits making the track.
Definition: CRTTrack.hh:28
float length
Track length.
Definition: CRTTrack.hh:45
double ts0_ns
Average T0 (nanosecond) of the two hits making the track.
Definition: CRTTrack.hh:26
float peshit
Total photoelectrons for this track (sum of PEs from the two CRTHits)
Definition: CRTTrack.hh:23
double ts0_ns_err_h1
T0 time error of first CRTHit.
Definition: CRTTrack.hh:49
bool complete
Whether or not the track is complete.
Definition: CRTTrack.hh:53
float z1_pos
Z position of first CRTHit.
Definition: CRTTrack.hh:37
double ts0_ns_h2
T0 time of second CRTHit.
Definition: CRTTrack.hh:50
float y_pos
position in y-direction (cm).
Definition: CRTHit.hh:40
float phizy
Track angle on the Z-Y plane.
Definition: CRTTrack.hh:47
int plane1
Plane ID of first CRTHit.
Definition: CRTTrack.hh:30
float z2_pos
Z position of second CRTHit.
Definition: CRTTrack.hh:43
float x_pos
position in x-direction (cm).
Definition: CRTHit.hh:38
float y2_pos
Y position of second CRTHit.
Definition: CRTTrack.hh:41
float y1_pos
Y position of first CRTHit.
Definition: CRTTrack.hh:35
float x2_pos
X position of second CRTHit.
Definition: CRTTrack.hh:39
float x1_err
X position error of first CRTHit.
Definition: CRTTrack.hh:34
int plane2
Plane ID of second CRTHit.
Definition: CRTTrack.hh:31
float z2_err
Z position error of second CRTHit.
Definition: CRTTrack.hh:44
float x2_err
X position error of second CRTHit.
Definition: CRTTrack.hh:40
float thetaxy
Track angle on the X-Y plane.
Definition: CRTTrack.hh:46
sbn::crt::CRTTrack sbnd::CRTTrackRecoAlg::FillCrtTrack ( sbn::crt::CRTHit  hit1,
sbn::crt::CRTHit  hit2,
size_t  nhits 
)

Definition at line 120 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

121 {
122 
123  bool complete = true;
124  if(nhits == 2){
125  // Track is incomplete if just between 2 top planes
126  if((hit1.tagger == "volTaggerTopHigh_0" && hit2.tagger == "volTaggerTopLow_0")
127  || (hit2.tagger == "volTaggerTopHigh_0" && hit1.tagger == "volTaggerTopLow_0")) complete = false;
128  return FillCrtTrack(hit1, hit2, complete);
129  }
130  else{
131  // Project track on to the limits of the CRT volume TODO errors
132  std::vector<double> crtLimits = fCrtGeo.CRTLimits();
133  TVector3 min (crtLimits[0], crtLimits[1], crtLimits[2]);
134  TVector3 max (crtLimits[3], crtLimits[4], crtLimits[5]);
135  TVector3 start (hit1.x_pos, hit1.y_pos, hit1.z_pos);
136  TVector3 end (hit2.x_pos, hit2.y_pos, hit2.z_pos);
137 
138  std::pair<TVector3, TVector3> intersection = CRTCommonUtils::CubeIntersection(min, max, start, end);
139  if(intersection.first.X() == -99999) return FillCrtTrack(hit1, hit2, complete);
140 
141  hit1.x_pos = intersection.first.X();
142  hit1.y_pos = intersection.first.Y();
143  hit1.z_pos = intersection.first.Z();
144  hit1.x_err = 0.;
145  hit1.y_err = 0.;
146  hit1.z_err = 0.;
147  hit2.x_pos = intersection.second.X();
148  hit2.y_pos = intersection.second.Y();
149  hit2.z_pos = intersection.second.Z();
150  hit1.x_err = 0.;
151  hit1.y_err = 0.;
152  hit1.z_err = 0.;
153  }
154 
155  return FillCrtTrack(hit1, hit2, complete);
156 }
float z_err
position uncertainty in z-direction (cm).
Definition: CRTHit.hh:43
float x_err
position uncertainty in x-direction (cm).
Definition: CRTHit.hh:39
float y_err
position uncertainty in y-direction (cm).
Definition: CRTHit.hh:41
float z_pos
position in z-direction (cm).
Definition: CRTHit.hh:42
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
float y_pos
position in y-direction (cm).
Definition: CRTHit.hh:40
sbn::crt::CRTTrack FillCrtTrack(sbn::crt::CRTHit hit1, sbn::crt::CRTHit hit2, bool complete)
float x_pos
position in x-direction (cm).
Definition: CRTHit.hh:38
std::pair< TVector3, TVector3 > CubeIntersection(TVector3 min, TVector3 max, TVector3 start, TVector3 end)
std::string tagger
Name of the CRT wall (in the form of strings).
Definition: CRTHit.hh:45
void sbnd::CRTTrackRecoAlg::reconfigure ( const Config config)

Definition at line 26 of file sbndcode/sbndcode/CRT/CRTUtils/CRTTrackRecoAlg.cc.

26  {
27 
28  fTimeLimit = config.TimeLimit();
29  fAverageHitDistance = config.AverageHitDistance();
30  fDistanceLimit = config.DistanceLimit();
31 
32  return;
33 
34 }

Member Data Documentation

double sbnd::CRTTrackRecoAlg::fAverageHitDistance
private
CRTGeoAlg sbnd::CRTTrackRecoAlg::fCrtGeo
private
double sbnd::CRTTrackRecoAlg::fDistanceLimit
private
double sbnd::CRTTrackRecoAlg::fTimeLimit
private
CRTHitRecoAlg sbnd::CRTTrackRecoAlg::hitAlg
private

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