CMS 3D CMS Logo

JetHTAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TkAlTools/JetHTAnalyzer
4 // Class: JetHTAnalyzer
5 //
13 //
14 // Original Author: Marco Musich
15 // Created: Fri, 29 Mar 2019 14:54:59 GMT
16 //
17 // Updated by: Jussi Viinikainen
18 //
19 
20 // system include files
21 #include <algorithm> // std::sort
22 #include <memory>
23 #include <string>
24 #include <vector> // std::vector
25 
26 // CMSSW includes
47 #include "DataFormats/Common/interface/TriggerResults.h" // Classes needed to print trigger results
49 
50 //
51 // class declaration
52 //
54 
55 class JetHTAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
56 public:
57  explicit JetHTAnalyzer(const edm::ParameterSet&);
58  ~JetHTAnalyzer() override = default;
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
61  static bool mysorter(reco::Track i, reco::Track j) { return (i.pt() > j.pt()); }
62 
63 private:
64  void beginJob() override;
65  void analyze(const edm::Event&, const edm::EventSetup&) override;
66  void endJob() override;
67 
68  // ----------member data ---------------------------
71 
74 
77 
78  const int printTriggerTable_;
79  const double minVtxNdf_;
80  const double minVtxWgt_;
81 
82  const std::vector<double> profilePtBorders_;
83  const std::vector<int> iovList_;
84 
85  // output histograms
87  TH1F* h_ntrks;
88  TH1F* h_probePt;
89  TH1F* h_probeEta;
90  TH1F* h_probePhi;
91  TH1F* h_probeDxy;
92  TH1F* h_probeDz;
94  TH1F* h_probeDzErr;
95 
97 
98  // for the conversions
99  static constexpr double cmToum = 10000;
100 };
101 
102 //
103 // Constructor
104 //
106  : pvsTag_(iConfig.getParameter<edm::InputTag>("vtxCollection")),
107  pvsToken_(consumes<reco::VertexCollection>(pvsTag_)),
108  tracksTag_(iConfig.getParameter<edm::InputTag>("trackCollection")),
109  tracksToken_(consumes<reco::TrackCollection>(tracksTag_)),
110  triggerTag_(iConfig.getParameter<edm::InputTag>("triggerResults")),
111  triggerToken_(consumes<edm::TriggerResults>(triggerTag_)),
112  printTriggerTable_(iConfig.getUntrackedParameter<int>("printTriggerTable")),
113  minVtxNdf_(iConfig.getUntrackedParameter<double>("minVertexNdf")),
114  minVtxWgt_(iConfig.getUntrackedParameter<double>("minVertexMeanWeight")),
115  profilePtBorders_(iConfig.getUntrackedParameter<std::vector<double>>("profilePtBorders")),
116  iovList_(iConfig.getUntrackedParameter<std::vector<int>>("iovList")) {
117  // Specify that TFileService is used by the class
118  usesResource(TFileService::kSharedResource);
119 }
120 
121 //
122 // member functions
123 //
124 
125 // ------------ method called for each event ------------
127  using namespace edm;
128 
129  const auto& vertices = iEvent.get(pvsToken_);
130  const reco::VertexCollection& pvtx = vertices;
131 
133 
134  // Find the IOV of the current event so that we can tag a histogram with the IOV
135  const int runNumber = iEvent.id().run();
136  std::string iovString = "iovNotFound";
137 
138  // Find the current IOV from the IOV list
139  if (runNumber >= iovList_.at(0)) { // If run number is smaller than the first item in the list, it is not in any IOV
140  for (std::vector<int>::size_type i = 1; i < iovList_.size(); i++) {
141  if (iovList_.at(i) > runNumber) {
142  iovString = Form("iov%d-%d", iovList_.at(i - 1), iovList_.at(i) - 1);
143  break;
144  }
145  }
146  }
147 
148  // Print the triggers to console
149  if (printTriggerTable_) {
150  const auto& triggerResults = iEvent.get(triggerToken_);
151 
152  const edm::TriggerNames& triggerNames = iEvent.triggerNames(triggerResults);
153  for (unsigned i = 0; i < triggerNames.size(); i++) {
154  const std::string& hltName = triggerNames.triggerName(i);
155  bool decision = triggerResults.accept(triggerNames.triggerIndex(hltName));
156  std::cout << hltName << " " << decision << std::endl;
157  }
158  }
159 
160  int counter = 0;
161  for (reco::VertexCollection::const_iterator pvIt = pvtx.begin(); pvIt != pvtx.end(); pvIt++) {
162  const reco::Vertex& iPV = *pvIt;
163  counter++;
164 
165  if (iPV.isFake())
166  continue;
168 
169  const math::XYZPoint pos_(iPV.x(), iPV.y(), iPV.z());
170 
171  // vertex selection as in bs code
172  if (iPV.ndof() < minVtxNdf_ || (iPV.ndof() + 3.) / iPV.tracksSize() < 2 * minVtxWgt_)
173  continue;
174 
176  for (trki = iPV.tracks_begin(); trki != iPV.tracks_end(); ++trki) {
177  if (trki->isNonnull()) {
178  reco::TrackRef trk_now(tracks, (*trki).key());
179  allTracks.push_back(*trk_now);
180  }
181  }
182 
183  // order with decreasing pt
184  std::sort(allTracks.begin(), allTracks.end(), mysorter);
185  uint ntrks = allTracks.size();
186  h_ntrks->Fill(ntrks);
187 
188  for (uint tracksIt = 0; tracksIt < ntrks; tracksIt++) {
189  auto tk = allTracks.at(tracksIt);
190 
191  double dxyRes = tk.dxy(pos_) * cmToum;
192  double dzRes = tk.dz(pos_) * cmToum;
193 
194  double dxy_err = tk.dxyError() * cmToum;
195  double dz_err = tk.dzError() * cmToum;
196 
197  float trackphi = tk.phi();
198  float tracketa = tk.eta();
199  float trackpt = tk.pt();
200 
201  h_probePt->Fill(trackpt);
202  h_probeEta->Fill(tracketa);
203  h_probePhi->Fill(trackphi);
204 
205  h_probeDxy->Fill(dxyRes);
206  h_probeDz->Fill(dzRes);
207  h_probeDxyErr->Fill(dxy_err);
208  h_probeDzErr->Fill(dz_err);
209 
210  mon.fillHisto("dxy", "all", dxyRes, 1.);
211  mon.fillHisto("dz", "all", dzRes, 1.);
212  mon.fillHisto("dxyerr", "all", dxy_err, 1.);
213  mon.fillHisto("dzerr", "all", dz_err, 1.);
214 
215  mon.fillProfile("dxyErrVsPt", "all", trackpt, dxy_err, 1.);
216  mon.fillProfile("dzErrVsPt", "all", trackpt, dz_err, 1.);
217 
218  mon.fillProfile("dxyErrVsPhi", "all", trackphi, dxy_err, 1.);
219  mon.fillProfile("dzErrVsPhi", "all", trackphi, dz_err, 1.);
220 
221  mon.fillProfile("dxyErrVsEta", "all", tracketa, dxy_err, 1.);
222  mon.fillProfile("dzErrVsEta", "all", tracketa, dz_err, 1.);
223 
224  // Integrated pT bins
225  for (std::vector<double>::size_type i = 0; i < profilePtBorders_.size(); i++) {
226  if (trackpt < profilePtBorders_.at(i))
227  break;
228  mon.fillProfile("dxyErrVsPtWide", "all", i, dxy_err, 1.);
229  mon.fillProfile("dzErrVsPtWide", "all", i, dz_err, 1.);
230  }
231 
232  // Fill IOV specific histograms
233  mon.fillHisto("dxy", iovString, dxyRes, 1.);
234  mon.fillHisto("dz", iovString, dzRes, 1.);
235  mon.fillHisto("dxyerr", iovString, dxy_err, 1.);
236  mon.fillHisto("dzerr", iovString, dz_err, 1.);
237 
238  mon.fillProfile("dxyErrVsPt", iovString, trackpt, dxy_err, 1.);
239  mon.fillProfile("dzErrVsPt", iovString, trackpt, dz_err, 1.);
240 
241  mon.fillProfile("dxyErrVsPhi", iovString, trackphi, dxy_err, 1.);
242  mon.fillProfile("dzErrVsPhi", iovString, trackphi, dz_err, 1.);
243 
244  mon.fillProfile("dxyErrVsEta", iovString, tracketa, dxy_err, 1.);
245  mon.fillProfile("dzErrVsEta", iovString, tracketa, dz_err, 1.);
246 
247  // Integrated pT bins
248  for (std::vector<double>::size_type i = 0; i < profilePtBorders_.size(); i++) {
249  if (trackpt < profilePtBorders_.at(i))
250  break;
251  mon.fillProfile("dxyErrVsPtWide", iovString, i, dxy_err, 1.);
252  mon.fillProfile("dzErrVsPtWide", iovString, i, dz_err, 1.);
253  }
254 
255  if (std::abs(tracketa) < 1.) {
256  mon.fillHisto("dxy", "central", dxyRes, 1.);
257  mon.fillHisto("dz", "central", dzRes, 1.);
258  mon.fillHisto("dxyerr", "central", dxy_err, 1.);
259  mon.fillHisto("dzerr", "central", dz_err, 1.);
260 
261  mon.fillProfile("dxyErrVsPt", "central", trackpt, dxy_err, 1.);
262  mon.fillProfile("dzErrVsPt", "central", trackpt, dz_err, 1.);
263 
264  mon.fillProfile("dxyErrVsPhi", "central", trackphi, dxy_err, 1.);
265  mon.fillProfile("dzErrVsPhi", "central", trackphi, dz_err, 1.);
266 
267  // Integrated pT bins
268  for (std::vector<double>::size_type i = 0; i < profilePtBorders_.size(); i++) {
269  if (trackpt < profilePtBorders_.at(i))
270  break;
271  mon.fillProfile("dxyErrVsPtWide", "central", i, dxy_err, 1.);
272  mon.fillProfile("dzErrVsPtWide", "central", i, dz_err, 1.);
273  }
274  }
275 
276  } // loop on tracks in vertex
277  } // loop on vertices
278 
279  mon.fillHisto("nvtx", "all", counter, 1.);
280 }
281 
282 // ------------ method called once each job just before starting event loop ------------
284  h_ntrks = outfile_->make<TH1F>("h_ntrks", "n. trks;n. of tracks/vertex;n. vertices", 100, 0, 100);
285  h_probePt = outfile_->make<TH1F>("h_probePt", "p_{T} of probe track;track p_{T} (GeV); tracks", 100, 0., 500.);
286  h_probeEta = outfile_->make<TH1F>("h_probeEta", "#eta of the probe track;track #eta;tracks", 54, -2.8, 2.8);
287  h_probePhi = outfile_->make<TH1F>("h_probePhi", "#phi of probe track;track #phi (rad);tracks", 100, -3.15, 3.15);
288 
289  h_probeDxy =
290  outfile_->make<TH1F>("h_probeDxy", "d_{xy}(PV) of the probe track;track d_{xy}(PV);tracks", 200, -100, 100);
291  h_probeDz = outfile_->make<TH1F>("h_probeDz", "d_{z}(PV) of the probe track;track d_{z}(PV);tracks", 200, -100, 100);
292  h_probeDxyErr = outfile_->make<TH1F>(
293  "h_probeDxyErr", "error on d_{xy}(PV) of the probe track;track error on d_{xy}(PV);tracks", 100, 0., 100);
294  h_probeDzErr = outfile_->make<TH1F>(
295  "h_probeDzErr", "error on d_{z}(PV) of the probe track;track error on d_{z}(PV);tracks", 100, 0., 100);
296 
297  mon.addHistogram(new TH1F("nvtx", ";Vertices;Events", 50, 0, 50));
298  mon.addHistogram(new TH1F("dxy", ";d_{xy};tracks", 100, -100, 100));
299  mon.addHistogram(new TH1F("dz", ";d_{z};tracks", 100, -100, 100));
300  mon.addHistogram(new TH1F("dxyerr", ";d_{xy} error;tracks", 100, 0., 200));
301  mon.addHistogram(new TH1F("dzerr", ";d_{z} error;tracks", 100, 0., 200));
302  mon.addHistogram(new TProfile("dxyErrVsPt", ";track p_{T};d_{xy} error", 100, 0., 200, 0., 100.));
303  mon.addHistogram(new TProfile("dzErrVsPt", ";track p_{T};d_{z} error", 100, 0., 200, 0., 100.));
304  mon.addHistogram(new TProfile("dxyErrVsPhi", ";track #varphi;d_{xy} error", 100, -M_PI, M_PI, 0., 100.));
305  mon.addHistogram(new TProfile("dzErrVsPhi", ";track #varphi;d_{z} error", 100, -M_PI, M_PI, 0., 100.));
306  mon.addHistogram(new TProfile("dxyErrVsEta", ";track #eta;d_{xy} error", 100, -2.5, 2.5, 0., 100.));
307  mon.addHistogram(new TProfile("dzErrVsEta", ";track #eta;d_{z} error", 100, -2.5, 2.5, 0., 100.));
308 
309  // Variable size histogram depending on the given pT bin borders
310  int nBins = profilePtBorders_.size();
312  new TProfile("dxyErrVsPtWide", ";track p_{T} wide bin;d_{xy} error", nBins, -0.5, nBins - 0.5, 0.0, 100.0));
314  new TProfile("dzErrVsPtWide", ";track p_{T} wide bin;d_{z} error", nBins, -0.5, nBins - 0.5, 0.0, 100.0));
315 }
316 
317 // ------------ method called once each job just after ending the event loop ------------
319 
320 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
323  desc.setComment("JetHT validation analyzer plugin.");
324  desc.add<edm::InputTag>("vtxCollection", edm::InputTag("offlinePrimaryVerticesFromRefittedTrks"));
325  desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults", "", "HLT"));
326  desc.add<edm::InputTag>("trackCollection", edm::InputTag("TrackRefitter"));
327  desc.addUntracked<int>("printTriggerTable", false);
328  desc.addUntracked<double>("minVertexNdf", 10.);
329  desc.addUntracked<double>("minVertexMeanWeight", 0.5);
330  desc.addUntracked<std::vector<double>>("profilePtBorders", {3, 5, 10, 20, 50, 100});
331  desc.addUntracked<std::vector<int>>("iovList", {0, 500000});
332  descriptions.addWithDefaultLabel(desc);
333 }
334 
335 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
static bool mysorter(reco::Track i, reco::Track j)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
TH1F * h_probeDzErr
bool fillHisto(std::string name, std::string tag, double valx, double weight, bool useBinWidth=false)
const edm::EDGetTokenT< edm::TriggerResults > triggerToken_
const std::vector< int > iovList_
double z() const
z coordinate
Definition: Vertex.h:133
bool isFake() const
Definition: Vertex.h:76
TH1F * h_probeDxyErr
JetHTAnalyzer(const edm::ParameterSet &)
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const int printTriggerTable_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
double ndof() const
Definition: Vertex.h:123
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
uint16_t size_type
const edm::InputTag pvsTag_
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.h:110
const edm::InputTag tracksTag_
SmartSelectionMonitor mon
size_t tracksSize() const
number of tracks
Definition: Vertex.h:112
int iEvent
Definition: GenABIO.cc:224
TH1 * addHistogram(TH1 *h, std::string tag)
const double minVtxNdf_
const edm::InputTag triggerTag_
static constexpr double cmToum
~JetHTAnalyzer() override=default
void endJob() override
bool fillProfile(std::string name, std::string tag, double valx, double valy, double weight)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.h:108
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double x() const
x coordinate
Definition: Vertex.h:129
static std::string const triggerResults
Definition: EdmProvDump.cc:47
const std::vector< double > profilePtBorders_
edm::Service< TFileService > outfile_
double y() const
y coordinate
Definition: Vertex.h:131
#define M_PI
void analyze(const edm::Event &, const edm::EventSetup &) override
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const double minVtxWgt_
const edm::EDGetTokenT< reco::VertexCollection > pvsToken_
void beginJob() override
fixed size matrix
HLT enums.
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38