CMS 3D CMS Logo

HLTL1TMuonSelector.cc
Go to the documentation of this file.
1 //-------------------------------------------------
2 //
13 //
14 //--------------------------------------------------
15 
16 // Class Header
17 #include "HLTL1TMuonSelector.h"
18 
19 // Framework
25 
26 using namespace std;
27 using namespace edm;
28 using namespace l1t;
29 
30 // constructors
32  : theSource(iConfig.getParameter<InputTag>("InputObjects")),
33  theL1MinPt(iConfig.getParameter<double>("L1MinPt")),
34  theL1MaxEta(iConfig.getParameter<double>("L1MaxEta")),
35  theL1MinQuality(iConfig.getParameter<unsigned int>("L1MinQuality")),
36  centralBxOnly_(iConfig.getParameter<bool>("CentralBxOnly")) {
37  muCollToken_ = consumes<MuonBxCollection>(theSource);
38 
39  produces<MuonBxCollection>();
40 }
41 
42 // destructor
44 
47  desc.add<edm::InputTag>("InputObjects", edm::InputTag("hltGmtStage2Digis"));
48  desc.add<double>("L1MinPt", -1.);
49  desc.add<double>("L1MaxEta", 5.0);
50  desc.add<unsigned int>("L1MinQuality", 0);
51  desc.add<bool>("CentralBxOnly", true);
52  descriptions.add("hltL1TMuonSelector", desc);
53 }
54 
56  const std::string metname = "Muon|RecoMuon|HLTL1TMuonSelector";
57 
58  unique_ptr<MuonBxCollection> output(new MuonBxCollection());
59 
60  // Muon particles
62  iEvent.getByToken(muCollToken_, muColl);
63  LogTrace(metname) << "Number of muons " << muColl->size() << endl;
64 
65  for (int ibx = muColl->getFirstBX(); ibx <= muColl->getLastBX(); ++ibx) {
66  if (centralBxOnly_ && (ibx != 0))
67  continue;
68  for (auto it = muColl->begin(ibx); it != muColl->end(ibx); it++) {
69  unsigned int quality = it->hwQual();
70  int valid_charge = it->hwChargeValid();
71 
72  float pt = it->pt();
73  float eta = it->eta();
74  float theta = 2 * atan(exp(-eta));
75  float phi = it->phi();
76  int charge = it->charge();
77  // Set charge=0 for the time being if the valid charge bit is zero
78  if (!valid_charge)
79  charge = 0;
80 
81  if (pt < theL1MinPt || fabs(eta) > theL1MaxEta)
82  continue;
83 
84  LogTrace(metname) << "L1 Muon Found";
85  LogTrace(metname) << "Pt = " << pt << " GeV/c";
86  LogTrace(metname) << "eta = " << eta;
87  LogTrace(metname) << "theta = " << theta << " rad";
88  LogTrace(metname) << "phi = " << phi << " rad";
89  LogTrace(metname) << "charge = " << charge;
90 
91  if (quality <= theL1MinQuality)
92  continue;
93  LogTrace(metname) << "quality = " << quality;
94 
95  output->push_back(ibx, *it);
96  }
97  }
98 
99  iEvent.put(std::move(output));
100 }
101 
102 // declare this class as a framework plugin
int getLastBX() const
HLTL1TMuonSelector(const edm::ParameterSet &)
Constructor.
int getFirstBX() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const std::string metname
delete x;
Definition: CaloConfig.h:22
#define LogTrace(id)
const_iterator begin(int bx) const
unsigned size(int bx) const
const unsigned theL1MinQuality
string quality
edm::InputTag theSource
edm::EDGetTokenT< l1t::MuonBxCollection > muCollToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BXVector< Muon > MuonBxCollection
Definition: Muon.h:11
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~HLTL1TMuonSelector() override
Destructor.
bool centralBxOnly_
use central bx only muons
const_iterator end(int bx) const
HLT enums.
const double theL1MaxEta
Geom::Theta< T > theta() const
def move(src, dest)
Definition: eostools.py:511