CMS 3D CMS Logo

L1TPhase2OuterTrackerTkMET.cc
Go to the documentation of this file.
1 // Package: SiOuterTracker
2 // Class: SiOuterTracker
3 //
4 // Author: Emily MacDonald (emily.kaelyn.macdonald@cern.ch)
5 
6 // system include files
7 #include <memory>
8 #include <vector>
9 #include <numeric>
10 #include <iostream>
11 #include <fstream>
12 #include <string>
13 
14 // user include files
24 
34 
36 
37 // constructors and destructor
39  : conf_(iConfig), m_topoToken(esConsumes()) {
40  topFolderName_ = conf_.getParameter<std::string>("TopFolderName");
42  consumes<std::vector<TTTrack<Ref_Phase2TrackerDigi_> > >(conf_.getParameter<edm::InputTag>("TTTracksTag"));
43  pvToken = consumes<l1t::VertexWordCollection>(conf_.getParameter<edm::InputTag>("L1VertexInputTag"));
44 
45  maxZ0 = conf_.getParameter<double>("maxZ0");
46  DeltaZ = conf_.getParameter<double>("DeltaZ");
47  chi2dofMax = conf_.getParameter<double>("chi2dofMax");
48  bendchi2Max = conf_.getParameter<double>("bendchi2Max");
49  minPt = conf_.getParameter<double>("minPt");
50  nStubsmin = iConfig.getParameter<int>("nStubsmin");
51  nStubsPSmin = iConfig.getParameter<int>("nStubsPSmin");
52  maxPt = conf_.getParameter<double>("maxPt");
53  maxEta = conf_.getParameter<double>("maxEta");
54  HighPtTracks = iConfig.getParameter<int>("HighPtTracks");
55 }
56 
58  // do anything here that needs to be done at desctruction time
59  // (e.g. close files, deallocate resources etc.)
60 }
61 
62 // member functions
63 
64 // ------------ method called for each event ------------
66  // L1 Primaries
68  iEvent.getByToken(pvToken, L1VertexHandle);
69 
71  iEvent.getByToken(ttTrackToken_, TTTrackHandle);
72 
73  // for PS stubs
74  // Tracker Topology
75  const TrackerTopology* const tTopo = &iSetup.getData(m_topoToken);
76 
77  // Adding protection
78  if (!TTTrackHandle.isValid()) {
79  edm::LogWarning("L1TPhase2OuterTrackerTkMET") << "cant find tracks" << std::endl;
80  return;
81  }
82  if (!L1VertexHandle.isValid()) {
83  edm::LogWarning("L1TPhase2OuterTrackerTkMET") << "cant find vertex" << std::endl;
84  return;
85  }
86  float sumPx = 0;
87  float sumPy = 0;
88  double sumPx_PU = 0;
89  double sumPy_PU = 0;
90  int nTracks_counter = 0;
91 
92  float zVTX = L1VertexHandle->begin()->z0();
93  unsigned int tkCnt = 0;
94  for (const auto& trackIter : *TTTrackHandle) {
95  edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_> > tempTrackPtr(TTTrackHandle, tkCnt++);
96  float pt = tempTrackPtr->momentum().perp();
97  float phi = tempTrackPtr->momentum().phi();
98  float eta = tempTrackPtr->momentum().eta();
99  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_> >, TTStub<Ref_Phase2TrackerDigi_> > >
100  theStubs = trackIter.getStubRefs();
101  int nstubs = (int)theStubs.size();
102 
103  float chi2dof = tempTrackPtr->chi2Red();
104  float bendchi2 = tempTrackPtr->stubPtConsistency();
105  float z0 = tempTrackPtr->z0();
106 
107  if (pt < minPt || fabs(z0) > maxZ0 || fabs(eta) > maxEta || chi2dof > chi2dofMax || bendchi2 > bendchi2Max)
108  continue;
109  if (maxPt > 0 && pt > maxPt) {
110  if (HighPtTracks == 0)
111  continue; // ignore these very high PT tracks: truncate
112  if (HighPtTracks == 1)
113  pt = maxPt; // saturate
114  }
115 
116  int nPS = 0.; // number of stubs in PS modules
117  // loop over the stubs
118  for (unsigned int istub = 0; istub < (unsigned int)theStubs.size(); istub++) {
119  DetId detId(theStubs.at(istub)->getDetId());
120  if (detId.det() == DetId::Detector::Tracker) {
121  if ((detId.subdetId() == StripSubdetector::TOB && tTopo->tobLayer(detId) <= 3) ||
122  (detId.subdetId() == StripSubdetector::TID && tTopo->tidRing(detId) <= 9))
123  nPS++;
124  }
125  }
126 
127  if (nstubs < nStubsmin || nPS < nStubsPSmin)
128  continue;
129 
130  // construct deltaZ cut to be based on track eta
131  if (fabs(eta) >= 0 && fabs(eta) < 0.7)
132  DeltaZ = 0.4;
133  else if (fabs(eta) >= 0.7 && fabs(eta) < 1.0)
134  DeltaZ = 0.6;
135  else if (fabs(eta) >= 1.0 && fabs(eta) < 1.2)
136  DeltaZ = 0.76;
137  else if (fabs(eta) >= 1.2 && fabs(eta) < 1.6)
138  DeltaZ = 1.0;
139  else if (fabs(eta) >= 1.6 && fabs(eta) < 2.0)
140  DeltaZ = 1.7;
141  else if (fabs(eta) >= 2.0 && fabs(eta) <= 2.4)
142  DeltaZ = 2.2;
143 
144  if (fabs(z0 - zVTX) <= DeltaZ) {
145  nTracks_counter++;
146  Track_Pt->Fill(pt);
147  Track_NStubs->Fill(nstubs);
148  Track_NPSstubs->Fill(nPS);
149  Track_Eta->Fill(eta);
150  Track_VtxZ->Fill(z0);
151  Track_Chi2Dof->Fill(chi2dof);
152  Track_BendChi2->Fill(bendchi2);
153 
154  sumPx += pt * cos(phi);
155  sumPy += pt * sin(phi);
156  } else { // PU sums
157  sumPx_PU += pt * cos(phi);
158  sumPy_PU += pt * sin(phi);
159  }
160  } // end loop over tracks
161 
162  Track_N->Fill(nTracks_counter);
163  float et = sqrt(sumPx * sumPx + sumPy * sumPy);
164  double etmiss_PU = sqrt(sumPx_PU * sumPx_PU + sumPy_PU * sumPy_PU);
165 
166  math::XYZTLorentzVector missingEt(-sumPx, -sumPy, 0, et);
167 
168  TkMET_QualityCuts->Fill(missingEt.Pt());
169  TkMET_PU->Fill(etmiss_PU);
170 
171 } // end of method
172 
173 // ------------ method called once each job just before starting event loop ------------
174 //Creating all histograms for DQM file output
176  edm::Run const& run,
177  edm::EventSetup const& es) {
179  iBooker.setCurrentFolder(topFolderName_ + "/TkMET_Tracks/");
180 
181  // Num of L1Tracks in tkMET selection
182  HistoName = "Track_N";
183  edm::ParameterSet psTrack_N = conf_.getParameter<edm::ParameterSet>("TH1_NTracks");
184  Track_N = iBooker.book1D(HistoName,
185  HistoName,
186  psTrack_N.getParameter<int32_t>("Nbinsx"),
187  psTrack_N.getParameter<double>("xmin"),
188  psTrack_N.getParameter<double>("xmax"));
189  Track_N->setAxisTitle("# L1 Tracks", 1);
190  Track_N->setAxisTitle("# Events", 2);
191 
192  //Pt of the tracks
193  edm::ParameterSet psTrack_Pt = conf_.getParameter<edm::ParameterSet>("TH1_Track_Pt");
194  HistoName = "Track_Pt";
195  Track_Pt = iBooker.book1D(HistoName,
196  HistoName,
197  psTrack_Pt.getParameter<int32_t>("Nbinsx"),
198  psTrack_Pt.getParameter<double>("xmin"),
199  psTrack_Pt.getParameter<double>("xmax"));
200  Track_Pt->setAxisTitle("p_{T} [GeV]", 1);
201  Track_Pt->setAxisTitle("# L1 Tracks", 2);
202 
203  //Eta
204  edm::ParameterSet psTrack_Eta = conf_.getParameter<edm::ParameterSet>("TH1_Track_Eta");
205  HistoName = "Track_Eta";
206  Track_Eta = iBooker.book1D(HistoName,
207  HistoName,
208  psTrack_Eta.getParameter<int32_t>("Nbinsx"),
209  psTrack_Eta.getParameter<double>("xmin"),
210  psTrack_Eta.getParameter<double>("xmax"));
211  Track_Eta->setAxisTitle("#eta", 1);
212  Track_Eta->setAxisTitle("# L1 Tracks", 2);
213 
214  //VtxZ
215  edm::ParameterSet psTrack_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1_Track_VtxZ");
216  HistoName = "Track_VtxZ";
217  Track_VtxZ = iBooker.book1D(HistoName,
218  HistoName,
219  psTrack_VtxZ.getParameter<int32_t>("Nbinsx"),
220  psTrack_VtxZ.getParameter<double>("xmin"),
221  psTrack_VtxZ.getParameter<double>("xmax"));
222  Track_VtxZ->setAxisTitle("L1 Track vertex position z [cm]", 1);
223  Track_VtxZ->setAxisTitle("# L1 Tracks", 2);
224 
225  //chi2dof
226  edm::ParameterSet psTrack_Chi2Dof = conf_.getParameter<edm::ParameterSet>("TH1_Track_Chi2Dof");
227  HistoName = "Track_Chi2Dof";
228  Track_Chi2Dof = iBooker.book1D(HistoName,
229  HistoName,
230  psTrack_Chi2Dof.getParameter<int32_t>("Nbinsx"),
231  psTrack_Chi2Dof.getParameter<double>("xmin"),
232  psTrack_Chi2Dof.getParameter<double>("xmax"));
233  Track_Chi2Dof->setAxisTitle("L1 Track #chi^{2}/D.O.F.", 1);
234  Track_Chi2Dof->setAxisTitle("# L1 Tracks", 2);
235 
236  //bend chi2
237  edm::ParameterSet psTrack_BendChi2 = conf_.getParameter<edm::ParameterSet>("TH1_Track_BendChi2");
238  HistoName = "Track_BendChi2";
239  Track_BendChi2 = iBooker.book1D(HistoName,
240  HistoName,
241  psTrack_BendChi2.getParameter<int32_t>("Nbinsx"),
242  psTrack_BendChi2.getParameter<double>("xmin"),
243  psTrack_BendChi2.getParameter<double>("xmax"));
244  Track_BendChi2->setAxisTitle("L1 Track Bend #chi^{2}", 1);
245  Track_BendChi2->setAxisTitle("# L1 Tracks", 2);
246 
247  //nstubs
248  edm::ParameterSet psTrack_NStubs = conf_.getParameter<edm::ParameterSet>("TH1_Track_NStubs");
249  HistoName = "Track_NStubs";
250  Track_NStubs = iBooker.book1D(HistoName,
251  HistoName,
252  psTrack_NStubs.getParameter<int32_t>("Nbinsx"),
253  psTrack_NStubs.getParameter<double>("xmin"),
254  psTrack_NStubs.getParameter<double>("xmax"));
255  Track_NStubs->setAxisTitle("# L1 Stubs", 1);
256  Track_NStubs->setAxisTitle("# L1 Tracks", 2);
257 
258  //nPSstubs
259  edm::ParameterSet psTrack_NPSstubs = conf_.getParameter<edm::ParameterSet>("TH1_Track_NPSstubs");
260  HistoName = "Track_NPSstubs";
261  Track_NPSstubs = iBooker.book1D(HistoName,
262  HistoName,
263  psTrack_NPSstubs.getParameter<int32_t>("Nbinsx"),
264  psTrack_NPSstubs.getParameter<double>("xmin"),
265  psTrack_NPSstubs.getParameter<double>("xmax"));
266  Track_NPSstubs->setAxisTitle("# PS Stubs", 1);
267  Track_NPSstubs->setAxisTitle("# L1 Tracks", 2);
268 
269  iBooker.setCurrentFolder(topFolderName_ + "/TkMET");
270  //loose tkMET
271  edm::ParameterSet psTrack_TkMET = conf_.getParameter<edm::ParameterSet>("TH1_Track_TkMET");
272  HistoName = "TkMET_QualityCuts";
274  HistoName,
275  psTrack_TkMET.getParameter<int32_t>("Nbinsx"),
276  psTrack_TkMET.getParameter<double>("xmin"),
277  psTrack_TkMET.getParameter<double>("xmax"));
278  TkMET_QualityCuts->setAxisTitle("L1 Track MET [GeV]", 1);
279  TkMET_QualityCuts->setAxisTitle("# Events", 2);
280 
281  //tkMET -- PU only
282  HistoName = "TkMET_PU";
283  TkMET_PU = iBooker.book1D(HistoName,
284  HistoName,
285  psTrack_TkMET.getParameter<int32_t>("Nbinsx"),
286  psTrack_TkMET.getParameter<double>("xmin"),
287  psTrack_TkMET.getParameter<double>("xmax"));
288  TkMET_PU->setAxisTitle("L1 Track MET (PU only) [GeV]", 1);
289  TkMET_PU->setAxisTitle("# Events", 2);
290 } //end of method
291 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int tobLayer(const DetId &id) const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void Fill(long long x)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > ttTrackToken_
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:23
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
L1TPhase2OuterTrackerTkMET(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static constexpr auto TOB
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
Definition: DetId.h:17
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< l1t::VertexWordCollection > pvToken
std::string HistoName
bool isValid() const
Definition: HandleBase.h:70
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > m_topoToken
unsigned int tidRing(const DetId &id) const
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
static constexpr auto TID
Definition: Run.h:45
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)