CMS 3D CMS Logo

RazorMonitor.cc
Go to the documentation of this file.
1 // -----------------------------
2 //
3 // Offline DQM for razor triggers. The razor inclusive analysis measures trigger efficiency
4 // in SingleElectron events (orthogonal to analysis), as a 2D function of the razor variables
5 // M_R and R^2. Also monitor dPhi_R, used offline for QCD and/or detector-related MET tail
6 // rejection.
7 // Based on DQMOffline/Trigger/plugins/METMonitor.cc
8 //
9 // -----------------------------
28 
29 #include <TVector3.h>
30 
31 class RazorMonitor : public DQMEDAnalyzer, public TriggerDQMBase {
32 public:
35 
37  ~RazorMonitor() throw() override;
38  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
39 
40  static double CalcMR(const math::XYZTLorentzVector& ja, const math::XYZTLorentzVector& jb);
41  static double CalcR(double MR,
45 
46 protected:
48  void analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) override;
49 
50 private:
52 
55 
59 
63 
64  ObjME MR_ME_;
65  ObjME Rsq_ME_;
66  ObjME dPhiR_ME_;
67  ObjME MRVsRsq_ME_;
68 
71 
74  unsigned int njets_;
75  float rsqCut_;
76  float mrCut_;
77 };
78 
79 // -----------------------------
80 //
81 // Offline DQM for razor triggers. The razor inclusive analysis measures trigger efficiency
82 // in SingleElectron events (orthogonal to analysis), as a 2D function of the razor variables
83 // M_R and R^2. Also monitor dPhi_R, used offline for QCD and/or detector-related MET tail
84 // rejection.
85 // Based on DQMOffline/Trigger/plugins/METMonitor.*
86 //
87 // -----------------------------
88 
90  : folderName_(iConfig.getParameter<std::string>("FolderName")),
91  requireValidHLTPaths_(iConfig.getParameter<bool>("requireValidHLTPaths")),
93  metToken_(consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("met"))),
94  jetToken_(mayConsume<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("jets"))),
96  consumes<std::vector<math::XYZTLorentzVector> >(iConfig.getParameter<edm::InputTag>("hemispheres"))),
97  rsq_binning_(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("rsqBins")),
98  mr_binning_(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("mrBins")),
100  iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("dphiRBins")),
102  iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet"), consumesCollector(), *this)),
104  iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet"), consumesCollector(), *this)),
105  metSelection_(iConfig.getParameter<std::string>("metSelection")),
106  jetSelection_(iConfig.getParameter<std::string>("jetSelection")),
107  njets_(iConfig.getParameter<unsigned int>("njets")),
108  rsqCut_(iConfig.getParameter<double>("rsqCut")),
109  mrCut_(iConfig.getParameter<double>("mrCut")) {}
110 
113  num_genTriggerEventFlag_.reset();
114  }
116  den_genTriggerEventFlag_.reset();
117  }
118 }
119 
120 void RazorMonitor::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) {
121  // Initialize the GenericTriggerEventFlag
123  num_genTriggerEventFlag_->initRun(iRun, iSetup);
125  den_genTriggerEventFlag_->initRun(iRun, iSetup);
126 
127  // check if every HLT path specified in numerator and denominator has a valid match in the HLT Menu
129  den_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->allHLTPathsAreValid() &&
130  den_genTriggerEventFlag_->allHLTPathsAreValid());
131 
132  // if valid HLT paths are required,
133  // create DQM outputs only if all paths are valid
135  return;
136  }
137 
138  std::string histname, histtitle;
139 
140  std::string currentFolder = folderName_;
141  ibooker.setCurrentFolder(currentFolder);
142 
143  // 1D hist, MR
144  histname = "MR";
145  histtitle = "PF MR";
146  bookME(ibooker, MR_ME_, histname, histtitle, mr_binning_);
147  setMETitle(MR_ME_, "PF M_{R} [GeV]", "events / [GeV]");
148 
149  // 1D hist, Rsq
150  histname = "Rsq";
151  histtitle = "PF Rsq";
152  bookME(ibooker, Rsq_ME_, histname, histtitle, rsq_binning_);
153  setMETitle(Rsq_ME_, "PF R^{2}", "events");
154 
155  // 1D hist, dPhiR
156  histname = "dPhiR";
157  histtitle = "dPhiR";
158  bookME(ibooker, dPhiR_ME_, histname, histtitle, dphiR_binning_);
159  setMETitle(dPhiR_ME_, "dPhi_{R}", "events");
160 
161  // 2D hist, MR & Rsq
162  histname = "MRVsRsq";
163  histtitle = "PF MR vs PF Rsq";
164  bookME(ibooker, MRVsRsq_ME_, histname, histtitle, mr_binning_, rsq_binning_);
165  setMETitle(MRVsRsq_ME_, "M_{R} [GeV]", "R^{2}");
166 }
167 
169  // if valid HLT paths are required,
170  // analyze event only if all paths are valid
172  return;
173  }
174 
175  // Filter out events if Trigger Filtering is requested
176  if (den_genTriggerEventFlag_->on() && !den_genTriggerEventFlag_->accept(iEvent, iSetup))
177  return;
178 
179  //met collection
181  iEvent.getByToken(metToken_, metHandle);
182  reco::PFMET pfmet = metHandle->front();
183  if (!metSelection_(pfmet))
184  return;
185 
186  //jet collection, track # of jets for two working points
188  iEvent.getByToken(jetToken_, jetHandle);
189  std::vector<reco::PFJet> jets;
190  if (jetHandle->size() < njets_)
191  return;
192  for (auto const& j : *jetHandle) {
193  if (jetSelection_(j))
194  jets.push_back(j);
195  }
196  if (jets.size() < njets_)
197  return;
198 
199  //razor hemisphere clustering from previous step
201  iEvent.getByToken(theHemispheres_, hemispheres);
202 
203  if (not hemispheres.isValid()) {
204  return;
205  }
206 
207  if (hemispheres
208  ->empty()) { // the Hemisphere Maker will produce an empty collection of hemispheres if # of jets is too high
209  edm::LogError("DQM_HLT_Razor") << "Cannot calculate M_R and R^2 because there are too many jets! (trigger passed "
210  "automatically without forming the hemispheres)"
211  << endl;
212  return;
213  }
214 
215  // should always have 2 hemispheres -- no muons included (c. 2017), if not return invalid hemisphere collection
216  // retaining check for hemisphere size 5 or 10 which correspond to the one or two muon case for possible future use
217  if (!hemispheres->empty() && hemispheres->size() != 2 && hemispheres->size() != 5 && hemispheres->size() != 10) {
218  edm::LogError("DQM_HLT_Razor") << "Invalid hemisphere collection! hemispheres->size() = " << hemispheres->size()
219  << endl;
220  return;
221  }
222 
223  //calculate razor variables, with hemispheres pT-ordered
224  double MR = 0, R = 0;
225  if (hemispheres->at(1).Pt() > hemispheres->at(0).Pt()) {
226  MR = CalcMR(hemispheres->at(1), hemispheres->at(0));
227  R = CalcR(MR, hemispheres->at(1), hemispheres->at(0), metHandle);
228  } else {
229  MR = CalcMR(hemispheres->at(0), hemispheres->at(1));
230  R = CalcR(MR, hemispheres->at(0), hemispheres->at(1), metHandle);
231  }
232 
233  double Rsq = R * R;
234  double dPhiR = abs(deltaPhi(hemispheres->at(0).Phi(), hemispheres->at(1).Phi()));
235 
236  //apply offline selection cuts
237  if (Rsq < rsqCut_ && MR < mrCut_)
238  return;
239 
240  // applying selection for denominator
241  if (den_genTriggerEventFlag_->on() && !den_genTriggerEventFlag_->accept(iEvent, iSetup))
242  return;
243 
244  // filling histograms (denominator)
245  if (Rsq >= rsqCut_) {
246  MR_ME_.denominator->Fill(MR);
247  }
248 
249  if (MR >= mrCut_) {
250  Rsq_ME_.denominator->Fill(Rsq);
251  }
252 
253  dPhiR_ME_.denominator->Fill(dPhiR);
254 
255  MRVsRsq_ME_.denominator->Fill(MR, Rsq);
256 
257  // applying selection for numerator
258  if (num_genTriggerEventFlag_->on() && !num_genTriggerEventFlag_->accept(iEvent, iSetup))
259  return;
260 
261  // filling histograms (numerator)
262  if (Rsq >= rsqCut_) {
263  MR_ME_.numerator->Fill(MR);
264  }
265 
266  if (MR >= mrCut_) {
267  Rsq_ME_.numerator->Fill(Rsq);
268  }
269 
270  dPhiR_ME_.numerator->Fill(dPhiR);
271 
272  MRVsRsq_ME_.numerator->Fill(MR, Rsq);
273 }
274 
277  desc.add<std::string>("FolderName", "HLT/SUSY/Razor");
278  desc.add<bool>("requireValidHLTPaths", true);
279 
280  desc.add<edm::InputTag>("met", edm::InputTag("pfMet"));
281  desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
282  desc.add<edm::InputTag>("hemispheres", edm::InputTag("hemispheresDQM"))
283  ->setComment("hemisphere jets used to compute razor variables");
284  desc.add<std::string>("metSelection", "pt > 0");
285 
286  // from 2016 offline selection
287  desc.add<std::string>("jetSelection", "pt > 80");
288  desc.add<unsigned int>("njets", 2);
289  desc.add<double>("mrCut", 300);
290  desc.add<double>("rsqCut", 0.15);
291 
294  desc.add<edm::ParameterSetDescription>("numGenericTriggerEventPSet", genericTriggerEventPSet);
295  desc.add<edm::ParameterSetDescription>("denGenericTriggerEventPSet", genericTriggerEventPSet);
296 
297  //binning from 2016 offline selection
299  std::vector<double> mrbins = {0., 100., 200., 300., 400., 500., 575., 650., 750., 900., 1200., 1600., 2500., 4000.};
300  histoPSet.add<std::vector<double> >("mrBins", mrbins);
301 
302  std::vector<double> rsqbins = {0., 0.05, 0.1, 0.15, 0.2, 0.25, 0.30, 0.41, 0.52, 0.64, 0.8, 1.5};
303  histoPSet.add<std::vector<double> >("rsqBins", rsqbins);
304 
305  std::vector<double> dphirbins = {0., 0.5, 1.0, 1.5, 2.0, 2.5, 2.8, 3.0, 3.2};
306  histoPSet.add<std::vector<double> >("dphiRBins", dphirbins);
307 
308  desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);
309 
310  descriptions.add("razorMonitoring", desc);
311 }
312 
313 //CalcMR and CalcR borrowed from HLTRFilter.cc
315  if (ja.Pt() <= 0.1)
316  return -1;
317 
318  double A = ja.P();
319  double B = jb.P();
320  double az = ja.Pz();
321  double bz = jb.Pz();
322  TVector3 jaT, jbT;
323  jaT.SetXYZ(ja.Px(), ja.Py(), 0.0);
324  jbT.SetXYZ(jb.Px(), jb.Py(), 0.0);
325  double ATBT = (jaT + jbT).Mag2();
326 
327  double MR = sqrt((A + B) * (A + B) - (az + bz) * (az + bz) -
328  (jbT.Dot(jbT) - jaT.Dot(jaT)) * (jbT.Dot(jbT) - jaT.Dot(jaT)) / (jaT + jbT).Mag2());
329  double mybeta = (jbT.Dot(jbT) - jaT.Dot(jaT)) / sqrt(ATBT * ((A + B) * (A + B) - (az + bz) * (az + bz)));
330 
331  double mygamma = 1. / sqrt(1. - mybeta * mybeta);
332 
333  //use gamma times MRstar
334  return MR * mygamma;
335 }
336 
337 double RazorMonitor::CalcR(double MR,
338  const math::XYZTLorentzVector& ja,
339  const math::XYZTLorentzVector& jb,
340  const edm::Handle<std::vector<reco::PFMET> >& inputMet) {
341  //now we can calculate MTR
342  const math::XYZVector met = (inputMet->front()).momentum();
343 
344  double MTR = sqrt(0.5 * (met.R() * (ja.Pt() + jb.Pt()) - met.Dot(ja.Vect() + jb.Vect())));
345 
346  //filter events
347  return float(MTR) / float(MR); //R
348 }
349 
350 // Define this as a plug-in
dqm::reco::MonitorElement MonitorElement
Definition: RazorMonitor.cc:33
edm::EDGetTokenT< reco::PFJetCollection > jetToken_
Definition: RazorMonitor.cc:57
edm::EDGetTokenT< reco::PFMETCollection > metToken_
Definition: RazorMonitor.cc:56
std::unique_ptr< GenericTriggerEventFlag > num_genTriggerEventFlag_
Definition: RazorMonitor.cc:69
Definition: APVGainStruct.h:7
static double CalcMR(const math::XYZTLorentzVector &ja, const math::XYZTLorentzVector &jb)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
bool hltPathsAreValid_
Definition: RazorMonitor.cc:54
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
Log< level::Error, false > LogError
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
void setMETitle(ObjME &me, const std::string &titleX, const std::string &titleY)
ObjME MRVsRsq_ME_
Definition: RazorMonitor.cc:67
void Fill(long long x)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
StringCutObjectSelector< reco::MET, true > metSelection_
Definition: RazorMonitor.cc:72
edm::EDGetTokenT< std::vector< math::XYZTLorentzVector > > theHemispheres_
Definition: RazorMonitor.cc:58
unsigned int njets_
Definition: RazorMonitor.cc:74
int iEvent
Definition: GenABIO.cc:224
std::vector< double > mr_binning_
Definition: RazorMonitor.cc:61
T sqrt(T t)
Definition: SSEVec.h:23
MonitorElement * denominator
MonitorElement * numerator
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< double > dphiR_binning_
Definition: RazorMonitor.cc:62
StringCutObjectSelector< reco::PFJet, true > jetSelection_
Definition: RazorMonitor.cc:73
std::vector< double > rsq_binning_
Definition: RazorMonitor.cc:60
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
std::unique_ptr< GenericTriggerEventFlag > den_genTriggerEventFlag_
Definition: RazorMonitor.cc:70
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void bookME(DQMStore::IBooker &, ObjME &me, const std::string &histname, const std::string &histtitle, const uint nbins, const double xmin, const double xmax, const bool bookDen=true)
dqm::reco::DQMStore DQMStore
Definition: RazorMonitor.cc:34
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
fixed size matrix
HLT enums.
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: APVGainStruct.h:7
~RazorMonitor() override
RazorMonitor(const edm::ParameterSet &)
Definition: RazorMonitor.cc:89
static double CalcR(double MR, const math::XYZTLorentzVector &ja, const math::XYZTLorentzVector &jb, const edm::Handle< std::vector< reco::PFMET > > &met)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
const std::string folderName_
Definition: RazorMonitor.cc:51
Definition: Run.h:45
const bool requireValidHLTPaths_
Definition: RazorMonitor.cc:53
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
Collection of PF MET.
Computes the MET from a collection of PFCandidates. HF missing!
Definition: PFMET.cc:24