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