CMS 3D CMS Logo

L1TkHTMissProducer.cc
Go to the documentation of this file.
1 // Original Author: Emmanuelle Perez,40 1-A28,+41227671915,
2 // Created: Tue Nov 12 17:03:19 CET 2013
3 
4 // system include files
5 #include <memory>
6 
7 // user include files
18 
19 using namespace l1t;
20 
22 public:
23  explicit L1TkHTMissProducer(const edm::ParameterSet&);
24  ~L1TkHTMissProducer() override;
25 
26 private:
27  virtual void beginJob();
28  void produce(edm::Event&, const edm::EventSetup&) override;
29  virtual void endJob();
30 
31  // ----------member data ---------------------------
32  const float jetMinPt_; // [GeV]
33  const float jetMaxEta_; // [rad]
34  const bool doVtxConstrain_; // require vertex constraint
35  const bool useCaloJets_; // Determines whether or not calo jets are used
36  const bool primaryVtxConstrain_; // use event primary vertex instead of leading jet (if doVtxConstrain)
37  const float deltaZ_; // for jets [cm] (if DoTvxConstrain)
38  const unsigned int minNtracksHighPt_;
39  const unsigned int minNtracksLowPt_;
40  const float minJetEtLowPt_; // for track jets, minimum et required, depending on number of low pT tracks
41  const float minJetEtHighPt_;
42  const bool displaced_; // Use prompt/displaced tracks
45 };
46 
48  : jetMinPt_((float)iConfig.getParameter<double>("jet_minPt")),
49  jetMaxEta_((float)iConfig.getParameter<double>("jet_maxEta")),
50  doVtxConstrain_(iConfig.getParameter<bool>("doVtxConstrain")),
51  useCaloJets_(iConfig.getParameter<bool>("useCaloJets")),
52  primaryVtxConstrain_(iConfig.getParameter<bool>("primaryVtxConstrain")),
53  deltaZ_((float)iConfig.getParameter<double>("deltaZ")),
54  minNtracksHighPt_(iConfig.getParameter<int>("jet_minNtracksHighPt")),
55  minNtracksLowPt_(iConfig.getParameter<int>("jet_minNtracksLowPt")),
56  minJetEtLowPt_(iConfig.getParameter<double>("jet_minJetEtLowPt")),
57  minJetEtHighPt_(iConfig.getParameter<double>("jet_minJetEtHighPt")),
58  displaced_(iConfig.getParameter<bool>("displaced")),
59  pvToken_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("L1VertexInputTag"))),
60  jetToken_(consumes<TkJetCollection>(iConfig.getParameter<edm::InputTag>("L1TkJetInputTag"))) {
61  if (useCaloJets_)
62  produces<TkHTMissCollection>("TkCaloHTMiss");
63  else if (displaced_)
64  produces<TkHTMissCollection>("L1TrackerHTMissExtended");
65  else
66  produces<TkHTMissCollection>("L1TrackerHTMiss");
67 }
68 
70 
72  using namespace edm;
73  std::unique_ptr<TkHTMissCollection> MHTCollection(new TkHTMissCollection);
74 
75  // L1 primary vertex
76  edm::Handle<VertexCollection> L1VertexHandle;
77  iEvent.getByToken(pvToken_, L1VertexHandle);
78 
79  // L1 track-trigger jets
80  edm::Handle<TkJetCollection> L1TkJetsHandle;
81  iEvent.getByToken(jetToken_, L1TkJetsHandle);
82  std::vector<TkJet>::const_iterator jetIter;
83 
84  if (!L1TkJetsHandle.isValid() && !displaced_) {
85  LogError("TkHTMissProducer") << "\nWarning: TkJetCollection not found in the event. Exit\n";
86  return;
87  }
88 
89  if (!L1TkJetsHandle.isValid() && displaced_) {
90  LogError("TkHTMissProducer") << "\nWarning: TkJetExtendedCollection not found in the event. Exit\n";
91  return;
92  }
93 
94  // ----------------------------------------------------------------------------------------------
95  // if primaryVtxConstrain_, use the primary vertex instead of z position from leading jet
96  // ----------------------------------------------------------------------------------------------
97  float evtZVtx = 999;
98  bool foundVtx = false;
99  edm::Ref<VertexCollection> L1VtxRef; // null reference
100 
101  if (useCaloJets_) {
103  if (!L1VertexHandle.isValid()) {
104  LogError("L1TkHTMissProducer") << "\nWarning: VertexCollection not found in the event. Exit\n";
105  return;
106  } else {
107  std::vector<Vertex>::const_iterator vtxIter = L1VertexHandle->begin();
108  // by convention, the first vertex in the collection is the one that should
109  // be used by default
110  evtZVtx = vtxIter->z0();
111  foundVtx = true;
112  int ivtx = 0;
113  edm::Ref<VertexCollection> vtxRef(L1VertexHandle, ivtx);
114  L1VtxRef = vtxRef;
115  }
116  } //endif primaryVtxConstrain_
117 
118  // ----------------------------------------------------------------------------------------------
119  // using z position of leading jet to define "event vertex"
120  // ----------------------------------------------------------------------------------------------
121  float zvtx_jetpt = -1.0; //pt of jet determining the event vertex
122  float jetVtxMax = 99.; //find z position of leading jet that has a z vertex!
123 
125  for (jetIter = L1TkJetsHandle->begin(); jetIter != L1TkJetsHandle->end(); ++jetIter) {
126  int ibx = jetIter->bx(); // only consider jets from the central BX
127  if (ibx != 0)
128  continue;
129 
130  float tmp_jet_vtx = jetIter->jetVtx();
131  float tmp_jet_pt = jetIter->pt();
132  float tmp_jet_eta = jetIter->eta();
133  if (tmp_jet_pt < jetMinPt_)
134  continue;
135  if (std::abs(tmp_jet_eta) > jetMaxEta_)
136  continue;
137  if (std::abs(tmp_jet_vtx) > jetVtxMax)
138  continue;
139 
140  // find vertex position of leading jet
141  if (tmp_jet_pt > zvtx_jetpt) {
142  evtZVtx = tmp_jet_vtx;
143  zvtx_jetpt = tmp_jet_pt;
144  foundVtx = true;
145  }
146  } //end loop over jets
147  } //endif z position from leading jet
148 
149  float sumPx_calo = 0;
150  float sumPy_calo = 0;
151  float HT_calo = 0;
152 
153  if (doVtxConstrain_ && !foundVtx)
154  LogWarning("L1TkHTMissProducer") << "Didn't find any z vertex (based on jet vertices) for this event!\n";
155 
156  // loop over jets
157  for (jetIter = L1TkJetsHandle->begin(); jetIter != L1TkJetsHandle->end(); ++jetIter) {
158  int ibx = jetIter->bx(); // only consider jets from the central BX
159  if (ibx != 0)
160  continue;
161 
162  float tmp_jet_px = jetIter->px();
163  float tmp_jet_py = jetIter->py();
164  float tmp_jet_et = jetIter->et();
165  float tmp_jet_vtx = jetIter->jetVtx();
166  if (jetIter->pt() < jetMinPt_)
167  continue;
168  if (std::abs(jetIter->eta()) > jetMaxEta_)
169  continue;
170 
171  // vertex consistency requirement
172  bool VtxRequirement = false;
173  if (foundVtx)
174  VtxRequirement = std::abs(tmp_jet_vtx - evtZVtx) < deltaZ_;
175 
176  if (!doVtxConstrain_ || VtxRequirement) {
177  sumPx_calo += tmp_jet_px;
178  sumPy_calo += tmp_jet_py;
179  HT_calo += tmp_jet_et;
180  }
181  } //end loop over jets
182 
183  // define missing HT
184  float et = sqrt(sumPx_calo * sumPx_calo + sumPy_calo * sumPy_calo);
185  math::XYZTLorentzVector missingEt(-sumPx_calo, -sumPy_calo, 0, et);
186  edm::RefProd<TkJetCollection> jetCollRef(L1TkJetsHandle);
187  TkHTMiss tkHTM(missingEt, HT_calo, jetCollRef, L1VtxRef);
188 
190  tkHTM.setVtx(evtZVtx);
191  }
192 
193  MHTCollection->push_back(tkHTM);
194  iEvent.put(std::move(MHTCollection), "L1TkCaloHTMiss");
195  }
196 
197  else { // Using standalone jets
198  float sumPx = 0;
199  float sumPy = 0;
200  float HT = 0;
201 
202  // loop over jets
203  for (jetIter = L1TkJetsHandle->begin(); jetIter != L1TkJetsHandle->end(); ++jetIter) {
204  float tmp_jet_px = jetIter->px();
205  float tmp_jet_py = jetIter->py();
206  float tmp_jet_et = jetIter->et();
207  float tmp_jet_pt = jetIter->pt();
208  if (tmp_jet_pt < jetMinPt_)
209  continue;
210  if (std::abs(jetIter->eta()) > jetMaxEta_)
211  continue;
212  if (jetIter->ntracks() < minNtracksLowPt_ && tmp_jet_et > minJetEtLowPt_)
213  continue;
214  if (jetIter->ntracks() < minNtracksHighPt_ && tmp_jet_et > minJetEtHighPt_)
215  continue;
216  sumPx += tmp_jet_px;
217  sumPy += tmp_jet_py;
218  HT += tmp_jet_pt;
219  } // end jet loop
220 
221  // define missing HT
222  float et = sqrt(sumPx * sumPx + sumPy * sumPy);
223  math::XYZTLorentzVector missingEt(-sumPx, -sumPy, 0, et);
224  edm::RefProd<TkJetCollection> jetCollRef(L1TkJetsHandle);
225  TkHTMiss tkHTM(missingEt, HT, jetCollRef, L1VtxRef);
226 
227  MHTCollection->push_back(tkHTM);
228  if (displaced_)
229  iEvent.put(std::move(MHTCollection), "L1TrackerHTMissExtended");
230  else
231  iEvent.put(std::move(MHTCollection), "L1TrackerHTMiss");
232  }
233 } //end producer
234 
236 
238 
~L1TkHTMissProducer() override
void produce(edm::Event &, const edm::EventSetup &) override
const edm::EDGetTokenT< TkJetCollection > jetToken_
const unsigned int minNtracksLowPt_
L1TkHTMissProducer(const edm::ParameterSet &)
std::vector< TkHTMiss > TkHTMissCollection
Definition: TkHTMissFwd.h:10
delete x;
Definition: CaloConfig.h:22
Log< level::Error, false > LogError
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
void beginJob()
Definition: Breakpoints.cc:14
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const unsigned int minNtracksHighPt_
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< pat::MHT > MHTCollection
Definition: MHT.h:39
void setVtx(const float &zvtx)
Definition: TkHTMiss.h:41
const bool primaryVtxConstrain_
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
const edm::EDGetTokenT< VertexCollection > pvToken_
Log< level::Warning, false > LogWarning
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
Definition: HT.h:21
def move(src, dest)
Definition: eostools.py:511