CMS 3D CMS Logo

HLTScoutingTrackProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTrigger/Muon
4 // Class: HLTScoutingTrackProducer
5 //
9 //
10 
11 #include <memory>
34 
36 public:
38  ~HLTScoutingTrackProducer() override;
39 
40  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
41 
42 private:
43  void produce(edm::StreamID sid, edm::Event& iEvent, edm::EventSetup const& setup) const final;
44 
47 
48  const int mantissaPrecision_;
49  const double vtxMinDist_;
50  const double ptMin_;
51 };
52 
53 //
54 // constructors and destructor
55 //
57  : otherTrackCollection_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("OtherTracks"))),
58  vertexCollection_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"))),
59  mantissaPrecision_(iConfig.getParameter<int>("mantissaPrecision")),
60  vtxMinDist_(iConfig.getParameter<double>("vtxMinDist")),
61  ptMin_(iConfig.getParameter<double>("ptMin")) {
62  //register products
63  produces<Run3ScoutingTrackCollection>();
64 }
65 
67 
68 // ------------ method called to produce the data ------------
70  using namespace edm;
71 
72  std::unique_ptr<Run3ScoutingTrackCollection> outTrack(new Run3ScoutingTrackCollection());
73 
74  Handle<reco::TrackCollection> otherTrackCollection;
76 
77  if (iEvent.getByToken(otherTrackCollection_, otherTrackCollection)) {
78  //match tracks to vertices
79  for (auto& trk : *otherTrackCollection) {
80  int vtxInd = -1;
81  double min_dist = vtxMinDist_;
82  int vtxIt = 0;
83 
84  if (trk.pt() < ptMin_)
85  continue;
86 
87  if (iEvent.getByToken(vertexCollection_, vertexCollection)) {
88  for (auto& vrt : *vertexCollection) {
89  double min_dist_tmp = pow(trk.dz(vrt.position()), 2); // hltPixelVertices only clustered in Z
90 
91  if (min_dist_tmp < min_dist) {
92  min_dist = min_dist_tmp;
93  vtxInd = vtxIt;
94  }
95 
96  vtxIt++;
97  }
98  }
99 
100  //fill track information
101  outTrack->emplace_back(
106  trk.ndof(),
107  trk.charge(),
110  trk.hitPattern().numberOfValidPixelHits(),
111  trk.hitPattern().trackerLayersWithMeasurement(),
112  trk.hitPattern().numberOfValidStripHits(),
132  vtxInd,
136  }
137  }
138 
139  iEvent.put(std::move(outTrack));
140 }
141 
142 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
145  desc.add<edm::InputTag>("OtherTracks", edm::InputTag("hltPixelTracksL3MuonNoVtx"));
146  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
147 
148  desc.add<int>("mantissaPrecision", 10)->setComment("default float16, change to 23 for float32");
149  desc.add<double>("vtxMinDist", 0.01);
150  desc.add<double>("ptMin", 0.3);
151  descriptions.add("hltScoutingTrackProducer", desc);
152 }
153 
154 // declare this class as a framework plugin
const edm::EDGetTokenT< reco::TrackCollection > otherTrackCollection_
const edm::EDGetTokenT< reco::VertexCollection > vertexCollection_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
void produce(edm::StreamID sid, edm::Event &iEvent, edm::EventSetup const &setup) const final
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
int iEvent
Definition: GenABIO.cc:224
~HLTScoutingTrackProducer() override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< Run3ScoutingTrack > Run3ScoutingTrackCollection
static float reduceMantissaToNbitsRounding(const float &f)
Definition: libminifloat.h:79
fixed size matrix
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLTScoutingTrackProducer(const edm::ParameterSet &)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511