File indexing completed on 2024-04-06 12:18:32
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "HLTrigger/JetMET/interface/HLTCaloTowerHtMhtProducer.h"
0011
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0015 #include "DataFormats/Common/interface/Handle.h"
0016
0017
0018 HLTCaloTowerHtMhtProducer::HLTCaloTowerHtMhtProducer(const edm::ParameterSet& iConfig)
0019 : usePt_(iConfig.getParameter<bool>("usePt")),
0020 minPtTowerHt_(iConfig.getParameter<double>("minPtTowerHt")),
0021 minPtTowerMht_(iConfig.getParameter<double>("minPtTowerMht")),
0022 maxEtaTowerHt_(iConfig.getParameter<double>("maxEtaTowerHt")),
0023 maxEtaTowerMht_(iConfig.getParameter<double>("maxEtaTowerMht")),
0024 towersLabel_(iConfig.getParameter<edm::InputTag>("towersLabel")) {
0025 m_theTowersToken = consumes<CaloTowerCollection>(towersLabel_);
0026
0027
0028 produces<reco::METCollection>();
0029 }
0030
0031
0032 HLTCaloTowerHtMhtProducer::~HLTCaloTowerHtMhtProducer() = default;
0033
0034
0035 void HLTCaloTowerHtMhtProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0036
0037 edm::ParameterSetDescription desc;
0038 desc.add<bool>("usePt", false);
0039 desc.add<double>("minPtTowerHt", 1.);
0040 desc.add<double>("minPtTowerMht", 1.);
0041 desc.add<double>("maxEtaTowerHt", 5.);
0042 desc.add<double>("maxEtaTowerMht", 5.);
0043 desc.add<edm::InputTag>("towersLabel", edm::InputTag("hltTowerMakerForAll"));
0044 descriptions.add("hltCaloTowerHtMhtProducer", desc);
0045 }
0046
0047
0048 void HLTCaloTowerHtMhtProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0049
0050 std::unique_ptr<reco::METCollection> result(new reco::METCollection());
0051
0052 edm::Handle<CaloTowerCollection> towers;
0053 iEvent.getByToken(m_theTowersToken, towers);
0054
0055 double ht = 0., mhx = 0., mhy = 0.;
0056
0057 if (!towers->empty()) {
0058 for (auto const& j : *towers) {
0059 double pt = usePt_ ? j.pt() : j.et();
0060 double eta = j.eta();
0061 double phi = j.phi();
0062 double px = usePt_ ? j.px() : j.et() * cos(phi);
0063 double py = usePt_ ? j.py() : j.et() * sin(phi);
0064
0065 if (pt > minPtTowerHt_ && std::abs(eta) < maxEtaTowerHt_) {
0066 ht += pt;
0067 }
0068
0069 if (pt > minPtTowerMht_ && std::abs(eta) < maxEtaTowerMht_) {
0070 mhx -= px;
0071 mhy -= py;
0072 }
0073 }
0074 }
0075
0076 reco::MET::LorentzVector p4(mhx, mhy, 0, sqrt(mhx * mhx + mhy * mhy));
0077 reco::MET::Point vtx(0, 0, 0);
0078 reco::MET htmht(ht, p4, vtx);
0079 result->push_back(htmht);
0080
0081
0082 iEvent.put(std::move(result));
0083 }