HoughLineFinderAna_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // HoughLineFinder class
4 //
5 // joshua.spitz@yale.edu
6 //
7 // This algorithm is designed to find lines (Houghclusters) from clusters found by DBSCAN after deconvolution and hit finding.
8 // The algorithm is based on:
9 // Queisser, A. "Computing the Hough Transform", C/C++ Users Journal 21, 12 (Dec. 2003).
10 // Niblack, W. and Petkovic, D. On Improving the Accuracy of the Hough Transform", Machine Vision and Applications 3, 87 (1990)
11 ////////////////////////////////////////////////////////////////////////
12 
13 #include <cmath>
14 #include <string>
15 
16 // ROOT includes
17 #include "TTree.h"
18 
19 // Framework includes
22 #include "canvas/Persistency/Common/FindManyP.h"
24 #include "fhiclcpp/ParameterSet.h"
29 #include "art_root_io/TFileService.h"
31 
32 // LArSoft includes
36 
37 namespace cluster {
38 
40 
41  public:
42 
43  explicit HoughLineFinderAna(fhicl::ParameterSet const& pset);
45 
46  private:
47 
48  void analyze(const art::Event&);
49  void beginJob();
50 
55  TTree* ftree;
56  int fm_run; // Run number
57  unsigned long int fm_run_timestamp; // Run number
58  int fm_event; // Event number
59  int fm_plane; // Plane number
60  int fm_dbsize;
61  int fm_clusterid; // Cluster ID
62  int fm_wirespan; // Wire spanned by track
63  int fm_sizeClusterZ; //Number of clusters
64  int fm_sizeHitZ; //Number of Hits
67  int *fm_wireZ;
68  int *fm_hitidZ;
69  float *fm_mipZ;
70  float *fm_drifttimeZ;
71  float *fm_widthZ;
72  float *fm_upadcZ;
73 
74  };
75 
76 
77 } // end namespace cluster
78 
79 //#endif // HoughLineFinderAna_H
80 
81 
82 namespace cluster {
83 
85  : EDAnalyzer(pset)
86  , fHoughModuleLabel (pset.get< std::string >("HoughModuleLabel"))
87  , fDigitModuleLabel (pset.get< std::string >("DigitModuleLabel"))
88  , fHitsModuleLabel (pset.get< std::string >("HitsModuleLabel"))
89  , fDBScanModuleLabel(pset.get< std::string >("DBScanModuleLabel"))
90  , fm_run(0)
91  , fm_event(0)
92  , fm_plane(0)
93  , fm_dbsize(0)
94  , fm_clusterid(0)
95  , fm_wirespan(0)
96  , fm_sizeClusterZ(10000)
97  , fm_sizeHitZ(10000)
98  , fm_clusterslope(0)
100  {
101  }
102 
103  //-------------------------------------------------
105  {
106  delete fm_hitidZ;
107  delete fm_mipZ;
108  delete fm_drifttimeZ;
109  delete fm_widthZ;
110  delete fm_upadcZ;
111  delete fm_wireZ;
112  }
113 
114  //-------------------------------------------------
116  {
117 
118 
119  // get access to the TFile service
121  ftree= tfs->make<TTree>("HoughTree","HoughTree");
122  fm_hitidZ = new int[fm_sizeHitZ];
123  fm_mipZ = new float[fm_sizeHitZ];
124  fm_drifttimeZ = new float[fm_sizeHitZ];
125  fm_widthZ = new float[fm_sizeHitZ];
126  fm_upadcZ = new float[fm_sizeHitZ];
127  fm_wireZ = new int[fm_sizeHitZ];
128  ftree->Branch("run", &fm_run, "run/I");
129  ftree->Branch("run_timestamp", &fm_run_timestamp, "run_timestamp/l"); //l is for ULong64_t
130  ftree->Branch("event", &fm_event, "event/I");
131  ftree->Branch("plane",&fm_plane,"plane/I");
132  ftree->Branch("dbsize",&fm_dbsize,"dbsize/I");
133  ftree->Branch("clusterid",&fm_clusterid,"clusterid/I");
134  ftree->Branch("clusterslope",&fm_clusterslope,"clusterslope/F");
135  ftree->Branch("clusterintercept",&fm_clusterintercept,"clusterintecept/F");
136  ftree->Branch("wirespan",&fm_wirespan,"wirespan/I");
137  ftree->Branch("numberHits",&fm_sizeHitZ,"numberHits/I");
138  ftree->Branch("numberClusters",&fm_sizeClusterZ,"numberClusters/I");
139  ftree->Branch("hitidZ",fm_hitidZ,"hitidZ[numberHits]/I");
140  ftree->Branch("wireZ",fm_wireZ,"wireZ[numberHits]/I");
141  ftree->Branch("mipZ",fm_mipZ,"mipZ[numberHits]/F");
142  ftree->Branch("drifttimeZ",fm_drifttimeZ,"drifttitmeZ[numberHits]/F");
143  ftree->Branch("widthZ",fm_widthZ,"widthZ[numberHits]/F");
144  }
145 
146 
148  {
149 
151  evt.getByLabel(fHoughModuleLabel,hlfListHandle);
152  art::Handle< std::vector<recob::Hit> > hitListHandle;
153  evt.getByLabel(fHitsModuleLabel,hitListHandle);
154  art::Handle< std::vector<recob::Cluster> > dbscanListHandle;
155  evt.getByLabel(fDBScanModuleLabel,dbscanListHandle);
156 
157  art::FindManyP<recob::Hit> fmh(dbscanListHandle, evt, fDBScanModuleLabel);
158  art::FindManyP<recob::Hit> fmhhl(hlfListHandle, evt, fHoughModuleLabel);
159 
162  // art::PtrVector<recob::Hit> hits;// unused, as yet. EC, 5-Oct-2010.
163 
164  for (size_t ii = 0; ii < hlfListHandle->size(); ++ii){
165  art::Ptr<recob::Cluster> cluster(hlfListHandle,ii);
166  clusters.push_back(cluster);
167  }
168 
169  for (size_t ii = 0; ii < dbscanListHandle->size(); ++ii){
170  art::Ptr<recob::Cluster> dbcluster(dbscanListHandle,ii);
171  dbclusters.push_back(dbcluster);
172  }
173 
174  MF_LOG_VERBATIM("HoughLineFinderAna") << "run : " << evt.id().run();
175  //std::cout << "subrun : " << evt.subRun() << std::endl;
176  MF_LOG_VERBATIM("HoughLineFinderAna") << "event : " << evt.id().event();
177  fm_run=evt.id().run();
178  fm_event=evt.id().event();
179  fm_run_timestamp=evt.time().value(); // won't cast, EC, 7-Oct-2010.
180  unsigned int firstwire=0;
181  unsigned int lastwire=0;
182  fm_sizeClusterZ=0;
183  fm_sizeHitZ=0;
184  fm_dbsize=0;
186 
187  for(auto view : geo->Views()){
188 
189  fm_dbsize = 0;
190  fm_sizeClusterZ = clusters.size();
191 
192  for(size_t j = 0; j < dbclusters.size(); ++j) {
193  if(dbclusters[j]->View() == view){
194  std::vector< art::Ptr<recob::Hit> > _dbhits = fmh.at(j);
195  fm_dbsize += _dbhits.size();
196  if(_dbhits.size() > 0) fm_plane = _dbhits.at(0)->WireID().Plane;
197  }
198  }
199 
200  for(size_t j = 0; j < clusters.size(); ++j) {
201  if(clusters[j]->View() == view){
202  fm_clusterid=clusters[j]->ID();
203  std::vector< art::Ptr<recob::Hit> > _hits = fmhhl.at(j);
204  fm_clusterslope=(double) std::tan(clusters[j]->StartAngle());
205  fm_clusterintercept=(double)clusters[j]->StartTick();
206  if(_hits.size()!=0){
207  fm_plane = _hits.at(0)->WireID().Plane;
208  firstwire = _hits[0]->WireID().Wire;
209  lastwire = _hits[_hits.size()-1]->WireID().Wire;
210  fm_wirespan = lastwire-firstwire;
211  fm_sizeHitZ = _hits.size();
212 
213  for(unsigned int i = 0; i < _hits.size(); ++i){
214 
215  fm_hitidZ[i] = i;
216  fm_wireZ[i] = _hits[i]->WireID().Wire;
217  fm_mipZ[i] = (double)_hits[i]->Integral();
218  fm_drifttimeZ[i] = (double)_hits[i]->PeakTime();
219  fm_widthZ[i] = (double) (2. * _hits[i]->RMS());
220  fm_upadcZ[i] = (double)_hits[i]->Integral();
221  }
222 
223  ftree->Fill();
224  }
225  }//end if in the correct view
226  }// end loop over clusters
227  }// end loop over views
228 
229  }
230 
231 }// end namespace
232 
233 
234 namespace cluster{
235 
237 
238 } // end namespace caldata
AdcChannelData::View View
std::string string
Definition: nybbler.cc:12
std::set< geo::View_t > const & Views() const
Returns a list of possible views in the detector.
STL namespace.
Cluster finding and building.
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
constexpr TimeValue_t value() const
Definition: Timestamp.h:23
RunNumber_t run() const
Definition: EventID.h:98
art framework interface to geometry description
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
Timestamp time() const
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
reference at(size_type n)
Definition: PtrVector.h:359
size_type size() const
Definition: PtrVector.h:302
HoughLineFinderAna(fhicl::ParameterSet const &pset)
Declaration of signal hit object.
#define MF_LOG_VERBATIM(category)
EventNumber_t event() const
Definition: EventID.h:116
TCEvent evt
Definition: DataStructs.cxx:7
auto const & get(AssnsNode< L, R, D > const &r)
Definition: AssnsNode.h:115
LArSoft geometry interface.
Definition: ChannelGeo.h:16
EventID id() const
Definition: Event.cc:34