CMS 3D CMS Logo

TrackMVAClassifier.h
Go to the documentation of this file.
1 #ifndef RecoTracker_FinalTrackSelectors_TrackMVAClassifierBase_h
2 #define RecoTracker_FinalTrackSelectors_TrackMVAClassifierBase_h
3 
4 
8 
9 
15 
17 
19 
20 #include <vector>
21 #include <memory>
22 
24 public:
25  explicit TrackMVAClassifierBase( const edm::ParameterSet & cfg );
26  ~TrackMVAClassifierBase() override;
27 protected:
28 
29  static void fill( edm::ParameterSetDescription& desc);
30 
31 
32  using MVACollection = std::vector<float>;
33  using QualityMaskCollection = std::vector<unsigned char>;
34 
35  virtual void initEvent(const edm::EventSetup& es) = 0;
36 
37  virtual void computeMVA(reco::TrackCollection const & tracks,
38  reco::BeamSpot const & beamSpot,
40  MVACollection & mvas) const = 0;
41 
42 private:
43  void produce(edm::Event& evt, const edm::EventSetup& es ) final;
44 
49 
51 
52  // MVA
53 
54  // qualitycuts (loose, tight, hp)
55  float qualityCuts[3];
56 
57 };
58 
59 template<typename MVA>
61 public:
64  mva(cfg.getParameter<edm::ParameterSet>("mva")){}
65 
66  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
68  fill(desc);
70  MVA::fillDescriptions(mvaDesc);
71  desc.add<edm::ParameterSetDescription>("mva",mvaDesc);
72  descriptions.add(MVA::name(), desc);
73  }
74 
75 
76 private:
77  void beginStream(edm::StreamID) final {
78  mva.beginStream();
79  }
80 
81  void initEvent(const edm::EventSetup& es) final {
82  mva.initEvent(es);
83  }
84 
86  reco::BeamSpot const & beamSpot,
88  MVACollection & mvas) const final {
89 
90  size_t current = 0;
91  for (auto const & trk : tracks) {
92  mvas[current++]= mva(trk,beamSpot,vertices);
93  }
94  }
95 
97 };
98 
99 
100 
101 #endif // RecoTracker_FinalTrackSelectors_TrackMVAClassifierBase_h
102 
TrackMVAClassifierBase(const edm::ParameterSet &cfg)
virtual void initEvent(const edm::EventSetup &es)=0
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
static void fill(edm::ParameterSetDescription &desc)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< float > MVACollection
virtual void computeMVA(reco::TrackCollection const &tracks, reco::BeamSpot const &beamSpot, reco::VertexCollection const &vertices, MVACollection &mvas) const =0
edm::EDGetTokenT< reco::BeamSpot > beamspot_
void produce(edm::Event &evt, const edm::EventSetup &es) final
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::VertexCollection > vertices_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void beginStream(edm::StreamID) final
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void computeMVA(reco::TrackCollection const &tracks, reco::BeamSpot const &beamSpot, reco::VertexCollection const &vertices, MVACollection &mvas) const final
HLT enums.
void initEvent(const edm::EventSetup &es) final
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< unsigned char > QualityMaskCollection
TrackMVAClassifier(const edm::ParameterSet &cfg)