LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
TrackContainmentAlg.cxx
Go to the documentation of this file.
1 #ifndef TRK_TRACKCONTAINMENTALG_CXX
2 #define TRK_TRACKCONTAINMENTALG_CXX
3 
4 #include "TrackContainmentAlg.hh"
5 
7 #include <iostream>
8 
10 
12  fTrackTree = tfs_tree_trk;
13  fTrackTree->SetObject(fTrackTree->GetName(),"Track Tree");
14  fTrackTree->Branch("run",&fRun);
15  fTrackTree->Branch("event",&fEvent);
16  fTrackTree->Branch("trk",&fTrackTreeObj,fTrackTreeObj.Leaflist().c_str());
17  fTrackTree->Branch("trk_col",&fCollection);
18  fTrackTree->Branch("trk_id",&fTrkID);
19  fTrackTree->Branch("trk_mindist",&fDistance);
20  fTrackTree->Branch("trk_containment",&fContainment);
21 }
22 
24 {
25  fZBuffer = p.get<double>("ZBuffer");
26  fYBuffer = p.get<double>("YBuffer");
27  fXBuffer = p.get<double>("XBuffer");
28  fIsolation = p.get<double>("Isolation");
29 
30  fDebug = p.get<bool>("Debug",false);
31  fMakeCosmicTags = p.get<bool>("MakeCosmicTags",true);
32  fFillOutputTree = p.get<bool>("FillOutputTree",true);
33 }
34 
36 {
37  if( track.Vertex().Z() < (0+fZBuffer) || track.Vertex().Z() > (geo.DetLength()-fZBuffer) )
38  return false;
39  if( track.End().Z() < (0+fZBuffer) || track.End().Z() > (geo.DetLength()-fZBuffer) )
40  return false;
41  if( track.Vertex().Y() < (-1*geo.DetHalfHeight()+fYBuffer) || track.Vertex().Y() > (geo.DetHalfHeight()-fYBuffer) )
42  return false;
43  if( track.End().Y() < (-1*geo.DetHalfHeight()+fYBuffer) || track.End().Y() > (geo.DetHalfHeight()-fYBuffer) )
44  return false;
45  if( track.Vertex().X() < (0+fXBuffer) || track.Vertex().X() > (2*geo.DetHalfWidth()-fXBuffer) )
46  return false;
47  if( track.End().X() < (0+fXBuffer) || track.End().X() > (2*geo.DetHalfWidth()-fXBuffer) )
48  return false;
49 
50  return true;
51 }
52 
54 {
55 
57 
58  if(track.Vertex().Z() < (0+fZBuffer) || track.Vertex().Z() > (geo.DetLength()-fZBuffer))
60  else if(track.Vertex().Y() < (-1*geo.DetHalfHeight()+fYBuffer) || track.Vertex().Y() > (geo.DetHalfHeight()-fYBuffer))
62  else if( (track.Vertex().X()>0 && track.Vertex().X() < (0+fXBuffer)) ||
63  (track.Vertex().X()<2*geo.DetHalfWidth() && track.Vertex().X() > (2*geo.DetHalfWidth()-fXBuffer)) )
65 
66  if( track.End().Z() < (0+fZBuffer) || track.End().Z() > (geo.DetLength()-fZBuffer) ){
75  }
76  else if( track.End().Y() < (-1*geo.DetHalfHeight()+fYBuffer) || track.End().Y() > (geo.DetHalfHeight()-fYBuffer) ){
85  }
86  else if( (track.End().X()>0 && track.End().X() < (0+fXBuffer)) ||
87  (track.End().X()<2*geo.DetHalfWidth() && track.End().X() > (2*geo.DetHalfWidth()-fXBuffer)) ){
96  }
97 
98  if(track.Vertex().X() < 0 || track.Vertex().X() > 2*geo.DetHalfWidth())
100  if(track.End().X() < 0 || track.End().X() > 2*geo.DetHalfWidth()){
103  else
105  }
106 
107  return id;
108 }
109 
111 {
112  double min_distance = 9e12;
113  double tmp;
114  for(size_t i_p=0; i_p<t_ref.NumberTrajectoryPoints(); ++i_p){
115  tmp =
116  (t_probe.Vertex().X()-t_ref.LocationAtPoint(i_p).X())*(t_probe.Vertex().X()-t_ref.LocationAtPoint(i_p).X()) +
117  (t_probe.Vertex().Y()-t_ref.LocationAtPoint(i_p).Y())*(t_probe.Vertex().Y()-t_ref.LocationAtPoint(i_p).Y()) +
118  (t_probe.Vertex().Z()-t_ref.LocationAtPoint(i_p).Z())*(t_probe.Vertex().Z()-t_ref.LocationAtPoint(i_p).Z());
119  //std::cout << "\t\t\ttmp=" << tmp << std::endl;
120  if(tmp<min_distance)
121  min_distance=tmp;
122  }
123  return std::sqrt(min_distance);
124 }
125 
127 {
128  double min_distance = 9e12;
129  double tmp;
130  for(size_t i_p=0; i_p<t_ref.NumberTrajectoryPoints(); ++i_p){
131  tmp =
132  (t_probe.End().X()-t_ref.LocationAtPoint(i_p).X())*(t_probe.End().X()-t_ref.LocationAtPoint(i_p).X()) +
133  (t_probe.End().Y()-t_ref.LocationAtPoint(i_p).Y())*(t_probe.End().Y()-t_ref.LocationAtPoint(i_p).Y()) +
134  (t_probe.End().Z()-t_ref.LocationAtPoint(i_p).Z())*(t_probe.End().Z()-t_ref.LocationAtPoint(i_p).Z());
135  if(tmp<min_distance)
136  min_distance=tmp;
137  }
138  return std::sqrt(min_distance);
139 }
140 
141 void trk::TrackContainmentAlg::SetRunEvent(unsigned int const& run, unsigned int const& event)
142 {
143  fRun = run;
144  fEvent = event;
145 }
146 
147 void trk::TrackContainmentAlg::ProcessTracks(std::vector< std::vector<recob::Track> > const& tracksVec,
148  geo::GeometryCore const& geo)
149 {
150 
151  if(fDebug){
152  std::cout << "Geometry:" << std::endl;
153  std::cout << "\t" << geo.DetHalfWidth() << " " << geo.DetHalfHeight() << " " << geo.DetLength() << std::endl;
154  std::cout << "\t z:(" << fZBuffer << "," << geo.DetLength()-fZBuffer << ")"
155  << "\t y:(" << -1.*geo.DetHalfHeight()+fYBuffer << "," << geo.DetHalfHeight()-fYBuffer << ")"
156  << "\t x:(" << 0+fXBuffer << "," << 2*geo.DetHalfWidth()-fXBuffer << ")" << std::endl;
157  }
158 
159 
160  int containment_level=0;
161  bool track_linked = false;
162  std::size_t n_tracks=0;
163 
164  fTrackContainmentLevel.clear();
165  fTrackContainmentLevel.resize(tracksVec.size());
166  fMinDistances.clear();
167  fMinDistances.resize(tracksVec.size());
168 
169  fTrackContainmentIndices.clear();
170  fTrackContainmentIndices.push_back( std::vector< std::pair<int,int> >() );
171 
172  fCosmicTags.clear();
173  fCosmicTags.resize(tracksVec.size());
174 
175  //first, loop through tracks and see what's not contained
176 
177  for(size_t i_tc=0; i_tc<tracksVec.size(); ++i_tc){
178  fTrackContainmentLevel[i_tc].resize(tracksVec[i_tc].size(),-1);
179  fMinDistances[i_tc].resize(tracksVec[i_tc].size(),9e12);
180  fCosmicTags[i_tc].resize(tracksVec[i_tc].size(),anab::CosmicTag(-1));
181  n_tracks += tracksVec[i_tc].size();
182  for(size_t i_t=0; i_t<tracksVec[i_tc].size(); ++i_t){
183 
184  if(!IsContained(tracksVec[i_tc][i_t],geo)){
185  if(!track_linked) track_linked=true;
186  fTrackContainmentLevel[i_tc][i_t] = 0;
187  fTrackContainmentIndices.back().emplace_back(i_tc,i_t);
188  if(fDebug){
189  std::cout << "\tTrack (" << i_tc << "," << i_t << ")"
190  << " " << containment_level << std::endl;
191  }
192 
193  }//end if contained
194  }//end loop over tracks
195 
196  }//end loop over track collections
197 
198 
199 
200  //now, while we are still linking tracks, loop over all tracks and note anything
201  //close to an uncontained (or linked) track
202  while(track_linked){
203 
204  track_linked=false;
205  ++containment_level;
206  fTrackContainmentIndices.push_back( std::vector< std::pair<int,int> >() );
207 
208  for(size_t i_tc=0; i_tc<tracksVec.size(); ++i_tc){
209  for(size_t i_t=0; i_t<tracksVec[i_tc].size(); ++i_t){
210  if(fTrackContainmentLevel[i_tc][i_t]>=0)
211  continue;
212  else
213  {
214  for(auto const& i_tr : fTrackContainmentIndices[containment_level-1]){
215 
216  /*
217  if(fDebug){
218  std::cout << "\t\t" << MinDistanceStartPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second]) << std::endl;
219  std::cout << "\t\t" << MinDistanceEndPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second]) << std::endl;
220  }
221  */
222 
223  if(MinDistanceStartPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second])<fMinDistances[i_tc][i_t])
224  fMinDistances[i_tc][i_t] = MinDistanceStartPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second]);
225  if(MinDistanceEndPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second])<fMinDistances[i_tc][i_t])
226  fMinDistances[i_tc][i_t] = MinDistanceEndPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second]);
227 
228  if(MinDistanceStartPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second])<fIsolation ||
229  MinDistanceEndPt(tracksVec[i_tc][i_t],tracksVec[i_tr.first][i_tr.second])<fIsolation){
230  if(!track_linked) track_linked=true;
231  fTrackContainmentLevel[i_tc][i_t] = containment_level;
232  fTrackContainmentIndices.back().emplace_back(i_tc,i_t);
233 
234  if(fDebug){
235  std::cout << "\tTrackPair (" << i_tc << "," << i_t << ") and (" << i_tr.first << "," << i_tr.second << ")"
236  << " " << containment_level << std::endl;
237  }
238 
239  }//end if track not isolated
240 
241  }//end loop over existing uncontained/linked tracks
242 
243  }//end if track not already linked
244 
245  }//end loops over tracks
246  }//end loop over track collections
247  }//end while linking tracks
248 
249 
250  if(fDebug)
251  std::cout << "All done! Now let's make the tree and tags!" << std::endl;
252 
253  //now we're going to will the tree and create tags if requested
254  for(size_t i_tc=0; i_tc<tracksVec.size(); ++i_tc){
255  for(size_t i_t=0; i_t<tracksVec[i_tc].size(); ++i_t){
256 
257  //fill ROOT Tree
258  if(fFillOutputTree){
259  fTrackTreeObj = TrackTree_t(tracksVec[i_tc][i_t]);
260  fDistance = fMinDistances[i_tc][i_t];
261  fCollection = i_tc;
262  fTrkID = i_t;
264  fTrackTree->Fill();
265  }
266 
267  if(fMakeCosmicTags){
268 
269  //default (if track looks contained and isolated)
270  float score=0;
272 
273  //overwrite if track is not contained or not isolated
274  if(fTrackContainmentLevel[i_tc][i_t]>=0){
275  score = 1./(1.+(float)fTrackContainmentLevel[i_tc][i_t]);
276  if(fTrackContainmentLevel[i_tc][i_t]==0)
277  id = GetCosmicTagID(tracksVec[i_tc][i_t],geo);
278  else
280  }
281 
282  fCosmicTags[i_tc][i_t] =
283  anab::CosmicTag(std::vector<float>{(float)tracksVec[i_tc][i_t].Vertex().X(),
284  (float)tracksVec[i_tc][i_t].Vertex().Y(),
285  (float)tracksVec[i_tc][i_t].Vertex().Z()},
286  std::vector<float>{(float)tracksVec[i_tc][i_t].End().X(),
287  (float)tracksVec[i_tc][i_t].End().Y(),
288  (float)tracksVec[i_tc][i_t].End().Z()},
289  score,id);
290  }//end cosmic tag making
291 
292 
293  //some debug work
294  if(fTrackContainmentLevel[i_tc][i_t]<0 && fDebug){
295  std::cout << "Track (" << i_tc << "," << i_t << ")"
296  << " " << fTrackContainmentLevel[i_tc][i_t]
297  << " " << fMinDistances[i_tc][i_t] << std::endl;
298  std::cout << "\tS_(X,Y,Z) = ("
299  << tracksVec[i_tc][i_t].Vertex().X() << ","
300  << tracksVec[i_tc][i_t].Vertex().Y() << ","
301  << tracksVec[i_tc][i_t].Vertex().Z() << ")" << std::endl;
302  std::cout << "\tNearest wire ..." << std::endl;
303  for(size_t i_p=0; i_p<geo.Nplanes(); ++i_p)
304  std::cout << "\t\tPlane " << i_p << " " << geo.NearestWireID(tracksVec[i_tc][i_t].Vertex(),i_p).Wire << std::endl;
305  std::cout << "\tE_(X,Y,Z) = ("
306  << tracksVec[i_tc][i_t].End().X() << ","
307  << tracksVec[i_tc][i_t].End().Y() << ","
308  << tracksVec[i_tc][i_t].End().Z() << ")" << std::endl;
309  std::cout << "\tNearest wire ..." << std::endl;
310  for(size_t i_p=0; i_p<geo.Nplanes(); ++i_p)
311  std::cout << "\t\tPlane " << i_p << " " << geo.NearestWireID(tracksVec[i_tc][i_t].End(),i_p).Wire << std::endl;
312  std::cout << "\tLength=" << tracksVec[i_tc][i_t].Length() << std::endl;
313  std::cout << "\tSimple_length=" << (tracksVec[i_tc][i_t].End()-tracksVec[i_tc][i_t].Vertex()).R() << std::endl;
314  }//end debug statements if track contained
315 
316  }//end loops over tracks
317  }//end loop over track collections
318 
319 }//end ProcessTracks
320 
321 
322 std::vector< std::vector<anab::CosmicTag> > const& trk::TrackContainmentAlg::GetTrackCosmicTags()
323 {
324  if(fMakeCosmicTags)
325  return fCosmicTags;
326  else
327  throw cet::exception("TrackContainmentAlg::GetTrackCosmicTags")
328  << "Cosmic tags not created. Set MakeCosmicTags to true in fcl paramters.";
329 }
330 
331 #endif
std::vector< std::vector< int > > fTrackContainmentLevel
double MinDistanceStartPt(recob::Track const &, recob::Track const &)
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
enum anab::cosmic_tag_id CosmicTagID_t
Point_t const & LocationAtPoint(size_t i) const
Access to track position at different points.
Definition: Track.h:129
std::vector< std::vector< double > > fMinDistances
size_t NumberTrajectoryPoints() const
Various functions related to the presence and the number of (valid) points.
Definition: Track.h:105
anab::CosmicTagID_t GetCosmicTagID(recob::Track const &, geo::GeometryCore const &)
Float_t tmp
Definition: plot.C:37
typename BeginEndPackage< L >::End End
void ProcessTracks(std::vector< std::vector< recob::Track > > const &, geo::GeometryCore const &)
std::vector< std::vector< std::pair< int, int > > > fTrackContainmentIndices
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
Access the description of detector geometry.
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
struct trk::TrackTree TrackTree_t
std::string Leaflist()
T get(std::string const &key) const
Definition: ParameterSet.h:231
Point_t const & Vertex() const
Access to track position at different points.
Definition: Track.h:127
Double_t R
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
void Configure(fhicl::ParameterSet const &)
void SetRunEvent(unsigned int const &, unsigned int const &)
Description of geometry of one entire detector.
double MinDistanceEndPt(recob::Track const &, recob::Track const &)
std::vector< std::vector< anab::CosmicTag > > fCosmicTags
std::vector< std::vector< anab::CosmicTag > > const & GetTrackCosmicTags()
TrackContainmentAlg()
Default constructor.
Point_t const & End() const
Access to track position at different points.
Definition: Track.h:128
geo::WireID NearestWireID(geo::Point_t const &point, geo::PlaneID const &planeid) const
Returns the ID of wire closest to position in the specified TPC.
Namespace collecting geometry-related classes utilities.
Float_t track
Definition: plot.C:34
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:52
bool IsContained(recob::Track const &, geo::GeometryCore const &)
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Event finding and building.