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::TkPrimaryVertexCollection>(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  float etTot = 0;
89  double sumPx_PU = 0;
90  double sumPy_PU = 0;
91  double etTot_PU = 0;
92  int nTracks_counter = 0;
93 
94  float zVTX = L1VertexHandle->begin()->zvertex();
95  unsigned int tkCnt = 0;
96  for (const auto& trackIter : *TTTrackHandle) {
97  edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_> > tempTrackPtr(TTTrackHandle, tkCnt++);
98  float pt = tempTrackPtr->momentum().perp();
99  float phi = tempTrackPtr->momentum().phi();
100  float eta = tempTrackPtr->momentum().eta();
101  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_> >, TTStub<Ref_Phase2TrackerDigi_> > >
102  theStubs = trackIter.getStubRefs();
103  int nstubs = (int)theStubs.size();
104 
105  float chi2dof = tempTrackPtr->chi2Red();
106  float bendchi2 = tempTrackPtr->stubPtConsistency();
107  float z0 = tempTrackPtr->z0();
108 
109  if (pt < minPt || fabs(z0) > maxZ0 || fabs(eta) > maxEta || chi2dof > chi2dofMax || bendchi2 > bendchi2Max)
110  continue;
111  if (maxPt > 0 && pt > maxPt) {
112  if (HighPtTracks == 0)
113  continue; // ignore these very high PT tracks: truncate
114  if (HighPtTracks == 1)
115  pt = maxPt; // saturate
116  }
117 
118  int nPS = 0.; // number of stubs in PS modules
119  // loop over the stubs
120  for (unsigned int istub = 0; istub < (unsigned int)theStubs.size(); istub++) {
121  DetId detId(theStubs.at(istub)->getDetId());
122  if (detId.det() == DetId::Detector::Tracker) {
123  if ((detId.subdetId() == StripSubdetector::TOB && tTopo->tobLayer(detId) <= 3) ||
124  (detId.subdetId() == StripSubdetector::TID && tTopo->tidRing(detId) <= 9))
125  nPS++;
126  }
127  }
128 
129  if (nstubs < nStubsmin || nPS < nStubsPSmin)
130  continue;
131 
132  // construct deltaZ cut to be based on track eta
133  if (fabs(eta) >= 0 && fabs(eta) < 0.7)
134  DeltaZ = 0.4;
135  else if (fabs(eta) >= 0.7 && fabs(eta) < 1.0)
136  DeltaZ = 0.6;
137  else if (fabs(eta) >= 1.0 && fabs(eta) < 1.2)
138  DeltaZ = 0.76;
139  else if (fabs(eta) >= 1.2 && fabs(eta) < 1.6)
140  DeltaZ = 1.0;
141  else if (fabs(eta) >= 1.6 && fabs(eta) < 2.0)
142  DeltaZ = 1.7;
143  else if (fabs(eta) >= 2.0 && fabs(eta) <= 2.4)
144  DeltaZ = 2.2;
145 
146  if (fabs(z0 - zVTX) <= DeltaZ) {
147  nTracks_counter++;
148  Track_Pt->Fill(pt);
149  Track_NStubs->Fill(nstubs);
150  Track_NPSstubs->Fill(nPS);
151  Track_Eta->Fill(eta);
152  Track_VtxZ->Fill(z0);
153  Track_Chi2Dof->Fill(chi2dof);
154  Track_BendChi2->Fill(bendchi2);
155 
156  sumPx += pt * cos(phi);
157  sumPy += pt * sin(phi);
158  etTot += pt;
159  } else { // PU sums
160  sumPx_PU += pt * cos(phi);
161  sumPy_PU += pt * sin(phi);
162  etTot_PU += pt;
163  }
164  } // end loop over tracks
165 
166  Track_N->Fill(nTracks_counter);
167  float et = sqrt(sumPx * sumPx + sumPy * sumPy);
168  double etmiss_PU = sqrt(sumPx_PU * sumPx_PU + sumPy_PU * sumPy_PU);
169 
170  math::XYZTLorentzVector missingEt(-sumPx, -sumPy, 0, et);
171 
172  TkMET_QualityCuts->Fill(missingEt.Pt());
173  TkMET_PU->Fill(etmiss_PU);
174 
175 } // end of method
176 
177 // ------------ method called once each job just before starting event loop ------------
178 //Creating all histograms for DQM file output
180  edm::Run const& run,
181  edm::EventSetup const& es) {
183  iBooker.setCurrentFolder(topFolderName_ + "/TkMET_Tracks/");
184 
185  // Num of L1Tracks in tkMET selection
186  HistoName = "Track_N";
187  edm::ParameterSet psTrack_N = conf_.getParameter<edm::ParameterSet>("TH1_NTracks");
188  Track_N = iBooker.book1D(HistoName,
189  HistoName,
190  psTrack_N.getParameter<int32_t>("Nbinsx"),
191  psTrack_N.getParameter<double>("xmin"),
192  psTrack_N.getParameter<double>("xmax"));
193  Track_N->setAxisTitle("# L1 Tracks", 1);
194  Track_N->setAxisTitle("# Events", 2);
195 
196  //Pt of the tracks
197  edm::ParameterSet psTrack_Pt = conf_.getParameter<edm::ParameterSet>("TH1_Track_Pt");
198  HistoName = "Track_Pt";
199  Track_Pt = iBooker.book1D(HistoName,
200  HistoName,
201  psTrack_Pt.getParameter<int32_t>("Nbinsx"),
202  psTrack_Pt.getParameter<double>("xmin"),
203  psTrack_Pt.getParameter<double>("xmax"));
204  Track_Pt->setAxisTitle("p_{T} [GeV]", 1);
205  Track_Pt->setAxisTitle("# L1 Tracks", 2);
206 
207  //Eta
208  edm::ParameterSet psTrack_Eta = conf_.getParameter<edm::ParameterSet>("TH1_Track_Eta");
209  HistoName = "Track_Eta";
210  Track_Eta = iBooker.book1D(HistoName,
211  HistoName,
212  psTrack_Eta.getParameter<int32_t>("Nbinsx"),
213  psTrack_Eta.getParameter<double>("xmin"),
214  psTrack_Eta.getParameter<double>("xmax"));
215  Track_Eta->setAxisTitle("#eta", 1);
216  Track_Eta->setAxisTitle("# L1 Tracks", 2);
217 
218  //VtxZ
219  edm::ParameterSet psTrack_VtxZ = conf_.getParameter<edm::ParameterSet>("TH1_Track_VtxZ");
220  HistoName = "Track_VtxZ";
221  Track_VtxZ = iBooker.book1D(HistoName,
222  HistoName,
223  psTrack_VtxZ.getParameter<int32_t>("Nbinsx"),
224  psTrack_VtxZ.getParameter<double>("xmin"),
225  psTrack_VtxZ.getParameter<double>("xmax"));
226  Track_VtxZ->setAxisTitle("L1 Track vertex position z [cm]", 1);
227  Track_VtxZ->setAxisTitle("# L1 Tracks", 2);
228 
229  //chi2dof
230  edm::ParameterSet psTrack_Chi2Dof = conf_.getParameter<edm::ParameterSet>("TH1_Track_Chi2Dof");
231  HistoName = "Track_Chi2Dof";
232  Track_Chi2Dof = iBooker.book1D(HistoName,
233  HistoName,
234  psTrack_Chi2Dof.getParameter<int32_t>("Nbinsx"),
235  psTrack_Chi2Dof.getParameter<double>("xmin"),
236  psTrack_Chi2Dof.getParameter<double>("xmax"));
237  Track_Chi2Dof->setAxisTitle("L1 Track #chi^{2}/D.O.F.", 1);
238  Track_Chi2Dof->setAxisTitle("# L1 Tracks", 2);
239 
240  //bend chi2
241  edm::ParameterSet psTrack_BendChi2 = conf_.getParameter<edm::ParameterSet>("TH1_Track_BendChi2");
242  HistoName = "Track_BendChi2";
243  Track_BendChi2 = iBooker.book1D(HistoName,
244  HistoName,
245  psTrack_BendChi2.getParameter<int32_t>("Nbinsx"),
246  psTrack_BendChi2.getParameter<double>("xmin"),
247  psTrack_BendChi2.getParameter<double>("xmax"));
248  Track_BendChi2->setAxisTitle("L1 Track Bend #chi^{2}", 1);
249  Track_BendChi2->setAxisTitle("# L1 Tracks", 2);
250 
251  //nstubs
252  edm::ParameterSet psTrack_NStubs = conf_.getParameter<edm::ParameterSet>("TH1_Track_NStubs");
253  HistoName = "Track_NStubs";
254  Track_NStubs = iBooker.book1D(HistoName,
255  HistoName,
256  psTrack_NStubs.getParameter<int32_t>("Nbinsx"),
257  psTrack_NStubs.getParameter<double>("xmin"),
258  psTrack_NStubs.getParameter<double>("xmax"));
259  Track_NStubs->setAxisTitle("# L1 Stubs", 1);
260  Track_NStubs->setAxisTitle("# L1 Tracks", 2);
261 
262  //nPSstubs
263  edm::ParameterSet psTrack_NPSstubs = conf_.getParameter<edm::ParameterSet>("TH1_Track_NPSstubs");
264  HistoName = "Track_NPSstubs";
265  Track_NPSstubs = iBooker.book1D(HistoName,
266  HistoName,
267  psTrack_NPSstubs.getParameter<int32_t>("Nbinsx"),
268  psTrack_NPSstubs.getParameter<double>("xmin"),
269  psTrack_NPSstubs.getParameter<double>("xmax"));
270  Track_NPSstubs->setAxisTitle("# PS Stubs", 1);
271  Track_NPSstubs->setAxisTitle("# L1 Tracks", 2);
272 
273  iBooker.setCurrentFolder(topFolderName_ + "/TkMET");
274  //loose tkMET
275  edm::ParameterSet psTrack_TkMET = conf_.getParameter<edm::ParameterSet>("TH1_Track_TkMET");
276  HistoName = "TkMET_QualityCuts";
278  HistoName,
279  psTrack_TkMET.getParameter<int32_t>("Nbinsx"),
280  psTrack_TkMET.getParameter<double>("xmin"),
281  psTrack_TkMET.getParameter<double>("xmax"));
282  TkMET_QualityCuts->setAxisTitle("L1 Track MET [GeV]", 1);
283  TkMET_QualityCuts->setAxisTitle("# Events", 2);
284 
285  //tkMET -- PU only
286  HistoName = "TkMET_PU";
287  TkMET_PU = iBooker.book1D(HistoName,
288  HistoName,
289  psTrack_TkMET.getParameter<int32_t>("Nbinsx"),
290  psTrack_TkMET.getParameter<double>("xmin"),
291  psTrack_TkMET.getParameter<double>("xmax"));
292  TkMET_PU->setAxisTitle("L1 Track MET (PU only) [GeV]", 1);
293  TkMET_PU->setAxisTitle("# Events", 2);
294 } //end of method
295 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
unsigned int tobLayer(const DetId &id) const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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:19
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 &)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static constexpr auto TOB
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
edm::EDGetTokenT< l1t::TkPrimaryVertexCollection > pvToken
Definition: DetId.h:17
void analyze(const edm::Event &, const edm::EventSetup &) override
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)