CMS 3D CMS Logo

Phase1L1TJetSumsProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1CaloTrigger
4 // Class: Phase1L1TJetSumsProducer
5 //
23 //
24 // Original Simone Bologna
25 //
26 
39 
40 #include <cmath>
41 
42 class Phase1L1TJetSumsProducer : public edm::one::EDProducer<edm::one::SharedResources> {
43 public:
45  ~Phase1L1TJetSumsProducer() override = default;
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48 
49 private:
50  void produce(edm::Event&, const edm::EventSetup&) override;
51 
52  // computes ht, adds jet pt to ht only if the pt of the jet is above the ht calculation threshold
53  l1t::EtSum computeHT(const edm::Handle<std::vector<reco::CaloJet> > inputJets) const;
54 
55  // computes MHT
56  // adds jet pt to mht only if the pt of the jet is above the mht calculation threshold
57  // performs some calculations with digitised/integer quantities to ensure agreement with firmware
58  l1t::EtSum computeMHT(const edm::Handle<std::vector<reco::CaloJet> > inputJets) const;
59 
61 
62  // holds the sin and cos for HLs LUT emulation
63  const std::vector<double> sinPhi_;
64  const std::vector<double> cosPhi_;
65  const unsigned int nBinsPhi_;
66 
67  // lower phi value
68  const double phiLow_;
69  // higher phi value
70  const double phiUp_;
71  // size of a phi bin
72  const double phiStep_;
73  // threshold for ht calculation
74  const double htPtThreshold_;
75  // threshold for ht calculation
76  const double mhtPtThreshold_;
77  // jet eta cut for ht calculation
78  const double htAbsEtaCut_;
79  // jet eta cut for mht calculation
80  const double mhtAbsEtaCut_;
81  // LSB of pt quantity
82  const double ptlsb_;
83  // label of sums
85 };
86 
87 // initialises plugin configuration and prepares ROOT file for saving the sums
89  : inputJetCollectionTag_{consumes<std::vector<reco::CaloJet> >(
90  iConfig.getParameter<edm::InputTag>("inputJetCollectionTag"))},
91  sinPhi_(iConfig.getParameter<std::vector<double> >("sinPhi")),
92  cosPhi_(iConfig.getParameter<std::vector<double> >("cosPhi")),
93  nBinsPhi_(iConfig.getParameter<unsigned int>("nBinsPhi")),
94  phiLow_(iConfig.getParameter<double>("phiLow")),
95  phiUp_(iConfig.getParameter<double>("phiUp")),
96  phiStep_((phiUp_ - phiLow_) / nBinsPhi_),
97  htPtThreshold_(iConfig.getParameter<double>("htPtThreshold")),
98  mhtPtThreshold_(iConfig.getParameter<double>("mhtPtThreshold")),
99  htAbsEtaCut_(iConfig.getParameter<double>("htAbsEtaCut")),
100  mhtAbsEtaCut_(iConfig.getParameter<double>("mhtAbsEtaCut")),
101  ptlsb_(iConfig.getParameter<double>("ptlsb")),
102  outputCollectionName_(iConfig.getParameter<std::string>("outputCollectionName")) {
103  usesResource("TFileService");
104  produces<std::vector<l1t::EtSum> >(outputCollectionName_).setBranchAlias(outputCollectionName_);
105 }
106 
108  const edm::Handle<std::vector<reco::CaloJet> >& jetCollectionHandle = iEvent.getHandle(inputJetCollectionTag_);
109 
110  // computing sums and storing them in sum object
111  l1t::EtSum lHT = computeHT(jetCollectionHandle);
112  l1t::EtSum lMHT = computeMHT(jetCollectionHandle);
113 
114  //packing sums in vector for event saving
115  std::unique_ptr<std::vector<l1t::EtSum> > lSumVectorPtr(new std::vector<l1t::EtSum>(0));
116  lSumVectorPtr->push_back(lHT);
117  lSumVectorPtr->push_back(lMHT);
118  iEvent.put(std::move(lSumVectorPtr), outputCollectionName_);
119 
120  return;
121 }
122 
123 l1t::EtSum Phase1L1TJetSumsProducer::computeHT(const edm::Handle<std::vector<reco::CaloJet> > inputJets) const {
124  double lHT = 0;
125  for (const auto& jet : *inputJets) {
126  double lJetPt = jet.pt();
127  double lJetPhi = jet.phi();
128  double lJetEta = jet.eta();
129  if ((lJetPhi < phiLow_) || (lJetPhi >= phiUp_))
130  continue;
131 
132  lHT += (lJetPt >= htPtThreshold_ && std::fabs(lJetEta) < htAbsEtaCut_) ? lJetPt : 0;
133  }
134 
136  lHTVector.SetPt(lHT);
137  lHTVector.SetEta(0);
138  lHTVector.SetPhi(0);
139  l1t::EtSum lHTSum(lHTVector, l1t::EtSum::EtSumType::kTotalHt, 0, 0, 0, 0);
140  return lHTSum;
141 }
142 
143 l1t::EtSum Phase1L1TJetSumsProducer::computeMHT(const edm::Handle<std::vector<reco::CaloJet> > inputJets) const {
144  int lTotalJetPx = 0;
145  int lTotalJetPy = 0;
146 
147  std::vector<unsigned int> jetPtInPhiBins(nBinsPhi_, 0);
148 
149  for (const auto& jet : *inputJets) {
150  double lJetPhi = jet.phi();
151 
152  if ((lJetPhi < phiLow_) || (lJetPhi >= phiUp_))
153  continue;
154 
155  unsigned int iPhi = (lJetPhi - phiLow_) / phiStep_;
156 
157  if (jet.pt() >= mhtPtThreshold_ && std::fabs(jet.eta()) < mhtAbsEtaCut_) {
158  unsigned int digiJetPt = floor(jet.pt() / ptlsb_);
159  jetPtInPhiBins[iPhi] += digiJetPt;
160  }
161  }
162 
163  for (unsigned int iPhi = 0; iPhi < jetPtInPhiBins.size(); ++iPhi) {
164  unsigned int digiJetPtSum = jetPtInPhiBins[iPhi];
165 
166  // retrieving sin cos from LUT emulator
167  double lSinPhi = sinPhi_[iPhi];
168  double lCosPhi = cosPhi_[iPhi];
169 
170  // checking if above threshold
171  lTotalJetPx += trunc(digiJetPtSum * lCosPhi);
172  lTotalJetPy += trunc(digiJetPtSum * lSinPhi);
173  }
174 
175  double lMHT = floor(sqrt(lTotalJetPx * lTotalJetPx + lTotalJetPy * lTotalJetPy)) * ptlsb_;
176  math::PtEtaPhiMLorentzVector lMHTVector(lMHT, 0, acos(lTotalJetPx / (lMHT / ptlsb_)), 0);
177  l1t::EtSum lMHTSum(lMHTVector, l1t::EtSum::EtSumType::kMissingHt, 0, 0, 0, 0);
178 
179  return lMHTSum;
180 }
181 
184  desc.add<edm::InputTag>("inputJetCollectionTag",
185  edm::InputTag("l1tPhase1JetCalibrator", "Phase1L1TJetFromPfCandidates"));
186  desc.add<std::vector<double> >("sinPhi");
187  desc.add<std::vector<double> >("cosPhi");
188  desc.add<unsigned int>("nBinsPhi", 72);
189  desc.add<double>("phiLow", -M_PI);
190  desc.add<double>("phiUp", M_PI);
191  desc.add<double>("htPtThreshold", 30);
192  desc.add<double>("mhtPtThreshold", 30);
193  desc.add<double>("htAbsEtaCut", 3);
194  desc.add<double>("mhtAbsEtaCut", 3);
195  desc.add<double>("ptlsb", 0.25), desc.add<string>("outputCollectionName", "Sums");
196  descriptions.add("Phase1L1TJetSumsProducer", desc);
197 }
198 
199 // creates the plugin for later use in python
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
l1t::EtSum computeHT(const edm::Handle< std::vector< reco::CaloJet > > inputJets) const
const std::vector< double > sinPhi_
l1t::EtSum computeMHT(const edm::Handle< std::vector< reco::CaloJet > > inputJets) const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
int iEvent
Definition: GenABIO.cc:224
const std::vector< double > cosPhi_
T sqrt(T t)
Definition: SSEVec.h:19
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Phase1L1TJetSumsProducer(const edm::ParameterSet &)
const edm::EDGetTokenT< std::vector< reco::CaloJet > > inputJetCollectionTag_
#define M_PI
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~Phase1L1TJetSumsProducer() override=default
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &, const edm::EventSetup &) override
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38