All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CRTTrackCosmicIdAna_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // Class: CRTTrackCosmicIdAna
3 // Module Type: analyzer
4 // File: CRTTrackCosmicIdAna_module.cc
5 //
6 // Tom Brooks (tbrooks@fnal.gov)
7 ////////////////////////////////////////////////////////////////////////
8 
9 // sbndcode includes
16 
17 // LArSoft includes
23 #include "nusimdata/SimulationBase/MCParticle.h"
24 
25 // Framework includes
26 #include "art/Framework/Core/EDAnalyzer.h"
27 #include "art/Framework/Principal/Event.h"
28 #include "art/Framework/Principal/Handle.h"
29 #include "art/Framework/Services/Registry/ServiceHandle.h"
30 #include "art_root_io/TFileService.h"
31 #include "canvas/Persistency/Common/FindManyP.h"
34 
35 // Utility libraries
36 #include "fhiclcpp/ParameterSet.h"
37 #include "fhiclcpp/types/Table.h"
38 #include "fhiclcpp/types/Atom.h"
39 
40 #include "Pandora/PdgTable.h"
41 
42 // ROOT includes. Note: To look up the properties of the ROOT classes,
43 // use the ROOT web site; e.g.,
44 // <https://root.cern.ch/doc/master/annotated.html>
45 #include "TH1.h"
46 #include "TVector3.h"
47 
48 // C++ includes
49 #include <map>
50 #include <vector>
51 #include <string>
52 
53 namespace sbnd {
54 
55  class CRTTrackCosmicIdAna : public art::EDAnalyzer {
56  public:
57 
58  // Describes configuration parameters of the module
59  struct Config {
60  using Name = fhicl::Name;
61  using Comment = fhicl::Comment;
62 
63  // One Atom for each parameter
64  fhicl::Atom<art::InputTag> SimModuleLabel {
65  Name("SimModuleLabel"),
66  Comment("tag of detector simulation data product")
67  };
68 
69  fhicl::Atom<art::InputTag> CRTTrackLabel {
70  Name("CRTTrackLabel"),
71  Comment("tag of CRT track producer data product")
72  };
73 
74  fhicl::Atom<art::InputTag> TPCTrackLabel {
75  Name("TPCTrackLabel"),
76  Comment("tag of tpc track producer data product")
77  };
78 
79  fhicl::Atom<art::InputTag> PandoraLabel {
80  Name("PandoraLabel"),
81  Comment("tag of pandora data product")
82  };
83 
84  fhicl::Atom<bool> Verbose {
85  Name("Verbose"),
86  Comment("Print information about what's going on")
87  };
88 
89  fhicl::Table<CRTTrackMatchAlg::Config> TrackMatchAlg {
90  Name("TrackMatchAlg"),
91  };
92 
93  fhicl::Table<CRTBackTracker::Config> CrtBackTrack {
94  Name("CrtBackTrack"),
95  };
96 
97  fhicl::Table<CrtTrackCosmicIdAlg::Config> CTTagAlg {
98  Name("CTTagAlg"),
99  };
100 
101  }; // Config
102 
103  using Parameters = art::EDAnalyzer::Table<Config>;
104 
105  // Constructor: configures module
106  explicit CRTTrackCosmicIdAna(Parameters const& config);
107 
108  // Called once, at start of the job
109  virtual void beginJob() override;
110 
111  // Called once per event
112  virtual void analyze(const art::Event& event) override;
113 
114  // Called once, at end of the job
115  virtual void endJob() override;
116 
117  typedef art::Handle< std::vector<recob::PFParticle> > PFParticleHandle;
118  typedef std::map< size_t, art::Ptr<recob::PFParticle> > PFParticleIdMap;
119 
120  private:
121 
122  void GetPFParticleIdMap(const PFParticleHandle &pfParticleHandle, PFParticleIdMap &pfParticleMap);
123 
124  // fcl file parameters
125  art::InputTag fSimModuleLabel; ///< name of detsim producer
126  art::InputTag fCRTTrackLabel; ///< name of CRT producer
127  art::InputTag fTPCTrackLabel; ///< name of CRT producer
128  art::InputTag fPandoraLabel;
129  bool fVerbose; ///< print information about what's going on
130 
132 
134 
137 
138  // Histograms
139  std::vector<std::string> types {"NuMuTrack", "CrTrack", "DirtTrack", "NuTrack", "NuMuPfp", "CrPfp", "DirtPfp", "NuPfp"};
140  std::map<std::string, TH1D*> hNumTrueMatches;
141  std::map<std::string, TH1D*> hMatchDCA;
142  std::map<std::string, TH1D*> hNoMatchDCA;
143  std::map<std::string, TH1D*> hMatchAngle;
144  std::map<std::string, TH1D*> hNoMatchAngle;
145  std::map<std::string, TH1D*> hDCATotal;
146  std::map<std::string, TH1D*> hDCATag;
147  std::map<std::string, TH1D*> hLengthTotal;
148  std::map<std::string, TH1D*> hLengthTag;
149 
150  }; // class CRTTrackCosmicIdAna
151 
152 
153  // Constructor
155  : EDAnalyzer(config)
156  , fSimModuleLabel (config().SimModuleLabel())
157  , fCRTTrackLabel (config().CRTTrackLabel())
158  , fTPCTrackLabel (config().TPCTrackLabel())
159  , fPandoraLabel (config().PandoraLabel())
160  , fVerbose (config().Verbose())
161  , trackAlg (config().TrackMatchAlg())
162  , fCrtBackTrack (config().CrtBackTrack())
163  , ctTag (config().CTTagAlg())
164  {
165 
166  } //CRTTrackCosmicIdAna()
167 
168 
170  {
171 
172  // Access tfileservice to handle creating and writing histograms
173  art::ServiceHandle<art::TFileService> tfs;
174  for(auto type : types){
175  hNumTrueMatches[type] = tfs->make<TH1D>(Form("%sNumTrueMatches", type.c_str()), "", 10, 0, 10);
176  hMatchDCA[type] = tfs->make<TH1D>(Form("%sMatchDCA", type.c_str()), "", 60, 0, 100);
177  hNoMatchDCA[type] = tfs->make<TH1D>(Form("%sNoMatchDCA", type.c_str()), "", 60, 0, 100);
178  hMatchAngle[type] = tfs->make<TH1D>(Form("%sMatchAngle", type.c_str()), "", 60, 0, 3.2);
179  hNoMatchAngle[type] = tfs->make<TH1D>(Form("%sNoMatchAngle", type.c_str()), "", 60, 0, 3.2);
180  hDCATotal[type] = tfs->make<TH1D>(Form("%sDCATotal", type.c_str()), "", 20, 0, 80);
181  hDCATag[type] = tfs->make<TH1D>(Form("%sDCATag", type.c_str()), "", 20, 0, 80);
182  hLengthTotal[type] = tfs->make<TH1D>(Form("%sLengthTotal", type.c_str()), "", 20, 0, 400);
183  hLengthTag[type] = tfs->make<TH1D>(Form("%sLengthTag", type.c_str()), "", 20, 0, 400);
184  }
185 
186  // Initial output
187  if(fVerbose) std::cout<<"----------------- CRT T0 Matching Ana Module -------------------"<<std::endl;
188 
189  } // CRTTrackCosmicIdAna::beginJob()
190 
191 
192  void CRTTrackCosmicIdAna::analyze(const art::Event& event)
193  {
194 
195  // Fetch basic event info
196  if(fVerbose){
197  std::cout<<"============================================"<<std::endl
198  <<"Run = "<<event.run()<<", SubRun = "<<event.subRun()<<", Event = "<<event.id().event()<<std::endl
199  <<"============================================"<<std::endl;
200  }
201 
202  //----------------------------------------------------------------------------------------------------------
203  // GETTING PRODUCTS
204  //----------------------------------------------------------------------------------------------------------
205 
206  // Get g4 particles
207  art::ServiceHandle<cheat::ParticleInventoryService> pi_serv;
208  auto particleHandle = event.getValidHandle<std::vector<simb::MCParticle>>(fSimModuleLabel);
209 
210  // Get CRT hits from the event
211  art::Handle< std::vector<sbn::crt::CRTTrack>> crtTrackHandle;
212  std::vector<art::Ptr<sbn::crt::CRTTrack> > crtTrackList;
213  if (event.getByLabel(fCRTTrackLabel, crtTrackHandle))
214  art::fill_ptr_vector(crtTrackList, crtTrackHandle);
215 
216  fCrtBackTrack.Initialize(event);
217  std::vector<sbn::crt::CRTTrack> crtTracks;
218  std::map<int, int> numCrtTrackMap;
219  int track_i = 0;
220  for(auto const& track : (crtTrackList)){
221  crtTracks.push_back(*track);
222  int trackTrueID = fCrtBackTrack.TrueIdFromTrackId(event, track_i);
223  track_i++;
224  double trackTime = track->ts1_ns * 1e-3;
225  if(trackTime > 0 && trackTime < 4) continue;
226  numCrtTrackMap[trackTrueID]++;
227  }
228 
229  // Get reconstructed tracks from the event
230  auto tpcTrackHandle = event.getValidHandle<std::vector<recob::Track>>(fTPCTrackLabel);
231  art::FindManyP<recob::Hit> findManyHits(tpcTrackHandle, event, fTPCTrackLabel);
232 
233  // Get PFParticles from pandora
234  PFParticleHandle pfParticleHandle;
235  event.getByLabel(fPandoraLabel, pfParticleHandle);
236  if( !pfParticleHandle.isValid() ){
237  if(fVerbose) std::cout<<"Failed to find the PFParticles."<<std::endl;
238  return;
239  }
240  PFParticleIdMap pfParticleMap;
241  this->GetPFParticleIdMap(pfParticleHandle, pfParticleMap);
242  // Get PFParticle to track associations
243  art::FindManyP< recob::Track > pfPartToTrackAssoc(pfParticleHandle, event, fTPCTrackLabel);
244 
245  //----------------------------------------------------------------------------------------------------------
246  // TRUTH MATCHING
247  //----------------------------------------------------------------------------------------------------------
248 
249  std::map<int, simb::MCParticle> particles;
250  std::vector<int> nuParticleIds;
251  std::vector<int> lepParticleIds;
252  std::vector<int> dirtParticleIds;
253  std::vector<int> crParticleIds;
254  // Loop over the true particles
255  for (auto const& particle: (*particleHandle)){
256 
257  // Make map with ID
258  int partID = particle.TrackId();
259  particles[partID] = particle;
260 
261  // Get MCTruth
262  art::Ptr<simb::MCTruth> truth = pi_serv->TrackIdToMCTruth_P(partID);
263  int pdg = std::abs(particle.PdgCode());
264 
265  // If origin is a neutrino
266  if(truth->Origin() == simb::kBeamNeutrino){
267  geo::Point_t vtx;
268  vtx.SetX(truth->GetNeutrino().Nu().Vx()); vtx.SetY(truth->GetNeutrino().Nu().Vy()); vtx.SetZ(truth->GetNeutrino().Nu().Vz());
269  // If neutrino vertex is not inside the TPC then call it a dirt particle
270  if(!fTpcGeo.InFiducial(vtx, 0, 0)){
271  dirtParticleIds.push_back(partID);
272  }
273  // If it's a primary muon
274  else if(pdg==13 && particle.Mother()==0){
275  lepParticleIds.push_back(partID);
276  }
277  // Other nu particles
278  else{
279  nuParticleIds.push_back(partID);
280  }
281  }
282 
283  // If origin is a cosmic ray
284  else if(truth->Origin() == simb::kCosmicRay){
285  crParticleIds.push_back(partID);
286  }
287 
288  }
289 
290  //----------------------------------------------------------------------------------------------------------
291  // DISTANCE OF CLOSEST APPROACH ANALYSIS
292  //----------------------------------------------------------------------------------------------------------
293  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(event);
294  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(event, clockData);
295 
296  // Loop over reconstructed tracks
297  for (auto const& tpcTrack : (*tpcTrackHandle)){
298  // Get the associated hits
299  std::vector<art::Ptr<recob::Hit>> hits = findManyHits.at(tpcTrack.ID());
300  int trackTrueID = RecoUtils::TrueParticleIDFromTotalRecoHits(clockData, hits, false);
301  std::string type = "none";
302  if(std::find(lepParticleIds.begin(), lepParticleIds.end(), trackTrueID) != lepParticleIds.end()) type = "NuMuTrack";
303  if(std::find(nuParticleIds.begin(), nuParticleIds.end(), trackTrueID) != nuParticleIds.end()) type = "NuTrack";
304  if(std::find(crParticleIds.begin(), crParticleIds.end(), trackTrueID) != crParticleIds.end()) type = "CrTrack";
305  if(std::find(dirtParticleIds.begin(), dirtParticleIds.end(), trackTrueID) != dirtParticleIds.end()) type = "DirtTrack";
306  if(type == "none") continue;
307 
308  if(numCrtTrackMap.find(trackTrueID) != numCrtTrackMap.end()){
309  hNumTrueMatches[type]->Fill(numCrtTrackMap[trackTrueID]);
310  }
311  else{
312  hNumTrueMatches[type]->Fill(0);
313  }
314 
315  // Calculate t0 from CRT track matching
316  std::pair<sbn::crt::CRTTrack, double> closestAngle = trackAlg.ClosestCRTTrackByAngle(detProp, tpcTrack, crtTracks, event);
317  std::pair<sbn::crt::CRTTrack, double> closestDCA = trackAlg.ClosestCRTTrackByDCA(detProp, tpcTrack, crtTracks, event);
318 
319  if(closestAngle.second != -99999){
320  int crtTrackTrueID = fCrtBackTrack.TrueIdFromTotalEnergy(event, closestAngle.first);
321  if(crtTrackTrueID == trackTrueID && crtTrackTrueID != -99999){
322  hMatchAngle[type]->Fill(closestAngle.second);
323  }
324  else{
325  hNoMatchAngle[type]->Fill(closestAngle.second);
326  }
327  }
328 
329  if(closestDCA.second != -99999){
330  int crtTrackTrueID = fCrtBackTrack.TrueIdFromTotalEnergy(event, closestDCA.first);
331  if(crtTrackTrueID == trackTrueID && crtTrackTrueID != -99999){
332  hMatchDCA[type]->Fill(closestDCA.second);
333  }
334  else{
335  hNoMatchDCA[type]->Fill(closestDCA.second);
336  }
337  }
338 
339 
340  int nbins = hDCATotal.begin()->second->GetNbinsX();
341  for(int i = 0; i < nbins; i++){
342  double DCAcut = hDCATotal.begin()->second->GetBinCenter(i);
343 
344  hDCATotal[type]->Fill(DCAcut);
345 
346  // If closest hit is below limit and track matches any hits then fill efficiency
347  if(closestDCA.second < DCAcut && closestDCA.second != -99999){
348  double trackTime = closestDCA.first.ts1_ns * 1e-3;
349  if(trackTime > 0 && trackTime < 4) continue;
350  hDCATag[type]->Fill(DCAcut);
351  }
352 
353  }
354 
355  hLengthTotal[type]->Fill(tpcTrack.Length());
356  if(ctTag.CrtTrackCosmicId(detProp, tpcTrack, crtTracks, event)){
357  hLengthTag[type]->Fill(tpcTrack.Length());
358  }
359  }
360 
361 
362  //Loop over the pfparticle map
363  for (PFParticleIdMap::const_iterator it = pfParticleMap.begin(); it != pfParticleMap.end(); ++it){
364 
365  const art::Ptr<recob::PFParticle> pParticle(it->second);
366  // Only look for primary particles
367  if (!pParticle->IsPrimary()) continue;
368  // Check if this particle is identified as the neutrino
369  const int pdg(pParticle->PdgCode());
370  const bool isNeutrino(std::abs(pdg) == pandora::NU_E || std::abs(pdg) == pandora::NU_MU || std::abs(pdg) == pandora::NU_TAU);
371  //Find neutrino pfparticle
372  if(!isNeutrino) continue;
373 
374  std::string type = "none";
375  std::vector<recob::Track> nuTracks;
376 
377  // Loop over daughters of pfparticle
378  for (const size_t daughterId : pParticle->Daughters()){
379 
380  // Get tracks associated with daughter
381  art::Ptr<recob::PFParticle> pDaughter = pfParticleMap.at(daughterId);
382  const std::vector< art::Ptr<recob::Track> > associatedTracks(pfPartToTrackAssoc.at(pDaughter.key()));
383  if(associatedTracks.size() != 1) continue;
384 
385  // Get the first associated track
386  recob::Track tpcTrack = *associatedTracks.front();
387  nuTracks.push_back(tpcTrack);
388 
389  // Truth match muon tracks and pfps
390  std::vector<art::Ptr<recob::Hit>> hits = findManyHits.at(tpcTrack.ID());
391  int trueId = RecoUtils::TrueParticleIDFromTotalRecoHits(clockData, hits, false);
392  if(std::find(lepParticleIds.begin(), lepParticleIds.end(), trueId) != lepParticleIds.end()){
393  type = "NuMuPfp";
394  }
395  else if(std::find(nuParticleIds.begin(), nuParticleIds.end(), trueId) != nuParticleIds.end()){
396  if(type != "NuMuPfp") type = "NuPfp";
397  }
398  else if(std::find(dirtParticleIds.begin(), dirtParticleIds.end(), trueId) != dirtParticleIds.end()){
399  if(type != "NuMuPfp" && type != "NuPfp") type = "DirtPfp";
400  }
401  else if(std::find(crParticleIds.begin(), crParticleIds.end(), trueId) != crParticleIds.end()){
402  if(type != "NuMuPfp" && type != "NuPfp" && type != "DirtPfp") type = "CrPfp";
403  }
404  }
405 
406  if(nuTracks.size() == 0) continue;
407 
408  std::sort(nuTracks.begin(), nuTracks.end(), [](auto& left, auto& right){
409  return left.Length() > right.Length();});
410 
411  recob::Track tpcTrack = nuTracks[0];
412  std::vector<art::Ptr<recob::Hit>> hits = findManyHits.at(tpcTrack.ID());
413  int trackTrueID = RecoUtils::TrueParticleIDFromTotalRecoHits(clockData, hits, false);
414 
415  if(numCrtTrackMap.find(trackTrueID) != numCrtTrackMap.end()){
416  hNumTrueMatches[type]->Fill(numCrtTrackMap[trackTrueID]);
417  }
418  else{
419  hNumTrueMatches[type]->Fill(0);
420  }
421 
422  std::pair<sbn::crt::CRTTrack, double> closestAngle = trackAlg.ClosestCRTTrackByAngle(detProp, tpcTrack, crtTracks, event);
423  std::pair<sbn::crt::CRTTrack, double> closestDCA = trackAlg.ClosestCRTTrackByDCA(detProp, tpcTrack, crtTracks, event);
424 
425  if(closestAngle.second != -99999){
426  int crtTrackTrueID = fCrtBackTrack.TrueIdFromTotalEnergy(event, closestAngle.first);
427  if(crtTrackTrueID == trackTrueID && crtTrackTrueID != -99999){
428  hMatchAngle[type]->Fill(closestAngle.second);
429  }
430  else{
431  hNoMatchAngle[type]->Fill(closestAngle.second);
432  }
433  }
434 
435  if(closestDCA.second != -99999){
436  int crtTrackTrueID = fCrtBackTrack.TrueIdFromTotalEnergy(event, closestDCA.first);
437  if(crtTrackTrueID == trackTrueID && crtTrackTrueID != -99999){
438  hMatchDCA[type]->Fill(closestDCA.second);
439  }
440  else{
441  hNoMatchDCA[type]->Fill(closestDCA.second);
442  }
443  }
444 
445 
446  int nbins = hDCATotal.begin()->second->GetNbinsX();
447  for(int i = 0; i < nbins; i++){
448  double DCAcut = hDCATotal.begin()->second->GetBinCenter(i);
449 
450  hDCATotal[type]->Fill(DCAcut);
451 
452  // If closest hit is below limit and track matches any hits then fill efficiency
453  if(closestDCA.second < DCAcut && closestDCA.second != -99999){
454  double trackTime = closestDCA.first.ts1_ns * 1e-3;
455  if(trackTime > 0 && trackTime < 4) continue;
456  hDCATag[type]->Fill(DCAcut);
457  }
458 
459  }
460 
461  hLengthTotal[type]->Fill(tpcTrack.Length());
462  if(ctTag.CrtTrackCosmicId(detProp, tpcTrack, crtTracks, event)){
463  hLengthTag[type]->Fill(tpcTrack.Length());
464  }
465  }
466 
467 
468  } // CRTTrackCosmicIdAna::analyze()
469 
470 
472 
473 
474  } // CRTTrackCosmicIdAna::endJob()
475 
476 
477  void CRTTrackCosmicIdAna::GetPFParticleIdMap(const PFParticleHandle &pfParticleHandle, PFParticleIdMap &pfParticleMap){
478  for (unsigned int i = 0; i < pfParticleHandle->size(); ++i){
479  const art::Ptr<recob::PFParticle> pParticle(pfParticleHandle, i);
480  if (!pfParticleMap.insert(PFParticleIdMap::value_type(pParticle->Self(), pParticle)).second){
481  std::cout << " Unable to get PFParticle ID map, the input PFParticle collection has repeat IDs!" <<"\n";
482  }
483  }
484  }
485 
486 
487  DEFINE_ART_MODULE(CRTTrackCosmicIdAna)
488 } // namespace sbnd
BEGIN_PROLOG Verbose
int TrueIdFromTrackId(const art::Event &event, int track_i)
art::Handle< std::vector< recob::PFParticle > > PFParticleHandle
var pdg
Definition: selectors.fcl:14
Declaration of signal hit object.
std::map< size_t, art::Ptr< recob::PFParticle > > PFParticleIdMap
walls no right
Definition: selectors.fcl:105
fhicl::Atom< art::InputTag > SimModuleLabel
art::InputTag fCRTTrackLabel
name of CRT producer
bool fVerbose
print information about what&#39;s going on
void GetPFParticleIdMap(const PFParticleHandle &pfParticleHandle, PFParticleIdMap &pfParticleMap)
process_name use argoneut_mc_hitfinder track
std::vector< std::string > types
int TrueIdFromTotalEnergy(const art::Event &event, const sbnd::crt::CRTData &data)
fhicl::Atom< art::InputTag > CRTTrackLabel
std::map< std::string, TH1D * > hDCATag
std::map< std::string, TH1D * > hDCATotal
T abs(T value)
art::InputTag fSimModuleLabel
name of detsim producer
art::EDAnalyzer::Table< Config > Parameters
std::pair< sbn::crt::CRTTrack, double > ClosestCRTTrackByDCA(detinfo::DetectorPropertiesData const &detProp, recob::Track tpcTrack, std::vector< sbn::crt::CRTTrack > crtTracks, const art::Event &event, double minAngle=0.)
fhicl::Atom< art::InputTag > PandoraLabel
bool CrtTrackCosmicId(detinfo::DetectorPropertiesData const &detProp, recob::Track track, std::vector< sbn::crt::CRTTrack > crtTracks, const art::Event &event)
const Cut kCosmicRay
BEGIN_PROLOG vertical distance to the surface Name
std::map< std::string, TH1D * > hNumTrueMatches
walls no left
Definition: selectors.fcl:105
int TrueParticleIDFromTotalRecoHits(detinfo::DetectorClocksData const &clockData, const std::vector< art::Ptr< recob::Hit > > &hits, bool rollup_unsaved_ids=1)
std::map< std::string, TH1D * > hNoMatchDCA
BEGIN_PROLOG don t mess with this TPCTrackLabel
Definition of data types for geometry description.
Provides recob::Track data product.
fhicl::Table< CRTBackTracker::Config > CrtBackTrack
std::map< std::string, TH1D * > hLengthTag
std::map< std::string, TH1D * > hMatchAngle
fhicl::Table< CrtTrackCosmicIdAlg::Config > CTTagAlg
virtual void analyze(const art::Event &event) override
std::pair< sbn::crt::CRTTrack, double > ClosestCRTTrackByAngle(detinfo::DetectorPropertiesData const &detProp, recob::Track tpcTrack, std::vector< sbn::crt::CRTTrack > crtTracks, const art::Event &event, double minDCA=0.)
do i e
fhicl::Table< CRTTrackMatchAlg::Config > TrackMatchAlg
stream1 can override from command line with o or output services user sbnd
art::ServiceHandle< art::TFileService > tfs
std::map< std::string, TH1D * > hMatchDCA
art::InputTag fTPCTrackLabel
name of CRT producer
std::map< std::string, TH1D * > hLengthTotal
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:184
CRTTrackCosmicIdAna(Parameters const &config)
fhicl::Atom< art::InputTag > TPCTrackLabel
BEGIN_PROLOG could also be cout
auto const detProp
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a &quot;fitted&quot; track:
bool InFiducial(geo::Point_t point, double fiducial)
std::map< std::string, TH1D * > hNoMatchAngle