File indexing completed on 2024-09-07 04:37:02
0001
0002
0003
0004
0005 #include <memory>
0006
0007
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/stream/EDProducer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "DataFormats/Math/interface/LorentzVector.h"
0015 #include "DataFormats/L1Trigger/interface/Vertex.h"
0016 #include "DataFormats/L1TCorrelator/interface/TkHTMiss.h"
0017 #include "DataFormats/L1TCorrelator/interface/TkHTMissFwd.h"
0018
0019 using namespace l1t;
0020
0021 class L1TkHTMissProducer : public edm::stream::EDProducer<> {
0022 public:
0023 explicit L1TkHTMissProducer(const edm::ParameterSet&);
0024 ~L1TkHTMissProducer() override;
0025
0026 private:
0027 virtual void beginJob();
0028 void produce(edm::Event&, const edm::EventSetup&) override;
0029 virtual void endJob();
0030
0031
0032 const float jetMinPt_;
0033 const float jetMaxEta_;
0034 const bool doVtxConstrain_;
0035 const bool useCaloJets_;
0036 const bool primaryVtxConstrain_;
0037 const float deltaZ_;
0038 const unsigned int minNtracksHighPt_;
0039 const unsigned int minNtracksLowPt_;
0040 const float minJetEtLowPt_;
0041 const float minJetEtHighPt_;
0042 const bool displaced_;
0043 const edm::EDGetTokenT<VertexCollection> pvToken_;
0044 const edm::EDGetTokenT<TkJetCollection> jetToken_;
0045 };
0046
0047 L1TkHTMissProducer::L1TkHTMissProducer(const edm::ParameterSet& iConfig)
0048 : jetMinPt_((float)iConfig.getParameter<double>("jet_minPt")),
0049 jetMaxEta_((float)iConfig.getParameter<double>("jet_maxEta")),
0050 doVtxConstrain_(iConfig.getParameter<bool>("doVtxConstrain")),
0051 useCaloJets_(iConfig.getParameter<bool>("useCaloJets")),
0052 primaryVtxConstrain_(iConfig.getParameter<bool>("primaryVtxConstrain")),
0053 deltaZ_((float)iConfig.getParameter<double>("deltaZ")),
0054 minNtracksHighPt_(iConfig.getParameter<int>("jet_minNtracksHighPt")),
0055 minNtracksLowPt_(iConfig.getParameter<int>("jet_minNtracksLowPt")),
0056 minJetEtLowPt_(iConfig.getParameter<double>("jet_minJetEtLowPt")),
0057 minJetEtHighPt_(iConfig.getParameter<double>("jet_minJetEtHighPt")),
0058 displaced_(iConfig.getParameter<bool>("displaced")),
0059 pvToken_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("L1VertexInputTag"))),
0060 jetToken_(consumes<TkJetCollection>(iConfig.getParameter<edm::InputTag>("L1TkJetInputTag"))) {
0061 if (useCaloJets_)
0062 produces<TkHTMissCollection>("TkCaloHTMiss");
0063 else if (displaced_)
0064 produces<TkHTMissCollection>("L1TrackerHTMissExtended");
0065 else
0066 produces<TkHTMissCollection>("L1TrackerHTMiss");
0067 }
0068
0069 L1TkHTMissProducer::~L1TkHTMissProducer() {}
0070
0071 void L1TkHTMissProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0072 using namespace edm;
0073 std::unique_ptr<TkHTMissCollection> MHTCollection(new TkHTMissCollection);
0074
0075
0076 edm::Handle<VertexCollection> L1VertexHandle;
0077 iEvent.getByToken(pvToken_, L1VertexHandle);
0078
0079
0080 edm::Handle<TkJetCollection> L1TkJetsHandle;
0081 iEvent.getByToken(jetToken_, L1TkJetsHandle);
0082 std::vector<TkJet>::const_iterator jetIter;
0083
0084 if (!L1TkJetsHandle.isValid() && !displaced_) {
0085 LogError("TkHTMissProducer") << "\nWarning: TkJetCollection not found in the event. Exit\n";
0086 return;
0087 }
0088
0089 if (!L1TkJetsHandle.isValid() && displaced_) {
0090 LogError("TkHTMissProducer") << "\nWarning: TkJetExtendedCollection not found in the event. Exit\n";
0091 return;
0092 }
0093
0094
0095
0096
0097 float evtZVtx = 999;
0098 bool foundVtx = false;
0099 edm::Ref<VertexCollection> L1VtxRef;
0100
0101 if (useCaloJets_) {
0102 if (doVtxConstrain_ && primaryVtxConstrain_) {
0103 if (!L1VertexHandle.isValid()) {
0104 LogError("L1TkHTMissProducer") << "\nWarning: VertexCollection not found in the event. Exit\n";
0105 return;
0106 } else {
0107 std::vector<Vertex>::const_iterator vtxIter = L1VertexHandle->begin();
0108
0109
0110 evtZVtx = vtxIter->z0();
0111 foundVtx = true;
0112 int ivtx = 0;
0113 edm::Ref<VertexCollection> vtxRef(L1VertexHandle, ivtx);
0114 L1VtxRef = vtxRef;
0115 }
0116 }
0117
0118
0119
0120
0121 float zvtx_jetpt = -1.0;
0122 float jetVtxMax = 99.;
0123
0124 if (doVtxConstrain_ && !primaryVtxConstrain_) {
0125 for (jetIter = L1TkJetsHandle->begin(); jetIter != L1TkJetsHandle->end(); ++jetIter) {
0126 int ibx = jetIter->bx();
0127 if (ibx != 0)
0128 continue;
0129
0130 float tmp_jet_vtx = jetIter->jetVtx();
0131 float tmp_jet_pt = jetIter->pt();
0132 float tmp_jet_eta = jetIter->eta();
0133 if (tmp_jet_pt < jetMinPt_)
0134 continue;
0135 if (std::abs(tmp_jet_eta) > jetMaxEta_)
0136 continue;
0137 if (std::abs(tmp_jet_vtx) > jetVtxMax)
0138 continue;
0139
0140
0141 if (tmp_jet_pt > zvtx_jetpt) {
0142 evtZVtx = tmp_jet_vtx;
0143 zvtx_jetpt = tmp_jet_pt;
0144 foundVtx = true;
0145 }
0146 }
0147 }
0148
0149 float sumPx_calo = 0;
0150 float sumPy_calo = 0;
0151 float HT_calo = 0;
0152
0153 if (doVtxConstrain_ && !foundVtx)
0154 LogWarning("L1TkHTMissProducer") << "Didn't find any z vertex (based on jet vertices) for this event!\n";
0155
0156
0157 for (jetIter = L1TkJetsHandle->begin(); jetIter != L1TkJetsHandle->end(); ++jetIter) {
0158 int ibx = jetIter->bx();
0159 if (ibx != 0)
0160 continue;
0161
0162 float tmp_jet_px = jetIter->px();
0163 float tmp_jet_py = jetIter->py();
0164 float tmp_jet_et = jetIter->et();
0165 float tmp_jet_vtx = jetIter->jetVtx();
0166 if (jetIter->pt() < jetMinPt_)
0167 continue;
0168 if (std::abs(jetIter->eta()) > jetMaxEta_)
0169 continue;
0170
0171
0172 bool VtxRequirement = false;
0173 if (foundVtx)
0174 VtxRequirement = std::abs(tmp_jet_vtx - evtZVtx) < deltaZ_;
0175
0176 if (!doVtxConstrain_ || VtxRequirement) {
0177 sumPx_calo += tmp_jet_px;
0178 sumPy_calo += tmp_jet_py;
0179 HT_calo += tmp_jet_et;
0180 }
0181 }
0182
0183
0184 float et = sqrt(sumPx_calo * sumPx_calo + sumPy_calo * sumPy_calo);
0185 math::XYZTLorentzVector missingEt(-sumPx_calo, -sumPy_calo, 0, et);
0186 edm::RefProd<TkJetCollection> jetCollRef(L1TkJetsHandle);
0187 TkHTMiss tkHTM(missingEt, HT_calo, jetCollRef, L1VtxRef);
0188
0189 if (doVtxConstrain_ && !primaryVtxConstrain_) {
0190 tkHTM.setVtx(evtZVtx);
0191 }
0192
0193 MHTCollection->push_back(tkHTM);
0194 iEvent.put(std::move(MHTCollection), "L1TkCaloHTMiss");
0195 }
0196
0197 else {
0198 float sumPx = 0;
0199 float sumPy = 0;
0200 float HT = 0;
0201
0202
0203 for (jetIter = L1TkJetsHandle->begin(); jetIter != L1TkJetsHandle->end(); ++jetIter) {
0204 float tmp_jet_px = jetIter->px();
0205 float tmp_jet_py = jetIter->py();
0206 float tmp_jet_et = jetIter->et();
0207 float tmp_jet_pt = jetIter->pt();
0208 if (tmp_jet_pt < jetMinPt_)
0209 continue;
0210 if (std::abs(jetIter->eta()) > jetMaxEta_)
0211 continue;
0212 if (jetIter->ntracks() < minNtracksLowPt_ && tmp_jet_et > minJetEtLowPt_)
0213 continue;
0214 if (jetIter->ntracks() < minNtracksHighPt_ && tmp_jet_et > minJetEtHighPt_)
0215 continue;
0216 sumPx += tmp_jet_px;
0217 sumPy += tmp_jet_py;
0218 HT += tmp_jet_pt;
0219 }
0220
0221
0222 float et = sqrt(sumPx * sumPx + sumPy * sumPy);
0223 math::XYZTLorentzVector missingEt(-sumPx, -sumPy, 0, et);
0224 edm::RefProd<TkJetCollection> jetCollRef(L1TkJetsHandle);
0225 TkHTMiss tkHTM(missingEt, HT, jetCollRef, L1VtxRef);
0226
0227 MHTCollection->push_back(tkHTM);
0228 if (displaced_)
0229 iEvent.put(std::move(MHTCollection), "L1TrackerHTMissExtended");
0230 else
0231 iEvent.put(std::move(MHTCollection), "L1TrackerHTMiss");
0232 }
0233 }
0234
0235 void L1TkHTMissProducer::beginJob() {}
0236
0237 void L1TkHTMissProducer::endJob() {}
0238
0239 DEFINE_FWK_MODULE(L1TkHTMissProducer);