CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 };
51 
52 //
53 // constructors and destructor
54 //
56  : otherTrackCollection_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("OtherTracks"))),
57  vertexCollection_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexCollection"))),
58  mantissaPrecision(iConfig.getParameter<int>("mantissaPrecision")),
59  vtxMinDist(iConfig.getParameter<double>("vtxMinDist")) {
60  //register products
61  produces<Run3ScoutingTrackCollection>();
62 }
63 
65 
66 // ------------ method called to produce the data ------------
68  using namespace edm;
69 
70  std::unique_ptr<Run3ScoutingTrackCollection> outTrack(new Run3ScoutingTrackCollection());
71 
72  Handle<reco::TrackCollection> otherTrackCollection;
74 
75  if (iEvent.getByToken(otherTrackCollection_, otherTrackCollection)) {
76  //match tracks to vertices
77  for (auto& trk : *otherTrackCollection) {
78  int vtxInd = -1;
79  double min_dist = vtxMinDist;
80  int vtxIt = 0;
81 
83  for (auto& vrt : *vertexCollection) {
84  double min_dist_tmp = pow(trk.dz(vrt.position()), 2); // hltPixelVertices only clustered in Z
85 
86  if (min_dist_tmp < min_dist) {
87  min_dist = min_dist_tmp;
88  vtxInd = vtxIt;
89  }
90 
91  vtxIt++;
92  }
93  }
94 
95  //fill track information
100  trk.ndof(),
101  trk.charge(),
104  trk.hitPattern().numberOfValidPixelHits(),
105  trk.hitPattern().trackerLayersWithMeasurement(),
106  trk.hitPattern().numberOfValidStripHits(),
126  vtxInd,
130  }
131  }
132 
133  iEvent.put(std::move(outTrack));
134 }
135 
136 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
139  desc.add<edm::InputTag>("OtherTracks", edm::InputTag("hltPixelTracksL3MuonNoVtx"));
140  desc.add<edm::InputTag>("vertexCollection", edm::InputTag("hltPixelVertices"));
141 
142  desc.add<int>("mantissaPrecision", 10)->setComment("default float16, change to 23 for float32");
143  desc.add<double>("vtxMinDist", 0.01);
144  descriptions.add("hltScoutingTrackProducer", desc);
145 }
146 
147 // declare this class as a framework plugin
const edm::EDGetTokenT< reco::TrackCollection > otherTrackCollection_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple mantissaPrecision
const edm::EDGetTokenT< reco::VertexCollection > vertexCollection_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#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:12
tuple vertexCollection
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
~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
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