CMS 3D CMS Logo

TrackMVAClassifierBase.cc
Go to the documentation of this file.
2 
4 
7 
8 #include <cassert>
9 
11  desc.add<edm::InputTag>("src", edm::InputTag());
12  desc.add<edm::InputTag>("beamspot", edm::InputTag("offlineBeamSpot"));
13  desc.add<edm::InputTag>("vertices", edm::InputTag("firstStepPrimaryVertices"));
14  desc.add<bool>("ignoreVertices", false);
15  // default cuts for "cut based classification"
16  std::vector<double> cuts = {-.7, 0.1, .7};
17  desc.add<std::vector<double>>("qualityCuts", cuts);
18 }
19 
21 
23  : src_(consumes<reco::TrackCollection>(cfg.getParameter<edm::InputTag>("src"))),
24  beamspot_(consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamspot"))),
25  vertices_(mayConsume<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("vertices"))),
26  ignoreVertices_(cfg.getParameter<bool>("ignoreVertices")) {
27  auto const& qv = cfg.getParameter<std::vector<double>>("qualityCuts");
28  assert(qv.size() == 3);
29  std::copy(std::begin(qv), std::end(qv), std::begin(qualityCuts));
30 
31  produces<MVACollection>("MVAValues");
32  produces<QualityMaskCollection>("QualityMasks");
33 }
34 
36  // Get tracks
38  evt.getByToken(src_, hSrcTrack);
39  auto const& tracks(*hSrcTrack);
40 
41  // looking for the beam spot
43  evt.getByToken(beamspot_, hBsp);
44 
45  // Select good primary vertices for use in subsequent track selection
47  evt.getByToken(vertices_, hVtx);
48 
49  initEvent(es);
50 
51  // products
52  auto mvaPairs = std::make_unique<MVAPairCollection>(tracks.size(), std::make_pair(-99.f, true));
53  auto mvas = std::make_unique<MVACollection>(tracks.size(), -99.f);
54  auto quals = std::make_unique<QualityMaskCollection>(tracks.size(), 0);
55 
56  if (hVtx.isValid() && !ignoreVertices_) {
57  computeMVA(tracks, *hBsp, *hVtx, *mvaPairs);
58  } else {
59  if (!ignoreVertices_)
60  edm::LogWarning("TrackMVAClassifierBase")
61  << "ignoreVertices is set to False in the configuration, but the vertex collection is not valid";
62  std::vector<reco::Vertex> vertices;
63  computeMVA(tracks, *hBsp, vertices, *mvaPairs);
64  }
65  assert((*mvaPairs).size() == tracks.size());
66 
67  unsigned int k = 0;
68  for (auto const& output : *mvaPairs) {
69  if (output.second) {
70  (*mvas)[k] = output.first;
71  } else {
72  // If the MVA value is known to be unreliable, force into generalTracks collection
73  (*mvas)[k] = std::max(output.first, float(qualityCuts[0] + 0.001));
74  }
75  float mva = (*mvas)[k];
76  (*quals)[k++] = (mva > qualityCuts[0]) << reco::TrackBase::loose |
79  }
80 
81  evt.put(std::move(mvas), "MVAValues");
82  evt.put(std::move(quals), "QualityMasks");
83 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
TrackMVAClassifierBase(const edm::ParameterSet &cfg)
virtual void initEvent(const edm::EventSetup &es)=0
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
virtual void computeMVA(reco::TrackCollection const &tracks, reco::BeamSpot const &beamSpot, reco::VertexCollection const &vertices, MVAPairCollection &mvas) const =0
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
static void fill(edm::ParameterSetDescription &desc)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
assert(be >=bs)
edm::EDGetTokenT< reco::BeamSpot > beamspot_
void produce(edm::Event &evt, const edm::EventSetup &es) final
double f[11][100]
edm::EDGetTokenT< reco::VertexCollection > vertices_
bool isValid() const
Definition: HandleBase.h:70
fixed size matrix
HLT enums.
Definition: output.py:1
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511