All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BeamFlashTrackMatchTaggerAlg.h
Go to the documentation of this file.
1 #ifndef BEAMFLASHTRACKMATCHTAGGERALG_H
2 #define BEAMFLASHTRACKMATCHTAGGERALG_H
3 /*!
4  * Title: Beam Flash<-->Track Match Algorithim Class
5  * Author: Wes Ketchum (wketchum@lanl.gov), based on code from Ben Jones
6  *
7  * Description: Algorithm that compares all tracks to the flash during the
8  * beam gate, and determines if that track is consistent with
9  * having produced that flash.
10  * Input: recob::OpFlash, recob::Track
11  * Output: anab::CosmicTag (and Assn<anab::CosmicTag,recob::Track>)
12 */
13 #include <iostream>
14 #include <string>
15 #include <vector>
16 
17 namespace fhicl { class ParameterSet; }
18 namespace geo { class GeometryCore; }
19 namespace detinfo { class LArProperties; }
20 namespace phot { class PhotonVisibilityService; }
21 namespace opdet { class OpDigiProperties; }
22 
27 
28 #include "nusimdata/SimulationBase/MCParticle.h"
29 
30 class TVector3;
31 class TH1F;
32 class TTree;
33 
34 namespace cosmic{
35  class BeamFlashTrackMatchTaggerAlg;
36 }
37 
38 
40  public:
41 
42  /// Pack of provider-interface supporting services we need
44 
45  BeamFlashTrackMatchTaggerAlg(fhicl::ParameterSet const& p);
46 
47  //how to run the algorithm
48  void RunCompatibilityCheck(std::vector<recob::OpFlash> const&,
49  std::vector<recob::Track> const&,
50  std::vector<anab::CosmicTag>&,
51  std::vector<size_t>&,
55 
56  void SetHypothesisComparisonTree(TTree*,TH1F*,TH1F*);
57 
58  void RunHypothesisComparison(unsigned int const,
59  unsigned int const,
60  std::vector<recob::OpFlash> const&,
61  std::vector<recob::Track> const&,
65 
66  void RunHypothesisComparison(unsigned int const,
67  unsigned int const,
68  std::vector<recob::OpFlash> const&,
69  std::vector<simb::MCParticle> const&,
73 
74  private:
75 
78  const bool DEBUG_FLAG;
79 
81  float fMinOpHitPE;
82  float fMIPdQdx;
86  unsigned int fCumulativeChannelCut;
87  float fIntegralCut;
88 
91 
92 
93  TTree* cTree;
94 
95  typedef struct FlashComparisonProperties{
96  unsigned int run;
97  unsigned int event;
98 
99  unsigned int flash_index;
101  float flash_y;
103  float flash_z;
105  unsigned int flash_nOpDet;
106 
107  unsigned int hyp_index;
108  float hyp_totalPE;
109  float hyp_y;
110  float hyp_sigmay;
111  float hyp_z;
112  float hyp_sigmaz;
113 
114  float trk_startx;
115  float trk_starty;
116  float trk_startz;
117 
118  float trk_endx;
119  float trk_endy;
120  float trk_endz;
121 
122  float chi2;
123 
124  std::string leaf_structure;
126  leaf_structure("run/i:event/i:flash_index/i:flash_totalPE/F:flash_y/F:flash_sigmay/F:flash_z/F:flash_sigmaz/F:flash_nOpDet/i:hyp_index/i:hyp_totalPE/F:hyp_y/F:hyp_sigmay/F:hyp_z/F:hyp_sigmaz/F:trk_startx/F:trk_starty/F:trk_startz/F:trk_endx/F:trk_endy/F:trk_endz/F:chi2/F"){}
127 
129 
131  std::vector<float> cOpDetVector_flash;
132  std::vector<float> cOpDetVector_hyp;
135 
136 
143 
144  //core functions
145  std::vector<float> GetMIPHypotheses(recob::Track const& track,
146  Providers_t providers,
149  float XOffset=0);
150 
151  std::vector<float> GetMIPHypotheses(simb::MCParticle const& particle,
152  size_t start_i,
153  size_t end_i,
154  Providers_t providers,
157  float XOffset=0);
158 
159  void AddLightFromSegment(TVector3 const& pt1,
160  TVector3 const& pt2,
161  std::vector<float> & lightHypothesis,
162  float & totalHypothesisPE,
163  geo::GeometryCore const& geom,
165  float const& PromptMIPScintYield,
166  float XOffset);
167 
168  void NormalizeLightHypothesis(std::vector<float> & lightHypothesis,
169  float const& totalHypothesisPE,
170  geo::GeometryCore const& geom);
171 
172  CompatibilityResultType CheckCompatibility(std::vector<float> const& lightHypothesis,
173  const recob::OpFlash* flashPointer,
174  geo::GeometryCore const& geom);
175 
176  bool InDetector(TVector3 const&, geo::GeometryCore const&);
177  bool InDriftWindow(double, double, geo::GeometryCore const&);
178 
179  void FillFlashProperties(std::vector<float> const& opdetVector,
180  float&,
181  float&, float&,
182  float&, float&,
183  geo::GeometryCore const& geom);
184 
185  float CalculateChi2(std::vector<float> const&,std::vector<float> const&);
186 
187  //debugging functions
188  void PrintTrackProperties(recob::Track const&, std::ostream* output=&std::cout);
189  void PrintFlashProperties(recob::OpFlash const&, std::ostream* output=&std::cout);
190  void PrintHypothesisFlashComparison(std::vector<float> const&,
191  const recob::OpFlash*,
192  geo::GeometryCore const& geom,
194  std::ostream* output=&std::cout);
195 
196 };
197 
198 #endif
CompatibilityResultType CheckCompatibility(std::vector< float > const &lightHypothesis, const recob::OpFlash *flashPointer, geo::GeometryCore const &geom)
BeamFlashTrackMatchTaggerAlg(fhicl::ParameterSet const &p)
enum anab::cosmic_tag_id CosmicTagID_t
pdgs p
Definition: selectors.fcl:22
void RunHypothesisComparison(unsigned int const, unsigned int const, std::vector< recob::OpFlash > const &, std::vector< recob::Track > const &, Providers_t, phot::PhotonVisibilityService const &, opdet::OpDigiProperties const &)
float CalculateChi2(std::vector< float > const &, std::vector< float > const &)
process_name use argoneut_mc_hitfinder track
struct cosmic::BeamFlashTrackMatchTaggerAlg::FlashComparisonProperties FlashComparisonProperties_t
void RunCompatibilityCheck(std::vector< recob::OpFlash > const &, std::vector< recob::Track > const &, std::vector< anab::CosmicTag > &, std::vector< size_t > &, Providers_t, phot::PhotonVisibilityService const &, opdet::OpDigiProperties const &)
lar::ProviderPack< geo::GeometryCore, detinfo::LArProperties > Providers_t
Pack of provider-interface supporting services we need.
std::vector< float > GetMIPHypotheses(recob::Track const &track, Providers_t providers, phot::PhotonVisibilityService const &pvs, opdet::OpDigiProperties const &, float XOffset=0)
void NormalizeLightHypothesis(std::vector< float > &lightHypothesis, float const &totalHypothesisPE, geo::GeometryCore const &geom)
void AddLightFromSegment(TVector3 const &pt1, TVector3 const &pt2, std::vector< float > &lightHypothesis, float &totalHypothesisPE, geo::GeometryCore const &geom, phot::PhotonVisibilityService const &pvs, float const &PromptMIPScintYield, float XOffset)
void PrintTrackProperties(recob::Track const &, std::ostream *output=&std::cout)
Description of geometry of one entire detector.
void SetHypothesisComparisonTree(TTree *, TH1F *, TH1F *)
Provides recob::Track data product.
void PrintHypothesisFlashComparison(std::vector< float > const &, const recob::OpFlash *, geo::GeometryCore const &geom, CompatibilityResultType, std::ostream *output=&std::cout)
Container for a list of pointers to providers.
Definition: ProviderPack.h:114
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG simSlidingORM6O6 effSlidingORW output
BEGIN_PROLOG cosmic
bool InDriftWindow(double, double, geo::GeometryCore const &)
void FillFlashProperties(std::vector< float > const &opdetVector, float &, float &, float &, float &, float &, geo::GeometryCore const &geom)
Data structure containing constant pointers to classes.
bool InDetector(TVector3 const &, geo::GeometryCore const &)
BEGIN_PROLOG could also be cout
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a &quot;fitted&quot; track:
void PrintFlashProperties(recob::OpFlash const &, std::ostream *output=&std::cout)